]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDQADataMakerRec.cxx
No need to liknk with lhapdf, pythia6 and microcern libraries
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDQADataMakerRec.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16
17 /* $Id$ */
18
19 // --- ROOT system ---
20 #include <TClonesArray.h>
21 #include <TFile.h> 
22 #include <TH1F.h> 
23 #include <TH2F.h>
24 #include <TProfile.h>
25 #include <Riostream.h>
26 // --- Standard library ---
27
28 // --- AliRoot header files ---
29 #include "AliESDCaloCluster.h"
30 #include "AliESDEvent.h"
31 #include "AliQAChecker.h"
32 #include "AliLog.h"
33 #include "AliHMPIDDigit.h"
34 #include "AliHMPIDHit.h"
35 #include "AliHMPIDCluster.h"
36 #include "AliHMPIDQADataMakerRec.h"
37 #include "AliHMPIDQAChecker.h"
38 #include "AliHMPIDParam.h"
39 #include "AliHMPIDRawStream.h"
40 #include "AliLog.h"
41
42 //.
43 // HMPID AliHMPIDQADataMakerRec base class
44 // for QA of reconstruction
45 // here also errors are calculated
46 //.
47
48 ClassImp(AliHMPIDQADataMakerRec)
49            
50 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
51   AliHMPIDQADataMakerRec::AliHMPIDQADataMakerRec() : 
52   AliQADataMakerRec(AliQA::GetDetName(AliQA::kHMPID), "HMPID Quality Assurance Data Maker"),fEvtRaw(0)
53 {
54   // ctor
55 }
56
57 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
58 AliHMPIDQADataMakerRec::AliHMPIDQADataMakerRec(const AliHMPIDQADataMakerRec& qadm) :
59   AliQADataMakerRec(),fEvtRaw(qadm.fEvtRaw)
60 {
61   //copy ctor 
62   SetName((const char*)qadm.GetName()) ; 
63   SetTitle((const char*)qadm.GetTitle()); 
64 }
65
66 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
67 AliHMPIDQADataMakerRec& AliHMPIDQADataMakerRec::operator = (const AliHMPIDQADataMakerRec& qadm )
68 {
69   // Equal operator.
70   this->~AliHMPIDQADataMakerRec();
71   new(this) AliHMPIDQADataMakerRec(qadm);
72   return *this;
73 }
74  
75 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
76 void AliHMPIDQADataMakerRec::InitRecPoints()
77 {
78   // create cluster histograms in RecPoint subdir
79
80   Bool_t expert = kTRUE;
81  
82   TProfile *hCluMult          = new TProfile("CluMult"   ,"Cluster multiplicity per chamber"    , 16, -1 , 7  , 0, 500);
83   Add2RecPointsList(hCluMult    , 0,expert);
84   
85   TH2F *hCluFlg          = new TH2F("CluFlg"      ,"Cluster flag  "                              ,  56  ,-1.5, 12.5, 70, -0.5, 6.5);
86   Add2RecPointsList(hCluFlg    , 1,expert);
87   
88   TH1F *hCluSizeMip[7], *hCluSizePho[7];
89   
90   TH1F *hCluQSect[42], *hCluQSectZoom[42];
91   
92   for(Int_t iCh =0; iCh <7; iCh++){
93    hCluSizeMip[iCh] = new TH1F(Form("CluSizeMipCh%i",iCh),Form("Cluster size  MIP  (cluster Q > 100 ADC) in Chamber %i",iCh),  50  , 0  , 50  );
94   Add2RecPointsList(hCluSizeMip[iCh], iCh+2,expert);
95
96    hCluSizePho[iCh]  = new TH1F(Form("CluSizePho%i",iCh ),Form("Cluster size  Phots(cluster Q < 100 ADC) in Chamber %i",iCh),  50  , 0  , 50  );
97   Add2RecPointsList(hCluSizePho[iCh], iCh+7+2,expert);
98
99     for(Int_t iSect =0; iSect < 6; iSect++){
100       hCluQSectZoom[iCh*6+iSect] = new TH1F(Form("QClusCh%iSect%iZoom",iCh,iSect) ,Form("Zoom on Cluster charge (ADC) in Chamber %i and sector %i",iCh,iSect),100,0,100);
101       Add2RecPointsList(hCluQSectZoom[iCh*6+iSect],2+14+iCh*6+iSect,expert);
102
103       hCluQSect[iCh*6+iSect] = new TH1F(Form("QClusCh%iSect%i",iCh,iSect) ,Form("Cluster charge (ADC) in Chamber %i and sector %i",iCh,iSect),250,0,5000);
104       Add2RecPointsList(hCluQSect[iCh*6+iSect],2+14+42+iCh*6+iSect);
105
106     }  
107   }
108 }
109 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
110 void AliHMPIDQADataMakerRec::InitRaws()
111 {
112 //
113 // Booking QA histo for Raw data
114 //
115 // All histograms implemented in InitRaws are used in AMORE. Any change here should be propagated to the amoreHMP-QA as well!!! (clm)
116 //  
117   
118   Bool_t expert = kTRUE;
119   const Int_t kNerr = (Int_t)AliHMPIDRawStream::kSumErr+1;
120   TH1F *hSumErr[14];
121   TH2F *hDilo[14];
122   TH2I *hPadMap[42]; //AMORE monitoring
123   TH1I *hPadQ[42]; //AMORE monitoring
124   for(Int_t iddl =0; iddl<AliHMPIDRawStream::kNDDL; iddl++) {
125     
126     hSumErr[iddl] = new TH1F(Form("hSumErrDDL%i",iddl), Form("Error summary for DDL %i",iddl), 2*kNerr,0,2*kNerr);
127     for(Int_t ilabel=0; ilabel< kNerr; ilabel++) {
128       hSumErr[iddl]->GetXaxis()->CenterLabels(kTRUE);
129       hSumErr[iddl]->GetXaxis()->SetBinLabel((2*ilabel+1),Form("%i  %s",ilabel+1,AliHMPIDRawStream::GetErrName(ilabel)));
130       }
131       
132     Add2RawsList(hSumErr[iddl],iddl,expert);
133     
134     hDilo[iddl] = new TH2F(Form("hDiloDDL%i",iddl),Form("Dilogic response at DDL;Row # ;Dilogic #",iddl),24,1,25,10,1,11);
135     Add2RawsList(hDilo[iddl],14+iddl,expert);
136   }//DDL loop
137  for(Int_t iCh = AliHMPIDParam::kMinCh; iCh <=AliHMPIDParam::kMaxCh ;iCh++) {
138    for(Int_t iPc = AliHMPIDParam::kMinPc; iPc <= AliHMPIDParam::kMaxPc ;iPc++) {
139       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);
140       Add2RawsList(hPadMap[iPc+6*iCh],28+iPc+6*iCh,expert); 
141       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);
142       Add2RawsList(hPadQ[iPc+6*iCh],70+iPc+6*iCh,expert); 
143     }//PC loop
144   }//Ch loop  
145
146     TH2I *hGeneralErrorSummary = new TH2I("GeneralErrorSummary"," DDL index vs Error type plot", 2*kNerr, 0, 2*kNerr, 2*AliHMPIDRawStream::kNDDL,0,2*AliHMPIDRawStream::kNDDL);
147   for(Int_t igenlabel =0 ; igenlabel< kNerr; igenlabel++) hGeneralErrorSummary->GetXaxis()->SetBinLabel((2*igenlabel+1),Form("%i  %s",igenlabel+1,AliHMPIDRawStream::GetErrName(igenlabel)));
148    Add2RawsList(hGeneralErrorSummary,14+14+42+42);
149   
150 }
151 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
152 void AliHMPIDQADataMakerRec::InitESDs()
153 {
154   //
155   //Booking ESDs histograms
156    TH2F*  hCkovP  = new TH2F("CkovP" , "#theta_{c}, [rad];P, [GeV]"   , 150,      0,  7  ,100, 0, 1)   ;
157    TH2F*  hSigP   = new TH2F("SigP"  ,"#sigma_{#theta_c} [mrad];[GeV]", 150,      0,  7  ,100, 0, 1)   ;
158    TH2F*  hDifXY  = new TH2F("DifXY" ,"diff"                          , 200,    -10, 10  ,200,-10,10)  ;
159    TH2F*  hMvsP = new TH2F("MvsP","Reconstructed Mass vs P",60,0,6,1000,0,1) ;
160    TH1F*  hPid[5];
161    hPid[0] = new TH1F("PidE" ,"electron response"              , 101, -0.005,1.005)             ;
162    hPid[1] = new TH1F("PidMu","#mu response"                   , 101, -0.005,1.005)             ;
163    hPid[2] = new TH1F("PidPi","#pi response"                   , 101, -0.005,1.005)             ;
164    hPid[3] = new TH1F("PidK" ,"K response"                     , 101, -0.005,1.005)             ;
165    hPid[4] = new TH1F("PidP" ,"p response"                     , 101, -0.005,1.005)             ;
166
167    Add2ESDsList(hCkovP,0);
168    Add2ESDsList(hSigP ,1,kTRUE);
169    Add2ESDsList(hDifXY,2);
170    Add2ESDsList(hMvsP,3,kTRUE);
171    for(Int_t i=0; i< 5; i++) Add2ESDsList(hPid[i],i+4,kTRUE);
172
173 }
174 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
175 void AliHMPIDQADataMakerRec::MakeRaws(AliRawReader *rawReader)
176 {
177 //
178 // Filling Raws QA histos
179 //
180           rawReader->Reset() ; 
181                 AliHMPIDRawStream stream(rawReader);
182
183     fEvtRaw++;
184
185     while(stream.Next())
186      {
187        UInt_t ddl=stream.GetDDLNumber(); //returns 0,1,2 ... 13
188        if(ddl > 13) continue;
189        for(Int_t iErr =1; iErr<(Int_t)AliHMPIDRawStream::kSumErr; iErr++){
190          Int_t numOfErr = stream.GetErrors(ddl,iErr);
191          GetRawsData(ddl)->Fill(iErr,numOfErr);
192          //Printf("err type %i   ddl number %i  Num of errors %i",iErr,ddl,numOfErr);
193         ((TH2I*)GetRawsData(14+14+42+42))->Fill(iErr,ddl,iErr); //
194        }
195        UInt_t word; Int_t Nddl, r, d, a;//pc,pcX,pcY;      
196        for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
197        AliHMPIDDigit dig(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);dig.Raw(word,Nddl,r,d,a);
198        GetRawsData(ddl+14)->Fill(r,d);
199        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));
200        GetRawsData(70+stream.Pc(Nddl,r,d,a)+6*AliHMPIDParam::DDL2C(ddl))->Fill(stream.GetChargeArray()[iPad]);
201        }      
202      }
203
204    stream.Delete();
205    
206 }
207 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
208 void AliHMPIDQADataMakerRec::MakeRecPoints(TTree * clustersTree)
209 {
210   //
211   //filling QA histos for clusters
212   //
213   AliHMPIDParam *pPar =AliHMPIDParam::Instance();
214   
215   static TClonesArray *clusters;
216   if(!clusters) clusters = new TClonesArray("AliHMPIDCluster");
217   
218   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
219     TBranch *branch = clustersTree->GetBranch(Form("HMPID%d",iCh));
220     branch->SetAddress(&clusters);
221     branch->GetEntry(0);
222     GetRecPointsData(0)->Fill(iCh,clusters->GetEntries());
223     TIter next(clusters);
224     AliHMPIDCluster *clu;
225     while ( (clu = dynamic_cast<AliHMPIDCluster *>(next())) ) {
226       GetRecPointsData(1)->Fill(clu->Status(),iCh);
227       Int_t sect =  pPar->InHVSector(clu->Y());
228       if(clu->Q()>100) GetRecPointsData(2+iCh)->Fill(clu->Size());
229       else {
230         GetRecPointsData(2+7+iCh)->Fill(clu->Size());
231         GetRecPointsData(2+14+iCh*6+sect)->Fill(clu->Q());
232       }    
233       GetRecPointsData(2+14+42+iCh*6+sect)->Fill(clu->Q());
234     }
235   }
236
237   clusters->Clear();
238   
239 }
240 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
241 void AliHMPIDQADataMakerRec::MakeESDs(AliESDEvent * esd)
242 {
243   //
244   //fills QA histos for ESD
245   //
246
247   for(Int_t iTrk = 0 ; iTrk < esd->GetNumberOfTracks() ; iTrk++){
248     AliESDtrack *pTrk = esd->GetTrack(iTrk) ;
249     GetESDsData(0)->Fill(pTrk->GetP(),pTrk->GetHMPIDsignal());
250     GetESDsData(1)->Fill( pTrk->GetP(),TMath::Sqrt(pTrk->GetHMPIDchi2()));
251     Float_t xm,ym; Int_t q,np;  
252     pTrk->GetHMPIDmip(xm,ym,q,np);                       //mip info
253     Float_t xRad,yRad,th,ph;        
254     pTrk->GetHMPIDtrk(xRad,yRad,th,ph);              //track info at the middle of the radiator
255     Float_t xPc = xRad+9.25*TMath::Tan(th)*TMath::Cos(ph); // temporar: linear extrapol (B=0!)
256     Float_t yPc = yRad+9.25*TMath::Tan(th)*TMath::Sin(ph); // temporar:          "
257     GetESDsData(2)->Fill(xm-xPc,ym-yPc); //track info
258     if(pTrk->GetHMPIDsignal()>0) {
259      Double_t a = 1.292*1.292*TMath::Cos(pTrk->GetHMPIDsignal())*TMath::Cos(pTrk->GetHMPIDsignal())-1.;
260      if(a > 0) {
261     Double_t mass = pTrk->P()*TMath::Sqrt(1.292*1.292*TMath::Cos(pTrk->GetHMPIDsignal())*TMath::Cos(pTrk->GetHMPIDsignal())-1.);
262     GetESDsData(3)->Fill( pTrk->GetP(),mass);
263      }
264     }
265    Double_t pid[5] ;      pTrk->GetHMPIDpid(pid) ;
266     for(Int_t i = 0 ; i < 5 ; i++) GetESDsData(4+i)->Fill(pid[i]) ;
267   }
268 }
269 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
270 void AliHMPIDQADataMakerRec::StartOfDetectorCycle()
271 {
272   //Detector specific actions at start of cycle
273   
274 }
275 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
276
277 void AliHMPIDQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray *histos)
278 {
279   //Detector specific actions at end of cycle
280   // do the QA checking
281   
282   if(task==AliQA::kRAWS) {
283     for(Int_t iddl=0;iddl<14;iddl++) {
284      TH1F *h = (TH1F*)histos->At(14+iddl); //ddl histos scaled by the number of events 
285      h->Scale(1./(Float_t)fEvtRaw);
286     }
287   }
288   
289    AliQAChecker::Instance()->Run(AliQA::kHMPID, task, histos);
290
291 }
292