1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //==============================================================================
17 // AliHMPIDTaskQA - Class representing a quality check tool of HMPID
18 // A set of histograms is created.
19 //==============================================================================
21 #ifndef AliHMPIDTaskQA_CXX
22 #define AliHMPIDTaskQA_CXX
29 #include "TGraphErrors.h"
30 #include "AliAnalysisManager.h"
31 #include "AliESDInputHandler.h"
32 #include "AliMCEventHandler.h"
33 #include "AliMCEvent.h"
34 #include "AliESDtrack.h"
37 #include "AliHMPIDTaskQA.h"
39 ClassImp(AliHMPIDTaskQA)
41 //__________________________________________________________________________
42 AliHMPIDTaskQA::AliHMPIDTaskQA() :
43 fESD(0x0),fMC(0x0),fUseMC(kTRUE),
45 fHmpPesdPhmp(0x0),fHmpCkovPesd(0x0),fHmpCkovPhmp(0x0),
46 fHmpMipCharge3cm(0x0),fHmpMipTrkDist(0x0),
66 for (Int_t i=0; i<7; i++){
69 fHmpPhotSin2th[i] = 0x0;
70 fHmpMipTrkDistPosX[i] = fHmpMipTrkDistNegX[i] = 0x0;
71 fHmpMipTrkDistPosY[i] = fHmpMipTrkDistNegY[i] = 0x0;
72 fHmpMipCharge[i] = 0x0;
76 //___________________________________________________________________________
77 AliHMPIDTaskQA::AliHMPIDTaskQA(const Char_t* name) :
78 AliAnalysisTaskSE(name),
79 fESD(0x0), fMC(0x0), fUseMC(kTRUE),
81 fHmpPesdPhmp(0x0),fHmpCkovPesd(0x0),fHmpCkovPhmp(0x0),
82 fHmpMipCharge3cm(0x0),fHmpMipTrkDist(0x0),
100 // Constructor. Initialization of Inputs and Outputs
102 for (Int_t i=0; i<7; i++){
103 fHmpPhotons[i] = 0x0;
105 fHmpPhotSin2th[i] = 0x0;
106 fHmpMipTrkDistPosX[i] = fHmpMipTrkDistNegX[i] = 0x0;
107 fHmpMipTrkDistPosY[i] = fHmpMipTrkDistNegY[i] = 0x0;
108 fHmpMipCharge[i] = 0x0;
111 DefineOutput(1,TList::Class());
114 //___________________________________________________________________________
115 AliHMPIDTaskQA& AliHMPIDTaskQA::operator=(const AliHMPIDTaskQA& c)
118 // Assignment operator
121 AliAnalysisTaskSE::operator=(c);
125 fHmpHistList = c.fHmpHistList;
126 fHmpPesdPhmp = c.fHmpPesdPhmp;
127 fHmpCkovPesd = c.fHmpCkovPesd;
128 fHmpCkovPhmp = c.fHmpCkovPhmp;
129 fHmpMipCharge3cm = c.fHmpMipCharge3cm;
130 fHmpMipTrkDist = c.fHmpMipTrkDist;
131 fHmpTrkFlags = c.fHmpTrkFlags;
134 fPionEff = c.fPionEff;
135 fKaonEff = c.fKaonEff;
136 fProtEff = c.fProtEff;
137 fPionTot = c.fPionTot;
138 fKaonTot = c.fKaonTot;
139 fProtTot = c.fProtTot;
140 fPionNot = c.fPionNot;
141 fKaonNot = c.fKaonNot;
142 fProtNot = c.fProtNot;
143 fPionCon = c.fPionCon;
144 fKaonCon = c.fKaonCon;
145 fProtCon = c.fProtCon;
146 for (Int_t i=0; i<7; i++){
147 fHmpPhotons[i] = c.fHmpPhotons[i];
148 fHmpPhotP[i] = c.fHmpPhotP[i];
149 fHmpPhotSin2th[i] = c.fHmpPhotSin2th[i];
150 fHmpMipTrkDistPosX[i] = c.fHmpMipTrkDistPosX[i];
151 fHmpMipTrkDistNegX[i] = c.fHmpMipTrkDistNegX[i];
152 fHmpMipTrkDistPosY[i] = c.fHmpMipTrkDistPosY[i];
153 fHmpMipTrkDistNegY[i] = c.fHmpMipTrkDistNegY[i];
154 fHmpMipCharge[i] = c.fHmpMipCharge[i];
160 //___________________________________________________________________________
161 AliHMPIDTaskQA::AliHMPIDTaskQA(const AliHMPIDTaskQA& c) :
162 AliAnalysisTaskSE(c),
163 fESD(c.fESD),fMC(c.fMC),fUseMC(c.fUseMC),
164 fHmpHistList(c.fHmpHistList),
165 fHmpPesdPhmp(c.fHmpPesdPhmp),fHmpCkovPesd(c.fHmpCkovPesd),fHmpCkovPhmp(c.fHmpCkovPhmp),
166 fHmpMipCharge3cm(c.fHmpMipCharge3cm),fHmpMipTrkDist(c.fHmpMipTrkDist),
167 fHmpTrkFlags(c.fHmpTrkFlags),
170 fPionEff(c.fPionEff),
171 fKaonEff(c.fKaonEff),
172 fProtEff(c.fProtEff),
173 fPionTot(c.fPionTot),
174 fKaonTot(c.fKaonTot),
175 fProtTot(c.fProtTot),
176 fPionNot(c.fPionNot),
177 fKaonNot(c.fKaonNot),
178 fProtNot(c.fProtNot),
179 fPionCon(c.fPionCon),
180 fKaonCon(c.fKaonCon),
186 for (Int_t i=0; i<7; i++){
187 fHmpPhotons[i] = c.fHmpPhotons[i];
188 fHmpPhotP[i] = c.fHmpPhotP[i];
189 fHmpPhotSin2th[i] = c.fHmpPhotSin2th[i];
190 fHmpMipTrkDistPosX[i] = c.fHmpMipTrkDistPosX[i];
191 fHmpMipTrkDistNegX[i] = c.fHmpMipTrkDistNegX[i];
192 fHmpMipTrkDistPosY[i] = c.fHmpMipTrkDistPosY[i];
193 fHmpMipTrkDistNegY[i] = c.fHmpMipTrkDistNegY[i];
194 fHmpMipCharge[i] = c.fHmpMipCharge[i];
198 //___________________________________________________________________________
199 AliHMPIDTaskQA::~AliHMPIDTaskQA() {
203 Info("~AliHMPIDTaskQA","Calling Destructor");
204 if (fHmpHistList && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) delete fHmpHistList;
207 //___________________________________________________________________________
208 void AliHMPIDTaskQA::ConnectInputData(Option_t *option)
210 AliAnalysisTaskSE::ConnectInputData(option);
213 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
215 AliDebug(2,Form("ERROR: Could not get ESDInputHandler"));
217 fESD = esdH->GetEvent();
221 AliMCEventHandler *mcH = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
223 AliDebug(2,Form("ERROR: Could not get MCEventHandler"));
229 //___________________________________________________________________________
230 void AliHMPIDTaskQA::UserExec(Option_t *)
232 Double_t priors[10]={1.,1.,1.,1.,1.,0.,0.,0.,0.,0.}; //{0.01,0.01,0.83,0.10,0.5};
234 AliPID *pPid = new AliPID();
235 pPid->SetPriors(priors);
236 AliESDtrack *track=0;
238 AliStack* pStack = 0;
241 AliMCEventHandler *mcH = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
243 fMC = mcH->MCEvent();
244 pStack = fMC->Stack();
249 // Main loop function, executed on Event basis
252 if( !fESD->GetPrimaryVertex()) return;
253 if( fESD->GetPrimaryVertex()->GetNContributors() < 1 ) return;
254 if( TMath::Abs(fESD->GetPrimaryVertex()->GetZ()) > 10.0 /* cm */) return;
257 for (Int_t iTrack = 0; iTrack < fESD->GetNumberOfTracks(); iTrack++) {
259 track = fESD->GetTrack(iTrack);
261 if(! track->IsOn(AliESDtrack::kITSrefit) || !track->IsOn(AliESDtrack::kITSrefit) ) continue;
263 Double_t rin[3], rout[3];
264 track->GetInnerXYZ(rin);
265 track->GetOuterXYZ(rout);
266 Double_t ktol = 0.001;
267 Double_t thetaCh = track->GetHMPIDsignal();
269 if(Equal(thetaCh,-20.,ktol)) fHmpTrkFlags->Fill(0);
270 else if(Equal(thetaCh,-9.,ktol)) fHmpTrkFlags->Fill(1);
271 else if(Equal(thetaCh,-5.,ktol)) fHmpTrkFlags->Fill(2);
272 else if(Equal(thetaCh,-11.,ktol)) fHmpTrkFlags->Fill(3);
273 else if(Equal(thetaCh,-22.,ktol)) fHmpTrkFlags->Fill(4);
274 else fHmpTrkFlags->Fill(4);
276 if(Equal(thetaCh,-20.,ktol)) continue;
277 if(track->GetHMPIDcluIdx() < 0) continue;
281 Float_t xpc, ypc, th, ph;
282 track->GetHMPIDmip(x,y,q,nph);
283 track->GetHMPIDtrk(xpc,ypc,th,ph);
285 if(Equal(x,0.,ktol) && Equal(y,0.,ktol) && Equal(xpc,0.,ktol) && Equal(ypc,0.,ktol)) continue;
287 Float_t sign = track->GetSign();
288 Int_t ch = track->GetHMPIDcluIdx()/1000000;
289 if( ch < 0 || ch > 6 ) continue; //Chamber numbering is [0,6]
291 if (sign > 0.) fHmpMipTrkDistPosX[ch]->Fill(xpc - x), fHmpMipTrkDistPosY[ch]->Fill(ypc - y);
292 if (sign < 0.) fHmpMipTrkDistNegX[ch]->Fill(xpc - x), fHmpMipTrkDistNegY[ch]->Fill(ypc - y);
293 fHmpMipCharge[ch]->Fill(q);
294 Double_t dist = TMath::Sqrt((xpc - x)*(xpc - x) + (ypc - y)*(ypc - y));
295 if (dist <= 3.0) fHmpMipCharge3cm->Fill(q);
296 fHmpMipTrkDist->Fill(dist);
298 track->GetHMPIDpid(probs);
299 pPid->SetProbabilities(probs);
301 if ((label = track->GetLabel()) < 0) continue;
302 pPart = pStack->Particle(label);
306 Double_t pHmp[3] = {0}, pHmp3 = 0;
307 if (track->GetOuterHmpPxPyPz(pHmp)) pHmp3 = TMath::Sqrt(pHmp[0]*pHmp[0]+pHmp[1]*pHmp[1]+pHmp[2]*pHmp[2]);
308 if (pHmp3) fHmpPesdPhmp->Fill(track->P(),pHmp3);
309 fHmpCkovPesd->Fill(track->P(),thetaCh);
310 if (pHmp3) fHmpCkovPhmp->Fill(pHmp3,thetaCh);
311 fHmpPhotons[ch]->Fill(nph);
312 fHmpPhotP[ch]->Fill(pHmp3,nph);
313 Double_t sin2 = TMath::Power(TMath::Sin(th),2);
314 fHmpPhotSin2th[ch]->Fill(sin2,nph);
316 if (fUseMC && dist<0.5 && TMath::Abs(th)<0.13){
317 if (!pStack->IsPhysicalPrimary(label)) continue;
318 if (track->Pt()<1. || track->Pt()>5.) continue;
319 Int_t pdgCode = TMath::Abs(pPart->GetPdgCode());
320 Int_t ptBin=(Int_t) (2*(track->Pt()-1));
321 if (pdgCode!=2212) fProtCon->Fill(ptBin);
323 fProtTot->Fill(ptBin);
324 fProtEff->Fill(ptBin,pPid->GetProbability(AliPID::kProton));
325 fPionNot->Fill(ptBin,pPid->GetProbability(AliPID::kPion));
326 fKaonNot->Fill(ptBin,pPid->GetProbability(AliPID::kKaon));
328 if (pdgCode!=211) fPionCon->Fill(ptBin);
329 if (pdgCode!=321) fKaonCon->Fill(ptBin);
332 Float_t weight=pPid->GetProbability(AliPID::kElectron)+
333 pPid->GetProbability(AliPID::kMuon)+
334 pPid->GetProbability(AliPID::kPion);
335 fPionTot->Fill(ptBin);
336 fPionEff->Fill(ptBin,weight);
337 fKaonNot->Fill(ptBin,pPid->GetProbability(AliPID::kKaon));
339 fProtNot->Fill(ptBin,pPid->GetProbability(AliPID::kProton));
343 fKaonTot->Fill(ptBin);
344 fKaonEff->Fill(ptBin,pPid->GetProbability(AliPID::kKaon));
345 fPionNot->Fill(ptBin,(pPid->GetProbability(AliPID::kPion)));
347 fProtNot->Fill(ptBin,(pPid->GetProbability(AliPID::kProton)));
354 /* PostData(0) is taken care of by AliAnalysisTaskSE */
355 PostData(1,fHmpHistList);
358 //___________________________________________________________________________
359 void AliHMPIDTaskQA::Terminate(Option_t*)
361 // The Terminate() function is the last function to be called during
362 // a query. It always runs on the client, it can be used to present
363 // the results graphically or save the results to file.
365 Info("Terminate"," ");
366 AliAnalysisTaskSE::Terminate();
370 //___________________________________________________________________________
371 void AliHMPIDTaskQA::UserCreateOutputObjects() {
373 //HERE ONE CAN CREATE OUTPUT OBJECTS
374 //TO BE SET BEFORE THE EXECUTION OF THE TASK
379 fHmpHistList = new TList();
380 fHmpHistList->SetOwner();
382 fHmpPesdPhmp = new TH2F("fHmpPesdPhmp","HMPID: ESD p - running p;HMP p (GeV/c);ESD p (Gev/c)",100,0,10,100,0,10);
383 fHmpHistList->Add(fHmpPesdPhmp);
385 fHmpCkovPesd = new TH2F("fHmpCkovPesd","HMPID: ThetaCherenkov vs P;p_esd (GeV/c);#Theta_C;Entries",500,0,10,500,0,1.0);
386 fHmpHistList->Add(fHmpCkovPesd);
388 fHmpCkovPhmp = new TH2F("fHmpCkovPhmp","HMPID: ThetaCherenkov vs P;p_hmp (GeV/c);#Theta_C;Entries",500,0,10,500,0,1.0);
389 fHmpHistList->Add(fHmpCkovPhmp);
391 for (Int_t i=0; i<7; i++){
392 TString title=Form("MIP-Track distance in local X, Chamber %d;distance (cm);Entries",i);
393 fHmpMipTrkDistPosX[i] = new TH1F(Form("fHmpMipTrkDistPosX%d",i),title.Data(),800,-20,20);
394 fHmpHistList->Add(fHmpMipTrkDistPosX[i]);
395 fHmpMipTrkDistNegX[i] = new TH1F(Form("fHmpMipTrkDistNegX%d",i),title.Data(),800,-20,20);
396 fHmpHistList->Add(fHmpMipTrkDistNegX[i]);
398 title=Form("MIP-Track distance in local Y, Chamber %d;distance (cm);Entries",i);
399 fHmpMipTrkDistPosY[i] = new TH1F(Form("fHmpMipTrkDistPosY%d",i),title.Data(),800,-20,20);
400 fHmpHistList->Add(fHmpMipTrkDistPosY[i]);
401 fHmpMipTrkDistNegY[i] = new TH1F(Form("fHmpMipTrkDistNegY%d",i),title.Data(),800,-20,20);
402 fHmpHistList->Add(fHmpMipTrkDistNegY[i]);
404 title=Form("Mip charge distribution, Chamber %d;Charge;Entries",i);
405 fHmpMipCharge[i] = new TH1F(Form("fHmpMipCharge%d",i),title.Data(),5001,-0.5,5000.5);
406 fHmpHistList->Add(fHmpMipCharge[i]);
408 title=Form("Photons, Chamber %d;N of Photons;Entries",i);
409 fHmpPhotons[i] = new TH1F(Form("fHmpPhotons%d",i),title.Data(),100,0,100);
410 fHmpHistList->Add(fHmpPhotons[i]);
412 title=Form("Photons versus HMP momentum, Chamber %d;P (GeV);N of Photons",i);
413 fHmpPhotP[i] = new TH2F(Form("fHmpPhotP%d",i),title.Data(),200,0.,20.,100,0,100);
414 fHmpHistList->Add(fHmpPhotP[i]);
416 title=Form("Photons versus sin(th)^2 (uncorrected), Chamber %d;P (GeV);N of Photons",i);
417 fHmpPhotSin2th[i] = new TH2F(Form("fHmpPhotSin2th%d",i),title.Data(),100,0.,1.,100,0,100);
418 fHmpHistList->Add(fHmpPhotSin2th[i]);
421 fHmpMipCharge3cm = new TH1F("fHmpMipCharge3cm","HMPID MIP Charge;MIP Charge (ADC);Entries",5001,-0.5,5000.5);
422 fHmpHistList->Add(fHmpMipCharge3cm);
424 fHmpMipTrkDist = new TH1F("fHmpMipTrkDist","Mip-track distance in all the chambers",100,0.,10.);
425 fHmpHistList->Add(fHmpMipTrkDist);
427 fHmpTrkFlags = new TH1F("fHmpTrkFlags","HMPID track flags",6,0,6);
428 TString summary[6] = {"NotPerformed","MipDistCut", "MipQdcCut", "NoPhotAccept", "kNoRad", "other"};
429 for(Int_t ibin = 0; ibin < 6; ibin++) fHmpTrkFlags->GetXaxis()->SetBinLabel(ibin+1,Form("%i %s",ibin+1,summary[ibin].Data()));
430 fHmpHistList->Add(fHmpTrkFlags);
432 fPionEff = new TH1F("PionEff","Identified pions",fN1,0,fN1);
433 fKaonEff = new TH1F("KaonEff","Identified kaons",fN1,0,fN1);
434 fProtEff = new TH1F("ProtEff","Identified protons",fN2,0,fN2);
435 fPionTot = new TH1I("PionTot","Total MC pions",fN1,0,fN1);
436 fKaonTot = new TH1I("KaonTot","Total MC kaons",fN1,0,fN1);
437 fProtTot = new TH1I("ProtTot","Total MC protons",fN2,0,fN2);
438 fPionNot = new TH1F("PionNot","Misidentified pions",fN1,0,fN1);
439 fKaonNot = new TH1F("KaonNot","Misidentified kaons",fN1,0,fN1);
440 fProtNot = new TH1F("ProtNot","Misidentified protons",fN2,0,fN2);
441 fPionCon = new TH1I("PionCon","Total not MC pions",fN1,0,fN1);
442 fKaonCon = new TH1I("KaonCon","Total not MC kaons",fN1,0,fN1);
443 fProtCon = new TH1I("ProtCon","Total not MC protons",fN2,0,fN2);
445 fHmpHistList->Add(fPionEff); fHmpHistList->Add(fKaonEff); fHmpHistList->Add(fProtEff);
446 fHmpHistList->Add(fPionTot); fHmpHistList->Add(fKaonTot); fHmpHistList->Add(fProtTot);
447 fHmpHistList->Add(fPionNot); fHmpHistList->Add(fKaonNot); fHmpHistList->Add(fProtNot);
448 fHmpHistList->Add(fPionCon); fHmpHistList->Add(fKaonCon); fHmpHistList->Add(fProtCon);
450 PostData(1,fHmpHistList);
453 //____________________________________________________________________________________________________________________________________
454 Bool_t AliHMPIDTaskQA::Equal(Double_t x, Double_t y, Double_t tolerance)
456 return abs(x - y) <= tolerance ;