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 Produces the data needed to calculate the TOF quality assurance.
17 QA objects are 1 & 2 Dimensional histograms.
21 #include <TClonesArray.h>
26 #include "AliESDEvent.h"
27 #include "AliESDtrack.h"
28 #include "AliTOFcluster.h"
29 #include "AliTOFdigit.h"
30 #include "AliTOFSDigit.h"
31 #include "AliTOFhitT0.h"
32 #include "AliTOFQADataMaker.h"
33 #include "AliRawReader.h"
34 #include "AliTOFRawStream.h"
35 #include "AliTOFrawData.h"
38 ClassImp(AliTOFQADataMaker)
40 //____________________________________________________________________________
41 AliTOFQADataMaker::AliTOFQADataMaker() :
42 AliQADataMaker(AliQA::GetDetName(AliQA::kTOF), "TOF Quality Assurance Data Maker")
49 //____________________________________________________________________________
50 AliTOFQADataMaker::AliTOFQADataMaker(const AliTOFQADataMaker& qadm) :
56 SetName((const char*)qadm.GetName()) ;
57 SetTitle((const char*)qadm.GetTitle());
60 //__________________________________________________________________
61 AliTOFQADataMaker& AliTOFQADataMaker::operator = (const AliTOFQADataMaker& qadm )
64 // assignment operator.
66 this->~AliTOFQADataMaker();
67 new(this) AliTOFQADataMaker(qadm);
71 //____________________________________________________________________________
72 void AliTOFQADataMaker::InitHits()
75 // create Hits histograms in Hits subdir
77 TH1F * h0 = new TH1F("hTOFHits", "Number of TOF Hits ",301, -1.02, 5.) ;
81 TH1F * h1 = new TH1F("hTOFHitsTime", "Hits Time Spectrum in TOF (ns)", 2000, 0., 200) ;
85 TH1F * h2 = new TH1F("hTOFHitsLength", "Length Spectrum in TOF (cm)", 500, 0., 500) ;
89 TH2F * h3 = new TH2F("hTOFHitsClusMap","Hits vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
94 //____________________________________________________________________________
95 void AliTOFQADataMaker::InitDigits()
98 // create Digits histograms in Digits subdir
100 TH1F * h0 = new TH1F("hTOFDigits", "Number of TOF Digits ",301, -1.02, 5.) ; h0->Sumw2() ;
101 Add2DigitsList(h0, 0) ;
103 TH1F * h1 = new TH1F("hTOFDigitsTime", "Digits Time Spectrum in TOF (ns)", 2000, 0., 200) ;
105 Add2DigitsList(h1, 1) ;
107 TH1F * h2 = new TH1F("hTOFDigitsToT", "Digits ToT Spectrum in TOF (ns)", 500, 0., 50) ;
109 Add2DigitsList(h2, 2) ;
111 TH2F * h3 = new TH2F("hTOFDigitsClusMap","Digits vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
113 Add2DigitsList(h3, 3) ;
117 //____________________________________________________________________________
118 void AliTOFQADataMaker::InitSDigits()
121 // create SDigits histograms in SDigits subdir
123 TH1F * h0 = new TH1F("hTOFSDigits", "Number of TOF SDigits ",301, -1.02, 5.) ; h0->Sumw2() ;
124 Add2SDigitsList(h0, 0) ;
126 TH1F * h1 = new TH1F("hTOFSDigitsTime", "SDigits Time Spectrum in TOF (ns)", 2000, 0., 200) ;
128 Add2SDigitsList(h1, 1) ;
130 TH2F * h2 = new TH2F("hTOFSDigitsClusMap","SDigits vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
132 Add2SDigitsList(h2, 2) ;
136 //____________________________________________________________________________
137 void AliTOFQADataMaker::InitRaws()
140 // create Raws histograms in Raws subdir
142 TH1F * h0 = new TH1F("hTOFRaws", "Number of TOF Raws ",301, -1.02, 5.) ; h0->Sumw2() ;
143 Add2RawsList(h0, 0) ;
145 TH1F * h1 = new TH1F("hTOFRawsTime", "Raws Time Spectrum in TOF (ns)", 2000, 0., 200) ;
147 Add2RawsList(h1, 1) ;
149 TH1F * h2 = new TH1F("hTOFRawsToT", "Raws ToT Spectrum in TOF (ns)", 500, 0., 50) ;
151 Add2RawsList(h2, 2) ;
153 TH2F * h3 = new TH2F("hTOFRawsClusMap","Raws vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
155 Add2RawsList(h3, 3) ;
159 //____________________________________________________________________________
160 void AliTOFQADataMaker::InitRecPoints()
163 // create RecPoints histograms in RecPoints subdir
165 TH1F * h0 = new TH1F("hTOFRecPoints", "Number of TOF RecPoints ",301, -1.02, 5.) ; h0->Sumw2() ;
166 Add2RecPointsList(h0, 0) ;
168 TH1F * h1 = new TH1F("hTOFRecPointsTime", "RecPoints Time Spectrum in TOF (ns)", 2000, 0., 200) ;
170 Add2RecPointsList(h1, 1) ;
172 TH1F * h2 = new TH1F("hTOFRecPointsRawTime", "RecPoints raw Time Spectrum in TOF (ns)", 2000, 0., 200) ;
174 Add2RecPointsList(h2, 2) ;
176 TH1F * h3 = new TH1F("hTOFRecPointsToT", "RecPoints ToT Spectrum in TOF (ns)", 500, 0., 50) ;
178 Add2RecPointsList(h3, 3) ;
180 TH2F * h4 = new TH2F("hTOFRecPointsClusMap","RecPoints vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ;
182 Add2RecPointsList(h4, 4) ;
186 //____________________________________________________________________________
187 void AliTOFQADataMaker::InitESDs()
190 //create ESDs histograms in ESDs subdir
192 TH1F * h0 = new TH1F("hTOFESDs", "Number of matched TOF tracks over ESDs", 250, -1., 4.) ;
194 Add2ESDsList(h0, 0) ;
196 TH1F * h1 = new TH1F("hTOFESDsTime", "Time Spectrum in TOF (ns)", 2000, 0., 200) ;
198 Add2ESDsList(h1, 1) ;
200 TH1F * h2 = new TH1F("hTOFESDsRawTime", "raw Time Spectrum in TOF (ns)", 2000, 0., 200) ;
202 Add2ESDsList(h2, 2) ;
204 TH1F * h3 = new TH1F("hTOFESDsToT", "ToT Spectrum in TOF (ns)", 500, 0., 50) ;
206 Add2ESDsList(h3, 3) ;
208 TH1F * h4 = new TH1F("hTOFESDsPID", "Fraction of matched TOF tracks with good PID glag", 100, 0., 1.) ;
210 Add2ESDsList(h4, 4) ;
214 //____________________________________________________________________________
215 void AliTOFQADataMaker::MakeHits(TClonesArray * hits)
218 //make QA data from Hits
224 Int_t nentries=hits->GetEntriesFast();
225 if(nentries==0)nentries=-1;
226 GetHitsData(0)->Fill(TMath::Log10(nentries)) ;
230 while ( (hit = dynamic_cast<AliTOFhitT0 *>(next())) ) {
232 GetHitsData(1)->Fill( hit->GetTof()*1.E9) ;//in ns
233 GetHitsData(2)->Fill( hit->GetLen()) ;//in cm
235 in[0] = hit->GetSector();
236 in[1] = hit->GetPlate();
237 in[2]= hit->GetStrip();
238 in[3]= hit->GetPadx();
239 in[4]= hit->GetPadz();
240 GetMapIndeces(in,out);
241 GetHitsData(3)->Fill( out[0],out[1]) ;//hit map
247 //____________________________________________________________________________
248 void AliTOFQADataMaker::MakeHits(TTree * hitTree)
251 // make QA data from Hit Tree
254 AliError("can't get the tree with TOF hits !");
258 TBranch * branch = hitTree->GetBranch("TOF") ;
261 AliError("TOF branch in Hit Tree not found") ;
265 TClonesArray * hits = new TClonesArray("AliTOFhitT0", 1000);
266 TClonesArray * dummy = new TClonesArray("AliTOFhitT0", 1000);
267 branch->SetAddress(&dummy);
269 for (Int_t ientry = 0 ; ientry < branch->GetEntries() ; ientry++) {
270 branch->GetEntry(ientry) ;
271 for (Int_t ihit = 0 ; ihit < dummy->GetEntries() ; ihit++) {
272 AliTOFhitT0 * hit = dynamic_cast<AliTOFhitT0 *> (dummy->At(ihit)) ;
273 new((*hits)[index]) AliTOFhitT0(*hit) ;
284 //____________________________________________________________________________
285 void AliTOFQADataMaker::MakeDigits(TClonesArray * digits)
288 // makes data from Digits
290 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
291 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
295 Int_t nentries=digits->GetEntriesFast();
296 if(nentries==0)nentries=-1;
297 GetDigitsData(0)->Fill(TMath::Log10(nentries)) ;
299 AliTOFdigit * digit ;
300 while ( (digit = dynamic_cast<AliTOFdigit *>(next())) ) {
302 GetDigitsData(1)->Fill( digit->GetTdc()*tdc2ns) ;//in ns
303 GetDigitsData(2)->Fill( digit->GetToT()*tot2ns) ;//in ns
305 in[0] = digit->GetSector();
306 in[1] = digit->GetPlate();
307 in[2] = digit->GetStrip();
308 in[3] = digit->GetPadx();
309 in[4]= digit->GetPadz();
310 GetMapIndeces(in,out);
311 GetDigitsData(3)->Fill( out[0],out[1]) ;//digit map
317 //____________________________________________________________________________
318 void AliTOFQADataMaker::MakeDigits(TTree * digitTree)
321 // makes data from Digit Tree
323 TClonesArray * digits = new TClonesArray("AliTOFdigit", 1000) ;
325 TBranch * branch = digitTree->GetBranch("TOF") ;
327 AliError("TOF branch in Digit Tree not found") ;
330 branch->SetAddress(&digits) ;
331 branch->GetEntry(0) ;
335 //____________________________________________________________________________
336 void AliTOFQADataMaker::MakeSDigits(TClonesArray * sdigits)
339 // makes data from SDigits
342 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
346 Int_t nentries=sdigits->GetEntriesFast();
347 if(nentries==0)nentries=-1;
348 GetSDigitsData(0)->Fill(TMath::Log10(nentries)) ;
350 TIter next(sdigits) ;
351 AliTOFSDigit * sdigit ;
352 while ( (sdigit = dynamic_cast<AliTOFSDigit *>(next())) ) {
354 for(Int_t i=0;i<sdigit->GetNDigits();i++){
355 GetSDigitsData(1)->Fill( sdigit->GetTdc(i)*tdc2ns) ;//in ns
358 in[0] = sdigit->GetSector();
359 in[1] = sdigit->GetPlate();
360 in[2] = sdigit->GetStrip();
361 in[3] = sdigit->GetPadx();
362 in[4]= sdigit->GetPadz();
363 GetMapIndeces(in,out);
364 GetSDigitsData(2)->Fill( out[0],out[1]) ;//sdigit map
368 //____________________________________________________________________________
369 void AliTOFQADataMaker::MakeSDigits(TTree * sdigitTree)
372 // makes data from SDigit Tree
374 TClonesArray * sdigits = new TClonesArray("AliTOFSDigit", 1000) ;
376 TBranch * branch = sdigitTree->GetBranch("TOF") ;
378 AliError("TOF branch in SDigit Tree not found") ;
381 branch->SetAddress(&sdigits) ;
382 branch->GetEntry(0) ;
383 MakeSDigits(sdigits) ;
386 //____________________________________________________________________________
387 void AliTOFQADataMaker::MakeRaws(AliRawReader* rawReader)
390 // makes data from Raws
393 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
394 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
401 TClonesArray * clonesRawData;
402 AliTOFRawStream tofInput(rawReader);
403 for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
405 tofInput.LoadRawData(iDDL);
406 clonesRawData = (TClonesArray*)tofInput.GetRawData();
407 for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
408 AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
409 if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
411 GetRawsData(1)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
412 GetRawsData(2)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
414 tofInput.EquipmentId2VolumeId(iDDL,
415 tofRawDatum->GetTRM(),
416 tofRawDatum->GetTRMchain(),
417 tofRawDatum->GetTDC(),
418 tofRawDatum->GetTDCchannel(),
421 GetMapIndeces(in,out);
422 GetRawsData(3)->Fill( out[0],out[1]) ;//raw map
426 clonesRawData->Clear();
431 if(nentries==0)nentries=-1;
432 GetRawsData(0)->Fill(TMath::Log10(nentries)) ;
436 //____________________________________________________________________________
437 void AliTOFQADataMaker::MakeRecPoints(TTree * clustersTree)
440 // Make data from Clusters
443 Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
444 Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
449 TBranch *branch=clustersTree->GetBranch("TOF");
451 AliError("can't get the branch with the TOF clusters !");
455 TObjArray * clusters = new TObjArray(100) ;
456 branch->SetAddress(&clusters);
457 clustersTree->GetEvent(0);
459 Int_t nentries=clusters->GetEntriesFast();
460 if(nentries==0)nentries=-1;
461 GetRecPointsData(0)->Fill(TMath::Log10(nentries)) ;
463 TIter next(clusters) ;
465 while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
466 GetRecPointsData(1)->Fill(c->GetTDC()*tdc2ns);
467 GetRecPointsData(2)->Fill(c->GetTDCRAW()*tdc2ns);
468 GetRecPointsData(3)->Fill(c->GetToT()*tot2ns);
470 in[0] = c->GetDetInd(0);
471 in[1] = c->GetDetInd(1);
472 in[2] = c->GetDetInd(2);
473 in[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
474 in[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
476 GetMapIndeces(in,out);
477 GetRecPointsData(4)->Fill(out[0],out[1]);
484 //____________________________________________________________________________
485 void AliTOFQADataMaker::MakeESDs(AliESDEvent * esd)
488 // make QA data from ESDs
490 Int_t ntrk = esd->GetNumberOfTracks() ;
494 AliESDtrack *track=esd->GetTrack(ntrk);
495 Double_t tofTime=track->GetTOFsignal()*1E-3;//in ns
496 Double_t tofTimeRaw=track->GetTOFsignalRaw()*1E-3;//in ns
497 Double_t tofToT=track->GetTOFsignalToT(); //in ns
498 if(!(tofTime>0))continue;
500 GetESDsData(1)->Fill(tofTime);
501 GetESDsData(2)->Fill(tofTimeRaw);
502 GetESDsData(3)->Fill(tofToT);
503 //check how many tracks where ESD PID is ok
504 UInt_t status=track->GetStatus();
505 if (((status&AliESDtrack::kESDpid)==0) ||
506 ((status&AliESDtrack::kTOFpid)==0)) continue;
511 if(nentries==0)nentries=-1;
512 GetESDsData(0)->Fill(TMath::Log10(nentries)) ;
513 if(ntof>0)GetESDsData(4)->Fill(ntofpid/ntof) ;
517 //____________________________________________________________________________
518 void AliTOFQADataMaker::StartOfDetectorCycle()
521 //Detector specific actions at start of cycle
525 //____________________________________________________________________________
526 void AliTOFQADataMaker::GetMapIndeces(Int_t* in , Int_t* out)
529 //return appropriate indeces for the theta-phi map
532 Int_t npadX = AliTOFGeometry::NpadX();
533 Int_t npadZ = AliTOFGeometry::NpadZ();
534 Int_t nStripA = AliTOFGeometry::NStripA();
535 Int_t nStripB = AliTOFGeometry::NStripB();
536 Int_t nStripC = AliTOFGeometry::NStripC();
538 Int_t isector = in[0];
539 Int_t iplate = in[1];
540 Int_t istrip = in[2];
544 Int_t stripOffset = 0;
550 stripOffset = nStripC;
553 stripOffset = nStripC+nStripB;
556 stripOffset = nStripC+nStripB+nStripA;
559 stripOffset = nStripC+nStripB+nStripA+nStripB;
562 AliError(Form("Wrong plate number in TOF (%d) !",iplate));
565 Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
566 Int_t phiindex=npadX*isector+ipadX+1;