]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG1/TRD/AliTRDcheckPID.cxx
prepare PID LQ ref maker for production
[u/mrichter/AliRoot.git] / PWG1 / TRD / AliTRDcheckPID.cxx
CommitLineData
1ee39b3a 1//////////////////////////////////////////////////////
2//
3// PID performance checker of the TRD
4//
5// Performs checks of ESD information for TRD-PID and recalculate PID based on OCDB information
6// Also provides performance plots on detector based on the PID information - for the moment only
7// MC source is used but V0 is also possible. Here is a list of detector performance checks
8// - Integrated dE/dx per chamber
9// - <PH> as function of time bin and local radial position
10// - number of clusters/tracklet
11// - number of tracklets/track
12//
13// Author : Alex Wilk <wilka@uni-muenster.de>
14// Alex Bercuci <A.Bercuci@gsi.de>
15// Markus Fasel <M.Fasel@gsi.de>
16//
17///////////////////////////////////////////////////////
18
19#include "TAxis.h"
20#include "TROOT.h"
21#include "TPDGCode.h"
22#include "TCanvas.h"
23#include "TF1.h"
24#include "TH1F.h"
25#include "TH1D.h"
26#include "TH2F.h"
27#include "TProfile.h"
28#include "TProfile2D.h"
29#include "TGraph.h"
30#include "TGraphErrors.h"
31#include "TLegend.h"
32
33#include <TClonesArray.h>
34#include <TObjArray.h>
35#include <TList.h>
36
37#include "AliESDEvent.h"
38#include "AliESDInputHandler.h"
39#include "AliTrackReference.h"
40
41#include "AliAnalysisTask.h"
42
43#include "AliTRDtrackerV1.h"
44#include "AliTRDtrackV1.h"
45#include "AliTRDcluster.h"
46#include "AliTRDReconstructor.h"
47#include "AliCDBManager.h"
48#include "AliTRDpidUtil.h"
49
2e6aa247 50#include "Cal/AliTRDCalPID.h"
51#include "Cal/AliTRDCalPIDNN.h"
1ee39b3a 52#include "AliTRDcheckPID.h"
53#include "info/AliTRDtrackInfo.h"
b91fdd71 54#include "info/AliTRDpidInfo.h"
1ee39b3a 55
2e6aa247 56Char_t const * AliTRDcheckPID::fgMethod[3] = {"LQ", "NN", "ESD"};
57
1ee39b3a 58ClassImp(AliTRDcheckPID)
59
60//________________________________________________________________________
61AliTRDcheckPID::AliTRDcheckPID()
62 :AliTRDrecoTask("checkPID", "TRD PID checker")
b91fdd71 63 ,fReconstructor(NULL)
64 ,fUtil(NULL)
65 ,fGraph(NULL)
66 ,fPID(NULL)
67 ,fMomentumAxis(NULL)
1ee39b3a 68 ,fMinNTracklets(AliTRDgeometry::kNlayer)
69 ,fMaxNTracklets(AliTRDgeometry::kNlayer)
70 {
71 //
72 // Default constructor
73 //
74
75 fReconstructor = new AliTRDReconstructor();
76 fReconstructor->SetRecoParam(AliTRDrecoParam::GetLowFluxParam());
77
78 // Initialize momentum axis with default values
79 Double_t defaultMomenta[AliTRDCalPID::kNMom+1];
80 for(Int_t imom = 0; imom < AliTRDCalPID::kNMom+1; imom++)
81 defaultMomenta[imom] = AliTRDCalPID::GetMomentumBinning(imom);
82 SetMomentumBinning(AliTRDCalPID::kNMom, defaultMomenta);
83
2e6aa247 84 memset(fEfficiency, 0, AliPID::kSPECIES*sizeof(TObjArray*));
85
1ee39b3a 86 fUtil = new AliTRDpidUtil();
87 InitFunctorList();
b91fdd71 88 DefineOutput(1, TObjArray::Class()); // PID for reference maker
1ee39b3a 89}
90
91
92//________________________________________________________________________
93AliTRDcheckPID::~AliTRDcheckPID()
94{
b91fdd71 95 if(fPID){fPID->Delete(); delete fPID;}
96 if(fGraph){fGraph->Delete(); delete fGraph;}
97 if(fReconstructor) delete fReconstructor;
98 if(fUtil) delete fUtil;
1ee39b3a 99}
100
101
102//________________________________________________________________________
103void AliTRDcheckPID::CreateOutputObjects()
104{
105 // Create histograms
106 // Called once
107
108 OpenFile(0, "RECREATE");
109 fContainer = Histos();
b91fdd71 110
111 fPID = new TObjArray();
112 fPID->SetOwner(kTRUE);
113}
114
115//________________________________________________________
116void AliTRDcheckPID::Exec(Option_t *opt)
117{
118 //
119 // Execution part
120 //
121
122 fPID->Delete();
123
124 AliTRDrecoTask::Exec(opt);
125
126 PostData(1, fPID);
1ee39b3a 127}
128
129
130//_______________________________________________________
131TObjArray * AliTRDcheckPID::Histos(){
132
133 //
134 // Create QA histograms
135 //
136 if(fContainer) return fContainer;
137
138 Int_t xBins = AliPID::kSPECIES*fMomentumAxis->GetNbins();
139 fContainer = new TObjArray(); fContainer->Expand(kNPlots);
140
141 const Float_t epsilon = 1./(2*(AliTRDpidUtil::kBins-1)); // get nice histos with bin center at 0 and 1
b91fdd71 142 TH1 *h = NULL;
1ee39b3a 143
2e6aa247 144 const Int_t kNmethodsPID=Int_t(sizeof(fgMethod)/sizeof(Char_t*));
145 // histos with posterior probabilities for all particle species
146 for(Int_t is=AliPID::kSPECIES; is--;){
147 fEfficiency[is] = new TObjArray(kNmethodsPID);
148 fEfficiency[is]->SetOwner();
149 fEfficiency[is]->SetName(Form("Eff_%s", AliPID::ParticleShortName(is)));
150 fContainer->AddAt(fEfficiency[is], is?kEfficiencyMu+is-1:kEfficiency);
151 for(Int_t im=kNmethodsPID; im--;){
152 if(!(h = (TH2F*)gROOT->FindObject(Form("PID_%s_%s", fgMethod[im], AliPID::ParticleShortName(is))))) {
153 h = new TH2F(Form("PID_%s_%s", fgMethod[im], AliPID::ParticleShortName(is)), "", xBins, -0.5, xBins - 0.5,
154 AliTRDpidUtil::kBins, 0.-epsilon, 1.+epsilon);
155 } else h->Reset();
156 fEfficiency[is]->AddAt(h, im);
157 }
158 }
1ee39b3a 159
160 // histos of the dE/dx distribution for all 5 particle species and 11 momenta
161 if(!(h = (TH2F*)gROOT->FindObject("dEdx"))){
162 h = new TH2F("dEdx", "",
163 xBins, -0.5, xBins - 0.5,
2e6aa247 164 200, 0, 15);
165// 200, 0, 10000);
1ee39b3a 166 } else h->Reset();
167 fContainer->AddAt(h, kdEdx);
168
169 // histos of the dE/dx slices for all 5 particle species and 11 momenta
170 if(!(h = (TH2F*)gROOT->FindObject("dEdxSlice"))){
171 h = new TH2F("dEdxSlice", "",
172 xBins*AliTRDpidUtil::kLQslices, -0.5, xBins*AliTRDpidUtil::kLQslices - 0.5,
173 200, 0, 5000);
174 } else h->Reset();
175 fContainer->AddAt(h, kdEdxSlice);
176
177 // histos of the pulse height distribution for all 5 particle species and 11 momenta
178 TObjArray *fPH = new TObjArray(2);
179 fPH->SetOwner(); fPH->SetName("PH");
180 fContainer->AddAt(fPH, kPH);
181 if(!(h = (TProfile2D*)gROOT->FindObject("PHT"))){
182 h = new TProfile2D("PHT", "",
183 xBins, -0.5, xBins - 0.5,
184 AliTRDtrackerV1::GetNTimeBins(), -0.5, AliTRDtrackerV1::GetNTimeBins() - 0.5);
185 } else h->Reset();
186 fPH->AddAt(h, 0);
187 if(!(h = (TProfile2D*)gROOT->FindObject("PHX"))){
188 h = new TProfile2D("PHX", "",
189 xBins, -0.5, xBins - 0.5,
190 AliTRDtrackerV1::GetNTimeBins(), 0., .5*AliTRDgeometry::CamHght()+AliTRDgeometry::CdrHght());
191 } else h->Reset();
192 fPH->AddAt(h, 1);
193
194 // histos of the number of clusters distribution for all 5 particle species and 11 momenta
195 if(!(h = (TH2F*)gROOT->FindObject("NClus"))){
196 h = new TH2F("NClus", "",
197 xBins, -0.5, xBins - 0.5,
198 50, -0.5, 49.5);
199 } else h->Reset();
200 fContainer->AddAt(h, kNClus);
201
202
203 // momentum distributions - absolute and in momentum bins
204 if(!(h = (TH1F*)gROOT->FindObject("hMom"))){
205 h = new TH1F("hMom", "momentum distribution", fMomentumAxis->GetNbins(), fMomentumAxis->GetXmin(), fMomentumAxis->GetXmax());
206 } else h->Reset();
207 fContainer->AddAt(h, kMomentum);
208
209 if(!(h = (TH1F*)gROOT->FindObject("hMomBin"))){
210 h = new TH1F("hMomBin", "momentum distribution in momentum bins", fMomentumAxis->GetNbins(), fMomentumAxis->GetXmin(), fMomentumAxis->GetXmax());
211 } else h->Reset();
212 fContainer->AddAt(h, kMomentumBin);
213
214 // Number of tracklets per track
215 if(!(h = (TH2F*)gROOT->FindObject("nTracklets"))){
216 h = new TH2F("nTracklets", "",
217 xBins, -0.5, xBins - 0.5,
218 AliTRDgeometry::kNlayer, 0.5, AliTRDgeometry::kNlayer+.5);
219 } else h->Reset();
220 fContainer->AddAt(h, kNTracklets);
221
222 return fContainer;
223}
224
225
226//________________________________________________________________________
227Bool_t AliTRDcheckPID::CheckTrackQuality(const AliTRDtrackV1* track) const
228{
229 //
230 // Check if the track is ok for PID
231 //
232
233 Int_t ntracklets = track->GetNumberOfTracklets();
234 if(ntracklets >= fMinNTracklets && ntracklets <= fMaxNTracklets) return 1;
235// if(!fESD)
236// return 0;
237
238 return 0;
239}
240
241//________________________________________________________________________
242Int_t AliTRDcheckPID::CalcPDG(AliTRDtrackV1* track)
243{
244// Documentation to come
245
246 track -> SetReconstructor(fReconstructor);
247 (const_cast<AliTRDrecoParam*>(fReconstructor->GetRecoParam()))->SetPIDNeuralNetwork();
248 track -> CookPID();
249
250 if(track -> GetPID(AliPID::kElectron) > track -> GetPID(AliPID::kMuon) + track -> GetPID(AliPID::kPion) + track -> GetPID(AliPID::kKaon) + track -> GetPID(AliPID::kProton)){
251 return kElectron;
252 }
253 else if(track -> GetPID(kProton) > track -> GetPID(AliPID::kPion) && track -> GetPID(AliPID::kProton) > track -> GetPID(AliPID::kKaon) && track -> GetPID(AliPID::kProton) > track -> GetPID(AliPID::kMuon)){
254 return kProton;
255 }
256 else if(track -> GetPID(AliPID::kKaon) > track -> GetPID(AliPID::kMuon) && track -> GetPID(AliPID::kKaon) > track -> GetPID(AliPID::kPion)){
257 return kKPlus;
258 }
259 else if(track -> GetPID(AliPID::kMuon) > track -> GetPID(AliPID::kPion)){
260 return kMuonPlus;
261 }
262 else{
263 return kPiPlus;
264 }
265}
266
267
268//_______________________________________________________
269TH1 *AliTRDcheckPID::PlotLQ(const AliTRDtrackV1 *track)
270{
271 //
272 // Plot the probabilities for electrons using 2-dim LQ
273 //
274
275 if(!fkESD){
276 AliWarning("No ESD info available.");
2e6aa247 277 return NULL;
1ee39b3a 278 }
279
280 if(track) fkTrack = track;
281 if(!fkTrack){
282 AliWarning("No Track defined.");
2e6aa247 283 return NULL;
1ee39b3a 284 }
285
286 ULong_t status;
287 status = fkESD -> GetStatus();
2e6aa247 288 if(!(status&AliESDtrack::kTRDin)) return NULL;
1ee39b3a 289
2e6aa247 290 if(!CheckTrackQuality(fkTrack)) return NULL;
1ee39b3a 291
292 AliTRDtrackV1 cTrack(*fkTrack);
293 cTrack.SetReconstructor(fReconstructor);
294
295 Int_t pdg = 0;
296 Float_t momentum = 0.;
297
298 if(fkMC){
299 if(fkMC->GetTrackRef()) momentum = fkMC->GetTrackRef()->P();
300 pdg = fkMC->GetPDG();
301 } else{
302 //AliWarning("No MC info available!");
303 momentum = cTrack.GetMomentum(0);
304 pdg = CalcPDG(&cTrack);
305 }
2e6aa247 306 if(!IsInRange(momentum)) return NULL;
1ee39b3a 307
308 (const_cast<AliTRDrecoParam*>(fReconstructor->GetRecoParam()))->SetPIDNeuralNetwork(kFALSE);
309 cTrack.CookPID();
2e6aa247 310 if(cTrack.GetNumberOfTrackletsPID() < fMinNTracklets) return NULL;
1ee39b3a 311 Int_t species = AliTRDpidUtil::Pdg2Pid(pdg);
2e6aa247 312
313 TH2F *hPIDLQ(NULL);
314 TObjArray *eff(NULL);
315 for(Int_t is=AliPID::kSPECIES; is--;){
316 if(!(eff = dynamic_cast<TObjArray *>(fContainer->At(is?kEfficiencyMu+is-1:kEfficiency)))){
317 AliWarning("No Histogram List defined.");
318 return NULL;
319 }
320 if(!(hPIDLQ = dynamic_cast<TH2F *>(eff->At(AliTRDpidUtil::kLQ)))){
321 AliWarning("No Histogram defined.");
322 return NULL;
323 }
324
325 hPIDLQ -> Fill(FindBin(species, momentum), cTrack.GetPID(is));
326 }
1ee39b3a 327 return hPIDLQ;
328}
329
330
2e6aa247 331
332
1ee39b3a 333//_______________________________________________________
334TH1 *AliTRDcheckPID::PlotNN(const AliTRDtrackV1 *track)
335{
336 //
337 // Plot the probabilities for electrons using 2-dim LQ
338 //
339
340 if(!fkESD){
341 AliWarning("No ESD info available.");
2e6aa247 342 return NULL;
1ee39b3a 343 }
344
345 if(track) fkTrack = track;
346 if(!fkTrack){
347 AliWarning("No Track defined.");
2e6aa247 348 return NULL;
1ee39b3a 349 }
350
351 ULong_t status;
352 status = fkESD -> GetStatus();
2e6aa247 353 if(!(status&AliESDtrack::kTRDin)) return NULL;
1ee39b3a 354
2e6aa247 355 if(!CheckTrackQuality(fkTrack)) return NULL;
1ee39b3a 356
1ee39b3a 357 AliTRDtrackV1 cTrack(*fkTrack);
358 cTrack.SetReconstructor(fReconstructor);
359
360 Int_t pdg = 0;
361 Float_t momentum = 0.;
362 if(fkMC){
363 if(fkMC->GetTrackRef()) momentum = fkMC->GetTrackRef()->P();
364 pdg = fkMC->GetPDG();
365 } else {
366 //AliWarning("No MC info available!");
367 momentum = cTrack.GetMomentum(0);
368 pdg = CalcPDG(&cTrack);
369 }
2e6aa247 370 if(!IsInRange(momentum)) return NULL;
1ee39b3a 371
372 (const_cast<AliTRDrecoParam*>(fReconstructor->GetRecoParam()))->SetPIDNeuralNetwork();
373 cTrack.CookPID();
2e6aa247 374 if(cTrack.GetNumberOfTrackletsPID() < fMinNTracklets) return NULL;
1ee39b3a 375
376 Int_t species = AliTRDpidUtil::Pdg2Pid(pdg);
2e6aa247 377
378
379 TH2F *hPIDNN(NULL);
380 TObjArray *eff(NULL);
381 for(Int_t is=AliPID::kSPECIES; is--;){
382 if(!(eff = dynamic_cast<TObjArray *>(fContainer->At(is?kEfficiencyMu+is-1:kEfficiency)))){
383 AliWarning("No Histogram List defined.");
384 return NULL;
385 }
386 if(!(hPIDNN = dynamic_cast<TH2F *>(eff->At(AliTRDpidUtil::kNN)))){
387 AliWarning("No Histogram defined.");
388 return NULL;
389 }
390
391 hPIDNN->Fill(FindBin(species, momentum), cTrack.GetPID(is));
392 }
1ee39b3a 393 return hPIDNN;
394}
395
396
2e6aa247 397
1ee39b3a 398//_______________________________________________________
399TH1 *AliTRDcheckPID::PlotESD(const AliTRDtrackV1 *track)
400{
401 //
402 // Plot the probabilities for electrons using 2-dim LQ
403 //
404
405 if(!fkESD){
406 AliWarning("No ESD info available.");
2e6aa247 407 return NULL;
1ee39b3a 408 }
409
410 if(track) fkTrack = track;
411 if(!fkTrack){
412 AliWarning("No Track defined.");
2e6aa247 413 return NULL;
1ee39b3a 414 }
415
416 ULong_t status;
417 status = fkESD -> GetStatus();
2e6aa247 418 if(!(status&AliESDtrack::kTRDin)) return NULL;
1ee39b3a 419
2e6aa247 420 if(!CheckTrackQuality(fkTrack)) return NULL;
421 if(fkTrack->GetNumberOfTrackletsPID() < fMinNTracklets) return NULL;
1ee39b3a 422
1ee39b3a 423 Int_t pdg = 0;
424 Float_t momentum = 0.;
425 if(fkMC){
426 if(fkMC->GetTrackRef()) momentum = fkMC->GetTrackRef()->P();
427 pdg = fkMC->GetPDG();
428 } else {
429 //AliWarning("No MC info available!");
430 AliTRDtrackV1 cTrack(*fkTrack);
431 cTrack.SetReconstructor(fReconstructor);
432 momentum = cTrack.GetMomentum(0);
433 pdg = CalcPDG(&cTrack);
434 }
2e6aa247 435 if(!IsInRange(momentum)) return NULL;
1ee39b3a 436
1ee39b3a 437 const Double32_t *pidESD = fkESD->GetResponseIter();
438 Int_t species = AliTRDpidUtil::Pdg2Pid(pdg);
2e6aa247 439
440 TH2F *hPID(NULL);
441 TObjArray *eff(NULL);
442 for(Int_t is=AliPID::kSPECIES; is--;){
443 if(!(eff = dynamic_cast<TObjArray *>(fContainer->At(is?kEfficiencyMu+is-1:kEfficiency)))){
444 AliWarning("No Histogram List defined.");
445 return NULL;
446 }
447 if(!(hPID = dynamic_cast<TH2F *>(eff->At(AliTRDpidUtil::kESD)))){
448 AliWarning("No Histogram defined.");
449 return NULL;
450 }
b91fdd71 451
452 AliTRDpidInfo *pid = new AliTRDpidInfo();
453 fPID->Add(pid);
2e6aa247 454 hPID->Fill(FindBin(species, momentum), pidESD[is]);
455 }
456 return hPID;
1ee39b3a 457}
458
459
2e6aa247 460
1ee39b3a 461//_______________________________________________________
462TH1 *AliTRDcheckPID::PlotdEdx(const AliTRDtrackV1 *track)
463{
464 //
465 // Plot the probabilities for electrons using 2-dim LQ
466 //
467
468 if(track) fkTrack = track;
469 if(!fkTrack){
470 AliWarning("No Track defined.");
2e6aa247 471 return NULL;
1ee39b3a 472 }
473
2e6aa247 474 if(!CheckTrackQuality(fkTrack)) return NULL;
1ee39b3a 475
476 TH2F *hdEdx;
477 if(!(hdEdx = dynamic_cast<TH2F *>(fContainer->At(kdEdx)))){
478 AliWarning("No Histogram defined.");
2e6aa247 479 return NULL;
1ee39b3a 480 }
481
482 AliTRDtrackV1 cTrack(*fkTrack);
483 cTrack.SetReconstructor(fReconstructor);
484 Int_t pdg = 0;
485 Float_t momentum = 0.;
486 if(fkMC){
487 if(fkMC->GetTrackRef()) momentum = fkMC->GetTrackRef()->P();
488 pdg = fkMC->GetPDG();
489 } else {
490 //AliWarning("No MC info available!");
491 momentum = cTrack.GetMomentum(0);
492 pdg = CalcPDG(&cTrack);
493 }
494 Int_t species = AliTRDpidUtil::Pdg2Pid(pdg);
2e6aa247 495 if(!IsInRange(momentum)) return NULL;
1ee39b3a 496
2e6aa247 497 (const_cast<AliTRDrecoParam*>(fReconstructor->GetRecoParam()))->SetPIDNeuralNetwork(kTRUE);
498// (const_cast<AliTRDrecoParam*>(fReconstructor->GetRecoParam()))->SetPIDNeuralNetwork(kFALSE);
1ee39b3a 499 Float_t sumdEdx = 0;
500 Int_t iBin = FindBin(species, momentum);
b91fdd71 501 AliTRDseedV1 *tracklet = NULL;
1ee39b3a 502 for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
503 sumdEdx = 0;
504 tracklet = cTrack.GetTracklet(iChamb);
505 if(!tracklet) continue;
2e6aa247 506 tracklet -> CookdEdx(AliTRDpidUtil::kNNslices);
507 for(Int_t i = AliTRDpidUtil::kNNslices; i--;) sumdEdx += tracklet->GetdEdx()[i];
508// tracklet -> CookdEdx(AliTRDpidUtil::kLQslices);
509// for(Int_t i = 3; i--;) sumdEdx += tracklet->GetdEdx()[i];
510 sumdEdx /= AliTRDCalPIDNN::kMLPscale;
1ee39b3a 511 hdEdx -> Fill(iBin, sumdEdx);
512 }
513 return hdEdx;
514}
515
516
517//_______________________________________________________
518TH1 *AliTRDcheckPID::PlotdEdxSlice(const AliTRDtrackV1 *track)
519{
520 //
521 // Plot the probabilities for electrons using 2-dim LQ
522 //
523
524 if(track) fkTrack = track;
525 if(!fkTrack){
526 AliWarning("No Track defined.");
2e6aa247 527 return NULL;
1ee39b3a 528 }
529
2e6aa247 530 if(!CheckTrackQuality(fkTrack)) return NULL;
1ee39b3a 531
532 TH2F *hdEdxSlice;
533 if(!(hdEdxSlice = dynamic_cast<TH2F *>(fContainer->At(kdEdxSlice)))){
534 AliWarning("No Histogram defined.");
2e6aa247 535 return NULL;
1ee39b3a 536 }
537
538 AliTRDtrackV1 cTrack(*fkTrack);
539 cTrack.SetReconstructor(fReconstructor);
540 Int_t pdg = 0;
541 Float_t momentum = 0.;
542 if(fkMC){
543 if(fkMC->GetTrackRef()) momentum = fkMC->GetTrackRef()->P();
544 pdg = fkMC->GetPDG();
545 } else {
546 //AliWarning("No MC info available!");
547 momentum = cTrack.GetMomentum(0);
548 pdg = CalcPDG(&cTrack);
549 }
2e6aa247 550 if(!IsInRange(momentum)) return NULL;
1ee39b3a 551
552 (const_cast<AliTRDrecoParam*>(fReconstructor->GetRecoParam()))->SetPIDNeuralNetwork(kFALSE);
553 Int_t iMomBin = fMomentumAxis->FindBin(momentum);
554 Int_t species = AliTRDpidUtil::Pdg2Pid(pdg);
555 Float_t *fdEdx;
b91fdd71 556 AliTRDseedV1 *tracklet = NULL;
1ee39b3a 557 for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
558 tracklet = cTrack.GetTracklet(iChamb);
559 if(!tracklet) continue;
560 tracklet -> CookdEdx(AliTRDpidUtil::kLQslices);
561 fdEdx = const_cast<Float_t *>(tracklet->GetdEdx());
562 for(Int_t iSlice = 0; iSlice < AliTRDpidUtil::kLQslices; iSlice++){
563 hdEdxSlice -> Fill(species * fMomentumAxis->GetNbins() * AliTRDpidUtil::kLQslices + (iMomBin-1) * AliTRDpidUtil::kLQslices + iSlice, fdEdx[iSlice]);
564 }
565 }
566
567 return hdEdxSlice;
568}
569
570
571//_______________________________________________________
572TH1 *AliTRDcheckPID::PlotPH(const AliTRDtrackV1 *track)
573{
574 //
575 // Plot the probabilities for electrons using 2-dim LQ
576 //
577
578 if(track) fkTrack = track;
579 if(!fkTrack){
580 AliWarning("No Track defined.");
2e6aa247 581 return NULL;
1ee39b3a 582 }
583
2e6aa247 584 if(!CheckTrackQuality(fkTrack)) return NULL;
1ee39b3a 585
b91fdd71 586 TObjArray *arr = NULL;
1ee39b3a 587 TProfile2D *hPHX, *hPHT;
588 if(!(arr = dynamic_cast<TObjArray *>(fContainer->At(kPH)))){
589 AliWarning("No Histogram defined.");
2e6aa247 590 return NULL;
1ee39b3a 591 }
592 hPHT = (TProfile2D*)arr->At(0);
593 hPHX = (TProfile2D*)arr->At(1);
594
595 Int_t pdg = 0;
596 Float_t momentum = 0.;
597 if(fkMC){
598 if(fkMC->GetTrackRef()) momentum = fkMC->GetTrackRef()->P();
599 pdg = fkMC->GetPDG();
600 } else {
601 //AliWarning("No MC info available!");
602 AliTRDtrackV1 cTrack(*fkTrack);
603 cTrack.SetReconstructor(fReconstructor);
604 momentum = cTrack.GetMomentum(0);
605 pdg = CalcPDG(&cTrack);
606 }
2e6aa247 607 if(!IsInRange(momentum)) return NULL;;
1ee39b3a 608
b91fdd71 609 AliTRDseedV1 *tracklet = NULL;
610 AliTRDcluster *cluster = NULL;
1ee39b3a 611 Int_t species = AliTRDpidUtil::Pdg2Pid(pdg);
612 Int_t iBin = FindBin(species, momentum);
613 for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
614 tracklet = fkTrack->GetTracklet(iChamb);
615 if(!tracklet) continue;
616 Float_t x0 = tracklet->GetX0();
617 for(Int_t ic = 0; ic < AliTRDtrackerV1::GetNTimeBins(); ic++){
618 if(!(cluster = tracklet->GetClusters(ic))) continue;
619 hPHT -> Fill(iBin, cluster->GetLocalTimeBin(), TMath::Abs(cluster->GetQ()));
620 hPHX -> Fill(iBin, x0 - cluster->GetX(), tracklet->GetdQdl(ic));
621 }
622 }
623 return hPHT;
624}
625
626
627//_______________________________________________________
628TH1 *AliTRDcheckPID::PlotNClus(const AliTRDtrackV1 *track)
629{
630 //
631 // Plot the probabilities for electrons using 2-dim LQ
632 //
633
634 if(track) fkTrack = track;
635 if(!fkTrack){
636 AliWarning("No Track defined.");
2e6aa247 637 return NULL;
1ee39b3a 638 }
639
2e6aa247 640 if(!CheckTrackQuality(fkTrack)) return NULL;
1ee39b3a 641
642 TH2F *hNClus;
643 if(!(hNClus = dynamic_cast<TH2F *>(fContainer->At(kNClus)))){
644 AliWarning("No Histogram defined.");
2e6aa247 645 return NULL;
1ee39b3a 646 }
647
648
649 Int_t pdg = 0;
650 Float_t momentum = 0.;
651 if(fkMC){
652 if(fkMC->GetTrackRef()) momentum = fkMC->GetTrackRef()->P();
653 pdg = fkMC->GetPDG();
654 } else {
655 //AliWarning("No MC info available!");
656 AliTRDtrackV1 cTrack(*fkTrack);
657 cTrack.SetReconstructor(fReconstructor);
658 momentum = cTrack.GetMomentum(0);
659 pdg = CalcPDG(&cTrack);
660 }
2e6aa247 661 if(!IsInRange(momentum)) return NULL;
1ee39b3a 662
663 Int_t species = AliTRDpidUtil::AliTRDpidUtil::Pdg2Pid(pdg);
664 Int_t iBin = FindBin(species, momentum);
b91fdd71 665 AliTRDseedV1 *tracklet = NULL;
1ee39b3a 666 for(Int_t iChamb = 0; iChamb < AliTRDgeometry::kNlayer; iChamb++){
667 tracklet = fkTrack->GetTracklet(iChamb);
668 if(!tracklet) continue;
669 hNClus -> Fill(iBin, tracklet->GetN());
670 }
671
672 return hNClus;
673}
674
675//_______________________________________________________
676TH1 *AliTRDcheckPID::PlotNTracklets(const AliTRDtrackV1 *track)
677{
678 //
679 // Plot the probabilities for electrons using 2-dim LQ
680 //
681
682 if(track) fkTrack = track;
683 if(!fkTrack){
684 AliWarning("No Track defined.");
2e6aa247 685 return NULL;
1ee39b3a 686 }
687
688 TH2F *hTracklets;
689 if(!(hTracklets = dynamic_cast<TH2F *>(fContainer->At(kNTracklets)))){
690 AliWarning("No Histogram defined.");
2e6aa247 691 return NULL;
1ee39b3a 692 }
693
694 AliTRDtrackV1 cTrack(*fkTrack);
695 cTrack.SetReconstructor(fReconstructor);
696 Int_t pdg = 0;
697 Float_t momentum = 0.;
698 if(fkMC){
699 if(fkMC->GetTrackRef()) momentum = fkMC->GetTrackRef()->P();
700 pdg = fkMC->GetPDG();
701 } else {
702 //AliWarning("No MC info available!");
703 momentum = cTrack.GetMomentum(0);
704 pdg = CalcPDG(&cTrack);
705 }
706 Int_t species = AliTRDpidUtil::Pdg2Pid(pdg);
2e6aa247 707 if(!IsInRange(momentum)) return NULL;
1ee39b3a 708
709 Int_t iBin = FindBin(species, momentum);
710 hTracklets -> Fill(iBin, cTrack.GetNumberOfTracklets());
711 return hTracklets;
712}
713
714//_______________________________________________________
715TH1 *AliTRDcheckPID::PlotMom(const AliTRDtrackV1 *track)
716{
717 //
718 // Plot the probabilities for electrons using 2-dim LQ
719 //
720
721 if(track) fkTrack = track;
722 if(!fkTrack){
723 AliWarning("No Track defined.");
2e6aa247 724 return NULL;
1ee39b3a 725 }
726
2e6aa247 727 if(!CheckTrackQuality(fkTrack)) return NULL;
1ee39b3a 728
729 TH1F *hMom;
730 if(!(hMom = dynamic_cast<TH1F *>(fContainer->At(kMomentum)))){
731 AliWarning("No Histogram defined.");
2e6aa247 732 return NULL;
1ee39b3a 733 }
734
735
736 Int_t pdg = 0;
737 Float_t momentum = 0.;
738 if(fkMC){
739 if(fkMC->GetTrackRef()) momentum = fkMC->GetTrackRef()->P();
740 pdg = fkMC->GetPDG();
741 } else {
742 //AliWarning("No MC info available!");
743 AliTRDtrackV1 cTrack(*fkTrack);
744 cTrack.SetReconstructor(fReconstructor);
745 momentum = cTrack.GetMomentum(0);
746 pdg = CalcPDG(&cTrack);
747 }
748 if(IsInRange(momentum)) hMom -> Fill(momentum);
749 return hMom;
750}
751
752
753//_______________________________________________________
754TH1 *AliTRDcheckPID::PlotMomBin(const AliTRDtrackV1 *track)
755{
756 //
757 // Plot the probabilities for electrons using 2-dim LQ
758 //
759
760 if(track) fkTrack = track;
761 if(!fkTrack){
762 AliWarning("No Track defined.");
2e6aa247 763 return NULL;
1ee39b3a 764 }
765
2e6aa247 766 if(!CheckTrackQuality(fkTrack)) return NULL;
1ee39b3a 767
768 TH1F *hMomBin;
769 if(!(hMomBin = dynamic_cast<TH1F *>(fContainer->At(kMomentumBin)))){
770 AliWarning("No Histogram defined.");
2e6aa247 771 return NULL;
1ee39b3a 772 }
773
774
775 Int_t pdg = 0;
776 Float_t momentum = 0.;
777
778 if(fkMC){
779 if(fkMC->GetTrackRef()) momentum = fkMC->GetTrackRef()->P();
780 pdg = fkMC->GetPDG();
781 } else {
782 //AliWarning("No MC info available!");
783 AliTRDtrackV1 cTrack(*fkTrack);
784 cTrack.SetReconstructor(fReconstructor);
785 momentum = cTrack.GetMomentum(0);
786 }
787 if(IsInRange(momentum)) hMomBin -> Fill(fMomentumAxis->FindBin(momentum));
788 return hMomBin;
789}
790
791
792//________________________________________________________
793Bool_t AliTRDcheckPID::GetRefFigure(Int_t ifig)
794{
795// Steering function to retrieve performance plots
796
797 Bool_t kFIRST = kTRUE;
b91fdd71 798 TGraphErrors *g = NULL;
799 TAxis *ax = NULL;
800 TObjArray *arr = NULL;
801 TH1 *h1 = NULL, *h=NULL;
802 TH2 *h2 = NULL;
803 TList *content = NULL;
1ee39b3a 804 switch(ifig){
805 case kEfficiency:{
806 gPad->Divide(2, 1, 1.e-5, 1.e-5);
807 TList *l=gPad->GetListOfPrimitives();
808 TVirtualPad *pad = ((TVirtualPad*)l->At(0));pad->cd();
809 pad->SetMargin(0.1, 0.01, 0.1, 0.01);
810
2e6aa247 811 TLegend *legEff = new TLegend(.64, .84, .98, .98);
812 legEff->SetBorderSize(1);legEff->SetTextSize(0.03255879);
1ee39b3a 813 legEff->SetFillColor(0);
814 h=new TH1S("hEff", "", 1, .5, 11.);
815 h->SetLineColor(1);h->SetLineWidth(1);
816 ax = h->GetXaxis();
817 ax->SetTitle("p [GeV/c]");
818 ax->SetRangeUser(.5, 11.);
819 ax->SetMoreLogLabels();
820 ax = h->GetYaxis();
821 ax->SetTitle("#epsilon_{#pi} [%]");
822 ax->CenterTitle();
823 ax->SetRangeUser(1.e-2, 10.);
824 h->Draw();
2e6aa247 825 content = (TList *)fGraph->FindObject(Form("Eff_%s", AliTRDCalPID::GetPartName(AliPID::kPion)));
1ee39b3a 826 if(!(g = (TGraphErrors*)content->At(AliTRDpidUtil::kLQ))) break;
827 if(!g->GetN()) break;
2e6aa247 828 legEff->SetHeader("PID Method [PION]");
1ee39b3a 829 g->Draw("pc"); legEff->AddEntry(g, "LQ 2D", "pl");
830 if(! (g = (TGraphErrors*)content->At(AliTRDpidUtil::kNN))) break;
831 g->Draw("pc"); legEff->AddEntry(g, "NN", "pl");
832 if(! (g = (TGraphErrors*)content->At(AliTRDpidUtil::kESD))) break;
833 g->Draw("p"); legEff->AddEntry(g, "ESD", "pl");
834 legEff->Draw();
835 gPad->SetLogy();
836 gPad->SetLogx();
837 gPad->SetGridy();
838 gPad->SetGridx();
839
2e6aa247 840
1ee39b3a 841 pad = ((TVirtualPad*)l->At(1));pad->cd();
842 pad->SetMargin(0.1, 0.01, 0.1, 0.01);
843 h=new TH1S("hThr", "", 1, .5, 11.);
844 h->SetLineColor(1);h->SetLineWidth(1);
845 ax = h->GetXaxis();
846 ax->SetTitle("p [GeV/c]");
847 ax->SetMoreLogLabels();
848 ax = h->GetYaxis();
849 ax->SetTitle("Threshold [%]");
850 ax->SetRangeUser(5.e-2, 1.);
851 h->Draw();
2e6aa247 852 content = (TList *)fGraph->FindObject("Thres");
1ee39b3a 853 if(!(g = (TGraphErrors*)content->At(AliTRDpidUtil::kLQ))) break;
854 if(!g->GetN()) break;
855 g->Draw("pc");
856 if(!(g = (TGraphErrors*)content->At(AliTRDpidUtil::kNN))) break;
857 g->Draw("pc");
858 if(!(g = (TGraphErrors*)content->At(AliTRDpidUtil::kESD))) break;
859 g->Draw("p");
860 gPad->SetLogx();
861 gPad->SetGridy();
862 gPad->SetGridx();
863 return kTRUE;
864 }
2e6aa247 865 case kEfficiencyKa:{
866 gPad->Divide(2, 1, 1.e-5, 1.e-5);
867 TList *l=gPad->GetListOfPrimitives();
868 TVirtualPad *pad = ((TVirtualPad*)l->At(0));pad->cd();
869 pad->SetMargin(0.1, 0.01, 0.1, 0.01);
870
871 TLegend *legEff = new TLegend(.64, .84, .98, .98);
872 legEff->SetBorderSize(1);legEff->SetTextSize(0.03255879);
873 legEff->SetFillColor(0);
874 h = (TH1S*)gROOT->FindObject("hEff");
875 h=(TH1S*)h->Clone("hEff_K");
876 h->SetYTitle("#epsilon_{K} [%]");
877 h->GetYaxis()->SetRangeUser(1.e-2, 1.e2);
878 h->Draw();
879 content = (TList *)fGraph->FindObject(Form("Eff_%s", AliTRDCalPID::GetPartName(AliPID::kKaon)));
880 if(!(g = (TGraphErrors*)content->At(AliTRDpidUtil::kLQ))) break;
881 if(!g->GetN()) break;
882 legEff->SetHeader("PID Method [KAON]");
883 g->Draw("pc"); legEff->AddEntry(g, "LQ 2D", "pl");
884 if(! (g = (TGraphErrors*)content->At(AliTRDpidUtil::kNN))) break;
885 g->Draw("pc"); legEff->AddEntry(g, "NN", "pl");
886 if(! (g = (TGraphErrors*)content->At(AliTRDpidUtil::kESD))) break;
887 g->Draw("p"); legEff->AddEntry(g, "ESD", "pl");
888 legEff->Draw();
889 gPad->SetLogy();
890 gPad->SetLogx();
891 gPad->SetGridy();
892 gPad->SetGridx();
893
894 TLegend *legEff2 = new TLegend(.64, .84, .98, .98);
895 legEff2->SetBorderSize(1);legEff2->SetTextSize(0.03255879);
896 legEff2->SetFillColor(0);
897 pad = ((TVirtualPad*)l->At(1));pad->cd();
898 pad->SetMargin(0.1, 0.01, 0.1, 0.01);
899 h=(TH1S*)h->Clone("hEff_p");
900 h->SetYTitle("#epsilon_{p} [%]");
901 h->GetYaxis()->SetRangeUser(1.e-2, 1.e2);
902 h->Draw();
903 content = (TList *)fGraph->FindObject(Form("Eff_%s", AliTRDCalPID::GetPartName(AliPID::kProton)));
904 if(!(g = (TGraphErrors*)content->At(AliTRDpidUtil::kLQ))) break;
905 if(!g->GetN()) break;
906 legEff2->SetHeader("PID Method [PROTON]");
907 g->Draw("pc"); legEff2->AddEntry(g, "LQ 2D", "pl");
908 if(! (g = (TGraphErrors*)content->At(AliTRDpidUtil::kNN))) break;
909 g->Draw("pc"); legEff2->AddEntry(g, "NN", "pl");
910 if(! (g = (TGraphErrors*)content->At(AliTRDpidUtil::kESD))) break;
911 g->Draw("p"); legEff2->AddEntry(g, "ESD", "pl");
912 legEff2->Draw();
913 gPad->SetLogy();
914 gPad->SetLogx();
915 gPad->SetGridy();
916 gPad->SetGridx();
917 return kTRUE;
918 }
1ee39b3a 919 case kdEdx:{
920 // save 2.0 GeV projection as reference
921 TLegend *legdEdx = new TLegend(.7, .7, .98, .98);
922 legdEdx->SetBorderSize(1);
923 kFIRST = kTRUE;
924 if(!(h2 = (TH2F*)(fContainer->At(kdEdx)))) break;
925 legdEdx->SetHeader("Particle Species");
926 gPad->SetMargin(0.1, 0.01, 0.1, 0.01);
927 for(Int_t is = AliPID::kSPECIES-1; is>=0; is--){
928 Int_t bin = FindBin(is, 2.);
929 h1 = h2->ProjectionY(Form("px%d", is), bin, bin);
930 if(!h1->GetEntries()) continue;
931 h1->Scale(1./h1->Integral());
932 h1->SetLineColor(AliTRDCalPID::GetPartColor(is));
933 if(kFIRST){
934 h1->GetXaxis()->SetTitle("dE/dx (a.u.)");
935 h1->GetYaxis()->SetTitle("<Entries>");
936 }
937 h = (TH1F*)h1->DrawClone(kFIRST ? "c" : "samec");
938 legdEdx->AddEntry(h, Form("%s", AliTRDCalPID::GetPartName(is)), "l");
939 kFIRST = kFALSE;
940 }
941 if(kFIRST) break;
942 legdEdx->Draw();
943 gPad->SetLogy();
944 gPad->SetLogx(0);
945 gPad->SetGridy();
946 gPad->SetGridx();
947 return kTRUE;
948 }
949 case kdEdxSlice:
950 break;
951 case kPH:{
952 gPad->Divide(2, 1, 1.e-5, 1.e-5);
953 TList *l=gPad->GetListOfPrimitives();
954
955 // save 2.0 GeV projection as reference
956 TLegend *legPH = new TLegend(.4, .7, .68, .98);
957 legPH->SetBorderSize(1);legPH->SetFillColor(0);
958 legPH->SetHeader("Particle Species");
959 if(!(arr = (TObjArray*)(fContainer->At(kPH)))) break;
960 if(!(h2 = (TProfile2D*)(arr->At(0)))) break;
961
962 TVirtualPad *pad = ((TVirtualPad*)l->At(0));pad->cd();
963 pad->SetMargin(0.1, 0.01, 0.1, 0.01);
964 kFIRST = kTRUE;
965 for(Int_t is=0; is<AliPID::kSPECIES; is++){
966 Int_t bin = FindBin(is, 2.);
967 h1 = h2->ProjectionY(Form("pyt%d", is), bin, bin);
968 if(!h1->GetEntries()) continue;
969 h1->SetMarkerStyle(24);
970 h1->SetMarkerColor(AliTRDCalPID::GetPartColor(is));
971 h1->SetLineColor(AliTRDCalPID::GetPartColor(is));
972 if(kFIRST){
973 h1->GetXaxis()->SetTitle("t_{drift} [100*ns]");
974 h1->GetYaxis()->SetTitle("<dQ/dt> [a.u.]");
975 }
976 h = (TH1F*)h1->DrawClone(kFIRST ? "c" : "samec");
977 legPH->AddEntry(h, Form("%s", AliTRDCalPID::GetPartName(is)), "pl");
978 kFIRST = kFALSE;
979 }
980
981 pad = ((TVirtualPad*)l->At(1));pad->cd();
982 pad->SetMargin(0.1, 0.01, 0.1, 0.01);
983 if(!(h2 = (TProfile2D*)(arr->At(1)))) break;
984 kFIRST = kTRUE;
985 for(Int_t is=0; is<AliPID::kSPECIES; is++){
986 Int_t bin = FindBin(is, 2.);
987 h1 = h2->ProjectionY(Form("pyx%d", is), bin, bin);
988 if(!h1->GetEntries()) continue;
989 h1->SetMarkerStyle(24);
990 h1->SetMarkerColor(AliTRDCalPID::GetPartColor(is));
991 h1->SetLineColor(AliTRDCalPID::GetPartColor(is));
992 if(kFIRST){
993 h1->GetXaxis()->SetTitle("x_{drift} [cm]");
994 h1->GetYaxis()->SetTitle("<dQ/dl> [a.u./cm]");
995 }
996 h = (TH1F*)h1->DrawClone(kFIRST ? "c" : "samec");
997 kFIRST = kFALSE;
998 }
999
1000 if(kFIRST) break;
1001 legPH->Draw();
1002 gPad->SetLogy(0);
1003 gPad->SetLogx(0);
1004 gPad->SetGridy();
1005 gPad->SetGridx();
1006 return kTRUE;
1007 }
1008 case kNClus:{
1009 // save 2.0 GeV projection as reference
1010 TLegend *legNClus = new TLegend(.13, .7, .4, .98);
1011 legNClus->SetBorderSize(1);
1012 legNClus->SetFillColor(0);
1013
1014 kFIRST = kTRUE;
1015 if(!(h2 = (TH2F*)(fContainer->At(kNClus)))) break;
1016 legNClus->SetHeader("Particle Species");
1017 for(Int_t is=0; is<AliPID::kSPECIES; is++){
1018 Int_t bin = FindBin(is, 2.);
1019 h1 = h2->ProjectionY(Form("pyNClus%d", is), bin, bin);
1020 if(!h1->GetEntries()) continue;
1021 h1->Scale(100./h1->Integral());
1022 //h1->SetMarkerStyle(24);
1023 //h1->SetMarkerColor(AliTRDCalPID::GetPartColor(is));
1024 h1->SetLineColor(AliTRDCalPID::GetPartColor(is));
1025 if(kFIRST){
1026 h1->GetXaxis()->SetTitle("N^{cl}/tracklet");
1027 h1->GetYaxis()->SetTitle("Prob. [%]");
1028 h = (TH1F*)h1->DrawClone("c");
1029 h->SetMaximum(55.);
1030 h->GetXaxis()->SetRangeUser(0., 35.);
1031 kFIRST = kFALSE;
1032 } else h = (TH1F*)h1->DrawClone("samec");
1033
1034 legNClus->AddEntry(h, Form("%s", AliTRDCalPID::GetPartName(is)), "l");
1035 }
1036 if(kFIRST) break;
1037 legNClus->Draw();
1038 gPad->SetLogy(0);
1039 gPad->SetLogx(0);
1040 gPad->SetGridy();
1041 gPad->SetGridx();
1042 return kTRUE;
1043 }
1044 case kMomentum:
1045 case kMomentumBin:
1046 break;
1047 case kNTracklets:{
1048 TLegend *legNClus = new TLegend(.4, .7, .68, .98);
1049 legNClus->SetBorderSize(1);
1050 kFIRST = kTRUE;
1051 if(!(h2 = (TH2F*)(fContainer->At(kNTracklets)))) break;
1052 legNClus->SetHeader("Particle Species");
1053 for(Int_t is=0; is<AliPID::kSPECIES; is++){
1054 Int_t bin = FindBin(is, 2.);
1055 h1 = h2->ProjectionY(Form("pyNTracklets%d", is), bin, bin);
1056 if(!h1->GetEntries()) continue;
1057 h1->Scale(100./h1->Integral());
1058 //h1->SetMarkerStyle(24);
1059 //h1->SetMarkerColor(AliTRDCalPID::GetPartColor(is));
1060 h1->SetLineColor(AliTRDCalPID::GetPartColor(is));
1061 if(kFIRST){
1062 h1->GetXaxis()->SetTitle("N^{trklt}/track");
1063 h1->GetXaxis()->SetRangeUser(1.,6.);
1064 h1->GetYaxis()->SetTitle("Prob. [%]");
1065 }
1066 h = (TH1F*)h1->DrawClone(kFIRST ? "c" : "samec");
1067 legNClus->AddEntry(h, Form("%s", AliTRDCalPID::GetPartName(is)), "l");
1068 kFIRST = kFALSE;
1069 }
1070 if(kFIRST) break;
1071 legNClus->Draw();
1072 gPad->SetLogy(0);
1073 gPad->SetLogx(0);
1074 gPad->SetGridy();
1075 gPad->SetGridx();
1076 return kTRUE;
1077 }
1078 }
1079 AliInfo(Form("Reference plot [%d] missing result", ifig));
1080 return kFALSE;
1081}
1082
1083//________________________________________________________________________
1084Bool_t AliTRDcheckPID::PostProcess()
1085{
1086 // Draw result to the screen
1087 // Called once at the end of the query
1088
1089 if (!fContainer) {
1090 Printf("ERROR: list not available");
1091 return kFALSE;
1092 }
2e6aa247 1093
1094 TObjArray *eff(NULL);
1ee39b3a 1095 if(!fGraph){
2e6aa247 1096 fGraph = new TObjArray(2*AliPID::kSPECIES);
1ee39b3a 1097 fGraph->SetOwner();
2e6aa247 1098
1099 if(!(eff = dynamic_cast<TObjArray *>(fContainer->At(kEfficiency)))){
1100 AliError("Efficiency container for Electrons missing.");
1101 return kFALSE;
1102 }
1103 EvaluateEfficiency(eff, fGraph, AliPID::kPion, 0.9);
1104 EvaluateEfficiency(eff, fGraph, AliPID::kKaon, 0.9);
1105 EvaluateEfficiency(eff, fGraph, AliPID::kProton, 0.9);
1ee39b3a 1106 }
2e6aa247 1107 fNRefFigures = 12;
1ee39b3a 1108 return kTRUE;
1109}
1110
1111//________________________________________________________________________
2e6aa247 1112void AliTRDcheckPID::EvaluateEfficiency(TObjArray * const histoContainer, TObjArray *results, Int_t species, Float_t electronEfficiency){
1ee39b3a 1113// Process PID information for pion efficiency
1114
1115 fUtil->SetElectronEfficiency(electronEfficiency);
1116
2e6aa247 1117
1118 const Int_t kNmethodsPID=Int_t(sizeof(fgMethod)/sizeof(Char_t*));
1119 Color_t colors[kNmethodsPID] = {kBlue, kGreen+2, kRed};
1120 Int_t markerStyle[kNmethodsPID] = {7, 7, 24};
1ee39b3a 1121 // efficiency graphs
2e6aa247 1122 TGraphErrors *g(NULL);
1123 TObjArray *eff = new TObjArray(kNmethodsPID); eff->SetOwner(); eff->SetName(Form("Eff_%s", AliTRDCalPID::GetPartName(species)));
1124 results->AddAt(eff, species);
1125 for(Int_t iMethod = 0; iMethod < kNmethodsPID; iMethod++){
1126 eff->AddAt(g = new TGraphErrors(), iMethod);
1127 g->SetName(Form("%s", fgMethod[iMethod]));
1ee39b3a 1128 g->SetLineColor(colors[iMethod]);
1129 g->SetMarkerColor(colors[iMethod]);
1130 g->SetMarkerStyle(markerStyle[iMethod]);
1131 }
1132
2e6aa247 1133 // Threshold graphs if not already
1134 TObjArray *thres(NULL);
1135 if(!(results->At(AliPID::kSPECIES))){
1136 thres = new TObjArray(kNmethodsPID); thres->SetOwner();
1137 thres->SetName("Thres");
1138 results->AddAt(thres, AliPID::kSPECIES);
1139 for(Int_t iMethod = 0; iMethod < kNmethodsPID; iMethod++){
1140 thres->AddAt(g = new TGraphErrors(), iMethod);
1141 g->SetName(Form("%s", fgMethod[iMethod]));
1142 g->SetLineColor(colors[iMethod]);
1143 g->SetMarkerColor(colors[iMethod]);
1144 g->SetMarkerStyle(markerStyle[iMethod]);
1145 }
1ee39b3a 1146 }
1ee39b3a 1147
2e6aa247 1148 TH2F *hPtr[kNmethodsPID]={
1149 (TH2F*)histoContainer->At(AliTRDpidUtil::kLQ),
1150 (TH2F*)histoContainer->At(AliTRDpidUtil::kNN),
1151 (TH2F*)histoContainer->At(AliTRDpidUtil::kESD)
1152 };
1ee39b3a 1153 for(Int_t iMom = 0; iMom < fMomentumAxis->GetNbins(); iMom++){
2e6aa247 1154 Float_t mom(fMomentumAxis->GetBinCenter(iMom+1));
1ee39b3a 1155
2e6aa247 1156 Int_t binEl(fMomentumAxis->GetNbins() * AliPID::kElectron + iMom + 1),
1157 binXX(fMomentumAxis->GetNbins() * species + iMom + 1);
1158
1159 for(Int_t iMethod = 0; iMethod < kNmethodsPID; iMethod++){
1160 // Calculate the Species Efficiency at electronEfficiency% electron efficiency for each Method
1161
1162 TH1D *histo1 = hPtr[iMethod] -> ProjectionY(Form("%s_el", fgMethod[iMethod]), binEl, binEl);
1163 TH1D *histo2 = hPtr[iMethod] -> ProjectionY(Form("%s_%s", fgMethod[iMethod], AliTRDCalPID::GetPartName(species)), binXX, binXX);
1ee39b3a 1164
1165 if(!fUtil->CalculatePionEffi(histo1, histo2)) continue;
1166
2e6aa247 1167 g=(TGraphErrors*)eff->At(iMethod);
1168 g->SetPoint(iMom, mom, 1.e2*fUtil->GetPionEfficiency());
1169 g->SetPointError(iMom, 0., 1.e2*fUtil->GetError());
1170 AliDebug(2, Form("%s Efficiency for %s is : %f +/- %f", AliTRDCalPID::GetPartName(species), fgMethod[iMethod], fUtil->GetPionEfficiency(), fUtil->GetError()));
1171
1172 if(!thres) continue;
1173 g=(TGraphErrors*)thres->At(iMethod);
1174 g->SetPoint(iMom, mom, fUtil->GetThreshold());
1175 g->SetPointError(iMom, 0., 0.);
1ee39b3a 1176 }
1177 }
1178}