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 /////////////////////////////////////////////////////
18 // Check basic detector results at ESD level
19 // - Geometrical efficiency
20 // - Tracking efficiency
25 // Alex Bercuci <A.Bercuci@gsi.de>
27 //////////////////////////////////////////////////////
29 #include <TClonesArray.h>
30 #include <TObjArray.h>
36 #include <TGraphErrors.h>
37 #include <TGraphAsymmErrors.h>
42 #include <TParticle.h>
45 #include "AliAnalysisManager.h"
46 #include "AliESDEvent.h"
47 #include "AliESDkink.h"
48 #include "AliMCEvent.h"
49 #include "AliESDInputHandler.h"
50 #include "AliMCEventHandler.h"
52 #include "AliESDtrack.h"
53 #include "AliMCParticle.h"
56 #include "AliTrackReference.h"
58 #include "AliTRDcheckESD.h"
60 ClassImp(AliTRDcheckESD)
62 const Float_t AliTRDcheckESD::fgkxTPC = 290.;
63 const Float_t AliTRDcheckESD::fgkxTOF = 365.;
64 const UChar_t AliTRDcheckESD::fgkNgraph[AliTRDcheckESD::kNrefs] ={
66 FILE* AliTRDcheckESD::fgFile = NULL;
68 const Float_t AliTRDcheckESD::fgkEvVertexZ = 15.;
69 const Int_t AliTRDcheckESD::fgkEvVertexN = 1;
70 const Float_t AliTRDcheckESD::fgkTrkDCAxy = 40.;
71 const Float_t AliTRDcheckESD::fgkTrkDCAz = 15.;
72 const Int_t AliTRDcheckESD::fgkNclTPC = 100;
73 const Float_t AliTRDcheckESD::fgkPt = 0.2;
74 const Float_t AliTRDcheckESD::fgkEta = 0.9;
76 //____________________________________________________________________
77 AliTRDcheckESD::AliTRDcheckESD():
87 // Default constructor
89 SetNameTitle("checkESD", "Check TRD @ ESD level");
93 //____________________________________________________________________
94 AliTRDcheckESD::AliTRDcheckESD(char* name):
95 AliAnalysisTaskSE(name)
104 // Default constructor
107 SetTitle("Check TRD @ ESD level");
108 DefineOutput(1, TObjArray::Class());
111 //____________________________________________________________________
112 AliTRDcheckESD::~AliTRDcheckESD()
125 //____________________________________________________________________
126 void AliTRDcheckESD::UserCreateOutputObjects()
129 // Create Output Containers (TObjectArray containing 1D histograms)
131 //OpenFile(0, "RECREATE");
137 //____________________________________________________________________
138 Bool_t AliTRDcheckESD::GetRefFigure(Int_t ifig)
140 if(ifig>=fNRefFigures){
141 AliWarning(Form("Ref plot %d not available. Valid only up to %d", ifig, fNRefFigures));
145 AliWarning("Please provide a canvas to draw results.");
148 gPad->SetLogx(0);gPad->SetLogy(0);
149 gPad->SetMargin(0.125, 0.015, 0.1, 0.015);
152 const Char_t *title[20];
154 if((hF=(TH1S*)gROOT->FindObject("hFcheckESD"))) delete hF;
156 TList *l(NULL); TVirtualPad *pad(NULL);
157 TGraphErrors *g(NULL);TGraphAsymmErrors *ga(NULL);
158 TObjArray *arr(NULL);
160 case kNCl: // number of clusters/track
161 if(!(arr = (TObjArray*)fResults->At(kNCl))) return kFALSE;
163 leg = new TLegend(.83, .7, .99, .96);
164 leg->SetHeader("Species");
165 leg->SetBorderSize(0); leg->SetFillStyle(0);
166 for(Int_t ig(0); ig<fgkNgraph[kNCl]; ig++){
167 if(!(g = (TGraphErrors*)arr->At(ig))) return kFALSE;
168 g->Draw(ig?"pc":"apc"); leg->AddEntry(g, g->GetTitle(), "pl");
170 hF=g->GetHistogram();
171 hF->SetXTitle("no of clusters");
172 hF->SetYTitle("entries");
173 hF->GetYaxis()->CenterTitle(1);
174 hF->GetYaxis()->SetTitleOffset(1.2);
177 leg->Draw(); gPad->SetLogy();
179 case kTRDstat: // Efficiency
180 if(!(arr = (TObjArray*)fResults->At(kTRDstat))) return kFALSE;
181 leg = new TLegend(.62, .77, .98, .98);
182 leg->SetHeader("TRD Efficiency");
183 leg->SetBorderSize(0); leg->SetFillStyle(0);
184 title[0] = "Geometrical (TRDin/TPCout)";
185 title[1] = "Tracking (TRDout/TRDin)";
186 title[2] = "PID (TRDpid/TRDin)";
187 title[3] = "Refit (TRDrefit/TRDin)";
188 hF = new TH1S("hFcheckESD", ";p [GeV/c];Efficiency", 10, 0.1, 10.);
190 hF->GetXaxis()->SetMoreLogLabels();
191 hF->GetYaxis()->CenterTitle(1);
193 for(Int_t ig(0); ig<fgkNgraph[kTRDstat]; ig++){
194 if(!(g = (TGraphErrors*)arr->At(ig))) return kFALSE;
195 g->Draw("pl"); leg->AddEntry(g, title[ig], "pl");
196 //PutTrendValue(name[id], g->GetMean(2));
197 //PutTrendValue(Form("%sRMS", name[id]), g->GetRMS(2));
199 leg->Draw(); gPad->SetLogx();
201 case kTRDmom: // Energy loss
202 if(!(arr = (TObjArray*)fResults->At(kTRDmom))) return kFALSE;
203 leg = new TLegend(.65, .7, .95, .99);
204 leg->SetHeader("Energy Loss");
205 leg->SetBorderSize(1); leg->SetFillColor(0);
206 title[0] = "Max & 90% quantile";
207 title[1] = "Mean & 60% quantile";
208 hF = new TH1S("hFcheckESD", ";layer;#Delta E", 6, -0.5, 5.5);
209 hF->SetMaximum(1.3);hF->SetMinimum(-.3);
211 for(Int_t ig(0); ig<fgkNgraph[kTRDmom]; ig++){
212 if(!(ga = (TGraphAsymmErrors*)arr->At(ig))) return kFALSE;
213 ga->Draw("pl"); leg->AddEntry(ga, title[ig], "pl");
214 //PutTrendValue(name[id], g->GetMean(2));
215 //PutTrendValue(Form("%sRMS", name[id]), g->GetRMS(2));
217 leg->Draw();gPad->SetLogx(kFALSE);
219 case kPtRes: // Pt resolution @ vertex
220 if(!(arr = (TObjArray*)fResults->At(kPtRes))) return kFALSE;
221 gPad->Divide(2, 1, 1.e-5, 1.e-5); l=gPad->GetListOfPrimitives();
222 pad = ((TVirtualPad*)l->At(0)); pad->cd(); pad->SetLogx();
223 pad->SetMargin(0.1, 0.022, 0.1, 0.023);
224 hF = new TH1S("hFcheckESD", "ITS+TPC+TRD;p_{t} [GeV/c];#Delta p_{t} / p_{t} [%]", 10, 0.2, 10.);
225 hF->SetMaximum(10.);hF->SetMinimum(-3.);
226 hF->GetXaxis()->SetMoreLogLabels();
227 hF->GetXaxis()->SetTitleOffset(1.2);
228 hF->GetYaxis()->CenterTitle();
230 //for(Int_t ig(0); ig<fgkNgraph[kPtRes]/2; ig++){
231 for(Int_t ig(2); ig<6; ig++){
232 if(!(g = (TGraphErrors*)arr->At(ig))) continue;
233 if(!g->GetN()) continue;
235 //PutTrendValue(name[id], g->GetMean(2));
236 //PutTrendValue(Form("%sRMS", name[id]), g->GetRMS(2));
238 pad = ((TVirtualPad*)l->At(1)); pad->cd(); pad->SetLogx();
239 pad->SetMargin(0.1, 0.22, 0.1, 0.023);
240 hF = (TH1*)hF->Clone("hFcheckESD1");
241 hF->SetTitle("ITS+TPC");
242 hF->SetMaximum(10.);hF->SetMinimum(-3.);
244 leg = new TLegend(.78, .1, .99, .98);
245 leg->SetHeader("P_{t} @ DCA");
246 leg->SetBorderSize(1); leg->SetFillColor(0);
247 leg->SetTextAlign(22);
248 leg->SetTextFont(12);
249 leg->SetTextSize(0.03813559);
252 //for(Int_t ig(fgkNgraph[kPtRes]/2); ig<fgkNgraph[kPtRes]; ig++){
253 for(Int_t ig(12); ig<16; ig++){
254 if(!(g = (TGraphErrors*)arr->At(ig))) continue;
255 if(!g->GetN()) continue;
257 g->Draw("pl"); leg->AddEntry(g, g->GetTitle(), "pl");
258 //PutTrendValue(name[id], g->GetMean(2));
259 //PutTrendValue(Form("%sRMS", name[id]), g->GetRMS(2));
261 if(nPlots) leg->Draw();
269 //____________________________________________________________________
270 void AliTRDcheckESD::UserExec(Option_t *){
274 fESD = dynamic_cast<AliESDEvent*>(InputEvent());
278 AliError("ESD event missing.");
282 // Get MC information if available
283 AliStack * fStack = NULL;
286 AliWarning("MC event missing");
289 if(!(fStack = fMC->Stack())){
290 AliWarning("MC stack missing");
297 AliESDtrack *esdTrack(NULL);
298 for(Int_t itrk = 0; itrk < fESD->GetNumberOfTracks(); itrk++){
299 esdTrack = fESD->GetTrack(itrk);
302 ULong_t status = esdTrack->GetStatus(); //PrintStatus(status);
305 Bool_t selected(kTRUE);
306 if(esdTrack->Pt() < fgkPt){
307 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Pt[%5.2f]", fESD->GetEventNumberInFile(), itrk, esdTrack->Pt()));
310 if(TMath::Abs(esdTrack->Eta()) > fgkEta){
311 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Eta[%5.2f]", fESD->GetEventNumberInFile(), itrk, TMath::Abs(esdTrack->Eta())));
314 if(!Bool_t(status & AliESDtrack::kTPCout)){
315 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] !TPCout", fESD->GetEventNumberInFile(), itrk));
318 if(esdTrack->GetKinkIndex(0) > 0){
319 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] Kink", fESD->GetEventNumberInFile(), itrk));
322 if(esdTrack->GetTPCNcls() < fgkNclTPC){
323 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] NclTPC[%d]", fESD->GetEventNumberInFile(), itrk, esdTrack->GetTPCNcls()));
326 Float_t par[2], cov[3];
327 esdTrack->GetImpactParameters(par, cov);
328 if(IsCollision()){ // cuts on DCA
329 if(TMath::Abs(par[0]) > fgkTrkDCAxy){
330 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] DCAxy[%f]", fESD->GetEventNumberInFile(), itrk, TMath::Abs(par[0])));
333 if(TMath::Abs(par[1]) > fgkTrkDCAz){
334 AliDebug(2, Form("Reject Ev[%4d] Trk[%3d] DCAz[%f]", fESD->GetEventNumberInFile(), itrk, TMath::Abs(par[1])));
338 if(!selected) continue;
340 //Int_t nTPC(esdTrack->GetNcls(1));
341 Int_t nTRD(esdTrack->GetNcls(2));
342 Double_t pt(esdTrack->Pt());
343 //Double_t eta(esdTrack->Eta());
344 //Double_t phi(esdTrack->Phi());
345 Double_t p[AliPID::kSPECIES]; esdTrack->GetTRDpid(p);
347 //esdTrack->GetTRDntrackletsPID();
348 Bool_t kBarrel = Bool_t(status & AliESDtrack::kTRDin);
350 // look at external track param
351 const AliExternalTrackParam *op = esdTrack->GetOuterParam();
352 const AliExternalTrackParam *ip = esdTrack->GetInnerParam();
354 Double_t pt0(0.), eta0(0.), phi0(0.), ptTRD(0.);
355 // read MC info if available
356 Bool_t kFOUND(kFALSE), kPhysPrim(kFALSE);
357 AliMCParticle *mcParticle(NULL);
359 AliTrackReference *ref(NULL);
360 Int_t fLabel(esdTrack->GetLabel());
361 Int_t fIdx(TMath::Abs(fLabel));
362 if(fIdx > fStack->GetNtrack()) continue;
365 if(!(mcParticle = (AliMCParticle*) fMC->GetTrack(fIdx))) {
366 AliWarning(Form("MC particle missing. Label[ %d].", fLabel));
369 pt0 = mcParticle->Pt();
370 eta0 = mcParticle->Eta();
371 phi0 = mcParticle->Phi();
372 kPhysPrim = fMC->IsPhysicalPrimary(fIdx);
374 // read track references
375 Int_t nRefs = mcParticle->GetNumberOfTrackReferences();
377 AliWarning(Form("No TR found for track @ Label[%d].", fLabel));
382 ref = mcParticle->GetTrackReference(iref);
383 if(ref->LocalX() > fgkxTPC) break;
387 if(ref->LocalX() > fgkxTOF){ // track skipping TRD fiducial volume
388 ref = mcParticle->GetTrackReference(TMath::Max(iref-1, 0));
390 } else { // track stopped in TPC
391 ref = mcParticle->GetTrackReference(TMath::Max(iref-1, 0));
393 ptTRD = ref->Pt();kFOUND=kTRUE;
394 } else { // use reconstructed values
396 Double_t x(op->GetX());
397 if(x<fgkxTOF && x>fgkxTPC){
410 h = (TH2I*)fHistos->At(kTRDstat);
411 if(status & AliESDtrack::kTPCout) h->Fill(ptTRD, kTPCout);
412 if(status & AliESDtrack::kTRDin) h->Fill(ptTRD, kTRDin);
413 if(kBarrel && (status & AliESDtrack::kTRDout)) h->Fill(ptTRD, kTRDout);
414 if(kBarrel && (status & AliESDtrack::kTRDpid)) h->Fill(ptTRD, kTRDpid);
415 if(kBarrel && (status & AliESDtrack::kTRDrefit)) h->Fill(ptTRD, kTRDref);
417 Int_t idx(HasMC() ? Pdg2Idx(TMath::Abs(mcParticle->PdgCode())): 0)
418 ,sgn(esdTrack->Charge()<0?0:1);
419 if(kBarrel && kPhysPrim) {
420 TH3 *h3 = (TH3S*)fHistos->At(kPtRes);
421 Int_t offset = (status & AliESDtrack::kTRDrefit) ? 0 : 10;
422 h3->Fill(pt0, 1.e2*(pt/pt0-1.),
423 offset + 2*idx + sgn);
425 ((TH1*)fHistos->At(kNCl))->Fill(nTRD, 2*idx + sgn);
427 h = (TH2I*)fHistos->At(kTRDmom);
429 for(Int_t ily=6; ily--;){
430 if((pTRD=esdTrack->GetTRDmomentum(ily))<0.) continue;
431 h->Fill(ip->GetP()-pTRD, ily);
435 PostData(1, fHistos);
438 //____________________________________________________________________
439 TObjArray* AliTRDcheckESD::Histos()
441 // Retrieve histograms array if already build or build it
443 if(fHistos) return fHistos;
445 fHistos = new TObjArray(kNhistos);
446 //fHistos->SetOwner(kTRUE);
450 // clusters per track
451 const Int_t kNpt(30);
453 Float_t binsPt[kNpt+1];
454 for(Int_t i=0;i<kNpt+1; i++,Pt+=(TMath::Exp(i*i*.001)-1.)) binsPt[i]=Pt;
455 if(!(h = (TH2S*)gROOT->FindObject("hNCl"))){
456 h = new TH2S("hNCl", "Clusters per TRD track;N_{cl}^{TRD};SPECIES;entries", 60, 0., 180., 10, -0.5, 9.5);
457 TAxis *ay(h->GetYaxis());
458 ay->SetLabelOffset(0.015);
459 for(Int_t i(0); i<AliPID::kSPECIES; i++){
460 ay->SetBinLabel(2*i+1, Form("%s^{-}", AliPID::ParticleLatexName(i)));
461 ay->SetBinLabel(2*i+2, Form("%s^{+}", AliPID::ParticleLatexName(i)));
464 fHistos->AddAt(h, kNCl); fNRefFigures++;
466 // status bits histogram
467 const Int_t kNbits(5);
469 Float_t binsBits[kNbits+1];
470 for(Int_t i=0; i<kNbits+1; i++,Bits+=1.) binsBits[i]=Bits;
471 if(!(h = (TH2I*)gROOT->FindObject("hTRDstat"))){
472 h = new TH2I("hTRDstat", "TRD status bits;p_{t} @ TRD [GeV/c];status;entries", kNpt, binsPt, kNbits, binsBits);
473 TAxis *ay(h->GetYaxis());
474 ay->SetBinLabel(1, "kTPCout");
475 ay->SetBinLabel(2, "kTRDin");
476 ay->SetBinLabel(3, "kTRDout");
477 ay->SetBinLabel(4, "kTRDpid");
478 ay->SetBinLabel(5, "kTRDrefit");
480 fHistos->AddAt(h, kTRDstat);
483 if(!(h = (TH2I*)gROOT->FindObject("hTRDmom"))){
484 h = new TH2I("hTRDmom", "TRD energy loss;p_{inner} - p_{ly} [GeV/c];ly;entries", 100, -1., 2., 6, -0.5, 5.5);
486 fHistos->AddAt(h, kTRDmom);
487 if(!HasMC()) return fHistos;
490 const Int_t kNdpt(100), kNspec(4*AliPID::kSPECIES);
491 Float_t DPt(-3.), Spec(-0.5);
492 Float_t binsDPt[kNdpt+1], binsSpec[kNspec+1];
493 for(Int_t i=0; i<kNdpt+1; i++,DPt+=6.e-2) binsDPt[i]=DPt;
494 for(Int_t i=0; i<kNspec+1; i++,Spec+=1.) binsSpec[i]=Spec;
495 if(!(h = (TH3S*)gROOT->FindObject("hPtRes"))){
496 h = new TH3S("hPtRes", "P_{t} resolution @ DCA;p_{t}^{MC} [GeV/c];#Delta p_{t}/p_{t}^{MC} [%];SPECIES", kNpt, binsPt, kNdpt, binsDPt, kNspec, binsSpec);
497 TAxis *az(h->GetZaxis());
498 az->SetLabelOffset(0.015);
499 for(Int_t i(0); i<AliPID::kSPECIES; i++){
500 az->SetBinLabel(2*i+1, Form("%s^{-}", AliPID::ParticleLatexName(i)));
501 az->SetBinLabel(2*i+2, Form("%s^{+}", AliPID::ParticleLatexName(i)));
502 az->SetBinLabel(10+2*i+1, Form("%s^{-}", AliPID::ParticleLatexName(i)));
503 az->SetBinLabel(10+2*i+2, Form("%s^{+}", AliPID::ParticleLatexName(i)));
506 fHistos->AddAt(h, kPtRes);
511 //____________________________________________________________________
512 Bool_t AliTRDcheckESD::Load(const Char_t *filename, const Char_t *name)
514 // Load data from performance file
516 if(!TFile::Open(filename)){
517 AliWarning(Form("Couldn't open file %s.", filename));
521 if(!(o = (TObjArray*)gFile->Get(name ? name : GetName()))){
522 AliWarning("Missing histogram container.");
525 fHistos = (TObjArray*)o->Clone(GetName());
527 SETBIT(fStatus, kLoad);
531 //_______________________________________________________
532 Bool_t AliTRDcheckESD::PutTrendValue(const Char_t *name, Double_t val)
534 // Dump trending value to default file
537 fgFile = fopen("TRD.Performance.txt", "at");
539 fprintf(fgFile, "%s_%s %f\n", GetName(), name, val);
543 //____________________________________________________________________
544 void AliTRDcheckESD::Terminate(Option_t *)
546 // Steer post-processing
548 fHistos = dynamic_cast<TObjArray *>(GetOutputData(1));
550 AliError("Histogram container not found in output");
555 const Char_t *name[kNrefs] = {
556 "Ncl", "Eff", "Eloss", "PtResDCA"
558 TObjArray *arr(NULL); TGraph *g(NULL);
560 fResults = new TObjArray(kNrefs);
561 fResults->SetOwner();
562 fResults->SetName("results");
563 for(Int_t iref(0); iref<kNrefs; iref++){
564 fResults->AddAt(arr = new TObjArray(fgkNgraph[iref]), iref);
565 arr->SetName(name[iref]); arr->SetOwner();
568 for(Int_t ig(0); ig<fgkNgraph[iref]; ig++){
569 arr->AddAt(g = new TGraphErrors(), ig);
570 g->SetLineColor(ig+1);
571 g->SetMarkerColor(ig+1);
572 g->SetMarkerStyle(ig+20);
573 g->SetNameTitle(Form("s%d", ig), ig? AliPID::ParticleLatexName(ig-1):"all");
577 for(Int_t ig(0); ig<fgkNgraph[iref]; ig++){
578 arr->AddAt(g = new TGraphAsymmErrors(), ig);
579 g->SetLineColor(ig+1);
580 g->SetMarkerColor(ig+1);
581 g->SetMarkerStyle(ig+20);
585 for(Int_t idx(0); idx<AliPID::kSPECIES; idx++){
587 arr->AddAt(g = new TGraphErrors(), ig);
588 g->SetLineColor(kRed-idx);
589 g->SetMarkerColor(kRed-idx);
590 g->SetMarkerStyle(20+idx);
591 g->SetNameTitle(Form("s%d", ig), Form("res %s", AliPID::ParticleLatexName(idx)));
592 arr->AddAt(g = new TGraphErrors(), ig+1);
593 g->SetLineColor(kBlue-idx);
594 g->SetMarkerColor(kBlue-idx);
595 g->SetMarkerStyle(20+idx);
596 g->SetNameTitle(Form("m%d", ig+1), Form("sys %s", AliPID::ParticleLatexName(idx)));
599 arr->AddAt(g = new TGraphErrors(), ig);
600 g->SetLineColor(kRed-idx);
601 g->SetMarkerColor(kRed-idx);
602 g->SetMarkerStyle(20+idx);
603 g->SetNameTitle(Form("s%d", ig), Form("sigma %s", AliPID::ParticleLatexName(idx)));
604 arr->AddAt(g = new TGraphErrors(), ig+1);
605 g->SetLineColor(kBlue-idx);
606 g->SetMarkerColor(kBlue-idx);
607 g->SetMarkerStyle(20+idx);
608 g->SetNameTitle(Form("m%d", ig+1), Form("mean %s", AliPID::ParticleLatexName(idx)));
612 for(Int_t ig(0); ig<fgkNgraph[iref]; ig++){
613 arr->AddAt(g = new TGraphErrors(), ig);
614 g->SetLineColor(ig+1);
615 g->SetMarkerColor(ig+1);
616 g->SetMarkerStyle(ig+20);
622 TH1 *h1[2] = {NULL, NULL};
627 if(!(h2 = (TH2I*)fHistos->At(kNCl))) return;
629 arr = (TObjArray*)fResults->At(kNCl);
630 h1[0] = h2->ProjectionX("Ncl_px");
631 TGraphErrors *ge=(TGraphErrors*)arr->At(0);
632 for(Int_t ib=2; ib<=ax->GetNbins(); ib++){
633 ge->SetPoint(ib-2, ax->GetBinCenter(ib), h1[0]->GetBinContent(ib));
635 for(Int_t is(1); is<=AliPID::kSPECIES; is++){
636 h1[0] = h2->ProjectionX("Ncl_px", 2*is-1, 2*is);
637 ge=(TGraphErrors*)arr->At(is);
638 for(Int_t ib=2; ib<=ax->GetNbins(); ib++){
639 ge->SetPoint(ib-2, ax->GetBinCenter(ib), h1[0]->GetBinContent(ib));
645 // geometrical efficiency
646 if(!(h2 = (TH2I*)fHistos->At(kTRDstat))) return;
647 arr = (TObjArray*)fResults->At(kTRDstat);
648 h1[0] = h2->ProjectionX("checkESDx0", kTPCout, kTPCout);
649 h1[1] = h2->ProjectionX("checkESDx1", kTRDin, kTRDin);
650 Process(h1, (TGraphErrors*)arr->At(0));
651 delete h1[0];delete h1[1];
652 // tracking efficiency
653 h1[0] = h2->ProjectionX("checkESDx0", kTRDin, kTRDin);
654 h1[1] = h2->ProjectionX("checkESDx1", kTRDout, kTRDout);
655 Process(h1, (TGraphErrors*)arr->At(1));
658 h1[1] = h2->ProjectionX("checkESDx1", kTRDpid, kTRDpid);
659 Process(h1, (TGraphErrors*)arr->At(2));
662 h1[1] = h2->ProjectionX("checkESDx1", kTRDref, kTRDref);
663 Process(h1, (TGraphErrors*)arr->At(3));
668 if(!(h2 = dynamic_cast<TH2I*>(fHistos->At(kTRDmom)))) return;
669 arr = (TObjArray*)fResults->At(kTRDmom);
670 TGraphAsymmErrors *g06 = (TGraphAsymmErrors*)arr->At(0), *g09 = (TGraphAsymmErrors*)arr->At(1);
673 const Double_t xq[nq] = {0.05, 0.2, 0.8, 0.95};
675 for(Int_t ily=6; ily--;){
676 h1[0] = h2->ProjectionX("checkESDp0", ily+1, ily+1);
677 h1[0]->GetQuantiles(nq,yq,xq);
678 g06->SetPoint(ily, Float_t(ily), ax->GetBinCenter(h1[0]->GetMaximumBin()));
679 g06->SetPointError(ily, 0., 0., TMath::Abs(yq[0]), yq[3]);
680 g09->SetPoint(ily, Float_t(ily), h1[0]->GetMean());
681 g09->SetPointError(ily, 0., 0., TMath::Abs(yq[1]), yq[2]);
683 //printf(" max[%f] mean[%f] q[%f %f %f %f]\n", ax->GetBinCenter(h1[0]->GetMaximumBin()), h1[0]->GetMean(), yq[0], yq[1], yq[2], yq[3]);
689 // Pt RESOLUTION @ DCA
690 TH3S* h3(NULL); TGraphErrors *gg[2] = {NULL,NULL};
691 if(!(h3 = dynamic_cast<TH3S*>(fHistos->At(kPtRes)))) return;
692 arr = (TObjArray*)fResults->At(kPtRes);
693 TAxis *az(h3->GetZaxis());
694 for(Int_t i(0); i<AliPID::kSPECIES; i++){
696 az->SetRange(idx+1, idx+2);
697 gg[1] = (TGraphErrors*)arr->At(idx);
698 gg[0] = (TGraphErrors*)arr->At(idx+1);
699 Process2D((TH2*)h3->Project3D("yx"), gg);
702 az->SetRange(idx+1, idx+2);
703 gg[1] = (TGraphErrors*)arr->At(idx);
704 gg[0] = (TGraphErrors*)arr->At(idx+1);
705 Process2D((TH2*)h3->Project3D("yx"), gg);
710 //____________________________________________________________________
711 Int_t AliTRDcheckESD::Pdg2Idx(Int_t pdg)
714 case kElectron: return AliPID::kElectron;
715 case kMuonMinus: return AliPID::kMuon;
716 case kPiPlus: return AliPID::kPion;
717 case kKPlus: return AliPID::kKaon;
718 case kProton: return AliPID::kProton;
723 //____________________________________________________________________
724 void AliTRDcheckESD::Process(TH1 **h1, TGraphErrors *g)
726 // Generic function to process one reference plot
728 Int_t n1 = 0, n2 = 0, ip=0;
731 TAxis *ax = h1[0]->GetXaxis();
732 for(Int_t ib=1; ib<=ax->GetNbins(); ib++){
733 if(!(n1 = (Int_t)h1[0]->GetBinContent(ib))) continue;
734 n2 = (Int_t)h1[1]->GetBinContent(ib);
735 eff = n2/Float_t(n1);
738 g->SetPoint(ip, ax->GetBinCenter(ib), eff);
739 g->SetPointError(ip, 0., n2 ? eff*TMath::Sqrt(1./n1+1./n2) : 0.);
742 //________________________________________________________
743 void AliTRDcheckESD::Process2D(TH2 * const h2, TGraphErrors **g)
750 if((n=g[0]->GetN())) for(;n--;) g[0]->RemovePoint(n);
751 if((n=g[1]->GetN())) for(;n--;) g[1]->RemovePoint(n);
752 TF1 f("fg", "gaus", -3.,3.);
753 for(Int_t ibin = 1; ibin <= h2->GetNbinsX(); ibin++){
754 Double_t x = h2->GetXaxis()->GetBinCenter(ibin);
755 TH1D *h = h2->ProjectionY("py", ibin, ibin);
756 if(h->GetEntries()<100) continue;
760 Int_t ip = g[0]->GetN();
761 g[0]->SetPoint(ip, x, f.GetParameter(1));
762 g[0]->SetPointError(ip, 0., f.GetParError(1));
763 g[1]->SetPoint(ip, x, f.GetParameter(2));
764 g[1]->SetPointError(ip, 0., f.GetParError(2));
768 //____________________________________________________________________
769 void AliTRDcheckESD::PrintStatus(ULong_t status)
771 // Dump track status to stdout
773 printf("ITS[i(%d) o(%d) r(%d)] TPC[i(%d) o(%d) r(%d) p(%d)] TRD[i(%d) o(%d) r(%d) p(%d) s(%d)] HMPID[o(%d) p(%d)]\n"
774 ,Bool_t(status & AliESDtrack::kITSin)
775 ,Bool_t(status & AliESDtrack::kITSout)
776 ,Bool_t(status & AliESDtrack::kITSrefit)
777 ,Bool_t(status & AliESDtrack::kTPCin)
778 ,Bool_t(status & AliESDtrack::kTPCout)
779 ,Bool_t(status & AliESDtrack::kTPCrefit)
780 ,Bool_t(status & AliESDtrack::kTPCpid)
781 ,Bool_t(status & AliESDtrack::kTRDin)
782 ,Bool_t(status & AliESDtrack::kTRDout)
783 ,Bool_t(status & AliESDtrack::kTRDrefit)
784 ,Bool_t(status & AliESDtrack::kTRDpid)
785 ,Bool_t(status & AliESDtrack::kTRDStop)
786 ,Bool_t(status & AliESDtrack::kHMPIDout)
787 ,Bool_t(status & AliESDtrack::kHMPIDpid)