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