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 **************************************************************************/
16 ///////////////////////////////////////////////////////////////////////
18 // Produces the data needed to calculate the TOF quality assurance. //
19 // QA objects are 1 & 2 Dimensional histograms. //
20 // author: S.Arcelli //
22 ///////////////////////////////////////////////////////////////////////
24 #include <TClonesArray.h>
31 #include "AliESDEvent.h"
32 #include "AliESDtrack.h"
33 #include "AliQAChecker.h"
34 #include "AliRawReader.h"
36 #include "AliTOFcluster.h"
37 #include "AliTOFQADataMakerRec.h"
38 #include "AliTOFRawStream.h"
39 #include "AliTOFrawData.h"
40 #include "AliTOFGeometry.h"
41 #include "AliTOFdigit.h"
43 ClassImp(AliTOFQADataMakerRec)
45 //____________________________________________________________________________
46 AliTOFQADataMakerRec::AliTOFQADataMakerRec() :
47 AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kTOF), "TOF Quality Assurance Data Maker")
54 //____________________________________________________________________________
55 AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
61 SetName((const char*)qadm.GetName()) ;
62 SetTitle((const char*)qadm.GetTitle());
65 //__________________________________________________________________
66 AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
69 // assignment operator.
71 this->~AliTOFQADataMakerRec();
72 new(this) AliTOFQADataMakerRec(qadm);
76 //____________________________________________________________________________
77 void AliTOFQADataMakerRec::InitRaws()
80 // create Raws histograms in Raws subdir
83 const Bool_t expert = kTRUE ;
84 const Bool_t saveCorr = kTRUE ;
85 const Bool_t image = kTRUE ;
87 TH1F * h0 = new TH1F("hTOFRaws", "Number of TOF Raws;TOF raw number [10 power];Counts ",301, -1.02, 5.) ;
89 Add2RawsList(h0, 0, !expert, image, !saveCorr) ;
91 TH1F * h1 = new TH1F("hTOFRawsTime", "Raws Time Spectrum in TOF (ns);Measured TOF time [ns];Counts", 2000, 0., 200) ;
93 Add2RawsList(h1, 1, !expert, image, !saveCorr) ;
95 TH1F * h2 = new TH1F("hTOFRawsToT", "Raws ToT Spectrum in TOF (ns);Measured TOT time [ns];Counts", 500, 0., 50) ;
97 Add2RawsList(h2, 2, !expert, image, !saveCorr) ;
99 TH2F * h3 = new TH2F("hTOFRawsClusMap","Raws vs TOF eta-phi;2*strip+padz (eta);48*sector+padx (phi)",183, -0.5, 182.5,865,-0.5,864.5) ;
101 h3->GetYaxis()->SetTitleOffset(1.15);
102 Add2RawsList(h3, 3, !expert, image, !saveCorr) ;
106 //____________________________________________________________________________
107 void AliTOFQADataMakerRec::InitDigits()
110 // create Digits histograms in Digits subdir
113 const Bool_t expert = kTRUE ;
114 const Bool_t image = kTRUE ;
116 TH1F * h0 = new TH1F("hTOFDigits", "Number of TOF Digits;TOF digit number [10 power];Counts ",301, -1.02, 5.) ;
118 Add2DigitsList(h0, 0, !expert, image) ;
120 TH1F * h1 = new TH1F("hTOFDigitsTime", "Digits Time Spectrum in TOF (ns);Digitized TOF time [ns];Counts", 2000, 0., 200) ;
122 Add2DigitsList(h1, 1, !expert, image) ;
124 TH1F * h2 = new TH1F("hTOFDigitsToT", "Digits ToT Spectrum in TOF (ns);Digitized TOF time [ns];Counts", 500, 0., 50) ;
126 Add2DigitsList(h2, 2, !expert, image) ;
128 TH2F * h3 = new TH2F("hTOFDigitsClusMap","Digits vs TOF eta-phi;2*strip+padz (eta);48*sector+padx (phi)",183, -0.5, 182.5,865,-0.5,864.5) ;
130 h3->GetYaxis()->SetTitleOffset(1.15);
131 Add2DigitsList(h3, 3, !expert, image) ;
135 //____________________________________________________________________________
136 void AliTOFQADataMakerRec::InitRecPoints()
139 // create RecPoints histograms in RecPoints subdir
142 const Bool_t expert = kTRUE ;
143 const Bool_t image = kTRUE ;
145 TH1F * h0 = new TH1F("hTOFRecPoints", "Number of TOF RecPoints;TOF recPoint number [10 power];Counts",301, -1.02, 5.) ;
147 Add2RecPointsList(h0, 0, !expert, image) ;
149 TH1F * h1 = new TH1F("hTOFRecPointsTime", "RecPoints Time Spectrum in TOF (ns);Calibrated TOF time [ns];Counts", 2000, 0., 200) ;
151 Add2RecPointsList(h1, 1, !expert, image) ;
153 TH1F * h2 = new TH1F("hTOFRecPointsRawTime", "RecPoints raw Time Spectrum in TOF (ns);Measured TOF time [ns];Counts", 2000, 0., 200) ;
155 Add2RecPointsList(h2, 2, !expert, image) ;
157 TH1F * h3 = new TH1F("hTOFRecPointsToT", "RecPoints ToT Spectrum in TOF (ns);Measured TOT [ns];Counts", 500, 0., 50) ;
159 Add2RecPointsList(h3, 3, !expert, image) ;
161 TH2F * h4 = new TH2F("hTOFRecPointsClusMap","RecPoints vs TOF phi-eta;2*strip+padz (eta);48*sector+padx (phi)",183, -0.5, 182.5,865,-0.5,864.5) ;
163 h4->GetYaxis()->SetTitleOffset(1.15);
164 Add2RecPointsList(h4, 4, !expert, image) ;
168 //____________________________________________________________________________
169 void AliTOFQADataMakerRec::InitESDs()
172 //create ESDs histograms in ESDs subdir
175 const Bool_t expert = kTRUE ;
176 const Bool_t image = kTRUE ;
177 TH1F * h0 = new TH1F("hTOFESDs", "Number of matched TOF tracks over ESDs;Number of TOF matched ESD tracks [10 power];Counts", 250, -1., 4.) ;
179 Add2ESDsList(h0, 0, !expert, image) ;
181 TH1F * h1 = new TH1F("hTOFESDsTime", "Time Spectrum in TOF (ns);Calibrated TOF time [ns];Counts", 2000, 0., 200) ;
183 Add2ESDsList(h1, 1, !expert, image) ;
185 TH1F * h2 = new TH1F("hTOFESDsRawTime", "raw Time Spectrum in TOF (ns);Measured TOF time [ns];Counts", 2000, 0., 200) ;
187 Add2ESDsList(h2, 2, !expert, image) ;
189 TH1F * h3 = new TH1F("hTOFESDsToT", "ToT Spectrum in TOF (ns);Measured TOF time [ns];Counts", 500, 0., 50) ;
191 Add2ESDsList(h3, 3, !expert, image) ;
193 TH1F * h4 = new TH1F("hTOFESDsPID", "Fraction of matched TOF tracks with good PID flag (%);Fraction of TOF matched ESD tracks with good flag [%];Counts", 101, 0., 101.) ;
195 Add2ESDsList(h4, 4, !expert, image) ;
198 //____________________________________________________________________________
199 void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
202 // makes data from Raws
205 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
206 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
213 TClonesArray * clonesRawData;
214 AliTOFRawStream tofInput(rawReader);
215 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
217 tofInput.LoadRawData(iDDL);
218 clonesRawData = (TClonesArray*)tofInput.GetRawData();
219 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
220 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
221 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
223 GetRawsData(1)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
224 GetRawsData(2)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
226 tofInput.EquipmentId2VolumeId(iDDL,
227 tofRawDatum->GetTRM(),
228 tofRawDatum->GetTRMchain(),
229 tofRawDatum->GetTDC(),
230 tofRawDatum->GetTDCchannel(),
233 GetMapIndeces(in,out);
234 GetRawsData(3)->Fill( out[0],out[1]) ;//raw map
238 clonesRawData->Clear();
244 GetRawsData(0)->Fill(-1.) ;
246 GetRawsData(0)->Fill(TMath::Log10(nentries)) ;
250 //____________________________________________________________________________
251 void AliTOFQADataMakerRec::MakeDigits(TClonesArray * digits)
254 // makes data from Digits
256 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
257 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
261 Int_t nentries=digits->GetEntriesFast();
263 GetDigitsData(0)->Fill(-1.) ;
265 GetDigitsData(0)->Fill(TMath::Log10(nentries)) ;
269 AliTOFdigit * digit ;
270 while ( (digit = dynamic_cast<AliTOFdigit *>(next())) ) {
272 GetDigitsData(1)->Fill( digit->GetTdc()*tdc2ns) ;//in ns
273 GetDigitsData(2)->Fill( digit->GetToT()*tot2ns) ;//in ns
275 in[0] = digit->GetSector();
276 in[1] = digit->GetPlate();
277 in[2] = digit->GetStrip();
278 in[3] = digit->GetPadx();
279 in[4]= digit->GetPadz();
280 GetMapIndeces(in,out);
281 GetDigitsData(3)->Fill( out[0],out[1]) ;//digit map
286 //____________________________________________________________________________
287 void AliTOFQADataMakerRec::MakeDigits(TTree * digitTree)
290 // makes data from Digit Tree
292 TClonesArray * digits = new TClonesArray("AliTOFdigit", 1000) ;
294 TBranch * branch = digitTree->GetBranch("TOF") ;
296 AliError("TOF branch in Digit Tree not found") ;
299 branch->SetAddress(&digits) ;
300 branch->GetEntry(0) ;
304 //____________________________________________________________________________
305 void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
308 // Make data from Clusters
311 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
312 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
317 TBranch *branch=clustersTree->GetBranch("TOF");
319 AliError("can't get the branch with the TOF clusters !");
323 static TClonesArray dummy("AliTOFcluster",10000);
325 TClonesArray *clusters=&dummy;
326 branch->SetAddress(&clusters);
329 clustersTree->GetEvent(0);
331 Int_t nentries=clusters->GetEntriesFast();
333 GetRecPointsData(0)->Fill(-1.) ;
335 GetRecPointsData(0)->Fill(TMath::Log10(nentries)) ;
338 TIter next(clusters) ;
340 while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
341 GetRecPointsData(1)->Fill(c->GetTDC()*tdc2ns);
342 GetRecPointsData(2)->Fill(c->GetTDCRAW()*tdc2ns);
343 GetRecPointsData(3)->Fill(c->GetToT()*tot2ns);
345 in[0] = c->GetDetInd(0);
346 in[1] = c->GetDetInd(1);
347 in[2] = c->GetDetInd(2);
348 in[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
349 in[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
351 GetMapIndeces(in,out);
352 GetRecPointsData(4)->Fill(out[0],out[1]);
357 //____________________________________________________________________________
358 void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
361 // make QA data from ESDs
363 Int_t ntrk = esd->GetNumberOfTracks() ;
367 AliESDtrack *track=esd->GetTrack(ntrk);
368 Double_t tofTime=track->GetTOFsignal()*1E-3;//in ns
369 Double_t tofTimeRaw=track->GetTOFsignalRaw()*1E-3;//in ns
370 Double_t tofToT=track->GetTOFsignalToT(); //in ns
371 if(!(tofTime>0))continue;
373 GetESDsData(1)->Fill(tofTime);
374 GetESDsData(2)->Fill(tofTimeRaw);
375 GetESDsData(3)->Fill(tofToT);
376 //check how many tracks where ESD PID is ok
377 UInt_t status=track->GetStatus();
378 if (((status&AliESDtrack::kESDpid)==0) ||
379 ((status&AliESDtrack::kTOFpid)==0)) continue;
385 GetESDsData(0)->Fill(-1.) ;
387 GetESDsData(0)->Fill(TMath::Log10(nentries)) ;
391 Float_t ratio = (Float_t)ntofpid/(Float_t)ntof*100.;
392 GetESDsData(4)->Fill(ratio) ;
397 //____________________________________________________________________________
398 void AliTOFQADataMakerRec::StartOfDetectorCycle()
401 //Detector specific actions at start of cycle
405 //____________________________________________________________________________
406 void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray ** list)
408 //Detector specific actions at end of cycle
409 // do the QA checking
411 AliQAChecker::Instance()->Run(AliQAv1::kTOF, task, list) ;
413 //____________________________________________________________________________
414 void AliTOFQADataMakerRec::GetMapIndeces(Int_t* in , Int_t* out)
417 //return appropriate indeces for the theta-phi map
420 Int_t npadX = AliTOFGeometry::NpadX();
421 Int_t npadZ = AliTOFGeometry::NpadZ();
422 Int_t nStripA = AliTOFGeometry::NStripA();
423 Int_t nStripB = AliTOFGeometry::NStripB();
424 Int_t nStripC = AliTOFGeometry::NStripC();
426 Int_t isector = in[0];
427 Int_t iplate = in[1];
428 Int_t istrip = in[2];
432 Int_t stripOffset = 0;
438 stripOffset = nStripC;
441 stripOffset = nStripC+nStripB;
444 stripOffset = nStripC+nStripB+nStripA;
447 stripOffset = nStripC+nStripB+nStripA+nStripB;
450 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
453 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
454 Int_t phiindex=npadX*isector+ipadX+1;