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