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