]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDQADataMakerRec.cxx
coding convention for type (_t) and access of reference data from alien OCDB
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDQADataMakerRec.cxx
CommitLineData
04236e67 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
04236e67 19// --- ROOT system ---
20#include <TClonesArray.h>
21#include <TFile.h>
22#include <TH1F.h>
23#include <TH2F.h>
cc1abc7f 24#include <TProfile.h>
04236e67 25#include <Riostream.h>
26// --- Standard library ---
27
28// --- AliRoot header files ---
29#include "AliESDCaloCluster.h"
30#include "AliESDEvent.h"
cc1abc7f 31#include "AliQAChecker.h"
04236e67 32#include "AliLog.h"
33#include "AliHMPIDDigit.h"
34#include "AliHMPIDHit.h"
35#include "AliHMPIDCluster.h"
36#include "AliHMPIDQADataMakerRec.h"
cc1abc7f 37#include "AliHMPIDParam.h"
38#include "AliHMPIDRawStream.h"
39#include "AliLog.h"
b38ac33a 40
606697a8 41//.
b38ac33a 42// HMPID AliHMPIDQADataMakerRec base class
43// for QA of reconstruction
44// here also errors are calculated
45//.
46
04236e67 47ClassImp(AliHMPIDQADataMakerRec)
48
49//____________________________________________________________________________
50 AliHMPIDQADataMakerRec::AliHMPIDQADataMakerRec() :
cc1abc7f 51 AliQADataMakerRec(AliQA::GetDetName(AliQA::kHMPID), "HMPID Quality Assurance Data Maker")
04236e67 52{
53 // ctor
04236e67 54}
55
56//____________________________________________________________________________
57AliHMPIDQADataMakerRec::AliHMPIDQADataMakerRec(const AliHMPIDQADataMakerRec& qadm) :
cc1abc7f 58 AliQADataMakerRec()
04236e67 59{
60 //copy ctor
04236e67 61 SetName((const char*)qadm.GetName()) ;
62 SetTitle((const char*)qadm.GetTitle());
63}
64
65//__________________________________________________________________
66AliHMPIDQADataMakerRec& AliHMPIDQADataMakerRec::operator = (const AliHMPIDQADataMakerRec& qadm )
67{
68 // Equal operator.
69 this->~AliHMPIDQADataMakerRec();
70 new(this) AliHMPIDQADataMakerRec(qadm);
71 return *this;
72}
73
74//____________________________________________________________________________
75
76void AliHMPIDQADataMakerRec::InitRecPoints()
77{
78 // create cluster histograms in RecPoint subdir
cc1abc7f 79
80 TH1F *hCluEvt=new TH1F("CluPerEvt","Cluster multiplicity" ,100,0,100);
81 TH1F *hCluChi2 =new TH1F("CluChi2" ,"Chi2 " ,1000,0,100);
82 TH1F *hCluFlg =new TH1F("CluFlg" ,"Cluster flag" ,14,-1.5,12.5); hCluFlg->SetFillColor(5);
83 TH1F *hCluSize =new TH1F("CluSize" ,"Cluster size ",100,0,100);
84 TH1F *hCluQ =new TH1F("CluQ" ,"Cluster charge (ADC)",1000,0,5000);
85
86 Add2RecPointsList(hCluEvt , 0);
87 Add2RecPointsList(hCluChi2, 1);
88 Add2RecPointsList(hCluFlg , 2);
89 Add2RecPointsList(hCluSize, 3);
90 Add2RecPointsList(hCluQ , 4);
91}
92//____________________________________________________________________________
93
94void AliHMPIDQADataMakerRec::InitRaws()
95{
96//
97// Booking QA histo for Raw data
98//
b38ac33a 99 const Int_t kNerr = (Int_t)AliHMPIDRawStream::kSumErr+1;
cc1abc7f 100 TH1F *hqPad[14], *hSumErr[14];
101
3c8e86a0 102 for(Int_t iddl =0; iddl<AliHMPIDRawStream::kNDDL; iddl++) {
103 hqPad[iddl] = new TH1F(Form("hqPadDDL%i",iddl), Form("Pad Q Entries at DDL %i",iddl), 500,0,5000);
104 Add2RawsList(hqPad[iddl],iddl);
b38ac33a 105 hSumErr[iddl] = new TH1F(Form("SumErrDDL%i",iddl), Form("Error summary for ddl %i",iddl), 2*kNerr,0,2*kNerr);
3c8e86a0 106 hSumErr[iddl]->SetYTitle("%");
107
b38ac33a 108 for(Int_t ilabel=0; ilabel< kNerr; ilabel++) {
3c8e86a0 109 hSumErr[iddl]->GetXaxis()->CenterLabels(kTRUE);
110 //hSumErr[iddl]->GetXaxis()->SetBinLabel((2*ilabel+1),Form("%i %s",ilabel+1,hnames[ilabel]));
111 hSumErr[iddl]->GetXaxis()->SetBinLabel((2*ilabel+1),Form("%i %s",ilabel+1,AliHMPIDRawStream::GetErrName(ilabel)));
112 }
113
114 Add2RawsList(hSumErr[iddl],iddl+14);
cc1abc7f 115 }
116 TH1F *hNevRaws = new TH1F("NevRaws","Events per DDL",15,0,15);
117 Add2RawsList(hNevRaws,28);
04236e67 118}
119//____________________________________________________________________________
120void AliHMPIDQADataMakerRec::InitESDs()
121{
cc1abc7f 122 //
123 //Booking ESDs histograms
124 TH2F* hCkovP = new TH2F("CkovP" , "#theta_{c}, [rad];P, [GeV]" , 150, 0, 7 ,100, 0, 1) ;
125 TH2F* hSigP = new TH2F("SigP" ,"#sigma_{#theta_c} [mrad];[GeV]", 150, 0, 7 ,100, 0, 1) ;
126 TH2F* hMipXY = new TH2F("MipXY" ,"mip position" , 260, 0,130 ,252, 0,126) ;
127 TH2F* hDifXY = new TH2F("DifXY" ,"diff" , 200, -10, 10 ,200,-10,10) ;
128 TH1F* hPid[5];
129 hPid[0] = new TH1F("PidE" ,"electron response" , 101, -0.005,1.005) ;
130 hPid[1] = new TH1F("PidMu","#mu response" , 101, -0.005,1.005) ;
131 hPid[2] = new TH1F("PidPi","#pi response" , 101, -0.005,1.005) ;
132 hPid[3] = new TH1F("PidK" ,"K response" , 101, -0.005,1.005) ;
133 hPid[4] = new TH1F("PidP" ,"p response" ,101, -0.005,1.005) ;
134
135Add2ESDsList(hCkovP,0);
136Add2ESDsList(hSigP ,1);
137Add2ESDsList(hMipXY,2);
138Add2ESDsList(hDifXY,3);
139for(Int_t i=0; i< 5; i++) Add2ESDsList(hPid[i],i+4);
04236e67 140}
04236e67 141//____________________________________________________________________________
cc1abc7f 142void AliHMPIDQADataMakerRec::MakeRaws(AliRawReader *rawReader)
04236e67 143{
cc1abc7f 144//
145// Filling Raws QA histos
146//
147 AliHMPIDRawStream stream(rawReader);
148
149 while(stream.Next())
150 {
151 UInt_t ddl=stream.GetDDLNumber(); //returns 0,1,2 ... 13
04236e67 152
cc1abc7f 153 for(Int_t iPad=0;iPad<stream.GetNPads();iPad++) {
154 GetRawsData(ddl)->Fill(stream.GetChargeArray()[iPad]);}
155
156 GetRawsData(28)->Fill(ddl);
157
158 for(Int_t iErr =1; iErr<(Int_t)AliHMPIDRawStream::kSumErr; iErr++){
159
b38ac33a 160 Int_t numOfErr = stream.GetErrors(ddl,iErr);
cc1abc7f 161
b38ac33a 162 GetRawsData(ddl+14)->Fill(iErr,numOfErr);
cc1abc7f 163 }
164 }
165 stream.Delete();
166}
167//___________________________________________________________________________
168void AliHMPIDQADataMakerRec::MakeRecPoints(TTree * clustersTree)
169{
170 //
171 //filling QA histos for clusters
172 //
04236e67 173 TClonesArray *clusters = new TClonesArray("AliHMPIDCluster");
174 for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){
175 TBranch *branch = clustersTree->GetBranch(Form("HMPID%d",i));
176 branch->SetAddress(&clusters);
177 branch->GetEntry(0);
178
cc1abc7f 179 GetRecPointsData(0)->Fill(i,clusters->GetEntries());
04236e67 180 TIter next(clusters);
181 AliHMPIDCluster *clu;
cc1abc7f 182 while ( (clu = dynamic_cast<AliHMPIDCluster *>(next())) ) {
183 GetRecPointsData(1)->Fill(clu->Chi2());
184 GetRecPointsData(2)->Fill(clu->Status());
185 GetRecPointsData(3)->Fill(clu->Size());
186 GetRecPointsData(4)->Fill(clu->Q());
04236e67 187 }
188 }
189
190 clusters->Delete();
191 delete clusters;
192}
193
194//____________________________________________________________________________
195void AliHMPIDQADataMakerRec::MakeESDs(AliESDEvent * esd)
196{
cc1abc7f 197 //
04236e67 198 //fills QA histos for ESD
cc1abc7f 199 //
04236e67 200 for(Int_t iTrk = 0 ; iTrk < esd->GetNumberOfTracks() ; iTrk++){
201 AliESDtrack *pTrk = esd->GetTrack(iTrk) ;
cc1abc7f 202 GetESDsData(0)->Fill(pTrk->GetP(),pTrk->GetHMPIDsignal());
203 GetESDsData(1)->Fill( pTrk->GetP(),TMath::Sqrt(pTrk->GetHMPIDchi2()));
04236e67 204 Float_t xm,ym; Int_t q,np;
205 pTrk->GetHMPIDmip(xm,ym,q,np); //mip info
cc1abc7f 206 GetESDsData(2)->Fill(xm,ym);
04236e67 207 Float_t xRad,yRad,th,ph;
208 pTrk->GetHMPIDtrk(xRad,yRad,th,ph); //track info at the middle of the radiator
209 Float_t xPc = xRad+9.25*TMath::Tan(th)*TMath::Cos(ph); // temporar: linear extrapol (B=0!)
210 Float_t yPc = yRad+9.25*TMath::Tan(th)*TMath::Sin(ph); // temporar: "
cc1abc7f 211 GetESDsData(3)->Fill(xm-xPc,ym-yPc); //track info
04236e67 212 Double_t pid[5] ; pTrk->GetHMPIDpid(pid) ;
cc1abc7f 213 for(Int_t i = 0 ; i < 5 ; i++) GetESDsData(4+i)->Fill(pid[i]) ;
214 }
215}
216//____________________________________________________________________________
217void AliHMPIDQADataMakerRec::StartOfDetectorCycle()
218{
219 //Detector specific actions at start of cycle
220
221}
222
223void AliHMPIDQADataMakerRec::EndOfDetectorCycle(AliQA::TASKINDEX, TObjArray *)
224{
225 //Detector specific actions at end of cycle
226 // do the QA checking
227 //AliQAChecker::Instance()->Run(AliQA::kHMPID, task, obj);
228 //if(task==AliQA::kRAWS){
229
71d0d954 230// for(Int_t iddl=0; iddl<14; iddl++)
231// {
232// if(GetRawsData(28)->GetBinContent(iddl)!=0) GetRawsData(iddl+14)->Scale(100./GetRawsData(28)->GetBinContent(iddl));
233// }
04236e67 234}
235