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 **************************************************************************/
18 /////////////////////////////////////////////////////////////////////////////////////////
20 // Implementation of the TPC pulser calibration //
22 // Origin: Jens Wiechula, Marian Ivanov J.Wiechula@gsi.de, Marian.Ivanov@cern.ch //
24 /////////////////////////////////////////////////////////////////////////////////////////
25 /***************************************************************************
27 ***************************************************************************
29 The AliTPCCalibPulser class is used to get calibration data concerning the FEE using
30 runs performed with the calibration pulser.
32 The information retrieved is
34 - Signal width differences
35 - Amplification variations
37 the seen differences arise from the manufacturing tolerances of the PASAs and are very small within
38 one chip and somewhat large between different chips.
41 For each ROC three TH2S histos 'Reference Histograms' (ROC channel vs. [Time0, signal width, Q sum]) is created when
42 it is filled for the first time (GetHisto[T0,RMS,Q](ROC,kTRUE)). The histos are stored in the
43 TObjArrays fHistoT0Array, fHistoRMSArray and fHistoQArray.
48 Raw calibration pulser data is processed by calling one of the ProcessEvent(...) functions
49 (see below). These in the end call the Update(...) function.
51 - the Update(...) function:
52 In this function the array fPadSignal is filled with the adc signals between the specified range
53 fFirstTimeBin and fLastTimeBin for the current pad.
54 before going to the next pad the ProcessPad() function is called, which analyses the data for one pad
57 - the ProcessPad() function:
58 Find Pedestal and Noise information
59 - use database information which has to be set by calling
60 SetPedestalDatabase(AliTPCCalPad *pedestalTPC, AliTPCCalPad *padNoiseTPC)
61 - if no information from the pedestal data base
62 is available the informaion is calculated on the fly ( see FindPedestal() function )
64 Find the Pulser signal information
65 - calculate mean = T0, RMS = signal width and Q sum in a range of -2+7 timebins around Q max
66 the Q sum is scaled by pad area
67 (see FindPulserSignal(...) function)
69 Fill a temprary array for the T0 information (GetPadTimesEvent(fCurrentSector,kTRUE)) (why see below)
70 Fill the Q sum and RMS values in the histograms (GetHisto[RMS,Q](ROC,kTRUE)),
72 At the end of each event the EndEvent() function is called
74 - the EndEvent() function:
75 calculate the mean T0 for each ROC and fill the Time0 histogram with Time0-<Time0 for ROC>
76 This is done to overcome syncronisation problems between the trigger and the fec clock.
78 After accumulating the desired statistics the Analyse() function has to be called.
79 - the Analyse() function
80 Whithin this function the mean values of T0, RMS, Q are calculated for each pad, using
81 the AliMathBase::GetCOG(...) function, and the calibration
82 storage classes (AliTPCCalROC) are filled for each ROC.
83 The calibration information is stored in the TObjArrays fCalRocArrayT0, fCalRocArrayRMS and
88 User interface for filling data:
89 --------------------------------
91 To Fill information one of the following functions can be used:
93 Bool_t ProcessEvent(eventHeaderStruct *event);
95 - use AliTPCRawReaderDate and call ProcessEvent(AliRawReader *rawReader)
97 Bool_t ProcessEvent(AliRawReader *rawReader);
98 - process AliRawReader event
99 - use AliTPCRawStream to loop over data and call ProcessEvent(AliTPCRawStream *rawStream)
101 Bool_t ProcessEvent(AliTPCRawStream *rawStream);
102 - process event from AliTPCRawStream
103 - call Update function for signal filling
105 Int_t Update(const Int_t isector, const Int_t iRow, const Int_t
106 iPad, const Int_t iTimeBin, const Float_t signal);
107 - directly fill signal information (sector, row, pad, time bin, pad)
108 to the reference histograms
110 It is also possible to merge two independently taken calibrations using the function
112 void Merge(AliTPCCalibPulser *sig)
113 - copy histograms in 'sig' if the do not exist in this instance
114 - Add histograms in 'sig' to the histograms in this instance if the allready exist
115 - After merging call Analyse again!
119 -- example: filling data using root raw data:
120 void fillSignal(Char_t *filename)
122 rawReader = new AliRawReaderRoot(fileName);
123 if ( !rawReader ) return;
124 AliTPCCalibPulser *calib = new AliTPCCalibPulser;
125 while (rawReader->NextEvent()){
126 calib->ProcessEvent(rawReader);
129 calib->DumpToFile("SignalData.root");
135 What kind of information is stored and how to retrieve them:
136 ------------------------------------------------------------
138 - Accessing the 'Reference Histograms' (Time0, signal width and Q sum information pad by pad):
140 TH2F *GetHistoT0(Int_t sector);
141 TH2F *GetHistoRMS(Int_t sector);
142 TH2F *GetHistoQ(Int_t sector);
144 - Accessing the calibration storage objects:
146 AliTPCCalROC *GetCalRocT0(Int_t sector); // for the Time0 values
147 AliTPCCalROC *GetCalRocRMS(Int_t sector); // for the signal width values
148 AliTPCCalROC *GetCalRocQ(Int_t sector); // for the Q sum values
150 example for visualisation:
151 if the file "SignalData.root" was created using the above example one could do the following:
153 TFile fileSignal("SignalData.root")
154 AliTPCCalibPulser *sig = (AliTPCCalibPulser*)fileSignal->Get("AliTPCCalibPulser");
155 sig->GetCalRocT0(0)->Draw("colz");
156 sig->GetCalRocRMS(0)->Draw("colz");
158 or use the AliTPCCalPad functionality:
159 AliTPCCalPad padT0(ped->GetCalPadT0());
160 AliTPCCalPad padSigWidth(ped->GetCalPadRMS());
161 padT0->MakeHisto2D()->Draw("colz"); //Draw A-Side Time0 Information
162 padSigWidth->MakeHisto2D()->Draw("colz"); //Draw A-Side signal width Information
167 #include <TObjArray.h>
171 #include <TVectorF.h>
174 #include <TDirectory.h>
180 #include "AliRawReader.h"
181 #include "AliRawReaderRoot.h"
182 #include "AliRawReaderDate.h"
183 #include "AliTPCRawStream.h"
184 #include "AliTPCRawStreamFast.h"
185 #include "AliTPCCalROC.h"
186 #include "AliTPCCalPad.h"
187 #include "AliTPCROC.h"
188 #include "AliTPCParam.h"
189 #include "AliTPCCalibPulser.h"
190 #include "AliTPCcalibDB.h"
191 #include "AliMathBase.h"
192 #include "TTreeStream.h"
200 ClassImp(AliTPCCalibPulser)
202 AliTPCCalibPulser::AliTPCCalibPulser() :
216 fOldRCUformat(kTRUE),
217 fROC(AliTPCROC::Instance()),
218 fParam(new AliTPCParam),
226 fCalRocArrayOutliers(72),
230 fPadTimesArrayEvent(72),
232 fPadRMSArrayEvent(72),
233 fPadPedestalArrayEvent(72),
243 fVTime0OffsetCounter(72),
249 // AliTPCSignal default constructor
253 //_____________________________________________________________________
254 AliTPCCalibPulser::AliTPCCalibPulser(const AliTPCCalibPulser &sig) :
256 fFirstTimeBin(sig.fFirstTimeBin),
257 fLastTimeBin(sig.fLastTimeBin),
258 fNbinsT0(sig.fNbinsT0),
259 fXminT0(sig.fXminT0),
260 fXmaxT0(sig.fXmaxT0),
261 fNbinsQ(sig.fNbinsQ),
264 fNbinsRMS(sig.fNbinsRMS),
265 fXminRMS(sig.fXminRMS),
266 fXmaxRMS(sig.fXmaxRMS),
268 fOldRCUformat(kTRUE),
269 fROC(AliTPCROC::Instance()),
270 fParam(new AliTPCParam),
278 fCalRocArrayOutliers(72),
282 fPadTimesArrayEvent(72),
284 fPadRMSArrayEvent(72),
285 fPadPedestalArrayEvent(72),
295 fVTime0OffsetCounter(72),
298 fDebugLevel(sig.fDebugLevel)
301 // AliTPCSignal default constructor
304 for (Int_t iSec = 0; iSec < 72; ++iSec){
305 const AliTPCCalROC *calQ = (AliTPCCalROC*)sig.fCalRocArrayQ.UncheckedAt(iSec);
306 const AliTPCCalROC *calT0 = (AliTPCCalROC*)sig.fCalRocArrayT0.UncheckedAt(iSec);
307 const AliTPCCalROC *calRMS = (AliTPCCalROC*)sig.fCalRocArrayRMS.UncheckedAt(iSec);
308 const AliTPCCalROC *calOut = (AliTPCCalROC*)sig.fCalRocArrayOutliers.UncheckedAt(iSec);
310 const TH2S *hQ = (TH2S*)sig.fHistoQArray.UncheckedAt(iSec);
311 const TH2S *hT0 = (TH2S*)sig.fHistoT0Array.UncheckedAt(iSec);
312 const TH2S *hRMS = (TH2S*)sig.fHistoRMSArray.UncheckedAt(iSec);
314 if ( calQ != 0x0 ) fCalRocArrayQ.AddAt(new AliTPCCalROC(*calQ), iSec);
315 if ( calT0 != 0x0 ) fCalRocArrayT0.AddAt(new AliTPCCalROC(*calT0), iSec);
316 if ( calRMS != 0x0 ) fCalRocArrayRMS.AddAt(new AliTPCCalROC(*calRMS), iSec);
317 if ( calOut != 0x0 ) fCalRocArrayOutliers.AddAt(new AliTPCCalROC(*calOut), iSec);
320 TH2S *hNew = new TH2S(*hQ);
321 hNew->SetDirectory(0);
322 fHistoQArray.AddAt(hNew,iSec);
325 TH2S *hNew = new TH2S(*hT0);
326 hNew->SetDirectory(0);
327 fHistoQArray.AddAt(hNew,iSec);
330 TH2S *hNew = new TH2S(*hRMS);
331 hNew->SetDirectory(0);
332 fHistoQArray.AddAt(hNew,iSec);
337 //_____________________________________________________________________
338 AliTPCCalibPulser& AliTPCCalibPulser::operator = (const AliTPCCalibPulser &source)
341 // assignment operator
343 if (&source == this) return *this;
344 new (this) AliTPCCalibPulser(source);
348 //_____________________________________________________________________
349 AliTPCCalibPulser::~AliTPCCalibPulser()
355 fCalRocArrayT0.Delete();
356 fCalRocArrayQ.Delete();
357 fCalRocArrayRMS.Delete();
358 fCalRocArrayOutliers.Delete();
360 fHistoQArray.Delete();
361 fHistoT0Array.Delete();
362 fHistoRMSArray.Delete();
364 fPadTimesArrayEvent.Delete();
365 fPadQArrayEvent.Delete();
366 fPadRMSArrayEvent.Delete();
367 fPadPedestalArrayEvent.Delete();
369 if ( fDebugStreamer) delete fDebugStreamer;
373 //_____________________________________________________________________
374 Int_t AliTPCCalibPulser::Update(const Int_t icsector,
377 const Int_t icTimeBin,
378 const Float_t csignal)
381 // Signal filling methode on the fly pedestal and time offset correction if necessary.
382 // no extra analysis necessary. Assumes knowledge of the signal shape!
383 // assumes that it is looped over consecutive time bins of one pad
385 if ( (icTimeBin>fLastTimeBin) || (icTimeBin<fFirstTimeBin) ) return 0;
387 Int_t iChannel = fROC->GetRowIndexes(icsector)[icRow]+icPad; // global pad position in sector
389 //init first pad and sector in this event
390 if ( fCurrentChannel == -1 ) {
391 fCurrentChannel = iChannel;
392 fCurrentSector = icsector;
396 //process last pad if we change to a new one
397 if ( iChannel != fCurrentChannel ){
399 fCurrentChannel = iChannel;
400 fCurrentSector = icsector;
404 //fill signals for current pad
405 fPadSignal[icTimeBin]=csignal;
406 if ( csignal > fMaxPadSignal ){
407 fMaxPadSignal = csignal;
408 fMaxTimeBin = icTimeBin;
412 //_____________________________________________________________________
413 void AliTPCCalibPulser::FindPedestal(Float_t part)
416 // find pedestal and noise for the current pad. Use either database or
417 // truncated mean with part*100%
419 Bool_t noPedestal = kTRUE;;
420 if (fPedestalTPC&&fPadNoiseTPC){
421 //use pedestal database
422 //only load new pedestals if the sector has changed
423 if ( fCurrentSector!=fLastSector ){
424 fPedestalROC = fPedestalTPC->GetCalROC(fCurrentSector);
425 fPadNoiseROC = fPadNoiseTPC->GetCalROC(fCurrentSector);
426 fLastSector=fCurrentSector;
429 if ( fPedestalROC&&fPadNoiseROC ){
430 fPadPedestal = fPedestalROC->GetValue(fCurrentChannel);
431 fPadNoise = fPadNoiseROC->GetValue(fCurrentChannel);
437 //if we are not running with pedestal database, or for the current sector there is no information
438 //available, calculate the pedestal and noise on the fly
440 const Int_t kPedMax = 100; //maximum pedestal value
449 UShort_t histo[kPedMax];
450 memset(histo,0,kPedMax*sizeof(UShort_t));
452 for (Int_t i=fFirstTimeBin; i<=fLastTimeBin; ++i){
453 padSignal = fPadSignal.GetMatrixArray()[i];
454 if (padSignal<=0) continue;
455 if (padSignal>max && i>10) {
459 if (padSignal>kPedMax-1) continue;
460 histo[Int_t(padSignal+0.5)]++;
464 for (Int_t i=1; i<kPedMax; ++i){
465 if (count1<count0*0.5) median=i;
470 Float_t count=histo[median] ,mean=histo[median]*median, rms=histo[median]*median*median ;
472 for (Int_t idelta=1; idelta<10; ++idelta){
473 if (median-idelta<=0) continue;
474 if (median+idelta>kPedMax) continue;
475 if (count<part*count1){
476 count+=histo[median-idelta];
477 mean +=histo[median-idelta]*(median-idelta);
478 rms +=histo[median-idelta]*(median-idelta)*(median-idelta);
479 count+=histo[median+idelta];
480 mean +=histo[median+idelta]*(median+idelta);
481 rms +=histo[median+idelta]*(median+idelta)*(median+idelta);
488 rms = TMath::Sqrt(TMath::Abs(rms/count-mean*mean));
494 //_____________________________________________________________________
495 void AliTPCCalibPulser::FindPulserSignal(TVectorD ¶m, Float_t &qSum)
498 // Find position, signal width and height of the CE signal (last signal)
499 // param[0] = Qmax, param[1] = mean time, param[2] = rms;
500 // maxima: array of local maxima of the pad signal use the one closest to the mean CE position
503 Float_t ceQmax =0, ceQsum=0, ceTime=0, ceRMS=0;
504 Int_t cemaxpos = fMaxTimeBin;
505 Float_t ceSumThreshold = 8.*fPadNoise; // threshold for the signal sum
506 const Int_t kCemin = 2; // range for the analysis of the ce signal +- channels from the peak
507 const Int_t kCemax = 7;
510 ceQmax = fPadSignal.GetMatrixArray()[cemaxpos]-fPadPedestal;
511 for (Int_t i=cemaxpos-kCemin; i<cemaxpos+kCemax; ++i){
512 Float_t signal = fPadSignal.GetMatrixArray()[i]-fPadPedestal;
513 if ( (i>fFirstTimeBin) && (i<fLastTimeBin) && (signal>0) ){
514 ceTime+=signal*(i+0.5);
515 ceRMS +=signal*(i+0.5)*(i+0.5);
520 if (ceQmax&&ceQsum>ceSumThreshold) {
522 ceRMS = TMath::Sqrt(TMath::Abs(ceRMS/ceQsum-ceTime*ceTime));
523 fVTime0Offset.GetMatrixArray()[fCurrentSector]+=ceTime; // mean time for each sector
524 fVTime0OffsetCounter.GetMatrixArray()[fCurrentSector]++;
526 //Normalise Q to the pad area
527 Float_t norm = fParam->GetPadPitchWidth(fCurrentSector)*fParam->GetPadPitchLength(fCurrentSector,fCurrentRow);
541 //_____________________________________________________________________
542 void AliTPCCalibPulser::ProcessPad()
545 // Process data of current pad
551 FindPulserSignal(param, qSum);
553 Double_t meanT = param[1];
554 Double_t sigmaT = param[2];
556 //Fill Event T0 counter
557 (*GetPadTimesEvent(fCurrentSector,kTRUE)).GetMatrixArray()[fCurrentChannel] = meanT;
560 GetHistoQ(fCurrentSector,kTRUE)->Fill( TMath::Sqrt(qSum), fCurrentChannel );
563 GetHistoRMS(fCurrentSector,kTRUE)->Fill( sigmaT, fCurrentChannel );
566 //Fill debugging info
567 if ( fDebugLevel>0 ){
568 (*GetPadPedestalEvent(fCurrentSector,kTRUE))[fCurrentChannel]=fPadPedestal;
569 (*GetPadRMSEvent(fCurrentSector,kTRUE))[fCurrentChannel]=sigmaT;
570 (*GetPadQEvent(fCurrentSector,kTRUE))[fCurrentChannel]=qSum;
574 //_____________________________________________________________________
575 void AliTPCCalibPulser::EndEvent()
578 // Process data of current event
581 //check if last pad has allready been processed, if not do so
582 if ( fMaxTimeBin>-1 ) ProcessPad();
584 //loop over all ROCs, fill Time0 histogram corrected for the mean Time0 of each ROC
585 for ( Int_t iSec = 0; iSec<72; ++iSec ){
586 TVectorF *vTimes = GetPadTimesEvent(iSec);
587 if ( !vTimes ) continue;
589 for ( UInt_t iChannel=0; iChannel<fROC->GetNChannels(iSec); ++iChannel ){
590 Float_t time0 = fVTime0Offset[iSec]/fVTime0OffsetCounter[iSec];
591 Float_t time = (*vTimes).GetMatrixArray()[iChannel];
593 GetHistoT0(iSec,kTRUE)->Fill( time-time0,iChannel );
597 if ( fDebugLevel>0 ){
598 if ( !fDebugStreamer ) {
600 TDirectory *backup = gDirectory;
601 fDebugStreamer = new TTreeSRedirector("deb2.root");
602 if ( backup ) backup->cd(); //we don't want to be cd'd to the debug streamer
609 Float_t q = (*GetPadQEvent(iSec)).GetMatrixArray()[iChannel];
610 Float_t rms = (*GetPadRMSEvent(iSec)).GetMatrixArray()[iChannel];
612 UInt_t channel=iChannel;
615 while ( channel > (fROC->GetRowIndexes(sector)[row]+fROC->GetNPads(sector,row)-1) ) row++;
616 pad = channel-fROC->GetRowIndexes(sector)[row];
617 padc = pad-(fROC->GetNPads(sector,row)/2);
619 TH1F *h1 = new TH1F(Form("hSignalD%d.%d.%d",sector,row,pad),
620 Form("hSignalD%d.%d.%d",sector,row,pad),
621 fLastTimeBin-fFirstTimeBin,
622 fFirstTimeBin,fLastTimeBin);
625 for (Int_t i=fFirstTimeBin; i<fLastTimeBin+1; ++i)
626 h1->Fill(i,fPadSignal(i));
628 (*fDebugStreamer) << "DataPad" <<
629 // "Event=" << fEvent <<
630 "Sector="<< sector <<
634 "PadSec="<< channel <<
648 //_____________________________________________________________________
649 Bool_t AliTPCCalibPulser::ProcessEventFast(AliTPCRawStreamFast *rawStreamFast)
652 // Event Processing loop - AliTPCRawStream
656 Bool_t withInput = kFALSE;
658 while ( rawStreamFast->NextDDL() ){
659 while ( rawStreamFast->NextChannel() ){
660 Int_t isector = rawStreamFast->GetSector(); // current sector
661 Int_t iRow = rawStreamFast->GetRow(); // current row
662 Int_t iPad = rawStreamFast->GetPad(); // current pad
663 Int_t startTbin = (Int_t)rawStreamFast->GetStartTimeBin();
664 Int_t endTbin = (Int_t)rawStreamFast->GetEndTimeBin();
666 while ( rawStreamFast->NextBunch() ){
667 for (Int_t iTimeBin = startTbin; iTimeBin < endTbin; iTimeBin++){
668 Float_t signal=(Float_t)rawStreamFast->GetSignals()[iTimeBin-startTbin];
669 Update(isector,iRow,iPad,iTimeBin+1,signal);
680 //_____________________________________________________________________
681 Bool_t AliTPCCalibPulser::ProcessEventFast(AliRawReader *rawReader)
684 // Event processing loop - AliRawReader
686 AliTPCRawStreamFast *rawStreamFast = new AliTPCRawStreamFast(rawReader);
687 Bool_t res=ProcessEventFast(rawStreamFast);
688 delete rawStreamFast;
691 //_____________________________________________________________________
692 Bool_t AliTPCCalibPulser::ProcessEvent(AliTPCRawStream *rawStream)
695 // Event Processing loop - AliTPCRawStream
698 rawStream->SetOldRCUFormat(fOldRCUformat);
702 Bool_t withInput = kFALSE;
704 while (rawStream->Next()) {
705 Int_t isector = rawStream->GetSector(); // current sector
706 Int_t iRow = rawStream->GetRow(); // current row
707 Int_t iPad = rawStream->GetPad(); // current pad
708 Int_t iTimeBin = rawStream->GetTime(); // current time bin
709 Float_t signal = rawStream->GetSignal(); // current ADC signal
711 Update(isector,iRow,iPad,iTimeBin,signal);
719 //_____________________________________________________________________
720 Bool_t AliTPCCalibPulser::ProcessEvent(AliRawReader *rawReader)
723 // Event processing loop - AliRawReader
727 AliTPCRawStream rawStream(rawReader);
729 rawReader->Select("TPC");
731 return ProcessEvent(&rawStream);
733 //_____________________________________________________________________
734 Bool_t AliTPCCalibPulser::ProcessEvent(eventHeaderStruct *event)
737 // Event processing loop - date event
739 AliRawReader *rawReader = new AliRawReaderDate((void*)event);
740 Bool_t result=ProcessEvent(rawReader);
745 //_____________________________________________________________________
746 TH2S* AliTPCCalibPulser::GetHisto(Int_t sector, TObjArray *arr,
747 Int_t nbinsY, Float_t ymin, Float_t ymax,
748 Char_t *type, Bool_t force)
751 // return pointer to Q histogram
752 // if force is true create a new histogram if it doesn't exist allready
754 if ( !force || arr->UncheckedAt(sector) )
755 return (TH2S*)arr->UncheckedAt(sector);
757 // if we are forced and histogram doesn't yes exist create it
758 Char_t name[255], title[255];
760 sprintf(name,"hCalib%s%.2d",type,sector);
761 sprintf(title,"%s calibration histogram sector %.2d",type,sector);
763 // new histogram with Q calib information. One value for each pad!
764 TH2S* hist = new TH2S(name,title,
766 fROC->GetNChannels(sector),0,fROC->GetNChannels(sector));
767 hist->SetDirectory(0);
768 arr->AddAt(hist,sector);
771 //_____________________________________________________________________
772 TH2S* AliTPCCalibPulser::GetHistoT0(Int_t sector, Bool_t force)
775 // return pointer to T0 histogram
776 // if force is true create a new histogram if it doesn't exist allready
778 TObjArray *arr = &fHistoT0Array;
779 return GetHisto(sector, arr, fNbinsT0, fXminT0, fXmaxT0, "T0", force);
781 //_____________________________________________________________________
782 TH2S* AliTPCCalibPulser::GetHistoQ(Int_t sector, Bool_t force)
785 // return pointer to Q histogram
786 // if force is true create a new histogram if it doesn't exist allready
788 TObjArray *arr = &fHistoQArray;
789 return GetHisto(sector, arr, fNbinsQ, fXminQ, fXmaxQ, "Q", force);
791 //_____________________________________________________________________
792 TH2S* AliTPCCalibPulser::GetHistoRMS(Int_t sector, Bool_t force)
795 // return pointer to Q histogram
796 // if force is true create a new histogram if it doesn't exist allready
798 TObjArray *arr = &fHistoRMSArray;
799 return GetHisto(sector, arr, fNbinsRMS, fXminRMS, fXmaxRMS, "RMS", force);
801 //_____________________________________________________________________
802 TVectorF* AliTPCCalibPulser::GetPadInfoEvent(Int_t sector, TObjArray *arr, Bool_t force)
805 // return pointer to Pad Info from 'arr' for the current event and sector
806 // if force is true create it if it doesn't exist allready
808 if ( !force || arr->UncheckedAt(sector) )
809 return (TVectorF*)arr->UncheckedAt(sector);
811 TVectorF *vect = new TVectorF(fROC->GetNChannels(sector));
812 arr->AddAt(vect,sector);
815 //_____________________________________________________________________
816 TVectorF* AliTPCCalibPulser::GetPadTimesEvent(Int_t sector, Bool_t force)
819 // return pointer to Pad Times Array for the current event and sector
820 // if force is true create it if it doesn't exist allready
822 TObjArray *arr = &fPadTimesArrayEvent;
823 return GetPadInfoEvent(sector,arr,force);
825 //_____________________________________________________________________
826 TVectorF* AliTPCCalibPulser::GetPadQEvent(Int_t sector, Bool_t force)
829 // return pointer to Pad Q Array for the current event and sector
830 // if force is true create it if it doesn't exist allready
831 // for debugging purposes only
834 TObjArray *arr = &fPadQArrayEvent;
835 return GetPadInfoEvent(sector,arr,force);
837 //_____________________________________________________________________
838 TVectorF* AliTPCCalibPulser::GetPadRMSEvent(Int_t sector, Bool_t force)
841 // return pointer to Pad RMS Array for the current event and sector
842 // if force is true create it if it doesn't exist allready
843 // for debugging purposes only
845 TObjArray *arr = &fPadRMSArrayEvent;
846 return GetPadInfoEvent(sector,arr,force);
848 //_____________________________________________________________________
849 TVectorF* AliTPCCalibPulser::GetPadPedestalEvent(Int_t sector, Bool_t force)
852 // return pointer to Pad RMS Array for the current event and sector
853 // if force is true create it if it doesn't exist allready
854 // for debugging purposes only
856 TObjArray *arr = &fPadPedestalArrayEvent;
857 return GetPadInfoEvent(sector,arr,force);
859 //_____________________________________________________________________
860 AliTPCCalROC* AliTPCCalibPulser::GetCalRoc(Int_t sector, TObjArray* arr, Bool_t force) const
863 // return pointer to ROC Calibration
864 // if force is true create a new histogram if it doesn't exist allready
866 if ( !force || arr->UncheckedAt(sector) )
867 return (AliTPCCalROC*)arr->UncheckedAt(sector);
869 // if we are forced and histogram doesn't yes exist create it
871 // new AliTPCCalROC for T0 information. One value for each pad!
872 AliTPCCalROC *croc = new AliTPCCalROC(sector);
873 arr->AddAt(croc,sector);
876 //_____________________________________________________________________
877 AliTPCCalROC* AliTPCCalibPulser::GetCalRocT0(Int_t sector, Bool_t force)
880 // return pointer to Carge ROC Calibration
881 // if force is true create a new histogram if it doesn't exist allready
883 TObjArray *arr = &fCalRocArrayT0;
884 return GetCalRoc(sector, arr, force);
886 //_____________________________________________________________________
887 AliTPCCalROC* AliTPCCalibPulser::GetCalRocQ(Int_t sector, Bool_t force)
890 // return pointer to T0 ROC Calibration
891 // if force is true create a new histogram if it doesn't exist allready
893 TObjArray *arr = &fCalRocArrayQ;
894 return GetCalRoc(sector, arr, force);
896 //_____________________________________________________________________
897 AliTPCCalROC* AliTPCCalibPulser::GetCalRocRMS(Int_t sector, Bool_t force)
900 // return pointer to signal width ROC Calibration
901 // if force is true create a new histogram if it doesn't exist allready
903 TObjArray *arr = &fCalRocArrayRMS;
904 return GetCalRoc(sector, arr, force);
906 //_____________________________________________________________________
907 AliTPCCalROC* AliTPCCalibPulser::GetCalRocOutliers(Int_t sector, Bool_t force)
910 // return pointer to Outliers
911 // if force is true create a new histogram if it doesn't exist allready
913 TObjArray *arr = &fCalRocArrayOutliers;
914 return GetCalRoc(sector, arr, force);
916 //_____________________________________________________________________
917 void AliTPCCalibPulser::ResetEvent()
920 // Reset global counters -- Should be called before each event is processed
929 fPadTimesArrayEvent.Delete();
930 fPadQArrayEvent.Delete();
931 fPadRMSArrayEvent.Delete();
932 fPadPedestalArrayEvent.Delete();
934 for ( Int_t i=0; i<72; ++i ){
936 fVTime0OffsetCounter[i]=0;
939 //_____________________________________________________________________
940 void AliTPCCalibPulser::ResetPad()
943 // Reset pad infos -- Should be called after a pad has been processed
945 for (Int_t i=fFirstTimeBin; i<fLastTimeBin+1; ++i)
952 //_____________________________________________________________________
953 void AliTPCCalibPulser::Merge(AliTPCCalibPulser *sig)
956 // Merge reference histograms of sig to the current AliTPCCalibPulser
960 for (Int_t iSec=0; iSec<72; ++iSec){
961 TH2S *hRefQmerge = sig->GetHistoQ(iSec);
962 TH2S *hRefT0merge = sig->GetHistoT0(iSec);
963 TH2S *hRefRMSmerge = sig->GetHistoRMS(iSec);
967 TDirectory *dir = hRefQmerge->GetDirectory(); hRefQmerge->SetDirectory(0);
968 TH2S *hRefQ = GetHistoQ(iSec);
969 if ( hRefQ ) hRefQ->Add(hRefQmerge);
971 TH2S *hist = new TH2S(*hRefQmerge);
972 hist->SetDirectory(0);
973 fHistoQArray.AddAt(hist, iSec);
975 hRefQmerge->SetDirectory(dir);
978 TDirectory *dir = hRefT0merge->GetDirectory(); hRefT0merge->SetDirectory(0);
979 TH2S *hRefT0 = GetHistoT0(iSec);
980 if ( hRefT0 ) hRefT0->Add(hRefT0merge);
982 TH2S *hist = new TH2S(*hRefT0merge);
983 hist->SetDirectory(0);
984 fHistoT0Array.AddAt(hist, iSec);
986 hRefT0merge->SetDirectory(dir);
989 TDirectory *dir = hRefRMSmerge->GetDirectory(); hRefRMSmerge->SetDirectory(0);
990 TH2S *hRefRMS = GetHistoRMS(iSec);
991 if ( hRefRMS ) hRefRMS->Add(hRefRMSmerge);
993 TH2S *hist = new TH2S(*hRefRMSmerge);
994 hist->SetDirectory(0);
995 fHistoRMSArray.AddAt(hist, iSec);
997 hRefRMSmerge->SetDirectory(dir);
1002 //_____________________________________________________________________
1003 void AliTPCCalibPulser::Analyse()
1006 // Calculate calibration constants
1010 TVectorD paramT0(3);
1011 TVectorD paramRMS(3);
1012 TMatrixD dummy(3,3);
1014 for (Int_t iSec=0; iSec<72; ++iSec){
1015 TH2S *hT0 = GetHistoT0(iSec);
1016 if (!hT0 ) continue;
1018 AliTPCCalROC *rocQ = GetCalRocQ (iSec,kTRUE);
1019 AliTPCCalROC *rocT0 = GetCalRocT0 (iSec,kTRUE);
1020 AliTPCCalROC *rocRMS = GetCalRocRMS(iSec,kTRUE);
1021 AliTPCCalROC *rocOut = GetCalRocOutliers(iSec,kTRUE);
1023 TH2S *hQ = GetHistoQ(iSec);
1024 TH2S *hRMS = GetHistoRMS(iSec);
1026 Short_t *arrayhQ = hQ->GetArray();
1027 Short_t *arrayhT0 = hT0->GetArray();
1028 Short_t *arrayhRMS = hRMS->GetArray();
1030 UInt_t nChannels = fROC->GetNChannels(iSec);
1038 for (UInt_t iChannel=0; iChannel<nChannels; ++iChannel){
1040 Float_t cogTime0 = -1000;
1041 Float_t cogQ = -1000;
1042 Float_t cogRMS = -1000;
1045 Int_t offsetQ = (fNbinsQ+2)*(iChannel+1)+1;
1046 Int_t offsetT0 = (fNbinsT0+2)*(iChannel+1)+1;
1047 Int_t offsetRMS = (fNbinsRMS+2)*(iChannel+1)+1;
1049 AliMathBase::FitGaus(arrayhQ+offsetQ,fNbinsQ,fXminQ,fXmaxQ,¶mQ,&dummy);
1050 AliMathBase::FitGaus(arrayhT0+offsetT0,fNbinsT0,fXminT0,fXmaxT0,¶mT0,&dummy);
1051 AliMathBase::FitGaus(arrayhRMS+offsetRMS,fNbinsRMS,fXminRMS,fXmaxRMS,¶mRMS,&dummy);
1053 cogTime0 = paramT0[1];
1054 cogRMS = paramRMS[1];
1056 cogQ = AliMathBase::GetCOG(arrayhQ+offsetQ,fNbinsQ,fXminQ,fXmaxQ);
1057 cogTime0 = AliMathBase::GetCOG(arrayhT0+offsetT0,fNbinsT0,fXminT0,fXmaxT0);
1058 cogRMS = AliMathBase::GetCOG(arrayhRMS+offsetRMS,fNbinsRMS,fXminRMS,fXmaxRMS);
1061 if ( (cogQ < ??) && (cogTime0 > ??) && (cogTime0<??) && ( cogRMS>??) ){
1068 rocQ->SetValue(iChannel, cogQ*cogQ);
1069 rocT0->SetValue(iChannel, cogTime0);
1070 rocRMS->SetValue(iChannel, cogRMS);
1071 rocOut->SetValue(iChannel, cogOut);
1075 if ( fDebugLevel > 0 ){
1076 if ( !fDebugStreamer ) {
1078 TDirectory *backup = gDirectory;
1079 fDebugStreamer = new TTreeSRedirector("deb2.root");
1080 if ( backup ) backup->cd(); //we don't want to be cd'd to the debug streamer
1083 while ( iChannel > (fROC->GetRowIndexes(iSec)[row]+fROC->GetNPads(iSec,row)-1) ) row++;
1084 pad = iChannel-fROC->GetRowIndexes(iSec)[row];
1085 padc = pad-(fROC->GetNPads(iSec,row)/2);
1087 (*fDebugStreamer) << "DataEnd" <<
1088 "Sector=" << iSec <<
1092 "PadSec=" << iChannel <<
1094 "T0=" << cogTime0 <<
1101 delete fDebugStreamer;
1102 fDebugStreamer = 0x0;
1104 //_____________________________________________________________________
1105 void AliTPCCalibPulser::DumpToFile(const Char_t *filename, const Char_t *dir, Bool_t append)
1108 // Write class to file
1117 option = "recreate";
1119 TDirectory *backup = gDirectory;
1120 TFile f(filename,option.Data());
1122 if ( !sDir.IsNull() ){
1123 f.mkdir(sDir.Data());
1129 if ( backup ) backup->cd();
1131 //_____________________________________________________________________
1132 //_________________________ Test Functions ___________________________
1133 //_____________________________________________________________________
1134 TObjArray* AliTPCCalibPulser::TestBinning()
1137 // Function to test the binning of the reference histograms
1138 // type: T0, Q or RMS
1139 // mode: 0 - number of filled bins per channel
1140 // 1 - number of empty bins between filled bins in one ROC
1141 // returns TObjArray with the test histograms type*2+mode:
1142 // position 0 = T0,0 ; 1 = T0,1 ; 2 = Q,0 ...
1145 TObjArray *histArray = new TObjArray(6);
1146 const Char_t *type[] = {"T0","Q","RMS"};
1147 Int_t fNbins[3] = {fNbinsT0,fNbinsQ,fNbinsRMS};
1149 for (Int_t itype = 0; itype<3; ++itype){
1150 for (Int_t imode=0; imode<2; ++imode){
1151 Int_t icount = itype*2+imode;
1152 histArray->AddAt(new TH1F(Form("hTestBinning%s%d",type[itype],imode),
1153 Form("Test Binning of '%s', mode - %d",type[itype],imode),
1162 for (Int_t itype = 0; itype<3; ++itype){
1163 for (Int_t iSec=0; iSec<72; ++iSec){
1164 if ( itype == 0 ) hRef = GetHistoT0(iSec);
1165 if ( itype == 1 ) hRef = GetHistoQ(iSec);
1166 if ( itype == 2 ) hRef = GetHistoRMS(iSec);
1167 if ( hRef == 0x0 ) continue;
1168 array = (hRef->GetArray());
1169 UInt_t nChannels = fROC->GetNChannels(iSec);
1172 for (UInt_t iChannel=0; iChannel<nChannels; ++iChannel){
1174 Int_t offset = (fNbins[itype]+2)*(iChannel+1)+1;
1177 for (Int_t iBin=0; iBin<fNbins[itype]; ++iBin){
1178 if ( array[offset+iBin]>0 ) {
1180 if ( c1 && c2 ) nempty++;
1183 else if ( c1 ) c2 = 1;
1186 ((TH1F*)histArray->At(itype*2))->Fill(nfilled);
1188 ((TH1F*)histArray->At(itype*2+1))->Fill(iSec,nempty);