Line data Source code
1 : /**************************************************************************
2 : * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 : * *
4 : * Author: The ALICE Off-line Project. *
5 : * Contributors are mentioned in the code where appropriate. *
6 : * *
7 : * Permission to use, copy, modify and distribute this software and its *
8 : * documentation strictly for non-commercial purposes is hereby granted *
9 : * without fee, provided that the above copyright notice appears in all *
10 : * copies and that both the copyright notice and this permission notice *
11 : * appear in the supporting documentation. The authors make no claims *
12 : * about the suitability of this software for any purpose. It is *
13 : * provided "as is" without express or implied warranty. *
14 : **************************************************************************/
15 :
16 : /// \class AliTPCdataQA
17 : ///
18 : /// July 2011:
19 : ///
20 : /// Changes to accomodate updates of general DQM/QA changes to have per trigger
21 : /// histograms (for a given event specie).
22 : ///
23 : /// AliTPCdataQA has a new flag for only keeping DQM info event by
24 : /// event!
25 : /// The expert/DA functionality has been kept exactly the same!
26 : ///
27 : ///
28 : /// June 2010
29 : ///
30 : /// This update should solve two problems mainly:
31 : /// * The vs event histograms have been limited to a fixed size for the
32 : /// DQM. The 500k seemed to be a big size but is no longer so, so we
33 : /// need to dynamically expand the range. The non-trivial point is that
34 : /// we also have to do it for the copy owned by AliTPCQADataMakerRec.
35 : /// * The amoreGui now remembers the x-range of the first visualization so
36 : /// the trick of setting the relevant event range as the histogram is
37 : /// filled no longer works.
38 : ///
39 : /// The fix is a bit crude but avoids creating a new histogram. Instead
40 : /// the range is expanded (max events and events per bin is doubled) but
41 : /// the number of bins is kept constant! In this way we can change just
42 : /// double the max of the X-axis of the hist and rebin the data. The
43 : /// same can easily be done for the copy owned by AliTPCQADataMakerRec.
44 : ///
45 : /// CAUTION:
46 : /// If we change the number of bins we could crash the whole system
47 : /// because ROOT does not create space for extra bins! (but we do not do
48 : /// this). In that way it is a crude solution.
49 : /// The rebinning in the code only works for an even number of bins.
50 : ///
51 : /// In addition to the above a bug in the reading of the config file was
52 : /// also found and corrected. fAdcMax was set instead of fEventsPerBin.
53 : ///
54 : /// Finally cout was changes to AliInfo.
55 : ///
56 : /// February 2008
57 : ///
58 : /// The code has been heavily modified so that now the RAW data is
59 : /// "expanded" for each sector and stored in a big signal array. Then a
60 : /// simple version of the code in AliTPCclusterer is used to identify
61 : /// the local maxima and these are then used for QA. This gives a better
62 : /// estimate of the charge (both max and total) and also limits the
63 : /// effect of noise.
64 : ///
65 : /// Implementation:
66 : ///
67 : /// In Update the RAW signals >= 3 ADC channels are stored in the arrays.
68 : /// There are 3 arrays:
69 : ///
70 : /// ~~~{.cpp}
71 : /// Float_t** fAllBins 2d array [row][bin(pad, time)] ADC signal
72 : /// Int_t** fAllSigBins 2d array [row][signal#] bin(with signal)
73 : /// Int_t* fAllNSigBins; 1d array [row] Nsignals
74 : /// ~~~
75 : ///
76 : /// This is done sector by sector.
77 : ///
78 : /// When data from a new sector is encountered, the method
79 : /// FindLocalMaxima is called on the data from the previous sector, and
80 : /// the calibration/data objects are updated with the "cluster"
81 : /// info. Finally the arrays are cleared.
82 : ///
83 : /// The requirements for a local maxima is:
84 : /// Charge in bin is >= 5 ADC channels.
85 : /// Charge in bin is larger than all the 8 neighboring bins.
86 : /// At least one of the two pad neighbors has a signal.
87 : /// At least one of the two time neighbors has a signal.
88 : ///
89 : /// Before accessing the data it is expected that the Analyse method is
90 : /// called. This normalizes some of the data objects to per event or per
91 : /// cluster.
92 : /// If more data is passed to the class after Analyse has been called
93 : /// the normalization is reversed and Analyse has to be called again.
94 :
95 : //Root includes
96 : #include <TH1F.h>
97 : #include <TString.h>
98 : #include <TMath.h>
99 : #include <TDirectory.h>
100 : #include <TFile.h>
101 : #include <TError.h>
102 : #include <TMap.h>
103 : #include <TProfile.h>
104 : #include <TObjArray.h>
105 : //AliRoot includes
106 : #include "AliRawReader.h"
107 : #include "AliRawReaderRoot.h"
108 : #include "AliRawReaderDate.h"
109 : #include "AliTPCRawStreamV3.h"
110 : #include "AliTPCCalROC.h"
111 : #include "AliTPCROC.h"
112 : #include "AliMathBase.h"
113 : #include "TTreeStream.h"
114 :
115 : //date
116 : #include "event.h"
117 : #include "AliTPCCalPad.h"
118 : #include "AliTPCPreprocessorOnline.h"
119 :
120 : //header file
121 : #include "AliTPCdataQA.h"
122 : #include "AliLog.h"
123 :
124 :
125 : /// \cond CLASSIMP
126 24 : ClassImp(AliTPCdataQA)
127 : /// \endcond
128 :
129 2 : AliTPCdataQA::AliTPCdataQA() : /*FOLD00*/
130 2 : fFirstTimeBin(60),
131 2 : fLastTimeBin(1000),
132 2 : fAdcMin(3),
133 2 : fAdcMax(1000),
134 2 : fMinQMax(5.f),
135 2 : fRequireNeighbouringPad(kTRUE),
136 2 : fMapping(NULL),
137 2 : fPedestal(0),
138 2 : fNoise(0),
139 2 : fNLocalMaxima(0),
140 2 : fMaxCharge(0),
141 2 : fMeanCharge(0),
142 2 : fNoThreshold(0),
143 2 : fNTimeBins(0),
144 2 : fNPads(0),
145 2 : fTimePosition(0),
146 2 : fOverThreshold10(0),
147 2 : fOverThreshold20(0),
148 2 : fOverThreshold30(0),
149 2 : fHistQVsTimeSideA(0),
150 2 : fHistQVsTimeSideC(0),
151 2 : fHistQMaxVsTimeSideA(0),
152 2 : fHistQMaxVsTimeSideC(0),
153 2 : fHistOccupancyVsEvent(0),
154 2 : fHistNclustersVsEvent(0),
155 2 : fEventCounter(0),
156 2 : fIsAnalysed(kFALSE),
157 2 : fMaxEvents(500000), // Max events for event histograms
158 2 : fEventsPerBin(1000), // Events per bin for event histograms
159 2 : fSignalCounter(0), // Signal counter
160 2 : fClusterCounter(0), // Cluster counter
161 2 : fActiveChambers(72),
162 2 : fAllBins(0),
163 2 : fAllSigBins(0),
164 2 : fAllNSigBins(0),
165 2 : fRowsMax(0),
166 2 : fPadsMax(0),
167 2 : fTimeBinsMax(0),
168 2 : fIsDQM(kFALSE),
169 2 : fHistOccVsSector(0x0),
170 2 : fHistOcc2dVsSector(0x0),
171 2 : fHistQVsSector(0x0),
172 2 : fHistQmaxVsSector(0x0),
173 2 : fOccVec(0x0),
174 2 : fOccMaxVec(0x0),
175 2 : fOccVecFine(0x0),
176 2 : fOccMaxVecFine(0x0)
177 10 : {
178 : //
179 : // default constructor
180 : //
181 :
182 436 : for (Int_t i=0; i<72; ++i) {fActiveChambers.SetBitNumber(i,kTRUE);}
183 4 : }
184 :
185 : //_____________________________________________________________________
186 : AliTPCdataQA::AliTPCdataQA(const AliTPCdataQA &ped) : /*FOLD00*/
187 0 : TH1F(ped),
188 0 : fFirstTimeBin(ped.GetFirstTimeBin()),
189 0 : fLastTimeBin(ped.GetLastTimeBin()),
190 0 : fAdcMin(ped.GetAdcMin()),
191 0 : fAdcMax(ped.GetAdcMax()),
192 0 : fMinQMax(ped.GetMinQMax()),
193 0 : fRequireNeighbouringPad(ped.GetRequireNeighbouringPad()),
194 0 : fMapping(NULL),
195 0 : fPedestal(0),
196 0 : fNoise(0),
197 0 : fNLocalMaxima(0),
198 0 : fMaxCharge(0),
199 0 : fMeanCharge(0),
200 0 : fNoThreshold(0),
201 0 : fNTimeBins(0),
202 0 : fNPads(0),
203 0 : fTimePosition(0),
204 0 : fOverThreshold10(0),
205 0 : fOverThreshold20(0),
206 0 : fOverThreshold30(0),
207 0 : fHistQVsTimeSideA(0),
208 0 : fHistQVsTimeSideC(0),
209 0 : fHistQMaxVsTimeSideA(0),
210 0 : fHistQMaxVsTimeSideC(0),
211 0 : fHistOccupancyVsEvent(0),
212 0 : fHistNclustersVsEvent(0),
213 0 : fEventCounter(ped.GetEventCounter()),
214 0 : fIsAnalysed(ped.GetIsAnalysed()),
215 0 : fMaxEvents(ped.GetMaxEvents()),
216 0 : fEventsPerBin(ped.GetEventsPerBin()),
217 0 : fSignalCounter(ped.GetSignalCounter()),
218 0 : fClusterCounter(ped.GetClusterCounter()),
219 0 : fActiveChambers(ped.fActiveChambers),
220 0 : fAllBins(0),
221 0 : fAllSigBins(0),
222 0 : fAllNSigBins(0),
223 0 : fRowsMax(0),
224 0 : fPadsMax(0),
225 0 : fTimeBinsMax(0),
226 0 : fIsDQM(ped.GetIsDQM()),
227 0 : fHistOccVsSector(0x0),
228 0 : fHistOcc2dVsSector(0x0),
229 0 : fHistQVsSector(0x0),
230 0 : fHistQmaxVsSector(0x0),
231 0 : fOccVec(0x0),
232 0 : fOccMaxVec(0x0),
233 0 : fOccVecFine(0x0),
234 0 : fOccMaxVecFine(0x0)
235 0 : {
236 : /// copy constructor
237 :
238 0 : if(ped.GetNLocalMaxima())
239 0 : fNLocalMaxima = new AliTPCCalPad(*ped.GetNLocalMaxima());
240 0 : if(ped.GetMaxCharge())
241 0 : fMaxCharge = new AliTPCCalPad(*ped.GetMaxCharge());
242 0 : if(ped.GetMeanCharge())
243 0 : fMeanCharge = new AliTPCCalPad(*ped.GetMeanCharge());
244 0 : if(ped.GetNoThreshold())
245 0 : fNoThreshold = new AliTPCCalPad(*ped.GetNoThreshold());
246 0 : if(ped.GetNTimeBins())
247 0 : fNTimeBins = new AliTPCCalPad(*ped.GetNTimeBins());
248 0 : if(ped.GetNPads())
249 0 : fNPads = new AliTPCCalPad(*ped.GetNPads());
250 0 : if(ped.GetTimePosition())
251 0 : fTimePosition = new AliTPCCalPad(*ped.GetTimePosition());
252 0 : if(ped.GetOverThreshold10())
253 0 : fOverThreshold10 = new AliTPCCalPad(*ped.GetOverThreshold10());
254 0 : if(ped.GetOverThreshold20())
255 0 : fOverThreshold20 = new AliTPCCalPad(*ped.GetOverThreshold20());
256 0 : if(ped.GetOverThreshold30())
257 0 : fOverThreshold30 = new AliTPCCalPad(*ped.GetOverThreshold30());
258 0 : if(ped.GetHistQVsTimeSideA()) {
259 0 : fHistQVsTimeSideA = new TProfile(*ped.GetHistQVsTimeSideA());
260 0 : fHistQVsTimeSideA->SetDirectory(0);
261 : }
262 0 : if(ped.GetHistQVsTimeSideC()) {
263 0 : fHistQVsTimeSideC = new TProfile(*ped.GetHistQVsTimeSideC());
264 0 : fHistQVsTimeSideC->SetDirectory(0);
265 : }
266 0 : if(ped.GetHistQMaxVsTimeSideA()) {
267 0 : fHistQMaxVsTimeSideA = new TProfile(*ped.GetHistQMaxVsTimeSideA());
268 0 : fHistQMaxVsTimeSideA->SetDirectory(0);
269 : }
270 0 : if(ped.GetHistQMaxVsTimeSideC()) {
271 0 : fHistQMaxVsTimeSideC = new TProfile(*ped.GetHistQMaxVsTimeSideC());
272 0 : fHistQMaxVsTimeSideC->SetDirectory(0);
273 : }
274 0 : if(ped.GetHistOccupancyVsEventConst()) {
275 0 : fHistOccupancyVsEvent = new TH1F(*ped.GetHistOccupancyVsEventConst());
276 0 : fHistOccupancyVsEvent->SetDirectory(0);
277 : }
278 0 : if(ped.GetHistNclustersVsEventConst()) {
279 0 : fHistNclustersVsEvent = new TH1F(*ped.GetHistNclustersVsEventConst());
280 0 : fHistNclustersVsEvent->SetDirectory(0);
281 : }
282 0 : }
283 :
284 : //_____________________________________________________________________
285 : AliTPCdataQA::AliTPCdataQA(const TMap *config) : /*FOLD00*/
286 0 : TH1F("TPCRAW","TPCRAW",100,0,100),
287 0 : fFirstTimeBin(60),
288 0 : fLastTimeBin(1000),
289 0 : fAdcMin(3),
290 0 : fAdcMax(1000),
291 0 : fMinQMax(5.f),
292 0 : fRequireNeighbouringPad(kTRUE),
293 0 : fMapping(NULL),
294 0 : fPedestal(0),
295 0 : fNoise(0),
296 0 : fNLocalMaxima(0),
297 0 : fMaxCharge(0),
298 0 : fMeanCharge(0),
299 0 : fNoThreshold(0),
300 0 : fNTimeBins(0),
301 0 : fNPads(0),
302 0 : fTimePosition(0),
303 0 : fOverThreshold10(0),
304 0 : fOverThreshold20(0),
305 0 : fOverThreshold30(0),
306 0 : fHistQVsTimeSideA(0),
307 0 : fHistQVsTimeSideC(0),
308 0 : fHistQMaxVsTimeSideA(0),
309 0 : fHistQMaxVsTimeSideC(0),
310 0 : fHistOccupancyVsEvent(0),
311 0 : fHistNclustersVsEvent(0),
312 0 : fEventCounter(0),
313 0 : fIsAnalysed(kFALSE),
314 0 : fMaxEvents(500000),
315 0 : fEventsPerBin(1000),
316 0 : fSignalCounter(0),
317 0 : fClusterCounter(0),
318 0 : fActiveChambers(72),
319 0 : fAllBins(0),
320 0 : fAllSigBins(0),
321 0 : fAllNSigBins(0),
322 0 : fRowsMax(0),
323 0 : fPadsMax(0),
324 0 : fTimeBinsMax(0),
325 0 : fIsDQM(kFALSE),
326 0 : fHistOccVsSector(0x0),
327 0 : fHistOcc2dVsSector(0x0),
328 0 : fHistQVsSector(0x0),
329 0 : fHistQmaxVsSector(0x0),
330 0 : fOccVec(0x0),
331 0 : fOccMaxVec(0x0),
332 0 : fOccVecFine(0x0),
333 0 : fOccMaxVecFine(0x0)
334 0 : {
335 : /// default constructor
336 :
337 0 : if (config->GetValue("FirstTimeBin")) fFirstTimeBin = ((TObjString*)config->GetValue("FirstTimeBin"))->GetString().Atoi();
338 0 : if (config->GetValue("LastTimeBin")) fLastTimeBin = ((TObjString*)config->GetValue("LastTimeBin"))->GetString().Atoi();
339 0 : if (config->GetValue("AdcMin")) fAdcMin = ((TObjString*)config->GetValue("AdcMin"))->GetString().Atoi();
340 0 : if (config->GetValue("AdcMax")) fAdcMax = ((TObjString*)config->GetValue("AdcMax"))->GetString().Atoi();
341 0 : if (config->GetValue("MinQMax")) fMinQMax = ((TObjString*)config->GetValue("MinQMax"))->GetString().Atof();
342 0 : if (config->GetValue("MaxEvents")) fMaxEvents = ((TObjString*)config->GetValue("MaxEvents"))->GetString().Atoi();
343 0 : if (config->GetValue("EventsPerBin")) fEventsPerBin = ((TObjString*)config->GetValue("EventsPerBin"))->GetString().Atoi();
344 0 : if (config->GetValue("RequireNeighbouringPad")) fRequireNeighbouringPad = ((TObjString*)config->GetValue("RequireNeighbouringPad"))->GetString().Atoi();
345 0 : for (Int_t i=0; i<72; ++i) {fActiveChambers.SetBitNumber(i,kTRUE);}
346 0 : }
347 :
348 : //_____________________________________________________________________
349 : AliTPCdataQA& AliTPCdataQA::operator = (const AliTPCdataQA &source)
350 : {
351 : /// assignment operator
352 :
353 0 : if (&source == this) return *this;
354 0 : new (this) AliTPCdataQA(source);
355 :
356 0 : return *this;
357 0 : }
358 :
359 :
360 : //_____________________________________________________________________
361 0 : AliTPCdataQA::~AliTPCdataQA() /*FOLD00*/
362 0 : {
363 : /// destructor
364 :
365 : // do not delete fMapping, because we do not own it.
366 : // do not delete fMapping, because we do not own it.
367 : // do not delete fNoise and fPedestal, because we do not own them.
368 :
369 0 : delete fNLocalMaxima;
370 0 : delete fMaxCharge;
371 0 : delete fMeanCharge;
372 0 : delete fNoThreshold;
373 0 : delete fNTimeBins;
374 0 : delete fNPads;
375 0 : delete fTimePosition;
376 0 : delete fOverThreshold10;
377 0 : delete fOverThreshold20;
378 0 : delete fOverThreshold30;
379 0 : delete fHistQVsTimeSideA;
380 0 : delete fHistQVsTimeSideC;
381 0 : delete fHistQMaxVsTimeSideA;
382 0 : delete fHistQMaxVsTimeSideC;
383 0 : delete fHistOccupancyVsEvent;
384 0 : delete fHistNclustersVsEvent;
385 :
386 : // DQM
387 0 : delete fHistOccVsSector;
388 0 : delete fHistOcc2dVsSector;
389 0 : delete fHistQVsSector;
390 0 : delete fHistQmaxVsSector;
391 0 : delete fOccVec;
392 0 : delete fOccMaxVec;
393 0 : delete fOccVecFine;
394 0 : delete fOccMaxVecFine;
395 :
396 0 : for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
397 0 : delete [] fAllBins[iRow];
398 0 : delete [] fAllSigBins[iRow];
399 : }
400 0 : delete [] fAllBins;
401 0 : delete [] fAllSigBins;
402 0 : delete [] fAllNSigBins;
403 0 : }
404 :
405 : //_____________________________________________________________________
406 : TH1F* AliTPCdataQA::GetHistOccupancyVsEvent()
407 : {
408 : /// Create Occupancy vs event histogram
409 : /// (we create this histogram differently then the other histograms
410 : /// because this we want to be able to access and copy
411 : /// from AliTPCQAMakerRec before it normally would be created)
412 :
413 0 : if(!fHistOccupancyVsEvent) {
414 :
415 0 : Int_t nBins = fMaxEvents/fEventsPerBin;
416 0 : fHistOccupancyVsEvent = new TH1F("hOccupancyVsEvent", "Occupancy vs event number (~time); Event number; Occupancy", nBins, 0, nBins*fEventsPerBin);
417 0 : fHistOccupancyVsEvent->SetDirectory(0);
418 0 : }
419 :
420 0 : return fHistOccupancyVsEvent;
421 0 : }
422 :
423 : //_____________________________________________________________________
424 : TH1F* AliTPCdataQA::GetHistNclustersVsEvent()
425 : {
426 : /// Create Nclusters vs event histogram
427 : /// (we create this histogram differently then the other histograms
428 : /// because this we want to be able to access and copy
429 : /// from AliTPCQAMakerRec before it normally would be created)
430 :
431 0 : if(!fHistNclustersVsEvent) {
432 :
433 0 : Int_t nBins = fMaxEvents/fEventsPerBin;
434 0 : fHistNclustersVsEvent = new TH1F("hNclustersVsEvent", "Nclusters vs event number (~time); Event number; Nclusters per event", nBins, 0, nBins*fEventsPerBin);
435 0 : fHistNclustersVsEvent->SetDirectory(0);
436 0 : }
437 :
438 0 : return fHistNclustersVsEvent;
439 0 : }
440 :
441 : //_____________________________________________________________________
442 : void AliTPCdataQA::UpdateEventHistograms()
443 : {
444 : /// Update histograms that display occupancy and
445 : /// number of clusters as a function of number of
446 : /// events
447 :
448 0 : if (!fHistOccupancyVsEvent)
449 0 : GetHistOccupancyVsEvent();
450 0 : if (!fHistNclustersVsEvent)
451 0 : GetHistNclustersVsEvent();
452 :
453 0 : if(fEventCounter > fMaxEvents) {
454 :
455 : // we have to expand the histogram to handle the larger number of
456 : // events. The way it is done now is to double the range and the
457 : // number of events per bin (so the number of histogram bins stays
458 : // constant)
459 0 : fEventsPerBin *= 2;
460 0 : fMaxEvents *= 2;
461 :
462 : // Change histogram limits
463 0 : const Int_t nBins = fHistOccupancyVsEvent->GetXaxis()->GetNbins();
464 0 : fHistOccupancyVsEvent->GetXaxis()->Set(nBins, fHistOccupancyVsEvent->GetXaxis()->GetNbins(), fMaxEvents);
465 0 : fHistNclustersVsEvent->GetXaxis()->Set(nBins, fHistNclustersVsEvent->GetXaxis()->GetNbins(), fMaxEvents);
466 :
467 : // Rebin the histogram
468 0 : for(Int_t bin = 1; bin <= nBins; bin+=2) {
469 :
470 0 : Int_t newBin = TMath::Nint(Float_t(bin+1)/2.0);
471 0 : Float_t newContent = (fHistOccupancyVsEvent->GetBinContent(bin)+
472 0 : fHistOccupancyVsEvent->GetBinContent(bin+1))/2.0;
473 0 : fHistOccupancyVsEvent->SetBinContent(newBin, newContent);
474 :
475 0 : newContent = (fHistNclustersVsEvent->GetBinContent(bin)+
476 0 : fHistNclustersVsEvent->GetBinContent(bin+1))/2.0;
477 0 : fHistNclustersVsEvent->SetBinContent(newBin, newContent);
478 : }
479 :
480 : // Set the remaining bins to 0
481 0 : Int_t lastHalf = nBins/2 +1;
482 0 : for(Int_t bin = lastHalf; bin <= nBins; bin++) {
483 :
484 0 : fHistOccupancyVsEvent->SetBinContent(bin, 0);
485 0 : fHistNclustersVsEvent->SetBinContent(bin, 0);
486 : }
487 :
488 : // In this case we should nut update but wait untill the new
489 : // number of events per bin is reached!
490 : return;
491 : }
492 :
493 0 : const Int_t bin = TMath::Nint(Float_t(fEventCounter)/fEventsPerBin);
494 :
495 : Float_t averageOccupancy =
496 0 : Float_t(fSignalCounter)/fEventsPerBin/(fLastTimeBin - fFirstTimeBin +1.0)
497 0 : / 570132.0; // 570,132 is number of pads
498 0 : fHistOccupancyVsEvent->SetBinContent(bin, averageOccupancy);
499 0 : fSignalCounter = 0;
500 :
501 : Float_t averageNclusters =
502 0 : Float_t(fClusterCounter)/fEventsPerBin;
503 0 : fHistNclustersVsEvent->SetBinContent(bin, averageNclusters);
504 0 : fClusterCounter = 0;
505 0 : }
506 :
507 : //_____________________________________________________________________
508 : Bool_t AliTPCdataQA::ProcessEvent(AliTPCRawStreamV3 *const rawStreamV3)
509 : {
510 : /// Event Processing loop - AliTPCRawStreamV3
511 :
512 : Bool_t withInput = kFALSE;
513 : Int_t nSignals = 0;
514 : Int_t lastSector = -1;
515 :
516 0 : Init();
517 :
518 0 : while ( rawStreamV3->NextDDL() ){
519 :
520 0 : while ( rawStreamV3->NextChannel() ){
521 :
522 0 : Int_t iSector = rawStreamV3->GetSector(); // current sector
523 0 : Int_t iRow = rawStreamV3->GetRow(); // current row
524 0 : Int_t iPad = rawStreamV3->GetPad(); // current pad
525 0 : Int_t iPatch = rawStreamV3->GetPatchIndex(); // current patch
526 0 : Int_t iBranch = rawStreamV3->GetBranch(); // current branch
527 0 : if (iRow<0 || iPad<0) continue;
528 : // Call local maxima finder if the data is in a new sector
529 0 : if(iSector != lastSector) {
530 :
531 0 : if(nSignals>0)
532 0 : FindLocalMaxima(lastSector);
533 :
534 0 : CleanArrays();
535 : lastSector = iSector;
536 : nSignals = 0;
537 0 : }
538 :
539 0 : while ( rawStreamV3->NextBunch() ){
540 :
541 0 : Int_t startTbin = (Int_t)rawStreamV3->GetStartTimeBin();
542 0 : Int_t bunchlength = (Int_t)rawStreamV3->GetBunchLength();
543 0 : const UShort_t *sig = rawStreamV3->GetSignals();
544 :
545 0 : for (Int_t iTimeBin = 0; iTimeBin<bunchlength; iTimeBin++){
546 0 : Float_t signal=(Float_t)sig[iTimeBin];
547 0 : nSignals += Update(iSector,iRow,iPad,startTbin--,signal, iPatch, iBranch);
548 : withInput = kTRUE;
549 : }
550 : }
551 0 : }
552 : }
553 :
554 0 : if (lastSector>=0&&nSignals>0)
555 0 : FindLocalMaxima(lastSector);
556 :
557 0 : CleanArrays();
558 :
559 0 : return withInput;
560 : }
561 :
562 : //_____________________________________________________________________
563 : Bool_t AliTPCdataQA::ProcessEvent(AliRawReader *const rawReader)
564 : {
565 : /// Event processing loop - AliRawReader
566 :
567 0 : AliTPCRawStreamV3 rawStreamV3(rawReader,(AliAltroMapping**)fMapping);
568 0 : Bool_t res=ProcessEvent(&rawStreamV3);
569 0 : if(res) {
570 0 : fEventCounter++; // only increment event counter if there is TPC data
571 :
572 0 : if(fEventCounter%fEventsPerBin==0)
573 0 : UpdateEventHistograms();
574 : }
575 : return res;
576 0 : }
577 :
578 : //_____________________________________________________________________
579 : Bool_t AliTPCdataQA::ProcessEvent(eventHeaderStruct *const event)
580 : {
581 : /// process date event
582 :
583 0 : AliRawReaderDate rawReader((void*)event);
584 0 : Bool_t result=ProcessEvent(&rawReader);
585 0 : return result;
586 0 : }
587 :
588 :
589 :
590 : //_____________________________________________________________________
591 : void AliTPCdataQA::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
592 : {
593 : /// Write class to file
594 :
595 0 : TString sDir(dir);
596 0 : TString option;
597 :
598 0 : if ( append )
599 0 : option = "update";
600 : else
601 0 : option = "recreate";
602 :
603 0 : TDirectory *backup = gDirectory;
604 0 : TFile f(filename,option.Data());
605 0 : f.cd();
606 0 : if ( !sDir.IsNull() ){
607 0 : f.mkdir(sDir.Data());
608 0 : f.cd(sDir);
609 : }
610 0 : this->Write();
611 0 : f.Close();
612 :
613 0 : if ( backup ) backup->cd();
614 0 : }
615 :
616 :
617 : //_____________________________________________________________________
618 : Int_t AliTPCdataQA::Update(const Int_t iSector, /*FOLD00*/
619 : const Int_t iRow,
620 : const Int_t iPad,
621 : const Int_t iTimeBin,
622 : Float_t signal,
623 : const Int_t iPatch,
624 : const Int_t iBranch)
625 : {
626 : /// Signal filling method
627 :
628 0 : if (!fActiveChambers[iSector]) return 0;
629 : //
630 : // TimeBin cut
631 : //
632 0 : if (iTimeBin<fFirstTimeBin) return 0;
633 0 : if (iTimeBin>fLastTimeBin) return 0;
634 :
635 : // if pedestal calibrations are loaded subtract pedestals
636 0 : if(fPedestal) {
637 :
638 0 : Float_t ped = fPedestal->GetCalROC(iSector)->GetValue(iRow, iPad);
639 : // Don't use data from pads where pedestals are abnormally small or large
640 0 : if(ped<10 || ped>90)
641 0 : return 0;
642 0 : signal -= ped;
643 0 : }
644 :
645 0 : if(fIsDQM) {
646 :
647 0 : fOccVec->GetArray()[iSector] += 1.0;
648 : // To change before committing
649 0 : if(iPatch>=0 && iBranch>=0 && iPatch<=5 && iBranch <= 1)
650 0 : fOccVecFine->GetArray()[(iSector%36)*12+iPatch*2+iBranch] += 1.0;
651 : } else {
652 : // In fNoThreshold we fill all data to estimate the ZS volume
653 0 : Float_t count = fNoThreshold->GetCalROC(iSector)->GetValue(iRow, iPad);
654 0 : fNoThreshold->GetCalROC(iSector)->SetValue(iRow, iPad,count+1);
655 : }
656 :
657 : // Require at least 3 ADC channels
658 0 : if (signal < fAdcMin)
659 0 : return 0;
660 :
661 : // if noise calibrations are loaded require at least 3*sigmaNoise
662 0 : if(fNoise) {
663 :
664 0 : Float_t noise = fNoise->GetCalROC(iSector)->GetValue(iRow, iPad);
665 :
666 0 : if(signal < noise*3.0)
667 0 : return 0;
668 0 : }
669 :
670 : //
671 : // This signal is ok and we store it in the cluster map
672 : //
673 :
674 0 : SetExpandDigit(iRow, iPad, iTimeBin, signal);
675 :
676 0 : fSignalCounter++;
677 :
678 0 : return 1; // signal was accepted
679 0 : }
680 :
681 : //_____________________________________________________________________
682 : void AliTPCdataQA::FindLocalMaxima(const Int_t iSector)
683 : {
684 : /// This method is called after the data from each sector has been
685 : /// exapanded into an array
686 : /// Loop over the signals and identify local maxima and fill the
687 : /// calibration objects with the information
688 :
689 0 : if (!fActiveChambers[iSector]) return;
690 0 : if (!fAllBins) return; // nothing has been expanded to the array
691 :
692 : Int_t nLocalMaxima = 0;
693 0 : const Int_t maxTimeBin = fTimeBinsMax+4; // Used to step between neighboring pads
694 : // Because we have tha pad-time data in a
695 : // 1d array
696 :
697 0 : for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
698 :
699 0 : Float_t* allBins = fAllBins[iRow];
700 0 : Int_t* sigBins = fAllSigBins[iRow];
701 0 : const Int_t nSigBins = fAllNSigBins[iRow];
702 :
703 0 : for (Int_t iSig = 0; iSig < nSigBins; iSig++) {
704 :
705 0 : Int_t bin = sigBins[iSig];
706 0 : Float_t *b = &allBins[bin];
707 :
708 : //
709 : // Now we check if this is a local maximum
710 : //
711 :
712 0 : Float_t qMax = b[0];
713 :
714 : // First check that the charge is bigger than the threshold
715 0 : if (qMax<fMinQMax)
716 0 : continue;
717 :
718 : // Require at least one neighboring pad with signal
719 0 : if (fRequireNeighbouringPad && (b[-maxTimeBin]+b[maxTimeBin]<=0) ) continue;
720 :
721 : // Require at least one neighboring time bin with signal
722 0 : if (b[-1]+b[1]<=0) continue;
723 :
724 : //
725 : // Check that this is a local maximum
726 : // Note that the checking is done so that if 2 charges has the same
727 : // qMax then only 1 cluster is generated
728 : // (that is why there is BOTH > and >=)
729 : //
730 0 : if (b[-maxTimeBin] >= qMax) continue;
731 0 : if (b[-1 ] >= qMax) continue;
732 0 : if (b[+maxTimeBin] > qMax) continue;
733 0 : if (b[+1 ] > qMax) continue;
734 0 : if (b[-maxTimeBin-1] >= qMax) continue;
735 0 : if (b[+maxTimeBin-1] >= qMax) continue;
736 0 : if (b[+maxTimeBin+1] > qMax) continue;
737 0 : if (b[-maxTimeBin+1] >= qMax) continue;
738 :
739 : //
740 : // Now we accept the local maximum and fill the calibration/data objects
741 : //
742 0 : ++nLocalMaxima;
743 :
744 0 : Int_t iPad, iTimeBin;
745 0 : GetPadAndTimeBin(bin, iPad, iTimeBin);
746 :
747 0 : if(!fIsDQM) {
748 0 : Float_t count = fNLocalMaxima->GetCalROC(iSector)->GetValue(iRow, iPad);
749 0 : fNLocalMaxima->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
750 :
751 0 : count = fTimePosition->GetCalROC(iSector)->GetValue(iRow, iPad);
752 0 : fTimePosition->GetCalROC(iSector)->SetValue(iRow, iPad, count+iTimeBin);
753 :
754 0 : Float_t charge = fMaxCharge->GetCalROC(iSector)->GetValue(iRow, iPad);
755 0 : fMaxCharge->GetCalROC(iSector)->SetValue(iRow, iPad, charge + qMax);
756 :
757 0 : if(qMax>=10) {
758 0 : count = fOverThreshold10->GetCalROC(iSector)->GetValue(iRow, iPad);
759 0 : fOverThreshold10->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
760 0 : }
761 0 : if(qMax>=20) {
762 0 : count = fOverThreshold20->GetCalROC(iSector)->GetValue(iRow, iPad);
763 0 : fOverThreshold20->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
764 0 : }
765 0 : if(qMax>=30) {
766 0 : count = fOverThreshold30->GetCalROC(iSector)->GetValue(iRow, iPad);
767 0 : fOverThreshold30->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
768 0 : }
769 0 : }
770 :
771 : //
772 : // Calculate the total charge as the sum over the region:
773 : //
774 : // o o o o o
775 : // o i i i o
776 : // o i C i o
777 : // o i i i o
778 : // o o o o o
779 : //
780 : // with qmax at the center C.
781 : //
782 : // The inner charge (i) we always add, but we only add the outer
783 : // charge (o) if the neighboring inner bin (i) has a signal.
784 : //
785 0 : Int_t minP = 0, maxP = 0, minT = 0, maxT = 0;
786 : Float_t qTot = qMax;
787 0 : for(Int_t i = -1; i<=1; i++) {
788 0 : for(Int_t j = -1; j<=1; j++) {
789 :
790 0 : if(i==0 && j==0)
791 : continue;
792 :
793 0 : Float_t charge1 = GetQ(b, i, j, maxTimeBin, minT, maxT, minP, maxP);
794 0 : qTot += charge1;
795 0 : if(charge1>0) {
796 : // see if the next neighbor is also above threshold
797 0 : if(i*j==0) {
798 0 : qTot += GetQ(b, 2*i, 2*j, maxTimeBin, minT, maxT, minP, maxP);
799 0 : } else {
800 : // we are in a diagonal corner
801 0 : qTot += GetQ(b, i, 2*j, maxTimeBin, minT, maxT, minP, maxP);
802 0 : qTot += GetQ(b, 2*i, j, maxTimeBin, minT, maxT, minP, maxP);
803 0 : qTot += GetQ(b, 2*i, 2*j, maxTimeBin, minT, maxT, minP, maxP);
804 : }
805 : }
806 0 : }
807 : }
808 :
809 0 : if(fIsDQM) {
810 0 : fHistQVsSector->Fill(iSector, qTot);
811 0 : fHistQmaxVsSector->Fill(iSector, qMax);
812 0 : } else {
813 0 : Float_t charge = fMeanCharge->GetCalROC(iSector)->GetValue(iRow, iPad);
814 0 : fMeanCharge->GetCalROC(iSector)->SetValue(iRow, iPad, charge + qTot);
815 :
816 0 : Float_t count = fNTimeBins->GetCalROC(iSector)->GetValue(iRow, iPad);
817 0 : fNTimeBins->GetCalROC(iSector)->SetValue(iRow, iPad, count + maxT-minT+1);
818 :
819 0 : count = fNPads->GetCalROC(iSector)->GetValue(iRow, iPad);
820 0 : fNPads->GetCalROC(iSector)->SetValue(iRow, iPad, count + maxP-minP+1);
821 :
822 0 : if((iSector%36)<18) { // A side
823 0 : fHistQVsTimeSideA->Fill(iTimeBin, qTot);
824 0 : fHistQMaxVsTimeSideA->Fill(iTimeBin, qMax);
825 0 : } else {
826 0 : fHistQVsTimeSideC->Fill(iTimeBin, qTot);
827 0 : fHistQMaxVsTimeSideC->Fill(iTimeBin, qMax);
828 : }
829 : }
830 0 : } // end loop over signals
831 : } // end loop over rows
832 :
833 0 : fClusterCounter += nLocalMaxima;
834 0 : }
835 :
836 : //_____________________________________________________________________
837 : void AliTPCdataQA::Analyse()
838 : {
839 : /// Calculate calibration constants
840 :
841 0 : AliInfo("Analyse called");
842 :
843 0 : if(fIsDQM == kTRUE) {
844 :
845 0 : AliInfo("DQM flas is set -> No 2d information to analyze");
846 0 : return;
847 : }
848 :
849 0 : if(fIsAnalysed == kTRUE) {
850 :
851 0 : AliInfo("No new data since Analyse was called last time");
852 0 : return;
853 : }
854 :
855 0 : if(fEventCounter==0) {
856 :
857 0 : AliInfo("EventCounter == 0, Cannot analyse");
858 0 : return;
859 : }
860 :
861 0 : Int_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
862 0 : AliInfo(Form("EventCounter: %d , TimeBins: %d", fEventCounter, nTimeBins));
863 :
864 0 : Float_t normalization = 1.0 / Float_t(fEventCounter * nTimeBins);
865 0 : fNoThreshold->Multiply(normalization);
866 :
867 0 : fMeanCharge->Divide(fNLocalMaxima);
868 0 : fMaxCharge->Divide(fNLocalMaxima);
869 0 : fNTimeBins->Divide(fNLocalMaxima);
870 0 : fNPads->Divide(fNLocalMaxima);
871 0 : fTimePosition->Divide(fNLocalMaxima);
872 :
873 0 : fIsAnalysed = kTRUE;
874 0 : }
875 :
876 :
877 : //_____________________________________________________________________
878 : void AliTPCdataQA::MakeTree(const char *fname) const {
879 : /// Export result to the tree -located in the file
880 : /// This file can be analyzed using AliTPCCalibViewer
881 :
882 0 : AliTPCPreprocessorOnline preprocesor;
883 :
884 0 : if (fNLocalMaxima) preprocesor.AddComponent(fNLocalMaxima);
885 0 : if (fMaxCharge) preprocesor.AddComponent(fMaxCharge);
886 0 : if (fMeanCharge) preprocesor.AddComponent(fMeanCharge);
887 0 : if (fNoThreshold) preprocesor.AddComponent(fNoThreshold);
888 0 : if (fNTimeBins) preprocesor.AddComponent(fNTimeBins);
889 0 : if (fNPads) preprocesor.AddComponent(fNPads);
890 0 : if (fTimePosition) preprocesor.AddComponent(fTimePosition);
891 0 : if (fOverThreshold10) preprocesor.AddComponent(fOverThreshold10);
892 0 : if (fOverThreshold20) preprocesor.AddComponent(fOverThreshold20);
893 0 : if (fOverThreshold30) preprocesor.AddComponent(fOverThreshold30);
894 :
895 0 : preprocesor.DumpToFile(fname);
896 0 : }
897 :
898 :
899 : //_____________________________________________________________________
900 : void AliTPCdataQA::MakeArrays(){
901 : /// The arrays for expanding the raw data are defined and
902 : /// som parameters are intialised
903 :
904 0 : AliTPCROC * roc = AliTPCROC::Instance();
905 : //
906 : // To make the array big enough for all sectors we take
907 : // the dimensions from the outer row of an OROC (the last sector)
908 : //
909 0 : fRowsMax = roc->GetNRows(roc->GetNSector()-1);
910 0 : fPadsMax = roc->GetNPads(roc->GetNSector()-1,fRowsMax-1);
911 0 : fTimeBinsMax = fLastTimeBin - fFirstTimeBin +1;
912 :
913 : //
914 : // Since we have added 2 pads (TimeBins) before and after the real pads (TimeBins)
915 : // to make sure that we can always query the exanded table even when the
916 : // max is on the edge
917 : //
918 :
919 :
920 0 : fAllBins = new Float_t*[fRowsMax];
921 0 : fAllSigBins = new Int_t*[fRowsMax];
922 0 : fAllNSigBins = new Int_t[fRowsMax];
923 :
924 0 : for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
925 : //
926 0 : Int_t maxBin = (fTimeBinsMax+4)*(fPadsMax+4);
927 0 : fAllBins[iRow] = new Float_t[maxBin];
928 0 : memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin); // set all values to 0
929 0 : fAllSigBins[iRow] = new Int_t[maxBin];
930 0 : fAllNSigBins[iRow] = 0;
931 : }
932 0 : }
933 :
934 :
935 : //_____________________________________________________________________
936 : void AliTPCdataQA::CleanArrays(){
937 : ///
938 :
939 0 : if (!fAllBins) return;
940 :
941 0 : for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
942 :
943 : // To speed up the performance by a factor 2 on cosmic data (and
944 : // presumably pp data as well) where the ocupancy is low, the
945 : // memset is only called if there is more than 1000 signals for a
946 : // row (of the order 1% occupancy)
947 0 : if(fAllNSigBins[iRow]<1000) {
948 :
949 0 : Float_t* allBins = fAllBins[iRow];
950 0 : Int_t* sigBins = fAllSigBins[iRow];
951 : const Int_t nSignals = fAllNSigBins[iRow];
952 0 : for(Int_t i = 0; i < nSignals; i++)
953 0 : allBins[sigBins[i]]=0;
954 0 : } else {
955 :
956 0 : Int_t maxBin = (fTimeBinsMax+4)*(fPadsMax+4);
957 0 : memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
958 : }
959 :
960 0 : fAllNSigBins[iRow]=0;
961 : }
962 0 : }
963 :
964 : //_____________________________________________________________________
965 : void AliTPCdataQA::GetPadAndTimeBin(Int_t bin, Int_t& iPad, Int_t& iTimeBin){
966 : /// Return pad and timebin for a given bin
967 :
968 : // Int_t bin = iPad*(fTimeBinsMax+4)+iTimeBin;
969 0 : iTimeBin = bin%(fTimeBinsMax+4);
970 0 : iPad = (bin-iTimeBin)/(fTimeBinsMax+4);
971 :
972 0 : iPad -= 2;
973 0 : iTimeBin -= 2;
974 0 : iTimeBin += fFirstTimeBin;
975 :
976 0 : R__ASSERT(iPad>=0 && iPad<=fPadsMax);
977 0 : R__ASSERT(iTimeBin>=fFirstTimeBin && iTimeBin<=fLastTimeBin);
978 0 : }
979 :
980 : //_____________________________________________________________________
981 : void AliTPCdataQA::SetExpandDigit(const Int_t iRow, Int_t iPad,
982 : Int_t iTimeBin, const Float_t signal)
983 : {
984 : ///
985 :
986 : // Make the arrays for expanding the data
987 0 : if (!fAllBins)
988 0 : MakeArrays();
989 :
990 0 : R__ASSERT(iRow>=0 && iRow<fRowsMax);
991 0 : R__ASSERT(iPad>=0 && iPad<=fPadsMax);
992 0 : R__ASSERT(iTimeBin>=fFirstTimeBin && iTimeBin<=fLastTimeBin);
993 :
994 0 : iTimeBin -= fFirstTimeBin;
995 0 : iPad += 2;
996 0 : iTimeBin += 2;
997 :
998 0 : Int_t bin = iPad*(fTimeBinsMax+4)+iTimeBin;
999 0 : fAllBins[iRow][bin] = signal;
1000 0 : fAllSigBins[iRow][fAllNSigBins[iRow]] = bin;
1001 0 : fAllNSigBins[iRow]++;
1002 0 : }
1003 :
1004 : //______________________________________________________________________________
1005 : Float_t AliTPCdataQA::GetQ(const Float_t* adcArray, const Int_t time,
1006 : const Int_t pad, const Int_t maxTimeBins,
1007 : Int_t& timeMin, Int_t& timeMax,
1008 : Int_t& padMin, Int_t& padMax) const
1009 : {
1010 : /// This methods return the charge in the bin time+pad*maxTimeBins
1011 : /// If the charge is above 0 it also updates the padMin, padMax, timeMin
1012 : /// and timeMax if necessary
1013 :
1014 0 : Float_t charge = adcArray[time + pad*maxTimeBins];
1015 0 : if(charge > 0) {
1016 0 : timeMin = TMath::Min(time, timeMin); timeMax = TMath::Max(time, timeMax);
1017 0 : padMin = TMath::Min(pad, padMin); padMax = TMath::Max(pad, padMax);
1018 0 : }
1019 0 : return charge;
1020 : }
1021 :
1022 : //______________________________________________________________________________
1023 0 : void AliTPCdataQA::Streamer(TBuffer &xRuub)
1024 : {
1025 : /// Automatic schema evolution was first used from revision 4
1026 : /// Code based on:
1027 : /// http://root.cern.ch/root/roottalk/roottalk02/3207.html
1028 :
1029 4 : UInt_t xRuus, xRuuc;
1030 2 : if (xRuub.IsReading()) {
1031 2 : Version_t xRuuv = xRuub.ReadVersion(&xRuus, &xRuuc);
1032 : //we use the automatic algorithm for class version > 3
1033 2 : if (xRuuv > 3) {
1034 0 : AliTPCdataQA::Class()->ReadBuffer(xRuub, this, xRuuv, xRuus,
1035 0 : xRuuc);
1036 0 : return;
1037 : }
1038 2 : TH1F::Streamer(xRuub);
1039 2 : xRuub >> fFirstTimeBin;
1040 2 : xRuub >> fLastTimeBin;
1041 2 : xRuub >> fAdcMin;
1042 2 : xRuub >> fAdcMax;
1043 2 : xRuub >> fNLocalMaxima;
1044 2 : xRuub >> fMaxCharge;
1045 2 : xRuub >> fMeanCharge;
1046 2 : xRuub >> fNoThreshold;
1047 2 : xRuub >> fNTimeBins;
1048 2 : xRuub >> fNPads;
1049 2 : xRuub >> fTimePosition;
1050 2 : xRuub >> fEventCounter;
1051 2 : xRuub >> fIsAnalysed;
1052 2 : xRuub.CheckByteCount(xRuus, xRuuc, AliTPCdataQA::IsA());
1053 2 : } else {
1054 0 : AliTPCdataQA::Class()->WriteBuffer(xRuub,this);
1055 : }
1056 4 : }
1057 :
1058 : //____________________________________________________________________________________________
1059 : void AliTPCdataQA::FillOccupancyProfile()
1060 : {
1061 : /// This has to be filled at the end of the loop over data
1062 :
1063 0 : if(!fIsDQM)
1064 0 : AliInfo("Method only meaningful for DQM");
1065 :
1066 0 : for(Int_t i = 0; i < 72; i++) {
1067 :
1068 0 : fOccVec->GetArray()[i] /= fOccMaxVec->GetArray()[i];
1069 0 : fHistOccVsSector->Fill(i, fOccVec->GetArray()[i]);
1070 : }
1071 :
1072 : const Int_t nBranches = 36*12;
1073 0 : for(Int_t i = 0; i < nBranches; i++) {
1074 :
1075 0 : fOccVecFine->GetArray()[i] /= fOccMaxVecFine->GetArray()[i];
1076 :
1077 0 : const Int_t fullSector = Int_t(i/12);
1078 :
1079 0 : Int_t branch = i - fullSector*12;
1080 0 : const Int_t patch = Int_t(branch/2);
1081 :
1082 0 : branch -= patch*2;
1083 :
1084 0 : fHistOcc2dVsSector->Fill(fullSector+0.5*branch+0.1, patch+0.5, fOccVecFine->GetArray()[i]);
1085 : }
1086 0 : }
1087 :
1088 : //____________________________________________________________________________________________
1089 : void AliTPCdataQA::ResetProfiles()
1090 : {
1091 0 : if(!fIsDQM)
1092 0 : AliInfo("Method only meaningful for DQM");
1093 :
1094 0 : if(fHistQVsSector)
1095 0 : fHistQVsSector->Reset();
1096 0 : if(fHistQmaxVsSector)
1097 0 : fHistQmaxVsSector->Reset();
1098 0 : if(fHistOccVsSector)
1099 0 : fHistOccVsSector->Reset();
1100 0 : if(fHistOcc2dVsSector)
1101 0 : fHistOcc2dVsSector->Reset();
1102 :
1103 0 : if(fOccVec)
1104 0 : for(Int_t i = 0; i < 72; i++)
1105 0 : fOccVec->GetArray()[i] = 0.0;
1106 0 : if(fOccVecFine)
1107 0 : for(Int_t i = 0; i < 36*12; i++)
1108 0 : fOccVecFine->GetArray()[i] = 0.0;
1109 0 : }
1110 :
1111 : //____________________________________________________________________________________________
1112 : void AliTPCdataQA::Init()
1113 : {
1114 : /// Define the calibration objects the first time Update is called
1115 : /// NB! This has to be done first even if the data is rejected by the time
1116 : /// cut to make sure that the objects are available in Analyse
1117 :
1118 : /// 2015-11-09 arrays for expending the data are not allocated in advance
1119 : /// but when they are needed. That saves about 100 MByte of memory.
1120 :
1121 0 : if(!fIsDQM) {
1122 :
1123 0 : if (!fNLocalMaxima){
1124 0 : TObjArray configArr(72);
1125 0 : fNLocalMaxima = new AliTPCCalPad(ConfigArrRocs(&configArr,"NLocalMaxima"));
1126 0 : fMaxCharge = new AliTPCCalPad(ConfigArrRocs(&configArr,"MaxCharge"));
1127 0 : fMeanCharge = new AliTPCCalPad(ConfigArrRocs(&configArr,"MeanCharge"));
1128 0 : fNoThreshold = new AliTPCCalPad(ConfigArrRocs(&configArr,"NoThreshold"));
1129 0 : fNTimeBins = new AliTPCCalPad(ConfigArrRocs(&configArr,"NTimeBins"));
1130 0 : fNPads = new AliTPCCalPad(ConfigArrRocs(&configArr,"NPads"));
1131 0 : fTimePosition = new AliTPCCalPad(ConfigArrRocs(&configArr,"TimePosition"));
1132 0 : fOverThreshold10 = new AliTPCCalPad(ConfigArrRocs(&configArr,"OverThreshold10"));
1133 0 : fOverThreshold20 = new AliTPCCalPad(ConfigArrRocs(&configArr,"OverThreshold20"));
1134 0 : fOverThreshold30 = new AliTPCCalPad(ConfigArrRocs(&configArr,"OverThreshold30"));
1135 :
1136 0 : fHistQVsTimeSideA = new TProfile("hQVsTimeSideA", "Q vs time (side A); Time [Timebin]; Q [ADC ch]", 100, 0, 1000);
1137 0 : fHistQVsTimeSideA->SetDirectory(0);
1138 0 : fHistQVsTimeSideC = new TProfile("hQVsTimeSideC", "Q vs time (side C); Time [Timebin]; Q [ADC ch]", 100, 0, 1000);
1139 0 : fHistQVsTimeSideC->SetDirectory(0);
1140 0 : fHistQMaxVsTimeSideA = new TProfile("hQMaxVsTimeSideA", "Q_{MAX} vs time (side A); Time [Timebin]; Q_{MAX} [ADC ch]", 100, 0, 1000);
1141 0 : fHistQMaxVsTimeSideA->SetDirectory(0);
1142 0 : fHistQMaxVsTimeSideC = new TProfile("hQMaxVsTimeSideC", "Q_{MAX} vs time (side C); Time [Timebin]; Q_{MAX} [ADC ch]", 100, 0, 1000);
1143 0 : fHistQMaxVsTimeSideC->SetDirectory(0);
1144 0 : }
1145 : } else { // DQM histograms and array
1146 :
1147 0 : if (!fHistOccVsSector) {
1148 0 : fHistOccVsSector = new TProfile("hOccVsSector", "Occupancy vs sector; Sector; Occupancy", 72, 0, 72);
1149 0 : fHistOccVsSector->SetDirectory(0);
1150 :
1151 0 : fHistOcc2dVsSector = new TProfile2D("hOcc2dVsSector", "Occupancy vs sector and patch; Sector; Patch", 72, 0, 36, 6, 0, 6);
1152 0 : fHistOcc2dVsSector->SetDirectory(0);
1153 :
1154 0 : fHistQVsSector = new TProfile("hQVsSector", "Q vs sector; Sector; Q [ADC ch]", 72, 0, 72);
1155 0 : fHistQVsSector->SetDirectory(0);
1156 :
1157 0 : fHistQmaxVsSector = new TProfile("hQmaxVsSector", "Qmax vs sector; Sector; Qmax [ADC ch]", 72, 0, 72);
1158 0 : fHistQmaxVsSector->SetDirectory(0);
1159 :
1160 0 : fOccVec = new TArrayD(72);
1161 0 : for(Int_t i = 0; i < 72; i++)
1162 0 : fOccVec->GetArray()[i] = 0;
1163 :
1164 0 : fOccMaxVec = new TArrayD(72);
1165 0 : const Double_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
1166 0 : for(Int_t i = 0; i < 72; i++)
1167 :
1168 0 : if(i<36) // IROCs (5504 pads)
1169 0 : fOccMaxVec->GetArray()[i] = nTimeBins*5504;
1170 : else // OROCs (9984 pads)
1171 0 : fOccMaxVec->GetArray()[i] = nTimeBins*9984;
1172 :
1173 : // 12 branches for each full sector
1174 : const Int_t nBranches = 36*12;
1175 0 : fOccVecFine = new TArrayD(nBranches);
1176 0 : for(Int_t i = 0; i < nBranches; i++)
1177 0 : fOccVecFine->GetArray()[i] = 0;
1178 :
1179 : // Pads per patch same for all sectors
1180 0 : Int_t nPads0[6] = {1152, 1536, 1152, 1280, 1280, 1280};
1181 0 : Int_t nPads1[6] = {1152, 1664, 1152, 1280, 1280, 1280};
1182 :
1183 0 : fOccMaxVecFine = new TArrayD(nBranches);
1184 0 : for(Int_t i = 0; i < nBranches; i++) {
1185 :
1186 0 : const Int_t fullSector = Int_t(i/12);
1187 0 : Int_t branch = i - fullSector*12;
1188 0 : R__ASSERT(branch>=0 && branch<12);
1189 :
1190 0 : const Int_t patch = Int_t(branch/2);
1191 0 : branch -= patch*2;
1192 :
1193 0 : R__ASSERT(branch>=0 && branch<2);
1194 0 : if(branch == 0)
1195 0 : fOccMaxVecFine->GetArray()[i] = nTimeBins*nPads0[patch];
1196 : else // OROCs (9984 pads)
1197 0 : fOccMaxVecFine->GetArray()[i] = nTimeBins*nPads1[patch];
1198 : }
1199 0 : }
1200 : }
1201 :
1202 : //
1203 : // If Analyse has been previously called we need now to denormalize the data
1204 : // as more data is coming
1205 : //
1206 0 : if(fIsAnalysed == kTRUE && !fIsDQM) {
1207 :
1208 0 : const Int_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
1209 0 : const Float_t denormalization = Float_t(fEventCounter * nTimeBins);
1210 0 : fNoThreshold->Multiply(denormalization);
1211 :
1212 0 : fMeanCharge->Multiply(fNLocalMaxima);
1213 0 : fMaxCharge->Multiply(fNLocalMaxima);
1214 0 : fNTimeBins->Multiply(fNLocalMaxima);
1215 0 : fNPads->Multiply(fNLocalMaxima);
1216 0 : fTimePosition->Multiply(fNLocalMaxima);
1217 0 : fIsAnalysed = kFALSE;
1218 0 : }
1219 0 : }
1220 :
1221 : //____________________________________________________________________________________________
1222 : void AliTPCdataQA::ResetData()
1223 : {
1224 : /// reset all data
1225 :
1226 0 : if(!fIsDQM) {
1227 :
1228 0 : if (fNLocalMaxima){
1229 0 : fNoThreshold->Reset();
1230 0 : fNLocalMaxima->Reset();
1231 0 : fMeanCharge->Reset();
1232 0 : fMaxCharge->Reset();
1233 0 : fNTimeBins->Reset();
1234 0 : fNPads->Reset();
1235 0 : fTimePosition->Reset();
1236 0 : fOverThreshold10->Reset();
1237 0 : fOverThreshold20->Reset();
1238 0 : fOverThreshold30->Reset();
1239 :
1240 0 : fHistQVsTimeSideA->Reset();
1241 0 : fHistQVsTimeSideC->Reset();
1242 0 : fHistQMaxVsTimeSideA->Reset();
1243 0 : fHistQMaxVsTimeSideC->Reset();
1244 :
1245 0 : fIsAnalysed = kFALSE;
1246 :
1247 0 : }
1248 : }
1249 :
1250 0 : fEventCounter=0;
1251 0 : fClusterCounter=0;
1252 0 : }
1253 :
1254 : TObjArray *AliTPCdataQA::ConfigArrRocs(TObjArray *arr, const Text_t* name)
1255 : {
1256 : /// GetArray with confiured ROCs
1257 :
1258 0 : arr->Clear();
1259 0 : arr->SetName(name);
1260 0 : for (Int_t i=0; i<72; ++i){
1261 0 : if (fActiveChambers[i]) arr->AddAt(new AliTPCCalROC(i),i);
1262 : }
1263 0 : return arr;
1264 0 : }
|