]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TOF/AliTOFQADataMakerRec.cxx
Splitting of the QA maker into simulation and reconstruction dependent parts (Yves)
[u/mrichter/AliRoot.git] / TOF / AliTOFQADataMakerRec.cxx
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 //                                                                   //
18 //  Produces the data needed to calculate the TOF quality assurance. //
19 //  QA objects are 1 & 2 Dimensional histograms.                     //
20 //  author: S.Arcelli                                                //
21 //                                                                   //
22 ///////////////////////////////////////////////////////////////////////
23
24 #include <TClonesArray.h>
25 #include <TFile.h> 
26 #include <TH1I.h> 
27 #include <TH1F.h> 
28 #include <TH2F.h> 
29 #include "AliESDEvent.h"
30 #include "AliESDtrack.h"
31 #include "AliTOFcluster.h"
32 #include "AliTOFQADataMakerRec.h"
33 #include "AliQAChecker.h"
34 #include "AliRawReader.h"
35 #include "AliTOFRawStream.h"
36 #include "AliTOFrawData.h"
37 #include "AliLog.h"
38
39 ClassImp(AliTOFQADataMakerRec)
40            
41 //____________________________________________________________________________ 
42   AliTOFQADataMakerRec::AliTOFQADataMakerRec() : 
43   AliQADataMakerRec(AliQA::GetDetName(AliQA::kTOF), "TOF Quality Assurance Data Maker")
44 {
45   //
46   // ctor
47   //
48 }
49
50 //____________________________________________________________________________ 
51 AliTOFQADataMakerRec::AliTOFQADataMakerRec(const AliTOFQADataMakerRec& qadm) :
52   AliQADataMakerRec()
53 {
54   //
55   //copy ctor 
56   //
57   SetName((const char*)qadm.GetName()) ; 
58   SetTitle((const char*)qadm.GetTitle()); 
59 }
60
61 //__________________________________________________________________
62 AliTOFQADataMakerRec& AliTOFQADataMakerRec::operator = (const AliTOFQADataMakerRec& qadm )
63 {
64   //
65   // assignment operator.
66   //
67   this->~AliTOFQADataMakerRec();
68   new(this) AliTOFQADataMakerRec(qadm);
69   return *this;
70 }
71  
72 //____________________________________________________________________________ 
73 void AliTOFQADataMakerRec::InitRaws()
74 {
75   //
76   // create Raws histograms in Raws subdir
77   //
78   TH1F * h0 = new TH1F("hTOFRaws",    "Number of TOF Raws ",301, -1.02, 5.) ;   h0->Sumw2() ;
79   Add2RawsList(h0, 0) ;
80
81   TH1F * h1  = new TH1F("hTOFRawsTime", "Raws Time Spectrum in TOF (ns)", 2000, 0., 200) ; 
82   h1->Sumw2() ;
83   Add2RawsList(h1, 1) ;
84
85   TH1F * h2  = new TH1F("hTOFRawsToT", "Raws ToT Spectrum in TOF (ns)", 500, 0., 50) ; 
86   h2->Sumw2() ;
87   Add2RawsList(h2, 2) ;
88
89   TH2F * h3  = new TH2F("hTOFRawsClusMap","Raws vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ; 
90   h3->Sumw2() ;
91   Add2RawsList(h3, 3) ;
92
93 }
94
95 //____________________________________________________________________________ 
96 void AliTOFQADataMakerRec::InitRecPoints()
97 {
98   //
99   // create RecPoints histograms in RecPoints subdir
100   //
101   TH1F * h0 = new TH1F("hTOFRecPoints",    "Number of TOF RecPoints ",301, -1.02, 5.) ;   h0->Sumw2() ;
102   Add2RecPointsList(h0, 0) ;
103
104   TH1F * h1  = new TH1F("hTOFRecPointsTime", "RecPoints Time Spectrum in TOF (ns)", 2000, 0., 200) ; 
105   h1->Sumw2() ;
106   Add2RecPointsList(h1, 1) ;
107
108   TH1F * h2  = new TH1F("hTOFRecPointsRawTime", "RecPoints raw Time Spectrum in TOF (ns)", 2000, 0., 200) ; 
109   h2->Sumw2() ;
110   Add2RecPointsList(h2, 2) ;
111
112   TH1F * h3  = new TH1F("hTOFRecPointsToT", "RecPoints ToT Spectrum in TOF (ns)", 500, 0., 50) ; 
113   h3->Sumw2() ;
114   Add2RecPointsList(h3, 3) ;
115
116   TH2F * h4  = new TH2F("hTOFRecPointsClusMap","RecPoints vs TOF eta-phi",183, -0.5, 182.5,865,-0.5,864.5) ; 
117   h4->Sumw2() ;
118   Add2RecPointsList(h4, 4) ;
119
120 }
121
122 //____________________________________________________________________________ 
123 void AliTOFQADataMakerRec::InitESDs()
124 {
125   //
126   //create ESDs histograms in ESDs subdir
127   //
128   TH1F * h0 = new TH1F("hTOFESDs",    "Number of matched TOF tracks over ESDs",       250, -1., 4.) ;  
129   h0->Sumw2() ; 
130   Add2ESDsList(h0, 0) ;
131
132   TH1F * h1  = new TH1F("hTOFESDsTime", "Time Spectrum in TOF (ns)", 2000, 0., 200) ; 
133   h1->Sumw2() ;
134   Add2ESDsList(h1, 1) ;
135
136   TH1F * h2  = new TH1F("hTOFESDsRawTime", "raw Time Spectrum in TOF (ns)", 2000, 0., 200) ; 
137   h2->Sumw2() ;
138   Add2ESDsList(h2, 2) ;
139
140   TH1F * h3  = new TH1F("hTOFESDsToT", "ToT Spectrum in TOF (ns)", 500, 0., 50) ; 
141   h3->Sumw2() ;
142   Add2ESDsList(h3, 3) ;
143
144   TH1F * h4 = new TH1F("hTOFESDsPID",    "Fraction of matched TOF tracks with good PID glag", 100, 0., 1.) ;  
145   h4->Sumw2() ; 
146   Add2ESDsList(h4, 4) ;
147 }
148
149 //____________________________________________________________________________
150 void AliTOFQADataMakerRec::MakeRaws(AliRawReader* rawReader)
151 {
152   //
153   // makes data from Raws
154   //
155
156   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
157   Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
158
159
160   Int_t ntof = 0 ; 
161   Int_t in[5];
162   Int_t out[5];
163
164   TClonesArray * clonesRawData;
165   AliTOFRawStream tofInput(rawReader);
166   for (Int_t iDDL = 0; iDDL < AliTOFGeometry::NDDL()*AliTOFGeometry::NSectors(); iDDL++){
167     rawReader->Reset();
168     tofInput.LoadRawData(iDDL);
169     clonesRawData = (TClonesArray*)tofInput.GetRawData();
170     for (Int_t iRawData = 0; iRawData<clonesRawData->GetEntriesFast(); iRawData++) {
171       AliTOFrawData *tofRawDatum = (AliTOFrawData*)clonesRawData->UncheckedAt(iRawData);
172       if (!tofRawDatum->GetTOT() || !tofRawDatum->GetTOF()) continue;
173       ntof++;
174       GetRawsData(1)->Fill( tofRawDatum->GetTOF()*tdc2ns) ;//in ns
175       GetRawsData(2)->Fill( tofRawDatum->GetTOT()*tot2ns) ;//in ns
176
177       tofInput.EquipmentId2VolumeId(iDDL, 
178                                     tofRawDatum->GetTRM(), 
179                                     tofRawDatum->GetTRMchain(),
180                                     tofRawDatum->GetTDC(), 
181                                     tofRawDatum->GetTDCchannel(), 
182                                     in);
183     
184       GetMapIndeces(in,out);
185       GetRawsData(3)->Fill( out[0],out[1]) ;//raw map
186       
187     } // while loop
188     
189     clonesRawData->Clear();
190     
191   } // DDL Loop
192   
193   Int_t nentries=ntof;
194   if(nentries<=0){
195     GetRawsData(0)->Fill(-1.) ; 
196   }else{
197     GetRawsData(0)->Fill(TMath::Log10(nentries)) ; 
198   }
199 }
200
201 //____________________________________________________________________________
202 void AliTOFQADataMakerRec::MakeRecPoints(TTree * clustersTree)
203 {
204   //
205   // Make data from Clusters
206   //
207
208   Double_t tdc2ns=AliTOFGeometry::TdcBinWidth()*1E-3;
209   Double_t tot2ns=AliTOFGeometry::ToTBinWidth()*1E-3;
210
211   Int_t in[5];
212   Int_t out[5];
213
214   TBranch *branch=clustersTree->GetBranch("TOF");
215   if (!branch) { 
216     AliError("can't get the branch with the TOF clusters !");
217     return;
218   }
219
220   TClonesArray dummy("AliTOFcluster",10000), *clusters=&dummy;
221   branch->SetAddress(&clusters);
222
223   // Import the tree
224   clustersTree->GetEvent(0);  
225   
226   Int_t nentries=clusters->GetEntriesFast();
227   if(nentries<=0){
228     GetRecPointsData(0)->Fill(-1.) ; 
229   }else{
230     GetRecPointsData(0)->Fill(TMath::Log10(nentries)) ; 
231   } 
232  
233   TIter next(clusters) ; 
234   AliTOFcluster * c ; 
235   while ( (c = dynamic_cast<AliTOFcluster *>(next())) ) {
236     GetRecPointsData(1)->Fill(c->GetTDC()*tdc2ns);
237     GetRecPointsData(2)->Fill(c->GetTDCRAW()*tdc2ns);
238     GetRecPointsData(3)->Fill(c->GetToT()*tot2ns);
239     
240     in[0] = c->GetDetInd(0);
241     in[1] = c->GetDetInd(1);
242     in[2] = c->GetDetInd(2);
243     in[3] = c->GetDetInd(4); //X and Z indeces inverted in RecPoints
244     in[4] = c->GetDetInd(3); //X and Z indeces inverted in RecPoints
245     
246     GetMapIndeces(in,out);
247     GetRecPointsData(4)->Fill(out[0],out[1]);
248     
249   }
250 }
251
252 //____________________________________________________________________________
253 void AliTOFQADataMakerRec::MakeESDs(AliESDEvent * esd)
254 {
255   //
256   // make QA data from ESDs
257   //  
258   Int_t ntrk = esd->GetNumberOfTracks() ; 
259   Int_t ntof=0;
260   Int_t ntofpid=0;
261   while (ntrk--) {
262     AliESDtrack *track=esd->GetTrack(ntrk);
263     Double_t tofTime=track->GetTOFsignal()*1E-3;//in ns
264     Double_t tofTimeRaw=track->GetTOFsignalRaw()*1E-3;//in ns
265     Double_t tofToT=track->GetTOFsignalToT(); //in ns
266     if(!(tofTime>0))continue;
267     ntof++;
268     GetESDsData(1)->Fill(tofTime);
269     GetESDsData(2)->Fill(tofTimeRaw); 
270     GetESDsData(3)->Fill(tofToT);
271     //check how many tracks where ESD PID is ok 
272     UInt_t status=track->GetStatus();
273     if (((status&AliESDtrack::kESDpid)==0) || 
274         ((status&AliESDtrack::kTOFpid)==0)) continue;
275     ntofpid++;
276   }
277   
278   Int_t nentries=ntof;
279   if(nentries<=0){
280     GetESDsData(0)->Fill(-1.) ;
281   }else{
282     GetESDsData(0)->Fill(TMath::Log10(nentries)) ;
283   }
284
285   if(ntof>0)GetESDsData(4)->Fill(ntofpid/ntof) ;
286
287 }
288
289 //____________________________________________________________________________ 
290 void AliTOFQADataMakerRec::StartOfDetectorCycle()
291 {
292   //
293   //Detector specific actions at start of cycle
294   //to be implemented  
295 }
296
297 //____________________________________________________________________________ 
298 void AliTOFQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX task, TObjArray * list)
299 {
300   //Detector specific actions at end of cycle
301   // do the QA checking
302
303   AliQAChecker::Instance()->Run(AliQA::kTOF, task, list) ;  
304 }
305 //____________________________________________________________________________
306 void AliTOFQADataMakerRec::GetMapIndeces(Int_t* in , Int_t* out)
307 {
308   //
309   //return appropriate indeces for the theta-phi map
310   //
311
312   Int_t npadX = AliTOFGeometry::NpadX();
313   Int_t npadZ = AliTOFGeometry::NpadZ();
314   Int_t nStripA = AliTOFGeometry::NStripA();
315   Int_t nStripB = AliTOFGeometry::NStripB();
316   Int_t nStripC = AliTOFGeometry::NStripC();
317
318   Int_t isector = in[0];
319   Int_t iplate = in[1];
320   Int_t istrip = in[2];
321   Int_t ipadX = in[3]; 
322   Int_t ipadZ = in[4]; 
323   
324   Int_t stripOffset = 0;
325   switch (iplate) {
326   case 0:
327     stripOffset = 0;
328       break;
329   case 1:
330     stripOffset = nStripC;
331     break;
332   case 2:
333     stripOffset = nStripC+nStripB;
334     break;
335   case 3:
336     stripOffset = nStripC+nStripB+nStripA;
337     break;
338   case 4:
339     stripOffset = nStripC+nStripB+nStripA+nStripB;
340     break;
341   default:
342     AliError(Form("Wrong plate number in TOF (%d) !",iplate));
343     break;
344   };
345   Int_t zindex=npadZ*(istrip+stripOffset)+(ipadZ+1);
346   Int_t phiindex=npadX*isector+ipadX+1;
347   out[0]=zindex;  
348   out[1]=phiindex;  
349   
350 }