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;
69 Bool_t AliPerformanceMatch::fgUseMergeTHnSparse = kFALSE;
71 //_____________________________________________________________________________
72 AliPerformanceMatch::AliPerformanceMatch():
73 AliPerformanceObject("AliPerformanceMatch"),
90 //_____________________________________________________________________________
91 AliPerformanceMatch::AliPerformanceMatch(Char_t* name="AliPerformanceMatch", Char_t* title="AliPerformanceMatch",Int_t analysisMode=0,Bool_t hptGenerator=kFALSE):
92 AliPerformanceObject(name,title),
107 SetAnalysisMode(analysisMode);
108 SetHptGenerator(hptGenerator);
113 //_____________________________________________________________________________
114 AliPerformanceMatch::~AliPerformanceMatch()
118 if(fResolHisto) delete fResolHisto; fResolHisto=0;
119 if(fPullHisto) delete fPullHisto; fPullHisto=0;
120 if(fTrackingEffHisto) delete fTrackingEffHisto; fTrackingEffHisto = 0x0;
122 if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
123 if(fFolderObj) delete fFolderObj; fFolderObj=0;
126 //_____________________________________________________________________________
127 void AliPerformanceMatch::Init(){
130 // Make performance histogrms
135 Double_t ptMin = 1.e-2, ptMax = 10.;
137 Double_t *binsPt = 0;
138 if (IsHptGenerator()) {
139 nPtBins = 100; ptMax = 100.;
140 binsPt = CreateLogAxis(nPtBins,ptMin,ptMax);
142 binsPt = CreateLogAxis(nPtBins,ptMin,ptMax);
145 Double_t yMin = -0.02, yMax = 0.02;
146 Double_t zMin = -12.0, zMax = 12.0;
147 if(GetAnalysisMode() == 0 || GetAnalysisMode() == 1) {
148 yMin = -100.; yMax = 100.;
149 zMin = -100.; zMax = 100.;
152 // res_y:res_z:res_phi,res_lambda:res_pt:y:z:eta:phi:pt:isRec
153 Int_t binsResolHisto[9]={100,100,100,100,100,90,30,nPtBins,2};
154 Double_t minResolHisto[9]={-1.,-1.,-0.03,-0.03,-0.2, 0., -1.5, ptMin,0};
155 Double_t maxResolHisto[9]={ 1., 1., 0.03, 0.03, 0.2, 2.*TMath::Pi(), 1.5, ptMax,2};
157 fResolHisto = new THnSparseF("fResolHisto","res_y:res_z:res_phi:res_lambda:res_pt:phi:eta:pt:isRec",9,binsResolHisto,minResolHisto,maxResolHisto);
158 fResolHisto->SetBinEdges(7,binsPt);
160 fResolHisto->GetAxis(0)->SetTitle("y-y_{ref} (cm)");
161 fResolHisto->GetAxis(1)->SetTitle("z-z_{ref} (cm)");
162 fResolHisto->GetAxis(2)->SetTitle("#phi-#phi_{ref} (rad)");
163 fResolHisto->GetAxis(3)->SetTitle("#lambda-#lambda_{ref} (rad)");
164 fResolHisto->GetAxis(4)->SetTitle("(p_{T}/p_{Tref}-1)");
165 fResolHisto->GetAxis(5)->SetTitle("#phi_{ref} (rad)");
166 fResolHisto->GetAxis(6)->SetTitle("#eta_{ref}");
167 fResolHisto->GetAxis(7)->SetTitle("p_{Tref} (GeV/c)");
168 fResolHisto->GetAxis(8)->SetTitle("isReconstructed");
169 fResolHisto->Sumw2();
171 //pull_y:pull_z:pull_snp:pull_tgl:pull_1pt:y:z:snp:tgl:1pt:isRec
172 Int_t binsPullHisto[9]={100,100,100,100,100,90,30,nPtBins,2};
173 Double_t minPullHisto[9]={-5.,-5.,-5.,-5.,-5., 0.,-1.5, ptMin,0};
174 Double_t maxPullHisto[9]={ 5., 5., 5., 5., 5., 2.*TMath::Pi(), 1.5, ptMax,2};
175 fPullHisto = new THnSparseF("fPullHisto","pull_y:pull_z:pull_snp:pull_tgl:pull_1pt:phi:eta:1pt:isRec",9,binsPullHisto,minPullHisto,maxPullHisto);
177 fPullHisto->GetAxis(0)->SetTitle("(y-y_{ref})/#sigma");
178 fPullHisto->GetAxis(1)->SetTitle("(z-z_{ref})/#sigma");
179 fPullHisto->GetAxis(2)->SetTitle("(sin#phi-sin#phi_{ref})/#sigma");
180 fPullHisto->GetAxis(3)->SetTitle("(tan#lambda-tan#lambda_{ref})/#sigma");
181 fPullHisto->GetAxis(4)->SetTitle("(p_{Tref}/p_{T}-1)/#sigma");
182 fPullHisto->GetAxis(5)->SetTitle("#phi_{ref} (rad)");
183 fPullHisto->GetAxis(6)->SetTitle("#eta_{ref}");
184 fPullHisto->GetAxis(7)->SetTitle("1/p_{Tref} (GeV/c)^{-1}");
185 fPullHisto->GetAxis(8)->SetTitle("isReconstructed");
188 // -> has match:y:z:snp:tgl:phi:pt:ITSclusters
189 Int_t binsTrackingEffHisto[5] = { 2, 90, 100, 30, 7 };
190 Double_t minTrackingEffHisto[5] = {-0.5, 0., 0, -1.5, -0.5 };
191 Double_t maxTrackingEffHisto[5] = { 1.5, 2*TMath::Pi(), 20, 1.5, 6.5 };
193 fTrackingEffHisto = new THnSparseF("fTrackingEffHisto","has match:phi:pt:eta:ITSclusters",5,binsTrackingEffHisto,minTrackingEffHisto,maxTrackingEffHisto);
194 fTrackingEffHisto->GetAxis(0)->SetTitle("IsMatching");
195 fTrackingEffHisto->GetAxis(1)->SetTitle("phi (rad)");
196 fTrackingEffHisto->GetAxis(2)->SetTitle("p_{T}");
197 fTrackingEffHisto->GetAxis(3)->SetTitle("eta");
198 fTrackingEffHisto->GetAxis(4)->SetTitle("number of ITS clusters");
199 fTrackingEffHisto->Sumw2();
203 AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
205 AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object");
208 fAnalysisFolder = CreateFolder("folderMatch","Analysis Matching Folder");
210 // save merge status in object
211 fMergeTHnSparseObj = fgMergeTHnSparse;
215 //_____________________________________________________________________________
216 void AliPerformanceMatch::ProcessITSTPC(Int_t iTrack, AliESDEvent *const esdEvent, AliStack* /*const stack*/, AliESDtrack *const esdTrack)
219 // addition to standard analysis - check if ITS stand-alone tracks have a match in the TPC
220 // Origin: A. Kalwait
221 // Modified: J. Otwinowski
222 if(!esdEvent) return;
223 if(!esdTrack) return;
224 // if(!esdFriendTrack) return;
226 // ITS stand alone tracks with SPD layers
227 if(!(esdTrack->GetStatus() & AliESDtrack::kITSpureSA)) return;
228 if(!(esdTrack->GetStatus() & AliESDtrack::kITSrefit)) return;
229 if(esdTrack->GetNcls(0)<4) return;
230 if(!esdTrack->HasPointOnITSLayer(0) && !esdTrack->HasPointOnITSLayer(1)) return;
232 const AliESDVertex* vtxESD = esdEvent->GetPrimaryVertexTracks();
233 // const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTPC();
234 AliESDtrack* tpcTrack2 = NULL;
235 Bool_t hasMatch = kFALSE;
236 for (Int_t jTrack = 0; jTrack < esdEvent->GetNumberOfTracks(); jTrack++) {
237 // loop for all those tracks and check if the corresponding TPC track is found
238 if (jTrack==iTrack) continue;
239 AliESDtrack *trackTPC = esdEvent->GetTrack(jTrack);
240 if (!trackTPC) continue;
241 if (!trackTPC->GetTPCInnerParam()) continue;
242 if(!(trackTPC->GetStatus() & AliESDtrack::kTPCrefit)) continue;
244 // TPC nClust/track after first tracking pass
245 // if(trackTPC->GetTPCNclsIter1()<fCutsRC->GetMinNClustersTPC()) continue;
246 tpcTrack2 = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent, jTrack);
247 if(!tpcTrack2) continue;
248 if(!tpcTrack2->RelateToVertex(vtxESD,esdEvent->GetMagneticField(),100.)) { delete tpcTrack2; continue; }
250 if(!fCutsRC->AcceptTrack(tpcTrack2)) { delete tpcTrack2; continue; }
252 if (TMath::Abs(esdTrack->GetY() - tpcTrack2->GetY()) > 3) { delete tpcTrack2; continue; }
253 if (TMath::Abs(esdTrack->GetSnp() - tpcTrack2->GetSnp()) > 0.2) { delete tpcTrack2; continue; }
254 if (TMath::Abs(esdTrack->GetTgl() - tpcTrack2->GetTgl()) > 0.2) { delete tpcTrack2; continue; }
260 FillHistograms(tpcTrack2,esdTrack,hasMatch);
267 //_____________________________________________________________________________
268 void AliPerformanceMatch::ProcessTPCITS(AliStack* /*const stack*/, AliESDtrack *const esdTrack)
271 // Match TPC and ITS min-bias tracks
272 // at radius between detectors
274 if(!esdTrack) return;
275 // if(!esdFriendTrack) return;
277 Bool_t isTPC = kFALSE;
278 Bool_t isMatch = kFALSE;
281 if(esdTrack->Charge()==0) return;
282 if(!esdTrack->GetTPCInnerParam()) return;
283 if(!(esdTrack->GetStatus()&AliESDtrack::kTPCrefit)) return;
285 if(!fCutsRC->AcceptTrack(esdTrack)) return;
289 if( (esdTrack->GetStatus()&AliESDtrack::kITSrefit))
293 Double_t vecTrackingEff[5] = { isMatch,esdTrack->Phi(), esdTrack->Pt(),esdTrack->Eta(),esdTrack->GetITSclusters(0) };
294 fTrackingEffHisto->Fill(vecTrackingEff);
298 //_____________________________________________________________________________
299 /*void AliPerformanceMatch::ProcessTPCTRD(AliStack* , AliESDtrack *const esdTrack, AliESDfriendTrack *const esdFriendTrack)
304 //_____________________________________________________________________________
305 void AliPerformanceMatch::FillHistograms(AliESDtrack *const refParam, AliESDtrack *const param, Bool_t isRec)
308 // fill performance histograms
309 // (TPC always as reference)
312 if(!refParam) return;
317 // Deltas (dy,dz,dphi,dtheta,dpt)
319 Float_t delta[5] = {0};
321 delta[0] = param->GetY()-refParam->GetY();
322 delta[1] = param->GetZ()-refParam->GetZ();
323 delta[2] = TMath::ATan2(param->Py(),param->Px())-TMath::ATan2(refParam->Py(),refParam->Px());
324 delta[3] = TMath::ATan2(param->Pz(),param->Pt())-TMath::ATan2(refParam->Pz(),refParam->Pt());
325 if(refParam->Pt()) delta[4] = (param->Pt()-refParam->Pt())/refParam->Pt();
328 // Pulls (y,z,snp,tanl,1/pt)
330 Float_t sigma[5] = {0};
331 Float_t pull[5] = {0};
333 sigma[0] = TMath::Sqrt(param->GetSigmaY2()+refParam->GetSigmaY2());
334 sigma[1] = TMath::Sqrt(param->GetSigmaZ2()+refParam->GetSigmaZ2());
335 sigma[2] = TMath::Sqrt(param->GetSigmaSnp2()+refParam->GetSigmaSnp2());
336 sigma[3] = TMath::Sqrt(param->GetSigmaTgl2()+refParam->GetSigmaTgl2());
337 sigma[4] = TMath::Sqrt(param->GetSigma1Pt2()+refParam->GetSigma1Pt2());
338 if(sigma[0]) pull[0] = delta[0] / sigma[0];
339 if(sigma[1]) pull[1] = delta[1] / sigma[1];
340 if(sigma[2]) pull[2] = (param->GetSnp()-refParam->GetSnp()) / sigma[2];
341 if(sigma[3]) pull[3] = (param->GetTgl()-refParam->GetTgl()) / sigma[3];
342 if(sigma[4]) pull[4] = (param->OneOverPt()-refParam->OneOverPt()) / sigma[4];
346 Double_t vResolHisto[9] = {delta[0],delta[1],delta[2],delta[3],delta[4],refParam->Phi(),refParam->Eta(),refParam->Pt(),isRec};
348 fResolHisto->Fill(vResolHisto);
350 Double_t vPullHisto[9] = {pull[0],pull[1],pull[2],pull[3],pull[4],refParam->Phi(),refParam->Eta(),refParam->OneOverPt(),isRec};
352 fPullHisto->Fill(vPullHisto);
355 //_____________________________________________________________________________
356 void AliPerformanceMatch::Exec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent, AliESDfriend */*const esdFriend*/, const Bool_t bUseMC, const Bool_t /*bUseESDfriend*/)
358 // Process comparison information
362 Error("Exec","esdEvent not available");
365 AliHeader* header = 0;
366 AliGenEventHeader* genHeader = 0;
373 Error("Exec","mcEvent not available");
376 // get MC event header
377 header = mcEvent->Header();
379 Error("Exec","Header not available");
383 stack = mcEvent->Stack();
385 Error("Exec","Stack not available");
389 genHeader = header->GenEventHeader();
391 Error("Exec","Could not retrieve genHeader from Header");
394 genHeader->PrimaryVertex(vtxMC);
398 /* if(bUseESDfriend) {
400 Error("Exec","esdFriend not available");
406 if(!bUseMC && GetTriggerClass()) {
407 Bool_t isEventTriggered = esdEvent->IsTriggerClassFired(GetTriggerClass());
408 if(!isEventTriggered) return;
411 // get TPC event vertex
412 const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTPC();
413 if(vtxESD && (vtxESD->GetStatus()<=0)) return;
416 for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++)
418 AliESDtrack *track = esdEvent->GetTrack(iTrack);
421 /*AliESDfriendTrack *friendTrack=0;
423 friendTrack=esdFriend->GetTrack(iTrack);
424 if(!friendTrack) continue;
427 if(GetAnalysisMode() == 0){
428 if(!IsUseTOFBunchCrossing())
429 ProcessTPCITS(stack,track);
431 if( track->GetTOFBunchCrossing(esdEvent->GetMagneticField())==0)
432 ProcessTPCITS(stack,track);
434 /* else if(GetAnalysisMode() == 2) ProcessTPCTRD(stack,track,friendTrack);*/
435 else if(GetAnalysisMode() == 1) {ProcessITSTPC(iTrack,esdEvent,stack,track);
438 printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
445 //_____________________________________________________________________________
446 TH1F* AliPerformanceMatch::MakeResol(TH2F * his, Int_t integ, Bool_t type, Int_t cut){
447 // Create resolution histograms
450 if (!gPad) new TCanvas;
451 hisr = AliTreeDraw::CreateResHistoII(his,&hism,integ,kTRUE,cut);
452 if (type) return hism;
457 //_____________________________________________________________________________
458 void AliPerformanceMatch::Analyse() {
459 // Analyse comparison information and store output histograms
460 // in the folder "folderMatch"
464 TH1::AddDirectory(kFALSE);
469 TObjArray *aFolderObj = new TObjArray;
471 // write results in the folder
472 // TCanvas * c = new TCanvas("Phi resol Tan","Phi resol Tan");
478 if(GetAnalysisMode()==1 || GetAnalysisMode()==2) {
480 fResolHisto->GetAxis(8)->SetRangeUser(1.0,2.0); // only reconstructed
481 fPullHisto->GetAxis(8)->SetRangeUser(1.0,2.0); // only reconstructed
482 for(Int_t i=0; i<5; i++)
484 for(Int_t j=5; j<8; j++)
486 //if(j!=8) fResolHisto->GetAxis(8)->SetRangeUser(-0.9,0.89); // eta window
487 if(j!=6) fResolHisto->GetAxis(6)->SetRangeUser(0.0,1.49); // eta window
488 else fResolHisto->GetAxis(6)->SetRangeUser(-1.5,1.49);
489 fResolHisto->GetAxis(7)->SetRangeUser(0.01,10.); // pt threshold
492 AddProjection(aFolderObj, "match", fResolHisto, i, j, &selString);
495 if(j!=6) fPullHisto->GetAxis(6)->SetRangeUser(0.0,1.49); // eta window
496 else fPullHisto->GetAxis(6)->SetRangeUser(-1.5,1.49); // eta window
497 fPullHisto->GetAxis(7)->SetRangeUser(0.01,10.); // pt threshold
499 AddProjection(aFolderObj, "match", fPullHisto, i, j, &selString);
505 // TPC efficiency wrt ITS
507 if(GetAnalysisMode()==0) {
508 selString = "trackingeff";
509 AddProjection(aFolderObj, "match", fTrackingEffHisto, 0, &selString);
511 for(Int_t i=1; i<5; i++)
513 // all ITS standalone tracks
514 fTrackingEffHisto->GetAxis(0)->SetRange(1,fTrackingEffHisto->GetAxis(0)->GetNbins());
515 selString = "trackingeff_all";
516 AddProjection(aFolderObj, "match", fTrackingEffHisto, i, 3,&selString);
518 // TPC tracks which has matching with TPC
519 fTrackingEffHisto->GetAxis(0)->SetRange(2,2);
520 selString = "trackingeff_tpc";
521 AddProjection(aFolderObj, "match", fTrackingEffHisto, i, 3,&selString);
525 printf("exportToFolder\n");
526 fAnalysisFolder = ExportToFolder(aFolderObj);
528 // delete only TObjArray
529 if(fFolderObj) delete fFolderObj;
530 fFolderObj = aFolderObj;
536 //_____________________________________________________________________________
537 TFolder* AliPerformanceMatch::ExportToFolder(TObjArray * array)
539 // recreate folder avery time and export objects to new one
541 AliPerformanceMatch * comp=this;
542 TFolder *folder = comp->GetAnalysisFolder();
545 TFolder *newFolder = 0;
547 Int_t size = array->GetSize();
550 // get name and title from old folder
551 name = folder->GetName();
552 title = folder->GetTitle();
558 newFolder = CreateFolder(name.Data(),title.Data());
559 newFolder->SetOwner();
561 // add objects to folder
563 newFolder->Add(array->At(i));
571 //_____________________________________________________________________________
572 TFolder* AliPerformanceMatch::CreateFolder(TString name,TString title) {
573 // create folder for analysed histograms
576 folder = new TFolder(name.Data(),title.Data());
581 //_____________________________________________________________________________
582 Long64_t AliPerformanceMatch::Merge(TCollection* const list)
584 // Merge list of objects (needed by PROOF)
592 Bool_t merge = ((fgUseMergeTHnSparse && fgMergeTHnSparse) || (!fgUseMergeTHnSparse && fMergeTHnSparseObj));
594 TIterator* iter = list->MakeIterator();
596 TObjArray* objArrayList = 0;
597 objArrayList = new TObjArray();
599 // collection of generated histograms
601 while((obj = iter->Next()) != 0)
603 AliPerformanceMatch* entry = dynamic_cast<AliPerformanceMatch*>(obj);
604 if (entry == 0) continue;
606 if ((fResolHisto) && (entry->fResolHisto)) { fResolHisto->Add(entry->fResolHisto); }
607 if ((fPullHisto) && (entry->fPullHisto)) { fPullHisto->Add(entry->fPullHisto); }
608 if ((fTrackingEffHisto) && (entry->fTrackingEffHisto)) { fTrackingEffHisto->Add(entry->fTrackingEffHisto); }
610 // the analysisfolder is only merged if present
611 if (entry->fFolderObj) { objArrayList->Add(entry->fFolderObj); }
615 if (fFolderObj) { fFolderObj->Merge(objArrayList); }
616 // to signal that track histos were not merged: reset
617 if (!merge) { fResolHisto->Reset(); fPullHisto->Reset(); fTrackingEffHisto->Reset(); }
619 if (objArrayList) delete objArrayList; objArrayList=0;