1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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
31 In Update the RAW signals >= 3 ADC channels are stored in the 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
38 This is done sector by sector.
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.
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.
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
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.
71 #include <TDirectory.h>
76 #include "AliRawReader.h"
77 #include "AliRawReaderRoot.h"
78 #include "AliRawReaderDate.h"
79 #include "AliTPCRawStream.h"
80 #include "AliTPCCalROC.h"
81 #include "AliTPCROC.h"
82 #include "AliMathBase.h"
83 #include "TTreeStream.h"
84 #include "AliTPCRawStreamFast.h"
88 #include "AliTPCCalPad.h"
89 #include "AliTPCPreprocessorOnline.h"
92 #include "AliTPCdataQA.h"
94 ClassImp(AliTPCdataQA)
96 AliTPCdataQA::AliTPCdataQA() : /*FOLD00*/
97 TH1F("TPCRAW","TPCRAW",100,0,100),
125 // default constructor
130 //_____________________________________________________________________
131 AliTPCdataQA::AliTPCdataQA(const AliTPCdataQA &ped) : /*FOLD00*/
133 fFirstTimeBin(ped.GetFirstTimeBin()),
134 fLastTimeBin(ped.GetLastTimeBin()),
135 fAdcMin(ped.GetAdcMin()),
136 fAdcMax(ped.GetAdcMax()),
150 fEventCounter(ped.GetEventCounter()),
151 fIsAnalysed(ped.GetIsAnalysed()),
162 if(ped.GetNLocalMaxima())
163 fNLocalMaxima = new AliTPCCalPad(*ped.GetNLocalMaxima());
164 if(ped.GetMaxCharge())
165 fMaxCharge = new AliTPCCalPad(*ped.GetMaxCharge());
166 if(ped.GetMeanCharge())
167 fMeanCharge = new AliTPCCalPad(*ped.GetMeanCharge());
168 if(ped.GetNoThreshold())
169 fNoThreshold = new AliTPCCalPad(*ped.GetNoThreshold());
170 if(ped.GetNTimeBins())
171 fNTimeBins = new AliTPCCalPad(*ped.GetNTimeBins());
173 fNPads = new AliTPCCalPad(*ped.GetNPads());
174 if(ped.GetTimePosition())
175 fTimePosition = new AliTPCCalPad(*ped.GetTimePosition());
176 if(ped.GetOverThreshold10())
177 fOverThreshold10 = new AliTPCCalPad(*ped.GetOverThreshold10());
178 if(ped.GetOverThreshold20())
179 fOverThreshold20 = new AliTPCCalPad(*ped.GetOverThreshold20());
180 if(ped.GetOverThreshold30())
181 fOverThreshold30 = new AliTPCCalPad(*ped.GetOverThreshold30());
184 //_____________________________________________________________________
185 AliTPCdataQA::AliTPCdataQA(TMap *config) : /*FOLD00*/
186 TH1F("TPCRAW","TPCRAW",100,0,100),
214 // default constructor
216 if (config->GetValue("FirstTimeBin")) fFirstTimeBin = ((TObjString*)config->GetValue("FirstTimeBin"))->GetString().Atoi();
217 if (config->GetValue("LastTimeBin")) fLastTimeBin = ((TObjString*)config->GetValue("LastTimeBin"))->GetString().Atoi();
218 if (config->GetValue("AdcMin")) fAdcMin = ((TObjString*)config->GetValue("AdcMin"))->GetString().Atoi();
219 if (config->GetValue("AdcMax")) fAdcMax = ((TObjString*)config->GetValue("AdcMax"))->GetString().Atoi();
222 //_____________________________________________________________________
223 AliTPCdataQA& AliTPCdataQA::operator = (const AliTPCdataQA &source)
226 // assignment operator
228 if (&source == this) return *this;
229 new (this) AliTPCdataQA(source);
235 //_____________________________________________________________________
236 AliTPCdataQA::~AliTPCdataQA() /*FOLD00*/
242 // do not delete fMapping, because we do not own it.
243 // do not delete fMapping, because we do not own it.
244 // do not delete fNoise and fPedestal, because we do not own them.
246 delete fNLocalMaxima;
252 delete fTimePosition;
253 delete fOverThreshold10;
254 delete fOverThreshold20;
255 delete fOverThreshold30;
257 for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
258 delete [] fAllBins[iRow];
259 delete [] fAllSigBins[iRow];
262 delete [] fAllSigBins;
263 delete [] fAllNSigBins;
266 //_____________________________________________________________________
267 Bool_t AliTPCdataQA::ProcessEventFast(AliTPCRawStreamFast *rawStreamFast)
270 // Event Processing loop - AliTPCRawStream
272 Bool_t withInput = kFALSE;
274 Int_t lastSector = -1;
276 while ( rawStreamFast->NextDDL() ){
277 while ( rawStreamFast->NextChannel() ){
279 Int_t iSector = rawStreamFast->GetSector(); // current sector
280 Int_t iRow = rawStreamFast->GetRow(); // current row
281 Int_t iPad = rawStreamFast->GetPad(); // current pad
282 // Call local maxima finder if the data is in a new sector
283 if(iSector != lastSector) {
286 FindLocalMaxima(lastSector);
289 lastSector = iSector;
293 while ( rawStreamFast->NextBunch() ){
294 Int_t startTbin = (Int_t)rawStreamFast->GetStartTimeBin();
295 Int_t endTbin = (Int_t)rawStreamFast->GetEndTimeBin();
297 for (Int_t iTimeBin = startTbin; iTimeBin < endTbin; iTimeBin++){
298 Float_t signal = rawStreamFast->GetSignals()[iTimeBin-startTbin];
299 nSignals += Update(iSector,iRow,iPad,iTimeBin+1,signal);
308 //_____________________________________________________________________
309 Bool_t AliTPCdataQA::ProcessEventFast(AliRawReader *rawReader)
312 // Event processing loop - AliRawReader
314 AliTPCRawStreamFast *rawStreamFast = new AliTPCRawStreamFast(rawReader, (AliAltroMapping**)fMapping);
315 Bool_t res=ProcessEventFast(rawStreamFast);
317 fEventCounter++; // only increment event counter if there is TPC data
318 // otherwise Analyse (called in QA) fails
320 delete rawStreamFast;
324 //_____________________________________________________________________
325 Bool_t AliTPCdataQA::ProcessEvent(AliTPCRawStream *rawStream)
328 // Event Processing loop - AliTPCRawStream
331 Bool_t withInput = kFALSE;
333 Int_t lastSector = -1;
335 while (rawStream->Next()) {
337 Int_t iSector = rawStream->GetSector(); // current ROC
338 Int_t iRow = rawStream->GetRow(); // current row
339 Int_t iPad = rawStream->GetPad(); // current pad
340 Int_t iTimeBin = rawStream->GetTime(); // current time bin
341 Float_t signal = rawStream->GetSignal(); // current ADC signal
343 // Call local maxima finder if the data is in a new sector
344 if(iSector != lastSector) {
347 FindLocalMaxima(lastSector);
350 lastSector = iSector;
354 // Sometimes iRow==-1 if there is problems to read the data
356 nSignals += Update(iSector,iRow,iPad,iTimeBin,signal);
365 //_____________________________________________________________________
366 Bool_t AliTPCdataQA::ProcessEvent(AliRawReader *rawReader)
369 // Event processing loop - AliRawReader
372 // if fMapping is NULL the rawstream will crate its own mapping
373 AliTPCRawStream rawStream(rawReader, (AliAltroMapping**)fMapping);
374 rawReader->Select("TPC");
375 Bool_t res = ProcessEvent(&rawStream);
378 fEventCounter++; // only increment event counter if there is TPC data
379 // otherwise Analyse (called in QA) fails
384 //_____________________________________________________________________
385 Bool_t AliTPCdataQA::ProcessEvent(eventHeaderStruct *event)
388 // process date event
391 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
392 Bool_t result=ProcessEvent(rawReader);
399 //_____________________________________________________________________
400 void AliTPCdataQA::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append) /*FOLD00*/
403 // Write class to file
414 TDirectory *backup = gDirectory;
415 TFile f(filename,option.Data());
417 if ( !sDir.IsNull() ){
418 f.mkdir(sDir.Data());
424 if ( backup ) backup->cd();
428 //_____________________________________________________________________
429 Int_t AliTPCdataQA::Update(const Int_t iSector, /*FOLD00*/
432 const Int_t iTimeBin,
436 // Signal filling method
440 // Define the calibration objects the first time Update is called
441 // NB! This has to be done first even if the data is rejected by the time
442 // cut to make sure that the objects are available in Analyse
444 if (!fNLocalMaxima) fNLocalMaxima = new AliTPCCalPad("NLocalMaxima","NLocalMaxima");
445 if (!fMaxCharge) fMaxCharge = new AliTPCCalPad("MaxCharge","MaxCharge");
446 if (!fMeanCharge) fMeanCharge = new AliTPCCalPad("MeanCharge","MeanCharge");
447 if (!fNoThreshold) fNoThreshold = new AliTPCCalPad("NoThreshold","NoThreshold");
448 if (!fNTimeBins) fNTimeBins = new AliTPCCalPad("NTimeBins","NTimeBins");
449 if (!fNPads) fNPads = new AliTPCCalPad("NPads","NPads");
450 if (!fTimePosition) fTimePosition = new AliTPCCalPad("TimePosition","TimePosition");
451 if (!fOverThreshold10) fOverThreshold10 = new AliTPCCalPad("OverThreshold10","OverThreshold10");
452 if (!fOverThreshold20) fOverThreshold20 = new AliTPCCalPad("OverThreshold20","OverThreshold20");
453 if (!fOverThreshold30) fOverThreshold30 = new AliTPCCalPad("OverThreshold30","OverThreshold30");
454 // Make the arrays for expanding the data
460 // If Analyse has been previously called we need now to denormalize the data
461 // as more data is coming
463 if(fIsAnalysed == kTRUE) {
465 const Int_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
466 const Float_t denormalization = Float_t(fEventCounter * nTimeBins);
467 fNoThreshold->Multiply(denormalization);
469 fMeanCharge->Multiply(fNLocalMaxima);
470 fMaxCharge->Multiply(fNLocalMaxima);
471 fNTimeBins->Multiply(fNLocalMaxima);
472 fNPads->Multiply(fNLocalMaxima);
473 fTimePosition->Multiply(fNLocalMaxima);
474 fIsAnalysed = kFALSE;
480 if (iTimeBin<fFirstTimeBin) return 0;
481 if (iTimeBin>fLastTimeBin) return 0;
483 // if pedestal calibrations are loaded subtract pedestals
486 Float_t ped = fPedestal->GetCalROC(iSector)->GetValue(iRow, iPad);
487 // Don't use data from pads where pedestals are abnormally small or large
493 // In fNoThreshold we fill all data to estimate the ZS volume
494 Float_t count = fNoThreshold->GetCalROC(iSector)->GetValue(iRow, iPad);
495 fNoThreshold->GetCalROC(iSector)->SetValue(iRow, iPad,count+1);
497 // Require at least 3 ADC channels
501 // if noise calibrations are loaded require at least 3*sigmaNoise
504 Float_t noise = fNoise->GetCalROC(iSector)->GetValue(iRow, iPad);
506 if(signal < noise*3.0)
511 // This signal is ok and we store it in the cluster map
514 SetExpandDigit(iRow, iPad, iTimeBin, signal);
516 return 1; // signal was accepted
519 //_____________________________________________________________________
520 void AliTPCdataQA::FindLocalMaxima(const Int_t iSector)
523 // This method is called after the data from each sector has been
524 // exapanded into an array
525 // Loop over the signals and identify local maxima and fill the
526 // calibration objects with the information
529 Int_t nLocalMaxima = 0;
530 const Int_t maxTimeBin = fTimeBinsMax+4; // Used to step between neighboring pads
531 // Because we have tha pad-time data in a
534 for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
536 Float_t* allBins = fAllBins[iRow];
537 Int_t* sigBins = fAllSigBins[iRow];
538 const Int_t nSigBins = fAllNSigBins[iRow];
540 for (Int_t iSig = 0; iSig < nSigBins; iSig++) {
542 Int_t bin = sigBins[iSig];
543 Float_t *b = &allBins[bin];
546 // Now we check if this is a local maximum
551 // First check that the charge is bigger than the threshold
555 // Require at least one neighboring pad with signal
556 if (b[-maxTimeBin]+b[maxTimeBin]<=0) continue;
558 // Require at least one neighboring time bin with signal
559 if (b[-1]+b[1]<=0) continue;
562 // Check that this is a local maximum
563 // Note that the checking is done so that if 2 charges has the same
564 // qMax then only 1 cluster is generated
565 // (that is why there is BOTH > and >=)
567 if (b[-maxTimeBin] >= qMax) continue;
568 if (b[-1 ] >= qMax) continue;
569 if (b[+maxTimeBin] > qMax) continue;
570 if (b[+1 ] > qMax) continue;
571 if (b[-maxTimeBin-1] >= qMax) continue;
572 if (b[+maxTimeBin-1] >= qMax) continue;
573 if (b[+maxTimeBin+1] > qMax) continue;
574 if (b[-maxTimeBin+1] >= qMax) continue;
577 // Now we accept the local maximum and fill the calibration/data objects
581 Int_t iPad, iTimeBin;
582 GetPadAndTimeBin(bin, iPad, iTimeBin);
584 Float_t count = fNLocalMaxima->GetCalROC(iSector)->GetValue(iRow, iPad);
585 fNLocalMaxima->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
587 count = fTimePosition->GetCalROC(iSector)->GetValue(iRow, iPad);
588 fTimePosition->GetCalROC(iSector)->SetValue(iRow, iPad, count+iTimeBin);
590 Float_t charge = fMaxCharge->GetCalROC(iSector)->GetValue(iRow, iPad);
591 fMaxCharge->GetCalROC(iSector)->SetValue(iRow, iPad, charge + qMax);
594 count = fOverThreshold10->GetCalROC(iSector)->GetValue(iRow, iPad);
595 fOverThreshold10->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
598 count = fOverThreshold20->GetCalROC(iSector)->GetValue(iRow, iPad);
599 fOverThreshold20->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
602 count = fOverThreshold30->GetCalROC(iSector)->GetValue(iRow, iPad);
603 fOverThreshold30->GetCalROC(iSector)->SetValue(iRow, iPad, count+1);
607 // Calculate the total charge as the sum over the region:
615 // with qmax at the center C.
617 // The inner charge (i) we always add, but we only add the outer
618 // charge (o) if the neighboring inner bin (i) has a signal.
620 Int_t minP = 0, maxP = 0, minT = 0, maxT = 0;
622 for(Int_t i = -1; i<=1; i++) {
623 for(Int_t j = -1; j<=1; j++) {
628 Float_t charge1 = GetQ(b, i, j, maxTimeBin, minT, maxT, minP, maxP);
631 // see if the next neighbor is also above threshold
633 qTot += GetQ(b, 2*i, 2*j, maxTimeBin, minT, maxT, minP, maxP);
635 // we are in a diagonal corner
636 qTot += GetQ(b, i, 2*j, maxTimeBin, minT, maxT, minP, maxP);
637 qTot += GetQ(b, 2*i, j, maxTimeBin, minT, maxT, minP, maxP);
638 qTot += GetQ(b, 2*i, 2*j, maxTimeBin, minT, maxT, minP, maxP);
644 charge = fMeanCharge->GetCalROC(iSector)->GetValue(iRow, iPad);
645 fMeanCharge->GetCalROC(iSector)->SetValue(iRow, iPad, charge + qTot);
647 count = fNTimeBins->GetCalROC(iSector)->GetValue(iRow, iPad);
648 fNTimeBins->GetCalROC(iSector)->SetValue(iRow, iPad, count + maxT-minT+1);
650 count = fNPads->GetCalROC(iSector)->GetValue(iRow, iPad);
651 fNPads->GetCalROC(iSector)->SetValue(iRow, iPad, count + maxP-minP+1);
653 } // end loop over signals
654 } // end loop over rows
656 // cout << "Number of local maximas found: " << nLocalMaxima << endl;
659 //_____________________________________________________________________
660 void AliTPCdataQA::Analyse()
663 // Calculate calibration constants
666 cout << "Analyse called" << endl;
668 if(fIsAnalysed == kTRUE) {
670 cout << "No new data since Analyse was called last time" << endl;
674 if(fEventCounter==0) {
676 cout << "EventCounter == 0, Cannot analyse" << endl;
680 Int_t nTimeBins = fLastTimeBin - fFirstTimeBin +1;
681 cout << "EventCounter: " << fEventCounter << endl
682 << "TimeBins: " << nTimeBins << endl;
684 Float_t normalization = 1.0 / Float_t(fEventCounter * nTimeBins);
685 fNoThreshold->Multiply(normalization);
687 fMeanCharge->Divide(fNLocalMaxima);
688 fMaxCharge->Divide(fNLocalMaxima);
689 fNTimeBins->Divide(fNLocalMaxima);
690 fNPads->Divide(fNLocalMaxima);
691 fTimePosition->Divide(fNLocalMaxima);
697 //_____________________________________________________________________
698 void AliTPCdataQA::MakeTree(const char *fname){
700 // Export result to the tree -located in the file
701 // This file can be analyzed using AliTPCCalibViewer
703 AliTPCPreprocessorOnline preprocesor;
705 if (fNLocalMaxima) preprocesor.AddComponent(fNLocalMaxima);
706 if (fMaxCharge) preprocesor.AddComponent(fMaxCharge);
707 if (fMeanCharge) preprocesor.AddComponent(fMeanCharge);
708 if (fNoThreshold) preprocesor.AddComponent(fNoThreshold);
709 if (fNTimeBins) preprocesor.AddComponent(fNTimeBins);
710 if (fNPads) preprocesor.AddComponent(fNPads);
711 if (fTimePosition) preprocesor.AddComponent(fTimePosition);
712 if (fOverThreshold10) preprocesor.AddComponent(fOverThreshold10);
713 if (fOverThreshold20) preprocesor.AddComponent(fOverThreshold20);
714 if (fOverThreshold30) preprocesor.AddComponent(fOverThreshold30);
716 preprocesor.DumpToFile(fname);
720 //_____________________________________________________________________
721 void AliTPCdataQA::MakeArrays(){
723 // The arrays for expanding the raw data are defined and
724 // som parameters are intialised
726 AliTPCROC * roc = AliTPCROC::Instance();
728 // To make the array big enough for all sectors we take
729 // the dimensions from the outer row of an OROC (the last sector)
731 fRowsMax = roc->GetNRows(roc->GetNSector()-1);
732 fPadsMax = roc->GetNPads(roc->GetNSector()-1,fRowsMax-1);
733 fTimeBinsMax = fLastTimeBin - fFirstTimeBin +1;
736 // Since we have added 2 pads (TimeBins) before and after the real pads (TimeBins)
737 // to make sure that we can always query the exanded table even when the
738 // max is on the edge
742 fAllBins = new Float_t*[fRowsMax];
743 fAllSigBins = new Int_t*[fRowsMax];
744 fAllNSigBins = new Int_t[fRowsMax];
746 for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
748 Int_t maxBin = (fTimeBinsMax+4)*(fPadsMax+4);
749 fAllBins[iRow] = new Float_t[maxBin];
750 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin); // set all values to 0
751 fAllSigBins[iRow] = new Int_t[maxBin];
752 fAllNSigBins[iRow] = 0;
757 //_____________________________________________________________________
758 void AliTPCdataQA::CleanArrays(){
763 for (Int_t iRow = 0; iRow < fRowsMax; iRow++) {
765 Int_t maxBin = (fTimeBinsMax+4)*(fPadsMax+4);
766 memset(fAllBins[iRow],0,sizeof(Float_t)*maxBin);
767 fAllNSigBins[iRow]=0;
771 //_____________________________________________________________________
772 void AliTPCdataQA::GetPadAndTimeBin(Int_t bin, Int_t& iPad, Int_t& iTimeBin){
774 // Return pad and timebin for a given bin
777 // Int_t bin = iPad*(fTimeBinsMax+4)+iTimeBin;
778 iTimeBin = bin%(fTimeBinsMax+4);
779 iPad = (bin-iTimeBin)/(fTimeBinsMax+4);
783 iTimeBin += fFirstTimeBin;
785 R__ASSERT(iPad>=0 && iPad<=fPadsMax);
786 R__ASSERT(iTimeBin>=fFirstTimeBin && iTimeBin<=fLastTimeBin);
789 //_____________________________________________________________________
790 void AliTPCdataQA::SetExpandDigit(const Int_t iRow, Int_t iPad,
791 Int_t iTimeBin, const Float_t signal)
796 R__ASSERT(iRow>=0 && iRow<fRowsMax);
797 R__ASSERT(iPad>=0 && iPad<=fPadsMax);
798 R__ASSERT(iTimeBin>=fFirstTimeBin && iTimeBin<=fLastTimeBin);
800 iTimeBin -= fFirstTimeBin;
804 Int_t bin = iPad*(fTimeBinsMax+4)+iTimeBin;
805 fAllBins[iRow][bin] = signal;
806 fAllSigBins[iRow][fAllNSigBins[iRow]] = bin;
807 fAllNSigBins[iRow]++;
810 Float_t AliTPCdataQA::GetQ(const Float_t* adcArray, const Int_t time,
811 const Int_t pad, const Int_t maxTimeBins,
812 Int_t& timeMin, Int_t& timeMax,
813 Int_t& padMin, Int_t& padMax)
816 // This methods return the charge in the bin time+pad*maxTimeBins
817 // If the charge is above 0 it also updates the padMin, padMax, timeMin
818 // and timeMax if necessary
820 Float_t charge = adcArray[time + pad*maxTimeBins];
822 timeMin = TMath::Min(time, timeMin); timeMax = TMath::Max(time, timeMax);
823 padMin = TMath::Min(pad, padMin); padMax = TMath::Max(pad, padMax);