]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TPC/AliTPCdataQA.cxx
Remove compilation warnings (Marian)
[u/mrichter/AliRoot.git] / TPC / AliTPCdataQA.cxx
CommitLineData
0ffacf98 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
17/* $Id$ */
18
266f8637 19/*
20 February 2008
21
22 The code has been heavily modified so that now the RAW data is
23 "expanded" for each sector and stored in a big signal array. Then a
24 simple version of the code in AliTPCclustererMI is used to identify
25 the local maxima and these are then used for QA. This gives a better
26 estimate of the charge (both max and total) and also limits the
27 effect of noise.
28
29 Implementation:
30
31 In Update the RAW signals >= 3 ADC channels are stored in the arrays.
32
33 There are 3 arrays:
34 Float_t** fAllBins 2d array [row][bin(pad, time)] ADC signal
35 Int_t** fAllSigBins 2d array [row][signal#] bin(with signal)
36 Int_t* fAllNSigBins; 1d array [row] Nsignals
37
38 This is done sector by sector.
39
40 When data from a new sector is encountered, the method
41 FindLocalMaxima is called on the data from the previous sector, and
42 the calibration/data objects are updated with the "cluster"
43 info. Finally the arrays are cleared.
44
45 The requirements for a local maxima is:
46 Charge in bin is >= 5 ADC channels.
47 Charge in bin is larger than all the 8 neighboring bins.
48 At least one of the two pad neighbors has a signal.
49 At least one of the two time neighbors has a signal.
50
51 Before accessing the data it is expected that the Analyse method is
52 called. This normalizes some of the data objects to per event or per
53 cluster.
54 If more data is passed to the class after Analyse has been called
55 the normalization is reversed and Analyse has to be called again.
56*/
57
0ffacf98 58
f11b3071 59// stl includes
60#include <iostream>
61
62using namespace std;
63
0ffacf98 64//Root includes
65#include <TH1F.h>
66#include <TH2F.h>
67#include <TString.h>
68#include <TMath.h>
69#include <TF1.h>
70#include <TRandom.h>
71#include <TDirectory.h>
72#include <TFile.h>
266f8637 73#include <TError.h>
0ffacf98 74//AliRoot includes
75#include "AliRawReader.h"
76#include "AliRawReaderRoot.h"
77#include "AliRawReaderDate.h"
78#include "AliTPCRawStream.h"
79#include "AliTPCCalROC.h"
80#include "AliTPCROC.h"
81#include "AliMathBase.h"
82#include "TTreeStream.h"
83#include "AliTPCRawStreamFast.h"
84
85//date
86#include "event.h"
87#include "AliTPCCalPad.h"
258cd111 88#include "AliTPCPreprocessorOnline.h"
0ffacf98 89
90//header file
91#include "AliTPCdataQA.h"
92
0ffacf98 93ClassImp(AliTPCdataQA)
94
336156cc 95AliTPCdataQA::AliTPCdataQA() : /*FOLD00*/
96 TH1F("TPCRAW","TPCRAW",100,0,100),
0ffacf98 97 fFirstTimeBin(60),
98 fLastTimeBin(1000),
99 fAdcMin(1),
100 fAdcMax(100),
101 fOldRCUformat(kTRUE),
0ffacf98 102 fMapping(NULL),
f11b3071 103 fPedestal(0),
104 fNoise(0),
266f8637 105 fNLocalMaxima(0),
0ffacf98 106 fMaxCharge(0),
f11b3071 107 fMeanCharge(0),
108 fNoThreshold(0),
266f8637 109 fNTimeBins(0),
110 fNPads(0),
111 fTimePosition(0),
c94a79e1 112 fOverThreshold10(0),
113 fOverThreshold20(0),
114 fOverThreshold30(0),
266f8637 115 fEventCounter(0),
116 fIsAnalysed(kFALSE),
117 fAllBins(0),
118 fAllSigBins(0),
119 fAllNSigBins(0),
120 fRowsMax(0),
121 fPadsMax(0),
122 fTimeBinsMax(0)
0ffacf98 123{
124 //
125 // default constructor
126 //
127}
128
129
130//_____________________________________________________________________
131AliTPCdataQA::AliTPCdataQA(const AliTPCdataQA &ped) : /*FOLD00*/
336156cc 132 TH1F(ped),
0ffacf98 133 fFirstTimeBin(ped.GetFirstTimeBin()),
134 fLastTimeBin(ped.GetLastTimeBin()),
135 fAdcMin(ped.GetAdcMin()),
136 fAdcMax(ped.GetAdcMax()),
266f8637 137 fOldRCUformat(ped.GetOldRCUformat()),
138 fMapping(NULL),
139 fPedestal(0),
140 fNoise(0),
141 fNLocalMaxima(0),
142 fMaxCharge(0),
143 fMeanCharge(0),
144 fNoThreshold(0),
266f8637 145 fNTimeBins(0),
146 fNPads(0),
147 fTimePosition(0),
c94a79e1 148 fOverThreshold10(0),
149 fOverThreshold20(0),
150 fOverThreshold30(0),
266f8637 151 fEventCounter(ped.GetEventCounter()),
152 fIsAnalysed(ped.GetIsAnalysed()),
153 fAllBins(0),
154 fAllSigBins(0),
155 fAllNSigBins(0),
156 fRowsMax(0),
157 fPadsMax(0),
158 fTimeBinsMax(0)
0ffacf98 159{
160 //
161 // copy constructor
162 //
266f8637 163 if(ped.GetNLocalMaxima())
164 fNLocalMaxima = new AliTPCCalPad(*ped.GetNLocalMaxima());
165 if(ped.GetMaxCharge())
166 fMaxCharge = new AliTPCCalPad(*ped.GetMaxCharge());
167 if(ped.GetMeanCharge())
168 fMeanCharge = new AliTPCCalPad(*ped.GetMeanCharge());
169 if(ped.GetNoThreshold())
170 fNoThreshold = new AliTPCCalPad(*ped.GetNoThreshold());
171 if(ped.GetNTimeBins())
172 fNTimeBins = new AliTPCCalPad(*ped.GetNTimeBins());
173 if(ped.GetNPads())
174 fNPads = new AliTPCCalPad(*ped.GetNPads());
175 if(ped.GetTimePosition())
176 fTimePosition = new AliTPCCalPad(*ped.GetTimePosition());
177 if(ped.GetOverThreshold10())
178 fOverThreshold10 = new AliTPCCalPad(*ped.GetOverThreshold10());
179 if(ped.GetOverThreshold20())
180 fOverThreshold20 = new AliTPCCalPad(*ped.GetOverThreshold20());
181 if(ped.GetOverThreshold30())
182 fOverThreshold30 = new AliTPCCalPad(*ped.GetOverThreshold30());
0ffacf98 183}
184
185
186//_____________________________________________________________________
187AliTPCdataQA& AliTPCdataQA::operator = (const AliTPCdataQA &source)
188{
189 //
190 // assignment operator
191 //
192 if (&source == this) return *this;
193 new (this) AliTPCdataQA(source);
194
195 return *this;
196}
197
198
199//_____________________________________________________________________
200AliTPCdataQA::~AliTPCdataQA() /*FOLD00*/
201{
202 //
203 // destructor
204 //
205
206 // do not delete fMapping, because we do not own it.
266f8637 207 // do not delete fMapping, because we do not own it.
208 // do not delete fNoise and fPedestal, because we do not own them.
209
210 delete fNLocalMaxima;
211 delete fMaxCharge;
212 delete fMeanCharge;
213 delete fNoThreshold;
214 delete fNTimeBins;
215 delete fNPads;
216 delete fTimePosition;
217 delete fOverThreshold10;
218 delete fOverThreshold20;
219 delete fOverThreshold30;
220
221 for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
222 delete [] fAllBins[iRow];
223 delete [] fAllSigBins[iRow];
224 }
225 delete [] fAllBins;
226 delete [] fAllSigBins;
227 delete [] fAllNSigBins;
0ffacf98 228}
229
0ffacf98 230//_____________________________________________________________________
231Bool_t AliTPCdataQA::ProcessEventFast(AliTPCRawStreamFast *rawStreamFast)
232{
233 //
234 // Event Processing loop - AliTPCRawStream
235 //
236 Bool_t withInput = kFALSE;
266f8637 237 Int_t nSignals = 0;
238 Int_t lastSector = -1;
0ffacf98 239
240 while ( rawStreamFast->NextDDL() ){
241 while ( rawStreamFast->NextChannel() ){
266f8637 242
243 Int_t iSector = rawStreamFast->GetSector(); // current sector
244 Int_t iRow = rawStreamFast->GetRow(); // current row
245 Int_t iPad = rawStreamFast->GetPad(); // current pad
246 // Call local maxima finder if the data is in a new sector
247 if(iSector != lastSector) {
248
249 if(nSignals>0)
250 FindLocalMaxima(lastSector);
251
252 CleanArrays();
253 lastSector = iSector;
254 nSignals = 0;
255 }
256
257 while ( rawStreamFast->NextBunch() ){
258 Int_t startTbin = (Int_t)rawStreamFast->GetStartTimeBin();
259 Int_t endTbin = (Int_t)rawStreamFast->GetEndTimeBin();
260
261 for (Int_t iTimeBin = startTbin; iTimeBin < endTbin; iTimeBin++){
262 Float_t signal = rawStreamFast->GetSignals()[iTimeBin-startTbin];
263 nSignals += Update(iSector,iRow,iPad,iTimeBin+1,signal);
264 withInput = kTRUE;
0ffacf98 265 }
266f8637 266 }
0ffacf98 267 }
268 }
266f8637 269
0ffacf98 270 return withInput;
271}
272//_____________________________________________________________________
273Bool_t AliTPCdataQA::ProcessEventFast(AliRawReader *rawReader)
274{
275 //
276 // Event processing loop - AliRawReader
277 //
278 AliTPCRawStreamFast *rawStreamFast = new AliTPCRawStreamFast(rawReader, (AliAltroMapping**)fMapping);
279 Bool_t res=ProcessEventFast(rawStreamFast);
c75bf2f1 280 if(res)
281 fEventCounter++; // only increment event counter if there is TPC data
282 // otherwise Analyse (called in QA) fails
283
0ffacf98 284 delete rawStreamFast;
285 return res;
286}
287
288//_____________________________________________________________________
289Bool_t AliTPCdataQA::ProcessEvent(AliTPCRawStream *rawStream)
290{
291 //
292 // Event Processing loop - AliTPCRawStream
293 //
294
295 rawStream->SetOldRCUFormat(fOldRCUformat);
296
297 Bool_t withInput = kFALSE;
266f8637 298 Int_t nSignals = 0;
299 Int_t lastSector = -1;
0ffacf98 300
301 while (rawStream->Next()) {
302
303 Int_t iSector = rawStream->GetSector(); // current ROC
304 Int_t iRow = rawStream->GetRow(); // current row
305 Int_t iPad = rawStream->GetPad(); // current pad
306 Int_t iTimeBin = rawStream->GetTime(); // current time bin
307 Float_t signal = rawStream->GetSignal(); // current ADC signal
266f8637 308
309 // Call local maxima finder if the data is in a new sector
310 if(iSector != lastSector) {
311
312 if(nSignals>0)
313 FindLocalMaxima(lastSector);
314
315 CleanArrays();
316 lastSector = iSector;
317 nSignals = 0;
318 }
0ffacf98 319
266f8637 320 // Sometimes iRow==-1 if there is problems to read the data
321 if(iRow>=0) {
322 nSignals += Update(iSector,iRow,iPad,iTimeBin,signal);
323 withInput = kTRUE;
324 }
0ffacf98 325 }
326
327 return withInput;
328}
329
330
331//_____________________________________________________________________
332Bool_t AliTPCdataQA::ProcessEvent(AliRawReader *rawReader)
333{
334 //
335 // Event processing loop - AliRawReader
336 //
337
338 // if fMapping is NULL the rawstream will crate its own mapping
339 AliTPCRawStream rawStream(rawReader, (AliAltroMapping**)fMapping);
340 rawReader->Select("TPC");
c75bf2f1 341 Bool_t res = ProcessEvent(&rawStream);
342
343 if(res)
344 fEventCounter++; // only increment event counter if there is TPC data
345 // otherwise Analyse (called in QA) fails
346 return res;
0ffacf98 347}
348
349
350//_____________________________________________________________________
351Bool_t AliTPCdataQA::ProcessEvent(eventHeaderStruct *event)
352{
353 //
354 // process date event
355 //
356
357 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
358 Bool_t result=ProcessEvent(rawReader);
359 delete rawReader;
360 return result;
361}
362
363
364
365//_____________________________________________________________________
366void AliTPCdataQA::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
367{
368 //
369 // Write class to file
370 //
371
372 TString sDir(dir);
373 TString option;
374
375 if ( append )
376 option = "update";
377 else
378 option = "recreate";
379
380 TDirectory *backup = gDirectory;
381 TFile f(filename,option.Data());
382 f.cd();
383 if ( !sDir.IsNull() ){
384 f.mkdir(sDir.Data());
385 f.cd(sDir);
386 }
387 this->Write();
388 f.Close();
389
390 if ( backup ) backup->cd();
391}
392
393
394//_____________________________________________________________________
266f8637 395Int_t AliTPCdataQA::Update(const Int_t iSector, /*FOLD00*/
396 const Int_t iRow,
397 const Int_t iPad,
398 const Int_t iTimeBin,
399 Float_t signal)
0ffacf98 400{
401 //
402 // Signal filling method
403 //
266f8637 404
405 //
406 // Define the calibration objects the first time Update is called
407 // NB! This has to be done first even if the data is rejected by the time
408 // cut to make sure that the objects are available in Analyse
409 //
410 if (!fNLocalMaxima) fNLocalMaxima = new AliTPCCalPad("NLocalMaxima","NLocalMaxima");
0ffacf98 411 if (!fMaxCharge) fMaxCharge = new AliTPCCalPad("MaxCharge","MaxCharge");
f11b3071 412 if (!fMeanCharge) fMeanCharge = new AliTPCCalPad("MeanCharge","MeanCharge");
413 if (!fNoThreshold) fNoThreshold = new AliTPCCalPad("NoThreshold","NoThreshold");
266f8637 414 if (!fNTimeBins) fNTimeBins = new AliTPCCalPad("NTimeBins","NTimeBins");
415 if (!fNPads) fNPads = new AliTPCCalPad("NPads","NPads");
416 if (!fTimePosition) fTimePosition = new AliTPCCalPad("TimePosition","TimePosition");
0ffacf98 417 if (!fOverThreshold10) fOverThreshold10 = new AliTPCCalPad("OverThreshold10","OverThreshold10");
418 if (!fOverThreshold20) fOverThreshold20 = new AliTPCCalPad("OverThreshold20","OverThreshold20");
419 if (!fOverThreshold30) fOverThreshold30 = new AliTPCCalPad("OverThreshold30","OverThreshold30");
266f8637 420 // Make the arrays for expanding the data
f11b3071 421
266f8637 422 if (!fAllBins)
423 MakeArrays();
424
425 //
426 // If Analyse has been previously called we need now to denormalize the data
427 // as more data is coming
428 //
429 if(fIsAnalysed == kTRUE) {
430
431 const Int_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
432 const Float_t denormalization = Float_t(fEventCounter * nTimeBins);
433 fNoThreshold->Multiply(denormalization);
434
435 fMeanCharge->Multiply(fNLocalMaxima);
436 fMaxCharge->Multiply(fNLocalMaxima);
437 fNTimeBins->Multiply(fNLocalMaxima);
438 fNPads->Multiply(fNLocalMaxima);
439 fTimePosition->Multiply(fNLocalMaxima);
440 fIsAnalysed = kFALSE;
441 }
f11b3071 442
266f8637 443 //
444 // TimeBin cut
445 //
446 if (iTimeBin<fFirstTimeBin) return 0;
447 if (iTimeBin>fLastTimeBin) return 0;
448
f11b3071 449 // if pedestal calibrations are loaded subtract pedestals
450 if(fPedestal) {
451
266f8637 452 Float_t ped = fPedestal->GetCalROC(iSector)->GetValue(iRow, iPad);
453 // Don't use data from pads where pedestals are abnormally small or large
454 if(ped<10 || ped>90)
f11b3071 455 return 0;
266f8637 456 signal -= ped;
f11b3071 457 }
266f8637 458
459 // In fNoThreshold we fill all data to estimate the ZS volume
460 Float_t count = fNoThreshold->GetCalROC(iSector)->GetValue(iRow, iPad);
461 fNoThreshold->GetCalROC(iSector)->SetValue(iRow, iPad,count+1);
462
f11b3071 463 // Require at least 3 ADC channels
266f8637 464 if (signal < 3.0)
f11b3071 465 return 0;
466
467 // if noise calibrations are loaded require at least 3*sigmaNoise
468 if(fNoise) {
266f8637 469
470 Float_t noise = fNoise->GetCalROC(iSector)->GetValue(iRow, iPad);
471
472 if(signal < noise*3.0)
f11b3071 473 return 0;
474 }
266f8637 475
f11b3071 476 //
266f8637 477 // This signal is ok and we store it in the cluster map
0ffacf98 478 //
f11b3071 479
266f8637 480 SetExpandDigit(iRow, iPad, iTimeBin, signal);
481
482 return 1; // signal was accepted
f11b3071 483}
266f8637 484
f11b3071 485//_____________________________________________________________________
266f8637 486void AliTPCdataQA::FindLocalMaxima(const Int_t iSector)
f11b3071 487{
488 //
266f8637 489 // This method is called after the data from each sector has been
490 // exapanded into an array
491 // Loop over the signals and identify local maxima and fill the
492 // calibration objects with the information
f11b3071 493 //
266f8637 494
495 Int_t nLocalMaxima = 0;
496 const Int_t maxTimeBin = fTimeBinsMax+4; // Used to step between neighboring pads
497 // Because we have tha pad-time data in a
498 // 1d array
499
500 for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
501
502 Float_t* allBins = fAllBins[iRow];
503 Int_t* sigBins = fAllSigBins[iRow];
504 const Int_t nSigBins = fAllNSigBins[iRow];
505
506 for (Int_t iSig = 0; iSig < nSigBins; iSig++) {
507
508 Int_t bin = sigBins[iSig];
509 Float_t *b = &allBins[bin];
510
511 //
512 // Now we check if this is a local maximum
513 //
514
515 Float_t qMax = b[0];
516
517 // First check that the charge is bigger than the threshold
518 if (qMax<5)
519 continue;
520
521 // Require at least one neighboring pad with signal
522 if (b[-maxTimeBin]+b[maxTimeBin]<=0) continue;
523
524 // Require at least one neighboring time bin with signal
525 if (b[-1]+b[1]<=0) continue;
526
527 //
528 // Check that this is a local maximum
529 // Note that the checking is done so that if 2 charges has the same
530 // qMax then only 1 cluster is generated
531 // (that is why there is BOTH > and >=)
532 //
533 if (b[-maxTimeBin] >= qMax) continue;
534 if (b[-1 ] >= qMax) continue;
535 if (b[+maxTimeBin] > qMax) continue;
536 if (b[+1 ] > qMax) continue;
537 if (b[-maxTimeBin-1] >= qMax) continue;
538 if (b[+maxTimeBin-1] >= qMax) continue;
539 if (b[+maxTimeBin+1] > qMax) continue;
540 if (b[-maxTimeBin+1] >= qMax) continue;
541
542 //
543 // Now we accept the local maximum and fill the calibration/data objects
544 //
545 nLocalMaxima++;
546
547 Int_t iPad, iTimeBin;
548 GetPadAndTimeBin(bin, iPad, iTimeBin);
549
550 Float_t count = fNLocalMaxima->GetCalROC(iSector)->GetValue(iRow, iPad);
551 fNLocalMaxima->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
552
553 count = fTimePosition->GetCalROC(iSector)->GetValue(iRow, iPad);
554 fTimePosition->GetCalROC(iSector)->SetValue(iRow, iPad, count+iTimeBin);
555
556 Float_t charge = fMaxCharge->GetCalROC(iSector)->GetValue(iRow, iPad);
557 fMaxCharge->GetCalROC(iSector)->SetValue(iRow, iPad, charge + qMax);
558
559 if(qMax>=10) {
560 count = fOverThreshold10->GetCalROC(iSector)->GetValue(iRow, iPad);
561 fOverThreshold10->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
562 }
563 if(qMax>=20) {
564 count = fOverThreshold20->GetCalROC(iSector)->GetValue(iRow, iPad);
565 fOverThreshold20->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
566 }
567 if(qMax>=30) {
568 count = fOverThreshold30->GetCalROC(iSector)->GetValue(iRow, iPad);
569 fOverThreshold30->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
570 }
571
572 //
573 // Calculate the total charge as the sum over the region:
574 //
575 // o o o o o
576 // o i i i o
577 // o i C i o
578 // o i i i o
579 // o o o o o
580 //
581 // with qmax at the center C.
582 //
583 // The inner charge (i) we always add, but we only add the outer
584 // charge (o) if the neighboring inner bin (i) has a signal.
585 //
586 Int_t minP = 0, maxP = 0, minT = 0, maxT = 0;
587 Float_t qTot = qMax;
588 for(Int_t i = -1; i<=1; i++) {
589 for(Int_t j = -1; j<=1; j++) {
590
591 if(i==0 && j==0)
592 continue;
593
594 Float_t charge = GetQ(b, i, j, maxTimeBin, minT, maxT, minP, maxP);
595 qTot += charge;
596 if(charge>0) {
597 // see if the next neighbor is also above threshold
598 if(i*j==0) {
599 qTot += GetQ(b, 2*i, 2*j, maxTimeBin, minT, maxT, minP, maxP);
600 } else {
601 // we are in a diagonal corner
602 qTot += GetQ(b, i, 2*j, maxTimeBin, minT, maxT, minP, maxP);
603 qTot += GetQ(b, 2*i, j, maxTimeBin, minT, maxT, minP, maxP);
604 qTot += GetQ(b, 2*i, 2*j, maxTimeBin, minT, maxT, minP, maxP);
605 }
606 }
607 }
608 }
609
610 charge = fMeanCharge->GetCalROC(iSector)->GetValue(iRow, iPad);
611 fMeanCharge->GetCalROC(iSector)->SetValue(iRow, iPad, charge + qTot);
612
613 count = fNTimeBins->GetCalROC(iSector)->GetValue(iRow, iPad);
614 fNTimeBins->GetCalROC(iSector)->SetValue(iRow, iPad, count + maxT-minT+1);
615
616 count = fNPads->GetCalROC(iSector)->GetValue(iRow, iPad);
617 fNPads->GetCalROC(iSector)->SetValue(iRow, iPad, count + maxP-minP+1);
618
619 } // end loop over signals
620 } // end loop over rows
f11b3071 621
266f8637 622 // cout << "Number of local maximas found: " << nLocalMaxima << endl;
0ffacf98 623}
11ccf1c1 624
f11b3071 625//_____________________________________________________________________
626void AliTPCdataQA::Analyse()
627{
11ccf1c1 628 //
f11b3071 629 // Calculate calibration constants
11ccf1c1 630 //
f11b3071 631
632 cout << "Analyse called" << endl;
633
266f8637 634 if(fIsAnalysed == kTRUE) {
635
636 cout << "No new data since Analyse was called last time" << endl;
637 return;
638 }
f11b3071 639
266f8637 640 if(fEventCounter==0) {
641
f11b3071 642 cout << "EventCounter == 0, Cannot analyse" << endl;
643 return;
644 }
266f8637 645
f11b3071 646 Int_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
c75bf2f1 647 cout << "EventCounter: " << fEventCounter << endl
f11b3071 648 << "TimeBins: " << nTimeBins << endl;
649
f11b3071 650 Float_t normalization = 1.0 / Float_t(fEventCounter * nTimeBins);
266f8637 651 fNoThreshold->Multiply(normalization);
652
653 fMeanCharge->Divide(fNLocalMaxima);
654 fMaxCharge->Divide(fNLocalMaxima);
655 fNTimeBins->Divide(fNLocalMaxima);
656 fNPads->Divide(fNLocalMaxima);
657 fTimePosition->Divide(fNLocalMaxima);
658
659 fIsAnalysed = kTRUE;
11ccf1c1 660}
258cd111 661
662
266f8637 663//_____________________________________________________________________
258cd111 664void AliTPCdataQA::MakeTree(const char *fname){
665 //
666 // Export result to the tree -located in the file
667 // This file can be analyzed using AliTPCCalibViewer
668 //
258cd111 669 AliTPCPreprocessorOnline preprocesor;
266f8637 670
671 if (fNLocalMaxima) preprocesor.AddComponent(fNLocalMaxima);
672 if (fMaxCharge) preprocesor.AddComponent(fMaxCharge);
673 if (fMeanCharge) preprocesor.AddComponent(fMeanCharge);
674 if (fNoThreshold) preprocesor.AddComponent(fNoThreshold);
675 if (fNTimeBins) preprocesor.AddComponent(fNTimeBins);
676 if (fNPads) preprocesor.AddComponent(fNPads);
677 if (fTimePosition) preprocesor.AddComponent(fTimePosition);
678 if (fOverThreshold10) preprocesor.AddComponent(fOverThreshold10);
679 if (fOverThreshold20) preprocesor.AddComponent(fOverThreshold20);
680 if (fOverThreshold30) preprocesor.AddComponent(fOverThreshold30);
681
258cd111 682 preprocesor.DumpToFile(fname);
683}
c322f08a 684
685
266f8637 686//_____________________________________________________________________
c322f08a 687void AliTPCdataQA::MakeArrays(){
688 //
266f8637 689 // The arrays for expanding the raw data are defined and
690 // som parameters are intialised
c322f08a 691 //
692 AliTPCROC * roc = AliTPCROC::Instance();
693 //
266f8637 694 // To make the array big enough for all sectors we take
695 // the dimensions from the outer row of an OROC (the last sector)
696 //
697 fRowsMax = roc->GetNRows(roc->GetNSector()-1);
698 fPadsMax = roc->GetNPads(roc->GetNSector()-1,fRowsMax-1);
699 fTimeBinsMax = fLastTimeBin - fFirstTimeBin +1;
700
701 //
702 // Since we have added 2 pads (TimeBins) before and after the real pads (TimeBins)
703 // to make sure that we can always query the exanded table even when the
704 // max is on the edge
705 //
706
c322f08a 707
266f8637 708 fAllBins = new Float_t*[fRowsMax];
709 fAllSigBins = new Int_t*[fRowsMax];
710 fAllNSigBins = new Int_t[fRowsMax];
711
712 for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
c322f08a 713 //
266f8637 714 Int_t maxBin = (fTimeBinsMax+4)*(fPadsMax+4);
c322f08a 715 fAllBins[iRow] = new Float_t[maxBin];
266f8637 716 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin); // set all values to 0
c322f08a 717 fAllSigBins[iRow] = new Int_t[maxBin];
266f8637 718 fAllNSigBins[iRow] = 0;
c322f08a 719 }
720}
721
722
266f8637 723//_____________________________________________________________________
c322f08a 724void AliTPCdataQA::CleanArrays(){
725 //
726 //
727 //
266f8637 728
729 for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
c322f08a 730 //
266f8637 731 Int_t maxBin = (fTimeBinsMax+4)*(fPadsMax+4);
c322f08a 732 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
733 fAllNSigBins[iRow]=0;
734 }
735}
736
266f8637 737//_____________________________________________________________________
738void AliTPCdataQA::GetPadAndTimeBin(Int_t bin, Int_t& iPad, Int_t& iTimeBin){
739 //
740 // Return pad and timebin for a given bin
c322f08a 741 //
266f8637 742
743 // Int_t bin = iPad*(fTimeBinsMax+4)+iTimeBin;
744 iTimeBin = bin%(fTimeBinsMax+4);
745 iPad = (bin-iTimeBin)/(fTimeBinsMax+4);
746
747 iPad -= 2;
748 iTimeBin -= 2;
749 iTimeBin += fFirstTimeBin;
750
751 R__ASSERT(iPad>=0 && iPad<=fPadsMax);
752 R__ASSERT(iTimeBin>=fFirstTimeBin && iTimeBin<=fLastTimeBin);
753}
754
755//_____________________________________________________________________
756void AliTPCdataQA::SetExpandDigit(const Int_t iRow, Int_t iPad,
757 Int_t iTimeBin, const Float_t signal)
758{
c322f08a 759 //
266f8637 760 //
c322f08a 761 //
266f8637 762 R__ASSERT(iRow>=0 && iRow<fRowsMax);
763 R__ASSERT(iPad>=0 && iPad<=fPadsMax);
764 R__ASSERT(iTimeBin>=fFirstTimeBin && iTimeBin<=fLastTimeBin);
765
766 iTimeBin -= fFirstTimeBin;
767 iPad += 2;
768 iTimeBin += 2;
c322f08a 769
266f8637 770 Int_t bin = iPad*(fTimeBinsMax+4)+iTimeBin;
771 fAllBins[iRow][bin] = signal;
772 fAllSigBins[iRow][fAllNSigBins[iRow]] = bin;
773 fAllNSigBins[iRow]++;
774}
775
776Float_t AliTPCdataQA::GetQ(const Float_t* adcArray, const Int_t time,
777 const Int_t pad, const Int_t maxTimeBins,
778 Int_t& timeMin, Int_t& timeMax,
779 Int_t& padMin, Int_t& padMax)
780{
781 //
782 // This methods return the charge in the bin time+pad*maxTimeBins
783 // If the charge is above 0 it also updates the padMin, padMax, timeMin
784 // and timeMax if necessary
785 //
786 Float_t charge = adcArray[time + pad*maxTimeBins];
787 if(charge > 0) {
788 timeMin = TMath::Min(time, timeMin); timeMax = TMath::Max(time, timeMax);
789 padMin = TMath::Min(pad, padMin); padMax = TMath::Max(pad, padMax);
790 }
791 return charge;
c322f08a 792}