]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDQADataMakerRec.cxx
QA local OCDB
[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   TProfile *hCluMult          = new TProfile("CluMult"   ,"Cluster multiplicity per chamber"    , 16, -1 , 7  , 0, 500);
81   Add2RecPointsList(hCluMult    , 0);
82   
83   TH2F *hCluFlg          = new TH2F("CluFlg"      ,"Cluster flag  "                              ,  56  ,-1.5, 12.5, 70, -0.5, 6.5);
84   Add2RecPointsList(hCluFlg    , 1);
85   
86   TH1F *hCluSizeMip[7], *hCluSizePho[7];
87   
88   TH1F *hCluQSect[42], *hCluQSectZoom[42];
89   
90   for(Int_t iCh =0; iCh <7; iCh++){
91    hCluSizeMip[iCh] = new TH1F(Form("CluSizeMipCh%i",iCh),Form("Cluster size  MIP  (cluster Q > 100 ADC) in Chamber %i",iCh),  50  , 0  , 50  );
92   Add2RecPointsList(hCluSizeMip[iCh], iCh+2);
93
94    hCluSizePho[iCh]  = new TH1F(Form("CluSizePho%i",iCh ),Form("Cluster size  Phots(cluster Q < 100 ADC) in Chamber %i",iCh),  50  , 0  , 50  );
95   Add2RecPointsList(hCluSizePho[iCh], iCh+7+2);
96
97     for(Int_t iSect =0; iSect < 6; iSect++){
98       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);
99       Add2RecPointsList(hCluQSectZoom[iCh*6+iSect],2+14+iCh*6+iSect);
100
101       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);
102       Add2RecPointsList(hCluQSect[iCh*6+iSect],2+14+42+iCh*6+iSect);
103
104     }  
105   }
106 }
107 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
108 void AliHMPIDQADataMakerRec::InitRaws()
109 {
110 //
111 // Booking QA histo for Raw data
112 //
113 // All histograms implemented in InitRaws are used in AMORE. Any change here should be propagated to the amoreHMP-QA as well!!! (clm)
114 //  
115   
116   const Int_t kNerr = (Int_t)AliHMPIDRawStream::kSumErr+1;
117   TH1F *hSumErr[14];
118   TH2F *hDilo[14];
119   TH2I *hPadMap[42]; //AMORE monitoring
120   TH1I *hPadQ[42]; //AMORE monitoring
121   for(Int_t iddl =0; iddl<AliHMPIDRawStream::kNDDL; iddl++) {
122     
123     hSumErr[iddl] = new TH1F(Form("hSumErrDDL%i",iddl), Form("Error summary for DDL %i",iddl), 2*kNerr,0,2*kNerr);
124     for(Int_t ilabel=0; ilabel< kNerr; ilabel++) {
125       hSumErr[iddl]->GetXaxis()->CenterLabels(kTRUE);
126       hSumErr[iddl]->GetXaxis()->SetBinLabel((2*ilabel+1),Form("%i  %s",ilabel+1,AliHMPIDRawStream::GetErrName(ilabel)));
127       }
128       
129     Add2RawsList(hSumErr[iddl],iddl);
130     
131     hDilo[iddl] = new TH2F(Form("hDiloDDL%i",iddl),Form("Dilogic response at DDL;Row # ;Dilogic #",iddl),24,1,25,10,1,11);
132     Add2RawsList(hDilo[iddl],14+iddl);
133   }//DDL loop
134  for(Int_t iCh = AliHMPIDParam::kMinCh; iCh <=AliHMPIDParam::kMaxCh ;iCh++) {
135    for(Int_t iPc = AliHMPIDParam::kMinPc; iPc <= AliHMPIDParam::kMaxPc ;iPc++) {
136       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);
137       Add2RawsList(hPadMap[iPc+6*iCh],28+iPc+6*iCh); 
138       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);
139       Add2RawsList(hPadQ[iPc+6*iCh],70+iPc+6*iCh); 
140     }//PC loop
141   }//Ch loop  
142   
143 }
144 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
145 void AliHMPIDQADataMakerRec::InitESDs()
146 {
147   //
148   //Booking ESDs histograms
149    TH2F*  hCkovP  = new TH2F("CkovP" , "#theta_{c}, [rad];P, [GeV]"   , 150,      0,  7  ,100, 0, 1)   ;
150    TH2F*  hSigP   = new TH2F("SigP"  ,"#sigma_{#theta_c} [mrad];[GeV]", 150,      0,  7  ,100, 0, 1)   ;
151    TH2F*  hDifXY  = new TH2F("DifXY" ,"diff"                          , 200,    -10, 10  ,200,-10,10)  ;
152    TH2F*  hMvsP = new TH2F("MvsP","Reconstructed Mass vs P",60,0,6,1000,0,1) ;
153    TH1F*  hPid[5];
154    hPid[0] = new TH1F("PidE" ,"electron response"              , 101, -0.005,1.005)             ;
155    hPid[1] = new TH1F("PidMu","#mu response"                   , 101, -0.005,1.005)             ;
156    hPid[2] = new TH1F("PidPi","#pi response"                   , 101, -0.005,1.005)             ;
157    hPid[3] = new TH1F("PidK" ,"K response"                     , 101, -0.005,1.005)             ;
158    hPid[4] = new TH1F("PidP" ,"p response"                     , 101, -0.005,1.005)             ;
159
160    Add2ESDsList(hCkovP,0);
161    Add2ESDsList(hSigP ,1);
162    Add2ESDsList(hDifXY,2);
163    Add2ESDsList(hMvsP,3);
164    for(Int_t i=0; i< 5; i++) Add2ESDsList(hPid[i],i+4);
165
166
167 }
168 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
169 void AliHMPIDQADataMakerRec::MakeRaws(AliRawReader *rawReader)
170 {
171 //
172 // Filling Raws QA histos
173 //
174     AliHMPIDRawStream stream(rawReader);
175
176     fEvtRaw++;
177
178     while(stream.Next())
179      {
180        UInt_t ddl=stream.GetDDLNumber(); //returns 0,1,2 ... 13
181        if(ddl > 13) continue;
182        for(Int_t iErr =1; iErr<(Int_t)AliHMPIDRawStream::kSumErr; iErr++){
183          Int_t numOfErr = stream.GetErrors(ddl,iErr);
184          GetRawsData(ddl)->Fill(iErr,numOfErr);
185        }
186        UInt_t word; Int_t Nddl, r, d, a;//pc,pcX,pcY;      
187        for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
188        AliHMPIDDigit dig(stream.GetPadArray()[iPad],stream.GetChargeArray()[iPad]);dig.Raw(word,Nddl,r,d,a);
189        GetRawsData(ddl+14)->Fill(r,d);
190        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));
191        GetRawsData(70+stream.Pc(Nddl,r,d,a)+6*AliHMPIDParam::DDL2C(ddl))->Fill(stream.GetChargeArray()[iPad]);
192        }      
193      }
194    stream.Delete();
195    
196 }
197 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
198 void AliHMPIDQADataMakerRec::MakeRecPoints(TTree * clustersTree)
199 {
200   //
201   //filling QA histos for clusters
202   //
203   AliHMPIDParam *pPar =AliHMPIDParam::Instance();
204   
205   static TClonesArray *clusters;
206   if(!clusters) clusters = new TClonesArray("AliHMPIDCluster");
207   
208   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
209     TBranch *branch = clustersTree->GetBranch(Form("HMPID%d",iCh));
210     branch->SetAddress(&clusters);
211     branch->GetEntry(0);
212     GetRecPointsData(0)->Fill(iCh,clusters->GetEntries());
213     TIter next(clusters);
214     AliHMPIDCluster *clu;
215     while ( (clu = dynamic_cast<AliHMPIDCluster *>(next())) ) {
216       GetRecPointsData(1)->Fill(clu->Status(),iCh);
217       Int_t sect =  pPar->InHVSector(clu->Y());
218       if(clu->Q()>100) GetRecPointsData(2+iCh)->Fill(clu->Size());
219       else {
220         GetRecPointsData(2+7+iCh)->Fill(clu->Size());
221         GetRecPointsData(2+14+iCh*6+sect)->Fill(clu->Q());
222       }    
223       GetRecPointsData(2+14+42+iCh*6+sect)->Fill(clu->Q());
224     }
225   }
226
227   clusters->Clear();
228   
229 }
230 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
231 void AliHMPIDQADataMakerRec::MakeESDs(AliESDEvent * esd)
232 {
233   //
234   //fills QA histos for ESD
235   //
236
237   for(Int_t iTrk = 0 ; iTrk < esd->GetNumberOfTracks() ; iTrk++){
238     AliESDtrack *pTrk = esd->GetTrack(iTrk) ;
239     GetESDsData(0)->Fill(pTrk->GetP(),pTrk->GetHMPIDsignal());
240     GetESDsData(1)->Fill( pTrk->GetP(),TMath::Sqrt(pTrk->GetHMPIDchi2()));
241     Float_t xm,ym; Int_t q,np;  
242     pTrk->GetHMPIDmip(xm,ym,q,np);                       //mip info
243     Float_t xRad,yRad,th,ph;        
244     pTrk->GetHMPIDtrk(xRad,yRad,th,ph);              //track info at the middle of the radiator
245     Float_t xPc = xRad+9.25*TMath::Tan(th)*TMath::Cos(ph); // temporar: linear extrapol (B=0!)
246     Float_t yPc = yRad+9.25*TMath::Tan(th)*TMath::Sin(ph); // temporar:          "
247     GetESDsData(2)->Fill(xm-xPc,ym-yPc); //track info
248     if(pTrk->GetHMPIDsignal()>0) {
249      Double_t a = 1.292*1.292*TMath::Cos(pTrk->GetHMPIDsignal())*TMath::Cos(pTrk->GetHMPIDsignal())-1.;
250      if(a > 0) {
251     Double_t mass = pTrk->P()*TMath::Sqrt(1.292*1.292*TMath::Cos(pTrk->GetHMPIDsignal())*TMath::Cos(pTrk->GetHMPIDsignal())-1.);
252     GetESDsData(3)->Fill( pTrk->GetP(),mass);
253      }
254     }
255    Double_t pid[5] ;      pTrk->GetHMPIDpid(pid) ;
256     for(Int_t i = 0 ; i < 5 ; i++) GetESDsData(4+i)->Fill(pid[i]) ;
257   }
258 }
259 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
260 void AliHMPIDQADataMakerRec::StartOfDetectorCycle()
261 {
262   //Detector specific actions at start of cycle
263   
264 }
265 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
266
267 void AliHMPIDQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX_t task, TObjArray *histos)
268 {
269   //Detector specific actions at end of cycle
270   // do the QA checking
271   
272   if(task==AliQA::kRAWS) {
273     for(Int_t iddl=0;iddl<14;iddl++) {
274      TH1F *h = (TH1F*)histos->At(14+iddl); //ddl histos scaled by the number of events 
275      h->Scale(1./(Float_t)fEvtRaw);
276     }
277   }
278   
279    AliQAChecker::Instance()->Run(AliQA::kHMPID, task, histos);
280
281 }
282