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