]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDQADataMaker.cxx
roll back of previous commit
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDQADataMaker.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 /* $Id$ */
18
19 // --- ROOT system ---
20 #include <TClonesArray.h>
21 #include <TFile.h> 
22 #include <TH1F.h> 
23 #include <TH2F.h>
24 #include <TProfile.h>
25 #include <Riostream.h>
26 // --- Standard library ---
27
28 // --- AliRoot header files ---
29 #include "AliESDCaloCluster.h"
30 #include "AliESDEvent.h"
31 #include "AliQAChecker.h"
32 #include "AliLog.h"
33 #include "AliHMPIDDigit.h"
34 #include "AliHMPIDHit.h"
35 #include "AliHMPIDCluster.h"
36 #include "AliHMPIDQADataMaker.h"
37 #include "AliHMPIDParam.h"
38 #include "AliHMPIDRawStream.h"
39 #include "AliLog.h"
40 ClassImp(AliHMPIDQADataMaker)
41            
42 //____________________________________________________________________________ 
43   AliHMPIDQADataMaker::AliHMPIDQADataMaker() : 
44   AliQADataMaker(AliQAv1::GetDetName(AliQAv1::kHMPID), "HMPID Quality Assurance Data Maker")
45 {
46   // ctor
47 }
48
49 //____________________________________________________________________________ 
50 AliHMPIDQADataMaker::AliHMPIDQADataMaker(const AliHMPIDQADataMaker& qadm) :
51   AliQADataMaker() 
52 {
53   //copy ctor 
54   SetName((const char*)qadm.GetName()) ; 
55   SetTitle((const char*)qadm.GetTitle()); 
56 }
57
58 //__________________________________________________________________
59 AliHMPIDQADataMaker& AliHMPIDQADataMaker::operator = (const AliHMPIDQADataMaker& qadm )
60 {
61   // Equal operator.
62   this->~AliHMPIDQADataMaker();
63   new(this) AliHMPIDQADataMaker(qadm);
64   return *this;
65 }
66  
67 //____________________________________________________________________________ 
68 void AliHMPIDQADataMaker::InitHits()
69 {
70   // create Hits histograms in Hits subdir
71   TH1F *hHitQdc=new TH1F("HitQdc","HMPID Hit Qdc all chamber;QDC",500,0,4000);
72   Add2HitsList(hHitQdc,0);
73   TH2F *hHitMap[7];
74   for(Int_t iCh=0;iCh<7;iCh++) {
75     hHitMap[iCh]=new TH2F(Form("HMPID HitMap%i",iCh),Form("Ch%i;x_{Hit};y_{Hit}",iCh),162,-1,161,146,-1,145);   
76     Add2HitsList(hHitMap[iCh],iCh+1);
77   }
78   //
79   ClonePerTrigClass(AliQAv1::kHITS); // this should be the last line  
80 }
81
82 //____________________________________________________________________________ 
83 void AliHMPIDQADataMaker::InitDigits()
84 {
85   // create Digits histograms in Digits subdir
86   TH1F *hDigPcEvt = new TH1F("hDigPcEvt","PC occupancy",156,-1,77);
87   TH1F *hDigQ     = new TH1F("Q        ","Charge of digits (ADC)     ",3000,0,3000);
88   TH1F *hDigChEvt = new TH1F("hDigChEvt","Chamber occupancy per event",AliHMPIDParam::kMaxCh+1,AliHMPIDParam::kMinCh,AliHMPIDParam::kMaxCh+1);
89   
90   TProfile *tDigHighQ = new TProfile("tDigHighQ","Highest charge in chamber  ",AliHMPIDParam::kMaxCh+1,AliHMPIDParam::kMinCh,AliHMPIDParam::kMaxCh+1);
91   TProfile *tDigChEvt = new TProfile("tDigChEvt","Chamber occupancy per event (profile)",AliHMPIDParam::kMaxCh+1,AliHMPIDParam::kMinCh,AliHMPIDParam::kMaxCh+1);
92   
93   Add2DigitsList(hDigPcEvt,0);
94   Add2DigitsList(hDigQ    ,1);
95   Add2DigitsList(hDigChEvt,2);
96   Add2DigitsList(tDigHighQ,3);
97   Add2DigitsList(tDigChEvt,4);
98   //
99   ClonePerTrigClass(AliQAv1::kDIGITS); // this should be the last line
100 }
101
102 //____________________________________________________________________________ 
103 void AliHMPIDQADataMaker::InitSDigits()
104 {
105   // create SDigits histograms in SDigits subdir
106   TH1F   *hSDigits     = new TH1F("hHmpidSDigits",    "SDigits Q  distribution in HMPID",  500, 0., 5000.) ; 
107   Add2SDigitsList(hSDigits,0);
108   //
109   ClonePerTrigClass(AliQAv1::kSDIGITS); // this should be the last line
110 }
111
112 //____________________________________________________________________________ 
113
114 void AliHMPIDQADataMaker::InitRecPoints()
115 {
116   // create cluster histograms in RecPoint subdir
117
118   TH1F *hCluEvt=new TH1F("CluPerEvt","Cluster multiplicity"   ,100,0,100);
119   TH1F *hCluChi2  =new TH1F("CluChi2"  ,"Chi2 "               ,1000,0,100);
120   TH1F *hCluFlg   =new TH1F("CluFlg"   ,"Cluster flag"        ,14,-1.5,12.5); hCluFlg->SetFillColor(5);
121   TH1F *hCluSize  =new TH1F("CluSize"  ,"Cluster size        ",100,0,100);
122   TH1F *hCluQ     =new TH1F("CluQ"     ,"Cluster charge (ADC)",1000,0,5000);
123
124   Add2RecPointsList(hCluEvt , 0);
125   Add2RecPointsList(hCluChi2, 1);
126   Add2RecPointsList(hCluFlg , 2);
127   Add2RecPointsList(hCluSize, 3);
128   Add2RecPointsList(hCluQ   , 4);
129   //
130   ClonePerTrigClass(AliQAv1::kRECPOINTS); // this should be the last line
131 }
132 //____________________________________________________________________________
133
134 void AliHMPIDQADataMaker::InitRaws()
135 {
136   //
137   // Booking QA histo for Raw data
138   //
139   TH1F *hqPad[14];
140   for(Int_t iddl =0; iddl<14; iddl++) {
141   hqPad[iddl] = new TH1F(Form("hqPadDDL%i",iddl), Form("Pad Q Entries at DDL %i",iddl), 500,0,5000);
142   Add2RawsList(hqPad[iddl],iddl);
143   }
144
145   const Int_t nerr = (Int_t)AliHMPIDRawStream::kSumErr+1;
146   const char *hnames[nerr]={"RawDataSize","RawMarkerSize","WrongRow","WrongDilogic","WrongPad","EoEFlag",
147                              "EoESize","EoEDILOGIC","EoERow","BadSegWord","WrongSeg","RowMarkerSize","NoErrors","Invalid"};
148
149   TH1F *hSumErr = new TH1F("SumErr","Summary of the returned errors",2*nerr,0,nerr);
150
151   for(Int_t ilabel=0; ilabel< nerr; ilabel++) {
152   hSumErr->GetXaxis()->CenterLabels(kTRUE);
153   hSumErr->GetXaxis()->SetBinLabel((2*ilabel+1),Form("%i  %s",ilabel+1,hnames[ilabel]));
154
155   }
156   Add2RawsList(hSumErr,14);
157   //
158   ClonePerTrigClass(AliQAv1::kRAWS); // this should be the last line
159 }
160
161 //____________________________________________________________________________
162 void AliHMPIDQADataMaker::InitESDs()
163 {
164   //
165   //Booking ESDs histograms
166   TH2F*  hCkovP  = new TH2F("CkovP" , "#theta_{c}, [rad];P, [GeV]"   , 150,      0,  7  ,100, 0, 1)   ;
167   TH2F*  hSigP   = new TH2F("SigP"  ,"#sigma_{#theta_c} [mrad];[GeV]", 150,      0,  7  ,100, 0, 1)   ;
168   TH2F*  hMipXY  = new TH2F("MipXY" ,"mip position"                  , 260,      0,130  ,252, 0,126)  ;
169   TH2F*  hDifXY  = new TH2F("DifXY" ,"diff"                          , 200,    -10, 10  ,200,-10,10)  ;
170   TH1F*  hPid[5];
171   hPid[0] = new TH1F("PidE" ,"electron response"              , 101, -0.005,1.005)             ;
172   hPid[1] = new TH1F("PidMu","#mu response"                   , 101, -0.005,1.005)             ;
173   hPid[2] = new TH1F("PidPi","#pi response"                   , 101, -0.005,1.005)             ;
174   hPid[3] = new TH1F("PidK" ,"K response"                     , 101, -0.005,1.005)             ;
175   hPid[4] = new TH1F("PidP" ,"p response"                         ,101, -0.005,1.005)             ;
176   
177   Add2ESDsList(hCkovP,0);
178   Add2ESDsList(hSigP ,1);
179   Add2ESDsList(hMipXY,2);
180   Add2ESDsList(hDifXY,3);
181   for(Int_t i=0; i< 5; i++) Add2ESDsList(hPid[i],i+4);
182   //
183   ClonePerTrigClass(AliQAv1::kESDS); // this should be the last line
184 }
185 //____________________________________________________________________________
186
187 void AliHMPIDQADataMaker::MakeHits(TClonesArray * data)
188 {
189  //
190  //filling QA histos for Hits
191  //
192   TClonesArray * hits = dynamic_cast<TClonesArray *>(data) ; 
193   if (!hits){
194     AliError("Wrong type of hits container") ; 
195   } else {
196     TIter next(hits); 
197     AliHMPIDHit * hit ; 
198     while ( (hit = dynamic_cast<AliHMPIDHit *>(next())) ) {
199       if(hit->Pid()<500000) FillHitsData(0,hit->Q()) ;
200       if(hit->Pid()<500000) FillHitsData(hit->Ch()+1,hit->LorsX(),hit->LorsY());
201     }
202   } 
203   //
204 }
205 //___________________________________________________________________________
206 void AliHMPIDQADataMaker::MakeHits(TTree * data)
207 {
208 //
209 //Opening of the Hit TTree 
210 //
211  TClonesArray *pHits=new TClonesArray("AliHMPIDHit");  data->SetBranchAddress("HMPID",&pHits);
212   for(Int_t iEnt=0;iEnt<data->GetEntriesFast();iEnt++){//entries loop
213     data->GetEntry(iEnt);
214     MakeHits(pHits);
215   }//entries loop
216   //
217   IncEvCountCycleHits();
218   IncEvCountTotalHits();
219   //
220 }
221
222 //____________________________________________________________________________
223 void AliHMPIDQADataMaker::MakeDigits(TClonesArray * data)
224 {
225  //
226  //filling QA histos for Digits
227  //
228   TObjArray *chamber = dynamic_cast<TObjArray*>(data);
229   if ( !chamber) {
230     AliError("Wrong type of digits container") ; 
231   } else {
232     for(Int_t i =0; i< chamber->GetEntries(); i++)
233       {
234         TClonesArray * digits = dynamic_cast<TClonesArray*>(chamber->At(i)); 
235         FillDigitsData(2,i,digits->GetEntriesFast()/(48.*80.*6.));
236         FillDigitsData(4,i,digits->GetEntriesFast()/(48.*80.*6.));
237         Double_t highQ=0;
238         TIter next(digits); 
239         AliHMPIDDigit * digit; 
240         while ( (digit = dynamic_cast<AliHMPIDDigit *>(next())) ) {
241           FillDigitsData(0,10.*i+digit->Pc(),1./(48.*80.));
242           FillDigitsData(1,digit->Q());
243           if(digit->Q()>highQ) highQ = digit->Q();
244         }  
245         FillDigitsData(3,i,highQ);      
246       }
247   }
248 }
249 //___________________________________________________________________________
250 void AliHMPIDQADataMaker::MakeDigits(TTree * data)
251 {
252 //
253 //Opening the Digit Tree
254 //
255  TObjArray *pObjDig=new TObjArray(AliHMPIDParam::kMaxCh+1);
256   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
257     TClonesArray *pCA=new TClonesArray("AliHMPIDDigit");
258     pObjDig->AddAt(pCA,iCh);
259   }
260
261   pObjDig->SetOwner(kTRUE);
262
263   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
264     data->SetBranchAddress(Form("HMPID%i",iCh),&(*pObjDig)[iCh]);
265   }
266   data->GetEntry(0);
267
268   MakeDigits((TClonesArray *)pObjDig);
269   //
270   IncEvCountCycleDigits();
271   IncEvCountTotalDigits();
272   //
273 }
274
275 //____________________________________________________________________________
276
277 void AliHMPIDQADataMaker::MakeRaws(AliRawReader *rawReader)
278 {
279 //
280 // Filling Raws QA histos
281 //
282   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++) {
283     AliHMPIDRawStream stream(rawReader);
284     while(stream.Next())
285       {
286       UInt_t ddl=stream.GetDDLNumber(); //returns 0,1,2 ... 13
287       if((UInt_t)(2*iCh)==ddl || (UInt_t)(2*iCh+1)==ddl) {
288        for(Int_t row = 1; row <= AliHMPIDRawStream::kNRows; row++){
289         for(Int_t dil = 1; dil <= AliHMPIDRawStream::kNDILOGICAdd; dil++){
290           for(Int_t pad = 0; pad < AliHMPIDRawStream::kNPadAdd; pad++){
291             if(stream.GetCharge(ddl,row,dil,pad) < 1) continue;
292               FillRawsData(ddl,stream.GetCharge(ddl,row,dil,pad));
293 //              Printf("charge %i",stream.GetCharge(ddl,row,dil,pad));
294             }//pad
295           }//dil
296         }//row
297       }//while
298     }
299     for(Int_t iErr =1; iErr<(Int_t)AliHMPIDRawStream::kSumErr; iErr++){
300       Int_t errflag = stream.GetErrors(iErr);
301
302        if(errflag < 0) FillRawsData(14,(Int_t)AliHMPIDRawStream::kSumErr+0.5);
303        else if(errflag == 0) FillRawsData(14,(Int_t)AliHMPIDRawStream::kSumErr-0.5);
304        else FillRawsData(14,iErr-0.5);
305      }
306     stream.Delete();
307   }//chamber loop
308   //
309   IncEvCountCycleRaws();
310   IncEvCountTotalRaws();
311   //
312 }
313
314 //___________________________________________________________________________
315
316 void AliHMPIDQADataMaker::MakeSDigits(TClonesArray * data)
317 {
318  //
319  //filling QA histos for SDigits
320  //
321   TClonesArray * sdigits = dynamic_cast<TClonesArray *>(data) ; 
322   if (!sdigits) {
323     AliError("Wrong type of sdigits container") ; 
324   } else {
325     TIter next(sdigits) ; 
326     AliHMPIDDigit * sdigit ; 
327     while ( (sdigit = dynamic_cast<AliHMPIDDigit *>(next())) ) {
328             FillSDigitsData(0,sdigit->Q());
329     } 
330   }
331 }
332 //___________________________________________________________________________
333 void AliHMPIDQADataMaker::MakeSDigits(TTree * data)
334 {
335   //
336   // Opening the SDigit Tree
337   //
338   TClonesArray * sdigits = new TClonesArray("AliHMPIDDigit", 1000) ;
339
340   TBranch * branch = data->GetBranch("HMPID") ;
341   if ( ! branch ) {
342     AliError("HMPID SDigit Tree not found") ;
343     return;
344   }
345   branch->SetAddress(&sdigits) ;
346   branch->GetEntry(0) ;
347   MakeSDigits(sdigits) ;
348   //
349   IncEvCountCycleSDigits();
350   IncEvCountTotalSDigits();
351   //
352 }
353 //____________________________________________________________________________
354 void AliHMPIDQADataMaker::MakeRecPoints(TTree * clustersTree)
355 {
356   //
357   //filling QA histos for clusters
358   //
359   TClonesArray *clusters = new TClonesArray("AliHMPIDCluster");
360   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){
361     TBranch *branch = clustersTree->GetBranch(Form("HMPID%d",i));
362     branch->SetAddress(&clusters);
363     branch->GetEntry(0);
364
365     FillRecPointsData(0,i,clusters->GetEntries());
366     TIter next(clusters);
367     AliHMPIDCluster *clu;
368     while ( (clu = dynamic_cast<AliHMPIDCluster *>(next())) ) {
369       FillRecPointsData(1,clu->Chi2());
370       FillRecPointsData(2,clu->Status());
371       FillRecPointsData(3,clu->Size());
372       FillRecPointsData(4,clu->Q()); 
373     }
374   }
375
376   clusters->Delete();
377   delete clusters;
378   IncEvCountCycleRecPoints();
379   IncEvCountTotalRecPoints();
380   //
381 }
382
383 //____________________________________________________________________________
384 void AliHMPIDQADataMaker::MakeESDs(AliESDEvent * esd)
385 {
386   //
387   //fills QA histos for ESD
388   //
389   for(Int_t iTrk = 0 ; iTrk < esd->GetNumberOfTracks() ; iTrk++){
390     AliESDtrack *pTrk = esd->GetTrack(iTrk) ;
391     Float_t thetaCkov = -999.;
392     if(pTrk->GetHMPIDsignal()<0.) thetaCkov = pTrk->GetHMPIDsignal();
393     else                          thetaCkov = pTrk->GetHMPIDsignal() - (Int_t)pTrk->GetHMPIDsignal();;
394     FillESDsData(0,pTrk->GetP(),thetaCkov);
395     FillESDsData(1, pTrk->GetP(),TMath::Sqrt(pTrk->GetHMPIDchi2()));
396     Float_t xm,ym; Int_t q,np;  
397     pTrk->GetHMPIDmip(xm,ym,q,np);                       //mip info
398     FillESDsData(2,xm,ym);
399     Float_t xRad,yRad,th,ph;        
400     pTrk->GetHMPIDtrk(xRad,yRad,th,ph);              //track info at the middle of the radiator
401     Float_t xPc = xRad+9.25*TMath::Tan(th)*TMath::Cos(ph); // temporar: linear extrapol (B=0!)
402     Float_t yPc = yRad+9.25*TMath::Tan(th)*TMath::Sin(ph); // temporar:          "
403     FillESDsData(3,xm-xPc,ym-yPc); //track info
404     Double_t pid[5] ;      pTrk->GetHMPIDpid(pid) ;
405     for(Int_t i = 0 ; i < 5 ; i++) FillESDsData(4+i,pid[i]) ;
406   }
407   IncEvCountCycleESDs();
408   IncEvCountTotalESDs();
409   //
410 }
411 //____________________________________________________________________________
412 void AliHMPIDQADataMaker::StartOfDetectorCycle()
413 {
414   //Detector specific actions at start of cycle
415   
416 }
417
418 void AliHMPIDQADataMaker::EndOfDetectorCycle(AliQAv1::TASKINDEX task, TObjArray * obj)
419 {
420   //Detector specific actions at end of cycle
421   // do the QA checking
422 //  AliQAChecker::Instance()->Run(AliQAv1::kHMPID, task, obj) ;  
423 }
424