]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDQADataMakerRec.cxx
a3345868ebece4015cbf5ba086137de74d1154a5
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDQADataMakerRec.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 "AliHMPIDDigit.h"
36 #include "AliHMPIDCluster.h"
37 #include "AliHMPIDQADataMakerRec.h"
38 #include "AliHMPIDQAChecker.h"
39 #include "AliHMPIDParam.h"
40 #include "AliHMPIDRawStream.h"
41 #include "AliLog.h"
42
43 //.
44 // HMPID AliHMPIDQADataMakerRec base class
45 // for QA of reconstruction
46 // here also errors are calculated
47 //.
48
49 ClassImp(AliHMPIDQADataMakerRec)
50            
51 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
52   AliHMPIDQADataMakerRec::AliHMPIDQADataMakerRec() : 
53   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kHMPID), "HMPID Quality Assurance Data Maker"),fEvtRaw(0)
54 {
55   // ctor
56 }
57
58 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
59 AliHMPIDQADataMakerRec::AliHMPIDQADataMakerRec(const AliHMPIDQADataMakerRec& qadm) :
60   AliQADataMakerRec(),fEvtRaw(qadm.fEvtRaw)
61 {
62   //copy ctor 
63   SetName((const char*)qadm.GetName()) ; 
64   SetTitle((const char*)qadm.GetTitle()); 
65 }
66
67 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
68 AliHMPIDQADataMakerRec& AliHMPIDQADataMakerRec::operator = (const AliHMPIDQADataMakerRec& qadm )
69 {
70   // Equal operator.
71   this->~AliHMPIDQADataMakerRec();
72   new(this) AliHMPIDQADataMakerRec(qadm);
73   return *this;
74 }
75  
76 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
77 void AliHMPIDQADataMakerRec::InitDigits()
78 {
79   // create Digits histograms in Digits subdir
80   const Bool_t expert   = kTRUE ; 
81   const Bool_t image    = kTRUE ; 
82   
83   TH1F *hDigChEvt = new TH1F("hDigChEvt","Chamber occupancy per event;Occupancy [%];Counts",AliHMPIDParam::kMaxCh+1,AliHMPIDParam::kMinCh,AliHMPIDParam::kMaxCh+1);
84   TH1F *hDigPcEvt = new TH1F("hDigPcEvt","PC occupancy",156,-1,77);
85   TH2F *hDigMap[7];
86   TH1F *hDigQ[42];
87   for(Int_t iCh =0; iCh < 7; iCh++){
88     hDigMap[iCh] = new TH2F(Form("MapCh%i",iCh),Form("Digit Map in Chamber %i",iCh),159,0,159,143,0,143);
89     for(Int_t iPc =0; iPc < 6; iPc++ ){
90       hDigQ[iCh*6+iPc] = new TH1F(Form("QCh%iPc%i        ",iCh,iPc),Form("Charge of digits (ADC) in Chamber %i and PC %i;Charge [ADC counts];Counts",iCh,iPc),4100,0,4100);
91     }
92   }
93   
94   Add2DigitsList(hDigChEvt,0, !expert, image);
95   Add2DigitsList(hDigPcEvt,1,expert, !image);
96   for(Int_t iMap=0; iMap < 7; iMap++) Add2DigitsList(hDigMap[iMap],2+iMap,expert, !image);
97   for(Int_t iH =0; iH < 42 ; iH++) Add2DigitsList(hDigQ[iH]    ,9+iH,expert,!image);
98 }
99
100 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
101 void AliHMPIDQADataMakerRec::InitRecPoints()
102 {
103   // create cluster histograms in RecPoint subdir
104   const Bool_t expert   = kTRUE ; 
105   const Bool_t image    = kTRUE ; 
106  
107   TProfile *hCluMult          = new TProfile("CluMult"   ,"Cluster multiplicity per chamber;Chamber Id;# of clusters"    , 16, -1 , 7  , 0, 500);
108   Add2RecPointsList(hCluMult    , 0,expert, !image);
109   
110   TH2F *hCluFlg          = new TH2F("CluFlg"      ,"Cluster flag;??;??"                              ,  56  ,-1.5, 12.5, 70, -0.5, 6.5);
111   Add2RecPointsList(hCluFlg    , 1,expert, !image);
112   
113   TH1F *hCluSizeMip[7], *hCluSizePho[7];
114   
115   TH1F *hCluQSect[42], *hCluQSectZoom[42];
116   
117   for(Int_t iCh =0; iCh <7; iCh++){
118     hCluSizeMip[iCh] = new TH1F(Form("CluSizeMipCh%i",iCh),Form("Cluster size  MIP  (cluster Q > 100 ADC) in Chamber %i;Size [MIP];Counts",iCh),  50  , 0  , 50  );
119     Add2RecPointsList(hCluSizeMip[iCh], iCh+2,expert,!image);
120
121     hCluSizePho[iCh]  = new TH1F(Form("CluSizePho%i",iCh ),Form("Cluster size  Phots(cluster Q < 100 ADC) in Chamber %i;Size [MIP];Counts",iCh),  50  , 0  , 50  );
122     Add2RecPointsList(hCluSizePho[iCh], iCh+7+2,expert,!image);
123     
124     for(Int_t iSect =0; iSect < 6; iSect++){
125       hCluQSectZoom[iCh*6+iSect] = new TH1F(Form("QClusCh%iSect%iZoom",iCh,iSect) ,Form("Zoom on Cluster charge (ADC) in Chamber %i and sector %i;Charge [ADC counts];Counts",iCh,iSect),100,0,100);
126       Add2RecPointsList(hCluQSectZoom[iCh*6+iSect],2+14+iCh*6+iSect,expert,!image);
127       
128       hCluQSect[iCh*6+iSect] = new TH1F(Form("QClusCh%iSect%i",iCh,iSect) ,Form("Cluster charge (ADC) in Chamber %i and sector %i;Charge [ADC counts];Counts",iCh,iSect),250,0,5000);
129       Add2RecPointsList(hCluQSect[iCh*6+iSect],2+14+42+iCh*6+iSect, !expert, image);
130     }  
131   }
132 }
133 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
134 void AliHMPIDQADataMakerRec::InitRaws()
135 {
136 //
137 // Booking QA histo for Raw data
138 //
139 // All histograms implemented in InitRaws are used in AMORE. Any change here should be propagated to the amoreHMP-QA as well!!! (clm)
140 //  
141   
142   const Bool_t expert   = kTRUE ; 
143   const Bool_t saveCorr = kTRUE ; 
144   const Bool_t image    = kTRUE ; 
145   
146   const Int_t kNerr = (Int_t)AliHMPIDRawStream::kSumErr+1;
147   TH1F *hSumErr[14];
148   TH2F *hDilo[14];
149   TH2I *hPadMap[42]; //AMORE monitoring
150   TH1I *hPadQ[42]; //AMORE monitoring
151   for(Int_t iddl =0; iddl<AliHMPIDRawStream::kNDDL; iddl++) {
152     
153     hSumErr[iddl] = new TH1F(Form("hSumErrDDL%i",iddl), Form("Error summary for DDL %i;??;??",iddl), 2*kNerr,0,2*kNerr);
154     for(Int_t ilabel=0; ilabel< kNerr; ilabel++) {
155       hSumErr[iddl]->GetXaxis()->CenterLabels(kTRUE);
156       hSumErr[iddl]->GetXaxis()->SetBinLabel((2*ilabel+1),Form("%i  %s",ilabel+1,AliHMPIDRawStream::GetErrName(ilabel)));
157     }
158     
159     Add2RawsList(hSumErr[iddl],iddl,expert,!image, !saveCorr);
160     
161     hDilo[iddl] = new TH2F(Form("hDiloDDL%i",iddl),Form("Dilogic response at DDL;Row # ;Dilogic #",iddl),24,1,25,10,1,11);
162     Add2RawsList(hDilo[iddl],14+iddl,expert,!image, !saveCorr);
163   }//DDL loop
164   for(Int_t iCh = AliHMPIDParam::kMinCh; iCh <=AliHMPIDParam::kMaxCh ;iCh++) {
165     for(Int_t iPc = AliHMPIDParam::kMinPc; iPc <= AliHMPIDParam::kMaxPc ;iPc++) {
166       hPadMap[iPc+6*iCh] = new TH2I(Form("hPadMap_Ch_%i_Pc%i",iCh,iPc),Form("Pad Map of Ch: %i Pc: %i;Pad X;Pad Y;",iCh,iPc),80,0,80,48,0,48);
167       Add2RawsList(hPadMap[iPc+6*iCh],28+iPc+6*iCh,expert,!image, !saveCorr); 
168       hPadQ[iPc+6*iCh]   = new TH1I(Form("hPadQ_Ch_%i_Pc%i",iCh,iPc),Form("Pad Charge of Ch: %i Pc: %i;Pad Q;Entries;",iCh,iPc),4100,0,4100);
169       Add2RawsList(hPadQ[iPc+6*iCh],70+iPc+6*iCh,expert,!image, !saveCorr); 
170     }//PC loop
171   }//Ch loop  
172   
173   TH2I *hGeneralErrorSummary = new TH2I("GeneralErrorSummary"," DDL index vs Error type plot", 2*kNerr, 0, 2*kNerr, 2*AliHMPIDRawStream::kNDDL,0,2*AliHMPIDRawStream::kNDDL);
174   for(Int_t igenlabel =0 ; igenlabel< kNerr; igenlabel++) hGeneralErrorSummary->GetXaxis()->SetBinLabel((2*igenlabel+1),Form("%i  %s",igenlabel+1,AliHMPIDRawStream::GetErrName(igenlabel)));
175   Add2RawsList(hGeneralErrorSummary,14+14+42+42, !expert, image, !saveCorr);
176 }
177 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
178 void AliHMPIDQADataMakerRec::InitESDs()
179 {
180   //
181   //Booking ESDs histograms
182   const Bool_t expert   = kTRUE ; 
183   const Bool_t image    = kTRUE ; 
184
185   TH2F*  hCkovP  = new TH2F("CkovP" , "#theta_{c}, [rad];P, [GeV]"   , 150,      0,  7  ,100, 0, 1)   ;
186   TH2F*  hSigP   = new TH2F("SigP"  ,"#sigma_{#theta_c} [mrad];[GeV]", 150,      0,  7  ,100, 0, 1)   ;
187   TH2F*  hDifXY  = new TH2F("DifXY" ,"diff"                          , 200,    -10, 10  ,200,-10,10)  ;
188   TH2F*  hMvsP = new TH2F("MvsP","Reconstructed Mass vs P",60,0,6,1000,0,1) ;
189   TH1F*  hPid[5];
190   hPid[0] = new TH1F("PidE" ,"electron response"              , 101, -0.005,1.005)             ;
191   hPid[1] = new TH1F("PidMu","#mu response"                   , 101, -0.005,1.005)             ;
192   hPid[2] = new TH1F("PidPi","#pi response"                   , 101, -0.005,1.005)             ;
193   hPid[3] = new TH1F("PidK" ,"K response"                     , 101, -0.005,1.005)             ;
194   hPid[4] = new TH1F("PidP" ,"p response"                     , 101, -0.005,1.005)             ;
195   
196   Add2ESDsList(hCkovP,0, !expert, image);
197   Add2ESDsList(hSigP ,1, expert, !image);
198   Add2ESDsList(hDifXY,2, !expert, image);
199   Add2ESDsList(hMvsP,3, expert, !image);
200   for(Int_t i=0; i< 5; i++) Add2ESDsList(hPid[i],i+4, expert, !image);
201 }
202 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
203 void AliHMPIDQADataMakerRec::MakeRaws(AliRawReader *rawReader)
204 {
205 //
206 // Filling Raws QA histos
207 //
208           rawReader->Reset() ; 
209                 AliHMPIDRawStream stream(rawReader);
210
211     fEvtRaw++;
212
213     while(stream.Next())
214      {
215        UInt_t ddl=stream.GetDDLNumber(); //returns 0,1,2 ... 13
216        if(ddl > 13) continue;
217        for(Int_t iErr =1; iErr<(Int_t)AliHMPIDRawStream::kSumErr; iErr++){
218          Int_t numOfErr = stream.GetErrors(ddl,iErr);
219          GetRawsData(ddl)->Fill(iErr,numOfErr);
220          //Printf("err type %i   ddl number %i  Num of errors %i",iErr,ddl,numOfErr);
221         ((TH2I*)GetRawsData(14+14+42+42))->Fill(iErr,ddl,iErr); //
222        }
223        UInt_t word; Int_t Nddl, r, d, a;//pc,pcX,pcY;      
224        for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
225        AliHMPIDDigit dig(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);dig.Raw(word,Nddl,r,d,a);
226        GetRawsData(ddl+14)->Fill(r,d);
227        GetRawsData(28+stream.Pc(Nddl,r,d,a)+6*AliHMPIDParam::DDL2C(ddl))->Fill(stream.PadPcX(Nddl,r,d,a),stream.PadPcY(Nddl,r,d,a));
228        GetRawsData(70+stream.Pc(Nddl,r,d,a)+6*AliHMPIDParam::DDL2C(ddl))->Fill(stream.GetChargeArray()[iPad]);
229        }      
230      }
231
232    stream.Delete();
233    
234 }
235
236 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
237 void AliHMPIDQADataMakerRec::MakeDigits(TClonesArray * data)
238 {
239   //
240   //filling QA histos for Digits
241   //
242   
243   TObjArray *chamber = dynamic_cast<TObjArray*>(data);
244   if ( !chamber) {
245     AliError("Wrong type of digits container") ; 
246   } else {
247     for(Int_t i =0; i< chamber->GetEntries(); i++)
248       {
249       TClonesArray * digits = dynamic_cast<TClonesArray*>(chamber->At(i)); 
250       GetDigitsData(0)->Fill(i,digits->GetEntriesFast()/(48.*80.*6.));
251       TIter next(digits); 
252       AliHMPIDDigit * digit; 
253       while ( (digit = dynamic_cast<AliHMPIDDigit *>(next())) ) {
254         GetDigitsData(1)->Fill(10.*i+digit->Pc(),1./(48.*80.));
255         GetDigitsData(2+i)->Fill(digit->PadChX(),digit->PadChY());
256         GetDigitsData(9+i*6+digit->Pc())->Fill(digit->Q());
257       }  
258       }
259   }
260 }
261
262 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
263 void AliHMPIDQADataMakerRec::MakeDigits(TTree * data)
264 {
265   //
266   //Opening the Digit Tree
267   //
268   TObjArray *pObjDig=new TObjArray(AliHMPIDParam::kMaxCh+1);
269   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
270     TClonesArray *pCA=new TClonesArray("AliHMPIDDigit");
271     pObjDig->AddAt(pCA,iCh);
272   }
273   
274   pObjDig->SetOwner(kTRUE);
275   
276   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
277     data->SetBranchAddress(Form("HMPID%i",iCh),&(*pObjDig)[iCh]);
278   }
279   data->GetEntry(0);
280   
281   MakeDigits((TClonesArray *)pObjDig);
282 }
283
284 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
285 void AliHMPIDQADataMakerRec::MakeRecPoints(TTree * clustersTree)
286 {
287   //
288   //filling QA histos for clusters
289   //
290   AliHMPIDParam *pPar =AliHMPIDParam::Instance();
291   
292   static TClonesArray *clusters;
293   if(!clusters) clusters = new TClonesArray("AliHMPIDCluster");
294   
295   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
296     TBranch *branch = clustersTree->GetBranch(Form("HMPID%d",iCh));
297     branch->SetAddress(&clusters);
298     branch->GetEntry(0);
299     GetRecPointsData(0)->Fill(iCh,clusters->GetEntries());
300     TIter next(clusters);
301     AliHMPIDCluster *clu;
302     while ( (clu = dynamic_cast<AliHMPIDCluster *>(next())) ) {
303       GetRecPointsData(1)->Fill(clu->Status(),iCh);
304       Int_t sect =  pPar->InHVSector(clu->Y());
305       if(clu->Q()>100) GetRecPointsData(2+iCh)->Fill(clu->Size());
306       else {
307         GetRecPointsData(2+7+iCh)->Fill(clu->Size());
308         GetRecPointsData(2+14+iCh*6+sect)->Fill(clu->Q());
309       }    
310       GetRecPointsData(2+14+42+iCh*6+sect)->Fill(clu->Q());
311     }
312   }
313
314   clusters->Clear();
315   
316 }
317 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
318 void AliHMPIDQADataMakerRec::MakeESDs(AliESDEvent * esd)
319 {
320   //
321   //fills QA histos for ESD
322   //
323
324   for(Int_t iTrk = 0 ; iTrk < esd->GetNumberOfTracks() ; iTrk++){
325     AliESDtrack *pTrk = esd->GetTrack(iTrk) ;
326     GetESDsData(0)->Fill(pTrk->GetP(),pTrk->GetHMPIDsignal());
327     GetESDsData(1)->Fill( pTrk->GetP(),TMath::Sqrt(pTrk->GetHMPIDchi2()));
328     Float_t xm,ym; Int_t q,np;  
329     pTrk->GetHMPIDmip(xm,ym,q,np);                       //mip info
330     Float_t xRad,yRad,th,ph;        
331     pTrk->GetHMPIDtrk(xRad,yRad,th,ph);              //track info at the middle of the radiator
332     Float_t xPc = xRad+9.25*TMath::Tan(th)*TMath::Cos(ph); // temporar: linear extrapol (B=0!)
333     Float_t yPc = yRad+9.25*TMath::Tan(th)*TMath::Sin(ph); // temporar:          "
334     GetESDsData(2)->Fill(xm-xPc,ym-yPc); //track info
335     if(pTrk->GetHMPIDsignal()>0) {
336      Double_t a = 1.292*1.292*TMath::Cos(pTrk->GetHMPIDsignal())*TMath::Cos(pTrk->GetHMPIDsignal())-1.;
337      if(a > 0) {
338     Double_t mass = pTrk->P()*TMath::Sqrt(1.292*1.292*TMath::Cos(pTrk->GetHMPIDsignal())*TMath::Cos(pTrk->GetHMPIDsignal())-1.);
339     GetESDsData(3)->Fill( pTrk->GetP(),mass);
340      }
341     }
342    Double_t pid[5] ;      pTrk->GetHMPIDpid(pid) ;
343     for(Int_t i = 0 ; i < 5 ; i++) GetESDsData(4+i)->Fill(pid[i]) ;
344   }
345 }
346 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
347 void AliHMPIDQADataMakerRec::StartOfDetectorCycle()
348 {
349   //Detector specific actions at start of cycle
350   
351 }
352 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
353
354 void AliHMPIDQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray **histos)
355 {
356   //Detector specific actions at end of cycle
357   // do the QA checking
358   
359   if(task==AliQAv1::kRAWS) {
360     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
361       for(Int_t iddl=0;iddl<14;iddl++) {
362         TH1F *h = (TH1F*)histos[specie]->At(14+iddl); //ddl histos scaled by the number of events 
363         h->Scale(1./(Float_t)fEvtRaw);
364       }
365     }
366   }
367   
368    AliQAChecker::Instance()->Run(AliQAv1::kHMPID, task, histos);
369
370 }
371