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