]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDQADataMakerRec.cxx
b31c627289d7395f7032855ec0a75a4cf34d740b
[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 <TLine.h> 
23 #include <TH1F.h> 
24 #include <TH2F.h>
25 #include <TProfile.h>
26 #include <TProfile2D.h>
27 #include <Riostream.h>
28 // --- Standard library ---
29
30 // --- AliRoot header files ---
31 #include "AliESDCaloCluster.h"
32 #include "AliESDEvent.h"
33 #include "AliQAChecker.h"
34 #include "AliLog.h"
35 #include "AliHMPIDDigit.h"
36 #include "AliHMPIDHit.h"
37 #include "AliHMPIDDigit.h"
38 #include "AliHMPIDCluster.h"
39 #include "AliHMPIDQADataMakerRec.h"
40 #include "AliHMPIDQAChecker.h"
41 #include "AliHMPIDParam.h"
42 #include "AliHMPIDRawStream.h"
43 #include "AliLog.h"
44
45 //.
46 // HMPID AliHMPIDQADataMakerRec base class
47 // for QA of reconstruction
48 // here also errors are calculated
49 //.
50
51 ClassImp(AliHMPIDQADataMakerRec)
52            
53 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
54   AliHMPIDQADataMakerRec::AliHMPIDQADataMakerRec() : 
55   AliQADataMakerRec(AliQAv1::GetDetName(AliQAv1::kHMPID), "HMPID Quality Assurance Data Maker"),fEvtRaw(0),  fChannel(0)
56 {
57   // ctor
58 }
59
60 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
61 AliHMPIDQADataMakerRec::AliHMPIDQADataMakerRec(const AliHMPIDQADataMakerRec& qadm) :
62   AliQADataMakerRec(),fEvtRaw(qadm.fEvtRaw), fChannel(qadm.fChannel)
63 {
64   //copy ctor 
65   SetName((const char*)qadm.GetName()) ; 
66   SetTitle((const char*)qadm.GetTitle()); 
67 }
68
69 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
70 AliHMPIDQADataMakerRec& AliHMPIDQADataMakerRec::operator = (const AliHMPIDQADataMakerRec& qadm )
71 {
72   // Equal operator.
73   this->~AliHMPIDQADataMakerRec();
74   new(this) AliHMPIDQADataMakerRec(qadm);
75   return *this;
76 }
77  
78 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
79 void AliHMPIDQADataMakerRec::InitDigits()
80 {
81   // create Digits histograms in Digits subdir
82   const Bool_t expert   = kTRUE ; 
83   const Bool_t image    = kTRUE ; 
84   
85   TH1F *hDigChEvt = new TH1F("hDigChEvt","Chamber occupancy per event;Occupancy [%];Counts",AliHMPIDParam::kMaxCh+1,AliHMPIDParam::kMinCh,AliHMPIDParam::kMaxCh+1);
86   TH1F *hDigPcEvt = new TH1F("hDigPcEvt","PC occupancy",156,-1,77);
87   TH2F *hDigMap[7];
88   TH1F *hDigQ[42];
89   for(Int_t iCh =0; iCh < 7; iCh++){
90     hDigMap[iCh] = new TH2F(Form("MapCh%i",iCh),Form("Digit Map in Chamber %i",iCh),159,0,159,143,0,143);
91     for(Int_t iPc =0; iPc < 6; iPc++ ){
92       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);
93     }
94   }
95   
96   Add2DigitsList(hDigChEvt,0, !expert, image);
97   Add2DigitsList(hDigPcEvt,1,expert, !image);
98   for(Int_t iMap=0; iMap < 7; iMap++) Add2DigitsList(hDigMap[iMap],2+iMap,expert, !image);
99   for(Int_t iH =0; iH < 42 ; iH++) Add2DigitsList(hDigQ[iH]    ,9+iH,expert,!image);
100 }
101
102 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
103 void AliHMPIDQADataMakerRec::InitRecPoints()
104 {
105   // create cluster histograms in RecPoint subdir
106   const Bool_t expert   = kTRUE ; 
107   const Bool_t image    = kTRUE ; 
108  
109   TProfile *hCluMult          = new TProfile("CluMult"   ,"Cluster multiplicity per chamber;Chamber Id;# of clusters"    , 16, -1 , 7  , 0, 500);
110   Add2RecPointsList(hCluMult    , 0,expert, !image);
111   
112   TH2F *hCluFlg          = new TH2F("CluFlg"      ,"Cluster flag;??;??"                              ,  56  ,-1.5, 12.5, 70, -0.5, 6.5);
113   Add2RecPointsList(hCluFlg    , 1,expert, !image);
114   
115   TH1F *hCluSizeMip[7], *hCluSizePho[7];
116   
117   TH1F *hCluQSect[42], *hCluQSectZoom[42];
118   
119   for(Int_t iCh =0; iCh <7; iCh++){
120     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  );
121     Add2RecPointsList(hCluSizeMip[iCh], iCh+2,expert,!image);
122
123     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  );
124     Add2RecPointsList(hCluSizePho[iCh], iCh+7+2,expert,!image);
125     
126     for(Int_t iSect =0; iSect < 6; iSect++){
127       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);
128       Add2RecPointsList(hCluQSectZoom[iCh*6+iSect],2+14+iCh*6+iSect,expert,!image);
129       
130       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);
131       Add2RecPointsList(hCluQSect[iCh*6+iSect],2+14+42+iCh*6+iSect, !expert, image);
132     }  
133   }
134 }
135 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
136 void AliHMPIDQADataMakerRec::InitRaws()
137 {
138 //
139 // Booking QA histo for Raw data
140 //
141 // All histograms implemented in InitRaws are used in AMORE. Any change here should be propagated to the amoreHMP-QA as well!!! (clm)
142 //  
143   
144   const Bool_t expert   = kTRUE ; 
145   const Bool_t saveCorr = kTRUE ; 
146   const Bool_t image    = kTRUE ; 
147   
148   const Int_t kNerr = (Int_t)AliHMPIDRawStream::kSumErr+1;
149   TH1F *hSumErr[14];
150   TH2F *hDilo[14];
151   TH2I *hPadMap[42]; //AMORE monitoring
152   TH1I *hPadQ[42]; //AMORE monitoring
153    
154   
155   for(Int_t iddl =0; iddl<AliHMPIDRawStream::kNDDL; iddl++) {
156     
157     hSumErr[iddl] = new TH1F(Form("hSumErrDDL%i",iddl), Form("Error summary for DDL %i;??;??",iddl), 2*kNerr,0,2*kNerr);
158     for(Int_t ilabel=0; ilabel< kNerr; ilabel++) {
159       hSumErr[iddl]->GetXaxis()->CenterLabels(kTRUE);
160       hSumErr[iddl]->GetXaxis()->SetBinLabel((2*ilabel+1),Form("%i  %s",ilabel+1,AliHMPIDRawStream::GetErrName(ilabel)));
161     }
162     
163    Add2RawsList(hSumErr[iddl],iddl,expert,!image, !saveCorr);
164     
165     hDilo[iddl] = new TH2F(Form("hDiloDDL%i",iddl),Form("Dilogic response at DDL;Row # ;Dilogic #",iddl),24,1,25,10,1,11);
166     Add2RawsList(hDilo[iddl],14+iddl,expert,!image, !saveCorr);
167   }//DDL loop
168   for(Int_t iCh = AliHMPIDParam::kMinCh; iCh <=AliHMPIDParam::kMaxCh ;iCh++) {
169     for(Int_t iPc = AliHMPIDParam::kMinPc; iPc <= AliHMPIDParam::kMaxPc ;iPc++) {
170       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);
171      Add2RawsList(hPadMap[iPc+6*iCh],28+iPc+6*iCh,expert,!image, !saveCorr); 
172       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);
173       Add2RawsList(hPadQ[iPc+6*iCh],70+iPc+6*iCh,expert,!image, !saveCorr); 
174     }//PC loop
175   }//Ch loop  
176   
177   TH2I *hGeneralErrorSummary = new TH2I("GeneralErrorSummary"," DDL index vs Error type plot", 2*kNerr, 0, 2*kNerr, 2*AliHMPIDRawStream::kNDDL,0,2*AliHMPIDRawStream::kNDDL);
178   for(Int_t igenlabel =0 ; igenlabel< kNerr; igenlabel++) hGeneralErrorSummary->GetXaxis()->SetBinLabel((2*igenlabel+1),Form("%i  %s",igenlabel+1,AliHMPIDRawStream::GetErrName(igenlabel)));
179   Add2RawsList(hGeneralErrorSummary,14+14+42+42, expert, !image, !saveCorr);
180   
181   //___ for DQM shifter and eLogBook ___ start
182   //___ Double booking of histograms since TProfile cannot be display in summary image
183   //___ hence TProfile plots will not be shown in QA and LogBook!
184  
185   TH1F* hHmpDdlDataSize = new TH1F("hHmpDdlDataSize","HMP Data Size per DDL;;Data Size (Bytes)",14,0.5,14.5);
186   hHmpDdlDataSize->Sumw2();
187   hHmpDdlDataSize->SetOption("P");
188   hHmpDdlDataSize->SetMinimum(0);
189   for(Int_t iddl=0;iddl<14;iddl++)  hHmpDdlDataSize->GetXaxis()->SetBinLabel(iddl+1,Form("DDL_%d",1535+iddl+1));
190   hHmpDdlDataSize->SetStats(0);hHmpDdlDataSize->SetMinimum(0);hHmpDdlDataSize->SetMarkerStyle(20);
191   
192   Add2RawsList(hHmpDdlDataSize,14+14+42+42+1,!expert,image,saveCorr);   //shifter, image
193   
194   TH1F *fHmpPadOcc = new TH1F("fHmpPadOcc","HMP Average pad occupancy per DDL;;Pad occupancy (%)",14,0.5,14.5);
195   fHmpPadOcc->Sumw2();fHmpPadOcc->SetMinimum(0);
196   fHmpPadOcc->SetMinimum(0);
197   for(Int_t iddl=0;iddl<14;iddl++)  fHmpPadOcc->GetXaxis()->SetBinLabel(iddl+1,Form("DDL_%d",1535+iddl+1));
198   fHmpPadOcc->SetStats(0);fHmpPadOcc->SetMinimum(0);fHmpPadOcc->SetMarkerStyle(20);
199   Add2RawsList(fHmpPadOcc,14+14+42+42+2,!expert,image,!saveCorr);       //shifter, image
200
201   TH2F* fHmpBigMap = new TH2F("hHmpBigMap","HMP Sum Q Maps Ch: 0-6;Ch 0-6: pad X;Ch0, Ch1, Ch2, Ch3, Ch4, Ch5, Ch6 pad Y ;Sum Q / Nevt",160,0,160,1008,0,1008);  
202   fHmpBigMap->SetStats(0);  fHmpBigMap->SetOption("COLZ");
203   Add2RawsList(fHmpBigMap,14+14+42+42+3,!expert,image,!saveCorr);       //shifter, image
204    
205   TH2F *fHmpHvSectorQ = new TH2F("fHmpHvSectorQ","HMP HV Sector vs Q; Q (ADC);HV Sector (Ch0-Sc0,Ch0-Sc1,...);Entries*Q/Nevt",410,1,4101,42,0,42);
206   fHmpHvSectorQ->SetStats(0); fHmpHvSectorQ->SetOption("colz");
207   Add2RawsList(fHmpHvSectorQ,14+14+42+42+4,!expert,image,!saveCorr);    //shifter, image
208   
209   // TProfiles
210   TProfile* hHmpDdlDataSizePrf = new TProfile("hHmpDdlDataSizePrf","HMP Data Size per DDL;;Data Size (Bytes)",14,0.5,14.5);
211   hHmpDdlDataSizePrf->Sumw2();
212   hHmpDdlDataSizePrf->SetOption("P");
213   hHmpDdlDataSizePrf->SetMinimum(0);
214   for(Int_t iddl=0;iddl<14;iddl++)  hHmpDdlDataSizePrf->GetXaxis()->SetBinLabel(iddl+1,Form("DDL_%d",1535+iddl+1));
215   hHmpDdlDataSizePrf->SetStats(0);hHmpDdlDataSizePrf->SetMinimum(0);hHmpDdlDataSizePrf->SetMarkerStyle(20);
216   Add2RawsList(hHmpDdlDataSizePrf,14+14+42+42+5,expert,!image,saveCorr);   //expert, no image
217   
218   TProfile *fHmpPadOccPrf = new TProfile("fHmpPadOccPrf","HMP Average pad occupancy per DDL;;Pad occupancy (%)",14,0.5,14.5);
219   fHmpPadOccPrf->Sumw2();fHmpPadOccPrf->SetMinimum(0);
220   fHmpPadOccPrf->SetMinimum(0);
221   for(Int_t iddl=0;iddl<14;iddl++)  fHmpPadOccPrf->GetXaxis()->SetBinLabel(iddl+1,Form("DDL_%d",1535+iddl+1));
222   fHmpPadOccPrf->SetStats(0);fHmpPadOccPrf->SetMinimum(0);fHmpPadOccPrf->SetMarkerStyle(20);
223   Add2RawsList(fHmpPadOccPrf,14+14+42+42+6,expert,!image,saveCorr);       //expert, no image
224
225   
226   //___ for DQM shifter and eLogBook ___ stop
227 }
228 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
229 void AliHMPIDQADataMakerRec::InitESDs()
230 {
231   //
232   //Booking ESDs histograms
233   const Bool_t expert   = kTRUE ; 
234   const Bool_t image    = kTRUE ; 
235
236   TH2F*  hCkovP  = new TH2F("CkovP" , "#theta_{c}, [rad];P, [GeV]"   , 150,      0,  7  ,100, 0, 1)   ;
237   TH2F*  hSigP   = new TH2F("SigP"  ,"#sigma_{#theta_c} [mrad];[GeV]", 150,      0,  7  ,100, 0, 1)   ;
238   TH2F*  hDifXY  = new TH2F("DifXY" ,"diff"                          , 200,    -10, 10  ,200,-10,10)  ;
239   TH2F*  hMvsP = new TH2F("MvsP","Reconstructed Mass vs P",60,0,6,1000,0,1) ;
240   TH1F*  hPid[5];
241   hPid[0] = new TH1F("PidE" ,"electron response"              , 101, -0.005,1.005)             ;
242   hPid[1] = new TH1F("PidMu","#mu response"                   , 101, -0.005,1.005)             ;
243   hPid[2] = new TH1F("PidPi","#pi response"                   , 101, -0.005,1.005)             ;
244   hPid[3] = new TH1F("PidK" ,"K response"                     , 101, -0.005,1.005)             ;
245   hPid[4] = new TH1F("PidP" ,"p response"                     , 101, -0.005,1.005)             ;
246   
247   Add2ESDsList(hCkovP,0, !expert, image);
248   Add2ESDsList(hSigP ,1, expert, !image);
249   Add2ESDsList(hDifXY,2, !expert, image);
250   Add2ESDsList(hMvsP,3, expert, !image);
251   for(Int_t i=0; i< 5; i++) Add2ESDsList(hPid[i],i+4, expert, !image);
252 }
253 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
254 void AliHMPIDQADataMakerRec::MakeRaws(AliRawReader *rawReader)
255 {
256 //
257 // Filling Raws QA histos
258 //
259     rawReader->Reset() ; 
260     Int_t hmpDaqId = AliDAQ::DetectorID("HMPID");                               // shoudl be number 6
261     const UInt_t *detPattern = rawReader->GetDetectorPattern(); 
262     UInt_t isHmpInRawData = ( ((1 << hmpDaqId) & detPattern[0]) >> hmpDaqId);   // check the 6th bit if HMP is there or not
263     if (! isHmpInRawData ) return;                                              // if HMP is not in the event then skip it
264     
265     AliHMPIDRawStream stream(rawReader);
266     Int_t ddlOcc[14]={0};  
267     Int_t isHMPin=0;
268     UInt_t word; Int_t Nddl, r, d, a;
269     Int_t numPadsInDdl;
270                
271     while(stream.Next())
272      {
273        UInt_t ddl=stream.GetDDLNumber(); //returns 0,1,2 ... 13   
274        if(ddl > 13) continue;
275  
276      //  GetRawsData(14+14+42+42+1)->Fill(ddl+1,stream.GetDdlDataSize());
277        GetRawsData(14+14+42+42+5)->Fill(ddl+1,stream.GetDdlDataSize());
278        if(stream.GetDdlDataSize() > 0) 
279          {
280            isHMPin++;
281            //___ fill error histo
282            for(Int_t iErr =1; iErr<(Int_t)AliHMPIDRawStream::kSumErr; iErr++){
283              Int_t numOfErr = stream.GetErrors(ddl,iErr);
284              GetRawsData(ddl)->Fill(iErr,numOfErr);
285              ((TH2I*)GetRawsData(14+14+42+42))->Fill(iErr,ddl,iErr); //
286            }
287            
288            numPadsInDdl= stream.GetNPads();
289            ddlOcc[ddl] = numPadsInDdl;
290            GetRawsData(14+14+42+42+6)->Fill(ddl+1,numPadsInDdl/11520.0*100.0);
291             
292            //___ loop on pads from raw data from a ddl
293            for(Int_t iPad=0;iPad<numPadsInDdl;iPad++) {
294              AliHMPIDDigit dig(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);dig.Raw(word,Nddl,r,d,a);    
295              //for DQM shifter 
296              ((TH2F*)GetRawsData(14+14+42+42+3))->Fill(dig.PadChX(), dig.Ch()*144+dig.PadChY(),dig.Q());
297              ((TH2F*)GetRawsData(14+14+42+42+4))->Fill(dig.Q(),(ddl/2*6)+dig.PadChY()/24,dig.Q());
298             
299              GetRawsData(ddl+14)->Fill(r,d);
300              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));
301              GetRawsData(70+stream.Pc(Nddl,r,d,a)+6*AliHMPIDParam::DDL2C(ddl))->Fill(stream.GetChargeArray()[iPad]);
302             // GetRawsData(14+14+42+42+6)->Fill(ddl+1,1);
303            }//pad loop
304          }     
305      }//next
306     
307     
308     if(isHMPin > 0) fEvtRaw++;
309      
310     
311 }//MakeRaws
312 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
313 void AliHMPIDQADataMakerRec::MakeDigits()
314 {
315   //
316   //filling QA histos for Digits
317   //
318
319   Int_t i = fChannel ; 
320   GetDigitsData(0)->Fill(i,fDigitsArray->GetEntriesFast()/(48.*80.*6.));
321   TIter next(fDigitsArray); 
322   AliHMPIDDigit * digit; 
323   while ( (digit = dynamic_cast<AliHMPIDDigit *>(next())) ) {
324     GetDigitsData(1)->Fill(10.*i+digit->Pc(),1./(48.*80.));
325     GetDigitsData(2+i)->Fill(digit->PadChX(),digit->PadChY());
326     GetDigitsData(9+i*6+digit->Pc())->Fill(digit->Q());
327   }  
328 }  
329   
330 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
331 void AliHMPIDQADataMakerRec::MakeDigits(TTree * digTree)
332 {
333   //
334   //Opening the Digit Tree
335   //
336
337   if(fDigitsArray) 
338     fDigitsArray->Clear() ; 
339   else
340     fDigitsArray=new TClonesArray("AliHMPIDDigit");
341   
342   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
343     fChannel = iCh ; 
344     fDigitsArray->Clear() ; 
345     TBranch *branch = digTree->GetBranch(Form("HMPID%d",iCh));
346     branch->SetAddress(&fDigitsArray);
347     branch->GetEntry(0); 
348     MakeDigits();
349   }
350 }
351
352 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
353 void AliHMPIDQADataMakerRec::MakeRecPoints(TTree * clustersTree)
354 {
355   //
356   //filling QA histos for clusters
357   //
358   AliHMPIDParam *pPar =AliHMPIDParam::Instance();
359  
360   if (fRecPointsArray) 
361     fRecPointsArray->Clear() ; 
362   else 
363     fRecPointsArray = new TClonesArray("AliHMPIDCluster");
364   
365   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
366     TBranch *branch = clustersTree->GetBranch(Form("HMPID%d",iCh));
367     branch->SetAddress(&fRecPointsArray);
368     branch->GetEntry(0);
369     GetRecPointsData(0)->Fill(iCh,fRecPointsArray->GetEntries());
370     TIter next(fRecPointsArray);
371     AliHMPIDCluster *clu;
372     while ( (clu = dynamic_cast<AliHMPIDCluster *>(next())) ) {
373       GetRecPointsData(1)->Fill(clu->Status(),iCh);
374       Int_t sect =  pPar->InHVSector(clu->Y());
375       if(clu->Q()>100) GetRecPointsData(2+iCh)->Fill(clu->Size());
376       else {
377         GetRecPointsData(2+7+iCh)->Fill(clu->Size());
378         GetRecPointsData(2+14+iCh*6+sect)->Fill(clu->Q());
379       }    
380       GetRecPointsData(2+14+42+iCh*6+sect)->Fill(clu->Q());
381     }
382   }
383 }
384 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
385 void AliHMPIDQADataMakerRec::MakeESDs(AliESDEvent * esd)
386 {
387   //
388   //fills QA histos for ESD
389   //
390  
391   for(Int_t iTrk = 0 ; iTrk < esd->GetNumberOfTracks() ; iTrk++){
392     AliESDtrack *pTrk = esd->GetTrack(iTrk) ;
393     GetESDsData(0)->Fill(pTrk->GetP(),pTrk->GetHMPIDsignal());
394     GetESDsData(1)->Fill( pTrk->GetP(),TMath::Sqrt(pTrk->GetHMPIDchi2()));
395     Float_t xm,ym; Int_t q,np;  
396     pTrk->GetHMPIDmip(xm,ym,q,np);                       //mip info
397     Float_t xRad,yRad,th,ph;        
398     pTrk->GetHMPIDtrk(xRad,yRad,th,ph);              //track info at the middle of the radiator
399     Float_t xPc = xRad+9.25*TMath::Tan(th)*TMath::Cos(ph); // temporar: linear extrapol (B=0!)
400     Float_t yPc = yRad+9.25*TMath::Tan(th)*TMath::Sin(ph); // temporar:          "
401     GetESDsData(2)->Fill(xm-xPc,ym-yPc); //track info
402     if(pTrk->GetHMPIDsignal()>0) {
403      Double_t a = 1.292*1.292*TMath::Cos(pTrk->GetHMPIDsignal())*TMath::Cos(pTrk->GetHMPIDsignal())-1.;
404      if(a > 0) {
405     Double_t mass = pTrk->P()*TMath::Sqrt(1.292*1.292*TMath::Cos(pTrk->GetHMPIDsignal())*TMath::Cos(pTrk->GetHMPIDsignal())-1.);
406     GetESDsData(3)->Fill( pTrk->GetP(),mass);
407      }
408     }
409    Double_t pid[5] ;      pTrk->GetHMPIDpid(pid) ;
410     for(Int_t i = 0 ; i < 5 ; i++) GetESDsData(4+i)->Fill(pid[i]) ;
411   }
412 }
413 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
414 void AliHMPIDQADataMakerRec::StartOfDetectorCycle()
415 {
416   //Detector specific actions at start of cycle
417   
418 }
419 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
420
421 void AliHMPIDQADataMakerRec::EndOfDetectorCycle(AliQAv1::TASKINDEX_t task, TObjArray **histos)
422 {
423   //Detector specific actions at end of cycle
424   // do the QA checking
425   if(task==AliQAv1::kRAWS) {
426     for (Int_t specie = 0 ; specie < AliRecoParam::kNSpecies ; specie++) {
427       if (! IsValidEventSpecie(specie, histos) )
428         continue ;
429       for(Int_t iddl=0;iddl<14;iddl++) {
430         TH1F *h = (TH1F*)histos[specie]->At(14+iddl); //ddl histos scaled by the number of events 
431         if(fEvtRaw>0) h->Scale(1.0/(Float_t)fEvtRaw);
432       }//ddl loop
433        
434       TH2F *h2 = (TH2F*)histos[specie]->At(14+14+42+42+3);
435       if(fEvtRaw>0) h2->Scale(1.0/(Float_t)fEvtRaw);
436       TLine *modline[6];
437       for(Int_t modcnt=0; modcnt < 6; modcnt++){ modline[modcnt] = new TLine(0,(1+modcnt)*144,160,(1+modcnt)*144); h2->GetListOfFunctions()->Add(modline[modcnt]); }
438      
439           
440       TH2F *h3 = (TH2F*)histos[specie]->At(14+14+42+42+4);
441       if(fEvtRaw>0) h3->Scale(1.0/(Float_t)fEvtRaw);
442       
443           
444       TLine *lineDdlDatSizeLow  = new TLine(0.5,932,14.5,932);   lineDdlDatSizeLow->SetLineColor(kGreen); lineDdlDatSizeLow->SetLineWidth(2);
445       TLine *lineDdlDatSizeUp   = new TLine(0.5,1500,14.5,1500); lineDdlDatSizeUp->SetLineColor(kGreen);  lineDdlDatSizeUp->SetLineWidth(2);
446       
447       Double_t binval=0,binerr=0;
448       
449       TH1F *h4 = (TH1F*)histos[specie]->At(14+14+42+42+1);
450       TProfile *h4prf = (TProfile*)histos[specie]->At(14+14+42+42+5);
451       for(Int_t iddl=1;iddl<=14;iddl++) 
452         {
453           binval=h4prf->GetBinContent(iddl);  binerr=h4prf->GetBinError(iddl);
454           h4->SetBinContent(iddl,binval);     h4->SetBinError(iddl,binerr);
455          }
456       //if(fEvtRaw>0) h4->Scale(1./(Float_t)fEvtRaw);
457       h4->GetListOfFunctions()->Add(lineDdlDatSizeLow);
458       h4->GetListOfFunctions()->Add(lineDdlDatSizeUp);
459
460       
461             
462       TLine *lineDdlPadOCcLow  = new TLine(0.5,0.086,14.5,0.086);   lineDdlPadOCcLow->SetLineColor(kGreen); lineDdlPadOCcLow->SetLineWidth(2);
463       TLine *lineDdlPadOCcUp   = new TLine(0.5,0.86,14.5,0.86); lineDdlPadOCcUp->SetLineColor(kGreen);  lineDdlPadOCcUp->SetLineWidth(2);
464       TH1F *h5 = (TH1F*)histos[specie]->At(14+14+42+42+2);
465       TProfile *h5prf = (TProfile*)histos[specie]->At(14+14+42+42+6);
466       for(Int_t iddl=1;iddl<=14;iddl++) 
467        {
468          binval=h5prf->GetBinContent(iddl);  binerr=h5prf->GetBinError(iddl);
469          h5->SetBinContent(iddl,binval);     h5->SetBinError(iddl,binerr);
470        }
471
472       //if(fEvtRaw>0) h5->Scale(1./(Float_t)fEvtRaw/11520.0*100.0);
473       h5->GetListOfFunctions()->Add(lineDdlPadOCcLow);
474       h5->GetListOfFunctions()->Add(lineDdlPadOCcUp);
475            
476      }//specie loop     
477   }//RAWS
478    
479   AliQAChecker::Instance()->Run(AliQAv1::kHMPID, task, histos);
480    
481    
482 }
483 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++