1 //------------------------------------------------------------------------------
2 // Implementation of AliPerformanceMatch class. It checks the track matching
3 // quality (residuals, pulls) bewteen tracking detectors TPC-ITS and TPC-TRD
4 // with TPC as the reference detector. In addition, the ITS and TRD track
5 // reconstruction efficiency with respect to TPC is calculated.
7 // The matchinig quality parameters are stored in the THnSparse histograms.
8 // The analysis of these histograms to extract reference information is done in
9 // the Analyse() function (post-processing).
11 // The histograms with reference information can be exported to the ROOT folder.
13 // Author: J.Otwinowski 17/10/2009
14 // Changes by M.Knichel 22/10/2010
15 //------------------------------------------------------------------------------
19 // after running comparison task, read the file, and get component
20 gROOT->LoadMacro("$ALICE_ROOT/PWG1/Macros/LoadMyLibs.C");
23 TFile f("Output.root");
24 AliPerformanceMatch * compObj = (AliPerformanceMatch*)coutput->FindObject("AliPerformanceMatchTPCTRD");
26 // analyse comparison data
29 // the output histograms/graphs will be stored in the folder "folderMatch"
30 compObj->GetAnalysisFolder()->ls("*");
32 // user can save whole comparison object (or only folder with anlysed histograms)
33 // in the seperate output file (e.g.)
34 TFile fout("Analysed_Match.root","recreate");
35 compObj->Write(); // compObj->GetAnalysisFolder()->Write();
46 #include "AliPerformanceMatch.h"
47 #include "AliESDEvent.h"
48 #include "AliESDVertex.h"
49 #include "AliESDtrack.h"
50 #include "AliESDfriendTrack.h"
51 #include "AliESDfriend.h"
53 #include "AliMCEvent.h"
54 #include "AliMCParticle.h"
55 #include "AliHeader.h"
56 #include "AliGenEventHeader.h"
58 #include "AliMCInfoCuts.h"
59 #include "AliRecInfoCuts.h"
60 #include "AliTracker.h"
61 #include "AliTRDtrackV1.h"
62 #include "AliTreeDraw.h"
66 ClassImp(AliPerformanceMatch)
68 Bool_t AliPerformanceMatch::fgMergeTHnSparse = kFALSE;
70 //_____________________________________________________________________________
71 AliPerformanceMatch::AliPerformanceMatch():
72 AliPerformanceObject("AliPerformanceMatch"),
89 //_____________________________________________________________________________
90 AliPerformanceMatch::AliPerformanceMatch(Char_t* name="AliPerformanceMatch", Char_t* title="AliPerformanceMatch",Int_t analysisMode=0,Bool_t hptGenerator=kFALSE):
91 AliPerformanceObject(name,title),
106 SetAnalysisMode(analysisMode);
107 SetHptGenerator(hptGenerator);
112 //_____________________________________________________________________________
113 AliPerformanceMatch::~AliPerformanceMatch()
117 if(fResolHisto) delete fResolHisto; fResolHisto=0;
118 if(fPullHisto) delete fPullHisto; fPullHisto=0;
119 if(fTrackingEffHisto) delete fTrackingEffHisto; fTrackingEffHisto = 0x0;
121 if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
122 if(fFolderObj) delete fFolderObj; fFolderObj=0;
125 //_____________________________________________________________________________
126 void AliPerformanceMatch::Init(){
129 // Make performance histogrms
134 Double_t ptMin = 1.e-2, ptMax = 10.;
136 Double_t *binsPt = 0;
137 if (IsHptGenerator()) {
138 nPtBins = 100; ptMax = 100.;
139 binsPt = CreateLogAxis(nPtBins,ptMin,ptMax);
141 binsPt = CreateLogAxis(nPtBins,ptMin,ptMax);
144 Double_t yMin = -0.02, yMax = 0.02;
145 Double_t zMin = -12.0, zMax = 12.0;
146 if(GetAnalysisMode() == 0 || GetAnalysisMode() == 1) {
147 yMin = -100.; yMax = 100.;
148 zMin = -100.; zMax = 100.;
151 // res_y:res_z:res_phi,res_lambda:res_pt:y:z:eta:phi:pt:isRec
152 Int_t binsResolHisto[11]={100,100,100,100,100,25,50,90,30,nPtBins,2};
153 Double_t minResolHisto[11]={-1.,-1.,-0.03,-0.03,-0.2, yMin, zMin, 0., -1.5, ptMin,0};
154 Double_t maxResolHisto[11]={ 1., 1., 0.03, 0.03, 0.2, yMax, zMax, 2.*TMath::Pi(), 1.5, ptMax,2};
156 fResolHisto = new THnSparseF("fResolHisto","res_y:res_z:res_phi:res_lambda:res_pt:y:z:phi:eta:pt:isRec",11,binsResolHisto,minResolHisto,maxResolHisto);
157 fResolHisto->SetBinEdges(9,binsPt);
159 fResolHisto->GetAxis(0)->SetTitle("y-y_{ref} (cm)");
160 fResolHisto->GetAxis(1)->SetTitle("z-z_{ref} (cm)");
161 fResolHisto->GetAxis(2)->SetTitle("#phi-#phi_{ref} (rad)");
162 fResolHisto->GetAxis(3)->SetTitle("#lambda-#lambda_{ref} (rad)");
163 fResolHisto->GetAxis(4)->SetTitle("(p_{T}/p_{Tref}-1)");
164 fResolHisto->GetAxis(5)->SetTitle("y_{ref} (cm)");
165 fResolHisto->GetAxis(6)->SetTitle("z_{ref} (cm)");
166 fResolHisto->GetAxis(7)->SetTitle("#phi_{ref} (rad)");
167 fResolHisto->GetAxis(8)->SetTitle("#eta_{ref}");
168 fResolHisto->GetAxis(9)->SetTitle("p_{Tref} (GeV/c)");
169 fResolHisto->GetAxis(10)->SetTitle("isReconstructed");
170 fResolHisto->Sumw2();
172 //pull_y:pull_z:pull_snp:pull_tgl:pull_1pt:y:z:snp:tgl:1pt:isRec
173 Int_t binsPullHisto[11]={100,100,100,100,100,50,50,50,50,nPtBins,2};
174 Double_t minPullHisto[11]={-5.,-5.,-5.,-5.,-5., yMin, zMin,-1.,-2.0, ptMin,0};
175 Double_t maxPullHisto[11]={ 5., 5., 5., 5., 5., yMax, zMax, 1., 2.0, ptMax,2};
176 fPullHisto = new THnSparseF("fPullHisto","pull_y:pull_z:pull_y:pull_z:pull_snp:pull_tgl:pull_1pt:y:z:snp:tgl:1pt:isRec",11,binsPullHisto,minPullHisto,maxPullHisto);
178 fPullHisto->GetAxis(0)->SetTitle("(y-y_{ref})/#sigma");
179 fPullHisto->GetAxis(1)->SetTitle("(z-z_{ref})/#sigma");
180 fPullHisto->GetAxis(2)->SetTitle("(sin#phi-sin#phi_{ref})/#sigma");
181 fPullHisto->GetAxis(3)->SetTitle("(tan#lambda-tan#lambda_{ref})/#sigma");
182 fPullHisto->GetAxis(4)->SetTitle("(p_{Tref}/p_{T}-1)/#sigma");
183 fPullHisto->GetAxis(5)->SetTitle("y_{ref} (cm)");
184 fPullHisto->GetAxis(6)->SetTitle("z_{ref} (cm)");
185 fPullHisto->GetAxis(7)->SetTitle("sin#phi_{ref}");
186 fPullHisto->GetAxis(8)->SetTitle("tan#lambda_{ref}");
187 fPullHisto->GetAxis(9)->SetTitle("1/p_{Tref} (GeV/c)^{-1}");
188 fPullHisto->GetAxis(10)->SetTitle("isReconstructed");
191 // -> has match:y:z:snp:tgl:phi:pt:ITSclusters
192 Int_t binsTrackingEffHisto[8] = { 2, 50, 100, 50, 50, 90, 100, 7 };
193 Double_t minTrackingEffHisto[8] = {-0.5, -25, -50, -1, -2, 0., 0, -0.5 };
194 Double_t maxTrackingEffHisto[8] = { 1.5, 25, 50, 1, 2, 2*TMath::Pi(), 20, 6.5 };
196 fTrackingEffHisto = new THnSparseF("fTrackingEffHisto","has match:y:z:snp:tgl:phi:pt:ITSclusters",8,binsTrackingEffHisto,minTrackingEffHisto,maxTrackingEffHisto);
197 fTrackingEffHisto->GetAxis(0)->SetTitle("IsMatching");
198 fTrackingEffHisto->GetAxis(1)->SetTitle("local y (cm)");
199 fTrackingEffHisto->GetAxis(2)->SetTitle("z (cm)");
200 fTrackingEffHisto->GetAxis(3)->SetTitle("sin(#phi)");
201 fTrackingEffHisto->GetAxis(4)->SetTitle("tan(#lambda)");
202 fTrackingEffHisto->GetAxis(5)->SetTitle("phi (rad)");
203 fTrackingEffHisto->GetAxis(6)->SetTitle("p_{T}");
204 fTrackingEffHisto->GetAxis(7)->SetTitle("number of ITS clusters");
205 fTrackingEffHisto->Sumw2();
209 AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
211 AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object");
214 fAnalysisFolder = CreateFolder("folderMatch","Analysis Matching Folder");
217 //_____________________________________________________________________________
218 void AliPerformanceMatch::ProcessITSTPC(Int_t iTrack, AliESDEvent *const esdEvent, AliStack* /*const stack*/, AliESDtrack *const esdTrack, AliESDfriendTrack *const esdFriendTrack)
221 // addition to standard analysis - check if ITS stand-alone tracks have a match in the TPC
222 // Origin: A. Kalwait
223 // Modified: J. Otwinowski
224 if(!esdEvent) return;
225 if(!esdTrack) return;
226 if(!esdFriendTrack) return;
228 // ITS stand alone tracks with SPD layers
229 if(!(esdTrack->GetStatus() & AliESDtrack::kITSpureSA)) return;
230 if(!(esdTrack->GetStatus() & AliESDtrack::kITSrefit)) return;
231 if(esdTrack->GetNcls(0)<fCutsRC->GetMinNClustersITS()) return;
232 if(!esdTrack->HasPointOnITSLayer(0) && !esdTrack->HasPointOnITSLayer(1)) return;
234 Bool_t hasMatch = kFALSE;
235 for (Int_t jTrack = 0; jTrack < esdEvent->GetNumberOfTracks(); jTrack++) {
236 // loop for all those tracks and check if the corresponding TPC track is found
237 if (jTrack==iTrack) continue;
238 AliESDtrack *trackTPC = esdEvent->GetTrack(jTrack);
239 if (!trackTPC) continue;
240 if (!trackTPC->GetTPCInnerParam()) continue;
242 // TPC nClust/track after first tracking pass
243 if(trackTPC->GetTPCNclsIter1()<fCutsRC->GetMinNClustersTPC()) continue;
245 AliExternalTrackParam *innerTPC = new AliExternalTrackParam(*(trackTPC->GetTPCInnerParam()));
246 if(!innerTPC) continue;
248 Double_t x[3]; trackTPC->GetXYZ(x);
249 Double_t b[3]; AliTracker::GetBxByBz(x,b);
250 Double_t dca[2],cov[3];
251 Bool_t isTPCOK = innerTPC->PropagateToDCABxByBz(esdEvent->GetPrimaryVertexSPD(),b,kVeryBig,dca,cov);
253 if(innerTPC) delete innerTPC;
260 Double_t dcaToVertex = -1;
261 if( fCutsRC->GetDCAToVertex2D() )
263 dcaToVertex = TMath::Sqrt(dca[0]*dca[0]/fCutsRC->GetMaxDCAToVertexXY()/fCutsRC->GetMaxDCAToVertexXY()+dca[1]*dca[1]/fCutsRC->GetMaxDCAToVertexZ()/fCutsRC->GetMaxDCAToVertexZ());
265 if(fCutsRC->GetDCAToVertex2D() && dcaToVertex > 1) {
266 delete innerTPC; continue; }
267 if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[0]) > fCutsRC->GetMaxDCAToVertexXY()) {
268 delete innerTPC; continue; }
269 if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[1]) > fCutsRC->GetMaxDCAToVertexZ()) {
270 delete innerTPC; continue; }
272 // rough checking if they match
274 printf("+++++++++++++++++++++++++++++++++++++++++++++++ ");
275 printf("esdTrack->GetY() %e, esdTrack->GetSnp() %e, esdTrack->GetTgl() %e \n",
276 esdTrack->GetY(), esdTrack->GetSnp(), esdTrack->GetTgl());
277 printf("innerTPC->GetY() %e, innerTPC->GetSnp() %e, innerTPC->GetTgl() %e \n",
278 innerTPC->GetY() , innerTPC->GetSnp() , innerTPC->GetTgl());
280 if (TMath::Abs(esdTrack->GetY() - innerTPC->GetY()) > 3) {
281 delete innerTPC; continue;
283 if (TMath::Abs(esdTrack->GetSnp() - innerTPC->GetSnp()) > 0.2) {
284 delete innerTPC; continue;
286 if (TMath::Abs(esdTrack->GetTgl() - innerTPC->GetTgl()) > 0.2) {
287 delete innerTPC; continue;
291 if(innerTPC) delete innerTPC;
293 //has match:y:z:snp:tgl:phi:pt:ITSclusters
294 Double_t vecTrackingEff[8] = { hasMatch,esdTrack->GetY(),esdTrack->GetZ(),esdTrack->GetSnp(),esdTrack->GetTgl(),esdTrack->Phi(), esdTrack->Pt(),esdTrack->GetITSclusters(0) };
295 fTrackingEffHisto->Fill(vecTrackingEff);
299 //_____________________________________________________________________________
300 void AliPerformanceMatch::ProcessTPCITS(AliStack* /*const stack*/, AliESDtrack *const esdTrack, AliESDfriendTrack *const esdFriendTrack)
303 // Match TPC and ITS min-bias tracks
304 // at radius between detectors
306 if(!esdTrack) return;
307 if(!esdFriendTrack) return;
310 // Propagate tracks to the radius between TPC-ITS
311 // using B-field and material budget
313 Double_t radius = fCutsRC->GetTPCITSMatchingRadius();
314 Double_t mass = esdTrack->GetMass();
315 Double_t step=1.0; // cm
318 // Propagate TPCinner (reference detector)
320 Bool_t isTPCOK=kFALSE;
321 AliExternalTrackParam *innerTPC=NULL;
323 Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
324 esdTrack->GetImpactParametersTPC(dca,cov);
329 Double_t dcaToVertex = -1;
330 if( fCutsRC->GetDCAToVertex2D() )
332 dcaToVertex = TMath::Sqrt(dca[0]*dca[0]/fCutsRC->GetMaxDCAToVertexXY()/fCutsRC->GetMaxDCAToVertexXY() + dca[1]*dca[1]/fCutsRC->GetMaxDCAToVertexZ()/fCutsRC->GetMaxDCAToVertexZ());
334 if(fCutsRC->GetDCAToVertex2D() && dcaToVertex > 1) return;
335 if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[0]) > fCutsRC->GetMaxDCAToVertexXY()) return;
336 if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[1]) > fCutsRC->GetMaxDCAToVertexZ()) return;
338 if( (esdTrack->GetTPCNclsIter1()>fCutsRC->GetMinNClustersTPC()) &&
339 (esdTrack->GetTPCInnerParam()) &&
340 (innerTPC=new AliExternalTrackParam(*(esdTrack->GetTPCInnerParam()))))
342 isTPCOK = AliTracker::PropagateTrackToBxByBz(innerTPC,radius,mass,step,kTRUE);
345 if(innerTPC) delete innerTPC;
350 // Propagate ITSouter
352 Bool_t isITSOK=kFALSE;
353 AliExternalTrackParam *outerITS=NULL;
355 if( (esdTrack->GetNcls(0)>fCutsRC->GetMinNClustersITS()) &&
356 (esdFriendTrack->GetITSOut()) &&
357 (outerITS=new AliExternalTrackParam(*(esdFriendTrack->GetITSOut()))))
359 isITSOK = AliTracker::PropagateTrackToBxByBz(outerITS,radius,mass,step,kTRUE);
363 // Fill histograms (TPC reference detector)
366 FillHistograms(innerTPC,outerITS,isITSOK);
368 if(outerITS) delete outerITS;
369 if(innerTPC) delete innerTPC;
372 //_____________________________________________________________________________
373 void AliPerformanceMatch::ProcessTPCTRD(AliStack* /*const stack*/, AliESDtrack *const esdTrack, AliESDfriendTrack *const esdFriendTrack)
376 // Match TPC and TRD min-bias tracks
377 // at radius between detectors. TPC is the reference detector.
379 if(!esdTrack) return;
380 if(!esdFriendTrack) return;
383 // Propagate tracks to the radius between TPC-TRD
384 // using B-field and material budget
386 Double_t radius = fCutsRC->GetTPCTRDMatchingRadius();
387 Double_t mass = esdTrack->GetMass();
388 Double_t step=1.0; // cm
391 // Propagate TPCouter (reference detector)
393 Bool_t isTPCOK=kFALSE;
394 AliExternalTrackParam *outerTPC=NULL;
396 Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
397 esdTrack->GetImpactParametersTPC(dca,cov);
402 Double_t dcaToVertex = -1;
403 if( fCutsRC->GetDCAToVertex2D() )
405 dcaToVertex = TMath::Sqrt(dca[0]*dca[0]/fCutsRC->GetMaxDCAToVertexXY()/fCutsRC->GetMaxDCAToVertexXY() + dca[1]*dca[1]/fCutsRC->GetMaxDCAToVertexZ()/fCutsRC->GetMaxDCAToVertexZ());
407 if(fCutsRC->GetDCAToVertex2D() && dcaToVertex > 1) return;
408 if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[0]) > fCutsRC->GetMaxDCAToVertexXY()) return;
409 if(!fCutsRC->GetDCAToVertex2D() && TMath::Abs(dca[1]) > fCutsRC->GetMaxDCAToVertexZ()) return;
411 if( (esdTrack->GetTPCNclsIter1()>fCutsRC->GetMinNClustersTPC()) &&
412 (esdFriendTrack->GetTPCOut()) &&
413 (outerTPC=new AliExternalTrackParam(*(esdFriendTrack->GetTPCOut()))))
415 isTPCOK = AliTracker::PropagateTrackToBxByBz(outerTPC,radius,mass,step,kTRUE);
418 if(outerTPC) delete outerTPC;
423 // Propagate TRDinner
425 Bool_t isTRDOK = kFALSE;
426 AliExternalTrackParam *innerTRD=NULL;
429 AliTRDtrackV1 *trdTrack=NULL; //esdFriendTrack = fESDfriend->GetTrack(itrk);
430 TObject *calObject=NULL;
432 while((calObject = esdFriendTrack->GetCalibObject(icalib++))) {
433 if(strcmp(calObject->IsA()->GetName(),"AliTRDtrackV1") != 0) continue; // Look for the TRDtrack
434 if(!(trdTrack = dynamic_cast<AliTRDtrackV1*>(calObject))) break;
438 (trdTrack->GetNumberOfTracklets()>fCutsRC->GetMinNTrackletsTRD()) &&
439 (trdTrack->GetTracklet(0)) &&
440 (esdFriendTrack->GetTRDIn()) &&
441 (innerTRD = new AliExternalTrackParam(*(esdFriendTrack->GetTRDIn()))))
443 isTRDOK = AliTracker::PropagateTrackToBxByBz(innerTRD,radius,mass,step,kTRUE);
447 // Fill histograms (TPC reference detector)
450 FillHistograms(outerTPC,innerTRD,isTRDOK);
452 if(outerTPC) delete outerTPC;
453 if(innerTRD) delete innerTRD;
456 //_____________________________________________________________________________
457 void AliPerformanceMatch::FillHistograms(AliExternalTrackParam *const refParam, AliExternalTrackParam *const param, Bool_t isRec)
460 // fill performance histograms
461 // (TPC always as reference)
463 if(!refParam) return;
466 // Deltas (dy,dz,dphi,dtheta,dpt)
468 Float_t delta[5] = {0};
470 delta[0] = param->GetY()-refParam->GetY();
471 delta[1] = param->GetZ()-refParam->GetZ();
472 delta[2] = TMath::ATan2(param->Py(),param->Px())-TMath::ATan2(refParam->Py(),refParam->Px());
473 delta[3] = TMath::ATan2(param->Pz(),param->Pt())-TMath::ATan2(refParam->Pz(),refParam->Pt());
474 if(refParam->Pt()) delta[4] = (param->Pt()-refParam->Pt())/refParam->Pt();
477 // Pulls (y,z,snp,tanl,1/pt)
479 Float_t sigma[5] = {0};
480 Float_t pull[5] = {0};
482 sigma[0] = TMath::Sqrt(param->GetSigmaY2()+refParam->GetSigmaY2());
483 sigma[1] = TMath::Sqrt(param->GetSigmaZ2()+refParam->GetSigmaZ2());
484 sigma[2] = TMath::Sqrt(param->GetSigmaSnp2()+refParam->GetSigmaSnp2());
485 sigma[3] = TMath::Sqrt(param->GetSigmaTgl2()+refParam->GetSigmaTgl2());
486 sigma[4] = TMath::Sqrt(param->GetSigma1Pt2()+refParam->GetSigma1Pt2());
487 if(sigma[0]) pull[0] = delta[0] / sigma[0];
488 if(sigma[1]) pull[1] = delta[1] / sigma[1];
489 if(sigma[2]) pull[2] = (param->GetSnp()-refParam->GetSnp()) / sigma[2];
490 if(sigma[3]) pull[3] = (param->GetTgl()-refParam->GetTgl()) / sigma[3];
491 if(sigma[4]) pull[4] = (param->OneOverPt()-refParam->OneOverPt()) / sigma[4];
495 Double_t vResolHisto[11] = {delta[0],delta[1],delta[2],delta[3],delta[4],refParam->GetY(),refParam->GetZ(),refParam->Phi(),refParam->Eta(),refParam->Pt(),isRec};
496 fResolHisto->Fill(vResolHisto);
498 Double_t vPullHisto[11] = {pull[0],pull[1],pull[2],pull[3],pull[4],refParam->GetY(),refParam->GetZ(),refParam->GetSnp(),refParam->GetTgl(),refParam->OneOverPt(),isRec};
499 fPullHisto->Fill(vPullHisto);
502 //_____________________________________________________________________________
503 void AliPerformanceMatch::Exec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent, AliESDfriend *const esdFriend, const Bool_t bUseMC, const Bool_t bUseESDfriend)
505 // Process comparison information
509 Error("Exec","esdEvent not available");
512 AliHeader* header = 0;
513 AliGenEventHeader* genHeader = 0;
520 Error("Exec","mcEvent not available");
523 // get MC event header
524 header = mcEvent->Header();
526 Error("Exec","Header not available");
530 stack = mcEvent->Stack();
532 Error("Exec","Stack not available");
536 genHeader = header->GenEventHeader();
538 Error("Exec","Could not retrieve genHeader from Header");
541 genHeader->PrimaryVertex(vtxMC);
547 Error("Exec","esdFriend not available");
553 if(!bUseMC && GetTriggerClass()) {
554 Bool_t isEventTriggered = esdEvent->IsTriggerClassFired(GetTriggerClass());
555 if(!isEventTriggered) return;
558 // get TPC event vertex
559 const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTPC();
560 if(vtxESD && (vtxESD->GetStatus()<=0)) return;
563 for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++)
565 AliESDtrack *track = esdEvent->GetTrack(iTrack);
568 AliESDfriendTrack *friendTrack=0;
570 friendTrack=esdFriend->GetTrack(iTrack);
571 if(!friendTrack) continue;
574 if(GetAnalysisMode() == 0) ProcessTPCITS(stack,track,friendTrack);
575 else if(GetAnalysisMode() == 1) ProcessTPCTRD(stack,track,friendTrack);
576 else if(GetAnalysisMode() == 2) ProcessITSTPC(iTrack,esdEvent,stack,track,friendTrack);
578 printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
584 //_____________________________________________________________________________
585 TH1F* AliPerformanceMatch::MakeResol(TH2F * his, Int_t integ, Bool_t type, Int_t cut){
586 // Create resolution histograms
589 if (!gPad) new TCanvas;
590 hisr = AliTreeDraw::CreateResHistoII(his,&hism,integ,kTRUE,cut);
591 if (type) return hism;
596 //_____________________________________________________________________________
597 void AliPerformanceMatch::Analyse() {
598 // Analyse comparison information and store output histograms
599 // in the folder "folderMatch"
603 TH1::AddDirectory(kFALSE);
608 TObjArray *aFolderObj = new TObjArray;
610 // write results in the folder
611 // TCanvas * c = new TCanvas("Phi resol Tan","Phi resol Tan");
617 if(GetAnalysisMode()==0 || GetAnalysisMode()==1) {
619 fResolHisto->GetAxis(10)->SetRangeUser(1.0,2.0); // only reconstructed
620 fPullHisto->GetAxis(10)->SetRangeUser(1.0,2.0); // only reconstructed
621 for(Int_t i=0; i<5; i++)
623 for(Int_t j=5; j<10; j++)
625 //if(j!=8) fResolHisto->GetAxis(8)->SetRangeUser(-0.9,0.89); // eta window
626 if(j!=8) fResolHisto->GetAxis(8)->SetRangeUser(0.0,0.89); // eta window
627 else fResolHisto->GetAxis(8)->SetRangeUser(-1.5,1.49);
628 fResolHisto->GetAxis(9)->SetRangeUser(0.1,100.); // pt threshold
631 AddProjection(aFolderObj, "match", fResolHisto, i, j, &selString);
634 h2D = (TH2F*)fResolHisto->Projection(i,j);
635 h = AliPerformanceMatch::MakeResol(h2D,1,0,100);
636 sprintf(name,"h_res_%d_vs_%d",i,j);
639 h->GetXaxis()->SetTitle(fResolHisto->GetAxis(j)->GetTitle());
640 sprintf(title,"%s %s",fResolHisto->GetAxis(i)->GetTitle(),"(resolution)");
641 h->GetYaxis()->SetTitle(title);
642 sprintf(title,"%s vs %s",title,fResolHisto->GetAxis(j)->GetTitle());
645 //if(j==9) h->SetBit(TH1::kLogX);
649 h = AliPerformanceMatch::MakeResol(h2D,1,1,100);
650 //h = (TH1F*)arr->At(1);
651 sprintf(name,"h_mean_res_%d_vs_%d",i,j);
654 h->GetXaxis()->SetTitle(fResolHisto->GetAxis(j)->GetTitle());
655 sprintf(title,"%s %s",fResolHisto->GetAxis(i)->GetTitle(),"(mean)");
656 h->GetYaxis()->SetTitle(title);
658 sprintf(title,"%s vs %s",title,fResolHisto->GetAxis(j)->GetTitle());
661 if(j==9) h->SetBit(TH1::kLogX);
665 if(j!=8) fPullHisto->GetAxis(8)->SetRangeUser(0.0,0.89); // eta window
666 else fPullHisto->GetAxis(8)->SetRangeUser(-1.5,1.49); // eta window
667 fPullHisto->GetAxis(9)->SetRangeUser(0.1,100.); // pt threshold
669 AddProjection(aFolderObj, "match", fPullHisto, i, j, &selString);
672 h2D = (TH2F*)fPullHisto->Projection(i,j);
674 h = AliPerformanceMatch::MakeResol(h2D,1,0,100);
675 sprintf(name,"h_pull_%d_vs_%d",i,j);
678 h->GetXaxis()->SetTitle(fPullHisto->GetAxis(j)->GetTitle());
679 sprintf(title,"%s %s",fPullHisto->GetAxis(i)->GetTitle(),"(resolution)");
680 h->GetYaxis()->SetTitle(title);
681 sprintf(title,"%s vs %s",title,fPullHisto->GetAxis(j)->GetTitle());
684 if(j==9) h->SetBit(TH1::kLogX);
687 h = AliPerformanceMatch::MakeResol(h2D,1,1,100);
688 sprintf(name,"h_mean_pull_%d_vs_%d",i,j);
691 h->GetXaxis()->SetTitle(fPullHisto->GetAxis(j)->GetTitle());
692 sprintf(title,"%s %s",fPullHisto->GetAxis(i)->GetTitle(),"(mean)");
693 h->GetYaxis()->SetTitle(title);
694 sprintf(title,"%s vs %s",title,fPullHisto->GetAxis(j)->GetTitle());
697 //if(j==9) h->SetBit(TH1::kLogX);
706 for(Int_t i=5; i<10; i++)
708 if(i!=8) fResolHisto->GetAxis(8)->SetRangeUser(-0.9,0.89); // eta window
709 else fResolHisto->GetAxis(8)->SetRangeUser(-1.5,1.49);
710 fResolHisto->GetAxis(9)->SetRangeUser(0.1,100.); // pt threshold
712 fResolHisto->GetAxis(10)->SetRange(1,fResolHisto->GetAxis(10)->GetNbins()); // all
713 selString = "eff_all";
714 AddProjection(aFolderObj, "match", fResolHisto, i, &selString);
715 // h = (TH1F*)fResolHisto->Projection(i);
717 fResolHisto->GetAxis(10)->SetRange(2,2); // only reconstructed
718 selString = "eff_rec";
719 AddProjection(aFolderObj, "match", fResolHisto, i, &selString);
720 //h2 = (TH1F*)fResolHisto->Projection(i);
723 TH1F* h2c = (TH1F*)h2->Clone();
724 h2c->Divide(h2,h,1,1,"B");
726 sprintf(name,"h_eff_%d",i);
729 h2c->GetXaxis()->SetTitle(h2c->GetXaxis()->GetTitle());
730 h2c->GetYaxis()->SetTitle("efficiency");
731 h2c->SetTitle("matching effciency");
733 aFolderObj->Add(h2c);
740 // TPC efficiency wrt ITS
742 if(GetAnalysisMode()==2) {
743 selString = "trackingeff";
744 AddProjection(aFolderObj, "match", fTrackingEffHisto, 0, &selString);
746 // h = (TH1F*)fTrackingEffHisto->Projection(0);
747 // aFolderObj->Add(h);
749 for(Int_t i=1; i<7; i++)
753 // calculate efficiency
756 // all ITS standalone tracks
757 fTrackingEffHisto->GetAxis(0)->SetRange(1,fTrackingEffHisto->GetAxis(0)->GetNbins());
758 //h = (TH1F*)fTrackingEffHisto->Projection(i);
759 selString = "trackingeff_all";
760 AddProjection(aFolderObj, "match", fTrackingEffHisto, i, &selString);
762 // TPC tracks which has matching with TPC
763 fTrackingEffHisto->GetAxis(0)->SetRange(2,2);
764 //h2 = (TH1F*)fTrackingEffHisto->Projection(i);
765 selString = "trackingeff_tpc";
766 AddProjection(aFolderObj, "match", fTrackingEffHisto, i, &selString);
769 TH1F* h2c = (TH1F*)h2->Clone();
770 h2c->Divide(h2,h,1,1,"B");
772 sprintf(name,"h_TPC_eff_%d",i);
775 h2c->GetXaxis()->SetTitle(h2c->GetXaxis()->GetTitle());
776 h2c->GetYaxis()->SetTitle("efficiency");
777 h2c->SetTitle("TPC effciency wrt ITS");
779 aFolderObj->Add(h2c);
784 printf("exportToFolder\n");
785 // export objects to analysis folder
786 fAnalysisFolder = ExportToFolder(aFolderObj);
788 // delete only TObjArray
789 if(fFolderObj) delete fFolderObj;
790 fFolderObj = aFolderObj;
795 //_____________________________________________________________________________
796 void AliPerformanceMatch::AnalyseFinal() {
798 printf("AliPerformanceMatch: no projections available to analyse\n");
801 TH1::AddDirectory(kFALSE);
805 TObjArray *aFolderObj = fFolderObj;
807 // write results in the folder
808 TCanvas * c = new TCanvas("Phi resol Tan","Phi resol Tan");
814 if(GetAnalysisMode()==0 || GetAnalysisMode()==1) {
816 for(Int_t i=0; i<5; i++)
818 for(Int_t j=5; j<10; j++)
820 sprintf(name,"h_tpc_match_resol_%d_%d",i,j);
821 h2D = dynamic_cast<TH2F*>(aFolderObj->FindObject(name));
823 h = AliPerformanceMatch::MakeResol(h2D,1,0,100);
824 sprintf(name,"h_res_%d_vs_%d",i,j);
826 h->GetXaxis()->SetTitle(fResolHisto->GetAxis(j)->GetTitle());
827 sprintf(title,"%s %s",fResolHisto->GetAxis(i)->GetTitle(),"(resolution)");
828 h->GetYaxis()->SetTitle(title);
829 sprintf(title,"%s vs %s",title,fResolHisto->GetAxis(j)->GetTitle());
832 //if(j==9) h->SetBit(TH1::kLogX);
835 h = AliPerformanceMatch::MakeResol(h2D,1,1,100);
836 //h = (TH1F*)arr->At(1);
837 sprintf(name,"h_mean_res_%d_vs_%d",i,j);
840 h->GetXaxis()->SetTitle(fResolHisto->GetAxis(j)->GetTitle());
841 sprintf(title,"%s %s",fResolHisto->GetAxis(i)->GetTitle(),"(mean)");
842 h->GetYaxis()->SetTitle(title);
844 sprintf(title,"%s vs %s",title,fResolHisto->GetAxis(j)->GetTitle());
847 if(j==9) h->SetBit(TH1::kLogX);
852 sprintf(name,"h_tpc_match_pull_%d_%d",i,j);
853 h2D = dynamic_cast<TH2F*>(aFolderObj->FindObject(name));
855 h = AliPerformanceMatch::MakeResol(h2D,1,0,100);
856 sprintf(name,"h_pull_%d_vs_%d",i,j);
859 h->GetXaxis()->SetTitle(fPullHisto->GetAxis(j)->GetTitle());
860 sprintf(title,"%s %s",fPullHisto->GetAxis(i)->GetTitle(),"(resolution)");
861 h->GetYaxis()->SetTitle(title);
862 sprintf(title,"%s vs %s",title,fPullHisto->GetAxis(j)->GetTitle());
865 if(j==9) h->SetBit(TH1::kLogX);
868 h = AliPerformanceMatch::MakeResol(h2D,1,1,100);
869 sprintf(name,"h_mean_pull_%d_vs_%d",i,j);
872 h->GetXaxis()->SetTitle(fPullHisto->GetAxis(j)->GetTitle());
873 sprintf(title,"%s %s",fPullHisto->GetAxis(i)->GetTitle(),"(mean)");
874 h->GetYaxis()->SetTitle(title);
875 sprintf(title,"%s vs %s",title,fPullHisto->GetAxis(j)->GetTitle());
878 //if(j==9) h->SetBit(TH1::kLogX);
888 for(Int_t i=5; i<10; i++)
890 sprintf(name,"h_tpc_match_eff_all_%d",i);
891 h = dynamic_cast<TH1F*>(aFolderObj->FindObject(name));
894 sprintf(name,"h_tpc_match_eff_rec_%d",i);
895 h2 = dynamic_cast<TH1F*>(aFolderObj->FindObject(name));
898 TH1F* h2c = (TH1F*)h2->Clone();
899 h2c->Divide(h2,h,1,1,"B");
901 sprintf(name,"h_eff_%d",i);
904 h2c->GetXaxis()->SetTitle(h2c->GetXaxis()->GetTitle());
905 h2c->GetYaxis()->SetTitle("efficiency");
906 h2c->SetTitle("matching effciency");
908 aFolderObj->Add(h2c);
917 // TPC efficiency wrt ITS
919 if(GetAnalysisMode()==2) {
921 for(Int_t i=1; i<7; i++)
925 // calculate efficiency
928 // all ITS standalone tracks
929 sprintf(name,"h_tpc_match_trackingeff_all_%d",i);
930 h = dynamic_cast<TH1F*>(aFolderObj->FindObject(name));
933 // TPC tracks which has matching with TPC
934 sprintf(name,"h_tpc_match_trackingeff_tpc_%d",i);
935 h2 = dynamic_cast<TH1F*>(aFolderObj->FindObject(name));
938 TH1F* h2c = (TH1F*)h2->Clone();
939 h2c->Divide(h2,h,1,1,"B");
941 sprintf(name,"h_TPC_eff_%d",i);
944 h2c->GetXaxis()->SetTitle(h2c->GetXaxis()->GetTitle());
945 h2c->GetYaxis()->SetTitle("efficiency");
946 h2c->SetTitle("TPC effciency wrt ITS");
948 aFolderObj->Add(h2c);
955 printf("exportToFolder\n");
956 // export objects to analysis folder
957 fAnalysisFolder = ExportToFolder(aFolderObj);
959 // delete only TObjArray
960 if(fFolderObj) delete fFolderObj;
961 fFolderObj = aFolderObj;
968 //_____________________________________________________________________________
969 TFolder* AliPerformanceMatch::ExportToFolder(TObjArray * array)
971 // recreate folder avery time and export objects to new one
973 AliPerformanceMatch * comp=this;
974 TFolder *folder = comp->GetAnalysisFolder();
977 TFolder *newFolder = 0;
979 Int_t size = array->GetSize();
982 // get name and title from old folder
983 name = folder->GetName();
984 title = folder->GetTitle();
990 newFolder = CreateFolder(name.Data(),title.Data());
991 newFolder->SetOwner();
993 // add objects to folder
995 newFolder->Add(array->At(i));
1003 //_____________________________________________________________________________
1004 TFolder* AliPerformanceMatch::CreateFolder(TString name,TString title) {
1005 // create folder for analysed histograms
1007 TFolder *folder = 0;
1008 folder = new TFolder(name.Data(),title.Data());
1013 //_____________________________________________________________________________
1014 Long64_t AliPerformanceMatch::Merge(TCollection* const list)
1016 // Merge list of objects (needed by PROOF)
1021 if (list->IsEmpty())
1024 TIterator* iter = list->MakeIterator();
1026 TObjArray* objArrayList = 0;
1027 objArrayList = new TObjArray();
1029 // collection of generated histograms
1031 while((obj = iter->Next()) != 0)
1033 AliPerformanceMatch* entry = dynamic_cast<AliPerformanceMatch*>(obj);
1034 if (entry == 0) continue;
1035 if (fgMergeTHnSparse) {
1036 if ((fResolHisto) && (entry->fResolHisto)) { fResolHisto->Add(entry->fResolHisto); }
1037 if ((fPullHisto) && (entry->fPullHisto)) { fPullHisto->Add(entry->fPullHisto); }
1038 if ((fTrackingEffHisto) && (entry->fTrackingEffHisto)) { fTrackingEffHisto->Add(entry->fTrackingEffHisto); }
1040 // the analysisfolder is only merged if present
1041 if (entry->fFolderObj) { objArrayList->Add(entry->fFolderObj); }
1045 if (fFolderObj) { fFolderObj->Merge(objArrayList); }
1046 // to signal that track histos were not merged: reset
1047 if (!fgMergeTHnSparse) { fResolHisto->Reset(); fPullHisto->Reset(); fTrackingEffHisto->Reset(); }
1049 if (objArrayList) delete objArrayList; objArrayList=0;