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