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/PWGPP/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(const Char_t* name, const Char_t* title, Int_t analysisMode, Bool_t hptGenerator):
73 AliPerformanceObject(name,title),
89 SetAnalysisMode(analysisMode);
90 SetHptGenerator(hptGenerator);
95 //_____________________________________________________________________________
96 AliPerformanceMatch::~AliPerformanceMatch()
100 if(fResolHisto) delete fResolHisto; fResolHisto=0;
101 if(fPullHisto) delete fPullHisto; fPullHisto=0;
102 if(fTrackingEffHisto) delete fTrackingEffHisto; fTrackingEffHisto = 0x0;
103 if(fTPCConstrain) delete fTPCConstrain; fTPCConstrain = 0x0;
105 if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
106 if(fFolderObj) delete fFolderObj; fFolderObj=0;
109 //_____________________________________________________________________________
110 void AliPerformanceMatch::Init(){
113 // Make performance histogrms
118 Double_t ptMin = 1.e-2, ptMax = 20.;
120 Double_t *binsPt = 0;
122 if (IsHptGenerator()) {
125 binsPt = CreateLogAxis(nPtBins,ptMin,ptMax);
128 // Double_t yMin = -0.02, yMax = 0.02;
129 // Double_t zMin = -12.0, zMax = 12.0;
130 // if(GetAnalysisMode() == 0 || GetAnalysisMode() == 1) {
131 // yMin = -100.; yMax = 100.;
132 // zMin = -100.; zMax = 100.;
136 //init ITS TPC Mactching
138 if(GetAnalysisMode()==1){
139 // res_y:res_z:res_phi,res_lambda:res_pt:y:z:eta:phi:pt:isRec
140 Int_t binsResolHisto[9]={100,100,100,100,100,90,30,nPtBins,2};
141 Double_t minResolHisto[9]={-1.,-1.,-0.03,-0.03,-0.2, 0., -1.5, ptMin,0};
142 Double_t maxResolHisto[9]={ 1., 1., 0.03, 0.03, 0.2, 2.*TMath::Pi(), 1.5, ptMax,2};
144 fResolHisto = new THnSparseF("fResolHisto","res_y:res_z:res_phi:res_lambda:res_pt:phi:eta:pt:isRec",9,binsResolHisto,minResolHisto,maxResolHisto);
145 fResolHisto->SetBinEdges(7,binsPt);
147 fResolHisto->GetAxis(0)->SetTitle("y-y_{ref} (cm)");
148 fResolHisto->GetAxis(1)->SetTitle("z-z_{ref} (cm)");
149 fResolHisto->GetAxis(2)->SetTitle("#phi-#phi_{ref} (rad)");
150 fResolHisto->GetAxis(3)->SetTitle("#lambda-#lambda_{ref} (rad)");
151 fResolHisto->GetAxis(4)->SetTitle("(p_{T}/p_{Tref}-1)");
152 fResolHisto->GetAxis(5)->SetTitle("#phi_{ref} (rad)");
153 fResolHisto->GetAxis(6)->SetTitle("#eta_{ref}");
154 fResolHisto->GetAxis(7)->SetTitle("p_{Tref} (GeV/c)");
155 fResolHisto->GetAxis(8)->SetTitle("isReconstructed");
156 fResolHisto->Sumw2();
158 //pull_y:pull_z:pull_snp:pull_tgl:pull_1pt:y:z:snp:tgl:1pt:isRec
159 Int_t binsPullHisto[9]={100,100,100,100,100,90,30,nPtBins,2};
160 Double_t minPullHisto[9]={-5.,-5.,-5.,-5.,-5., 0.,-1.5, ptMin,0};
161 Double_t maxPullHisto[9]={ 5., 5., 5., 5., 5., 2.*TMath::Pi(), 1.5, ptMax,2};
162 fPullHisto = new THnSparseF("fPullHisto","pull_y:pull_z:pull_snp:pull_tgl:pull_1pt:phi:eta:1pt:isRec",9,binsPullHisto,minPullHisto,maxPullHisto);
164 fPullHisto->GetAxis(0)->SetTitle("(y-y_{ref})/#sigma");
165 fPullHisto->GetAxis(1)->SetTitle("(z-z_{ref})/#sigma");
166 fPullHisto->GetAxis(2)->SetTitle("(sin#phi-sin#phi_{ref})/#sigma");
167 fPullHisto->GetAxis(3)->SetTitle("(tan#lambda-tan#lambda_{ref})/#sigma");
168 fPullHisto->GetAxis(4)->SetTitle("(p_{Tref}/p_{T}-1)/#sigma");
169 fPullHisto->GetAxis(5)->SetTitle("#phi_{ref} (rad)");
170 fPullHisto->GetAxis(6)->SetTitle("#eta_{ref}");
171 fPullHisto->GetAxis(7)->SetTitle("1/p_{Tref} (GeV/c)^{-1}");
172 fPullHisto->GetAxis(8)->SetTitle("isReconstructed");
179 if(GetAnalysisMode()==0){
180 // -> has match:y:z:snp:tgl:phi:pt:ITSclusters
181 Int_t binsTrackingEffHisto[5] = { 2, 90, nPtBins, 30, 7 };
182 Double_t minTrackingEffHisto[5] = {-0.5, 0., ptMin, -1.5, -0.5 };
183 Double_t maxTrackingEffHisto[5] = { 1.5, 2*TMath::Pi(), ptMax, 1.5, 6.5 };
185 fTrackingEffHisto = new THnSparseF("fTrackingEffHisto","has match:phi:pt:eta:ITSclusters",5,binsTrackingEffHisto,minTrackingEffHisto,maxTrackingEffHisto);
186 fTrackingEffHisto->SetBinEdges(2,binsPt);
187 fTrackingEffHisto->GetAxis(0)->SetTitle("IsMatching");
188 fTrackingEffHisto->GetAxis(1)->SetTitle("phi (rad)");
189 fTrackingEffHisto->GetAxis(2)->SetTitle("p_{T}");
190 fTrackingEffHisto->GetAxis(3)->SetTitle("eta");
191 fTrackingEffHisto->GetAxis(4)->SetTitle("number of ITS clusters");
192 fTrackingEffHisto->Sumw2();
196 //TPC constraining to vertex
198 if(GetAnalysisMode()==2){
199 //initilization of fTPCConstrain
200 Int_t binsTPCConstrain[4] = {100, 90, nPtBins, 30};
201 Double_t minTPCConstrain[4] = {-5, 0, ptMin, -1.5};
202 Double_t maxTPCConstrain[4] = {5, 2*TMath::Pi(), ptMax, 1.5};
204 fTPCConstrain = new THnSparseF("fTPCConstrain","pull_phi:phi:pt:eta",4, binsTPCConstrain,minTPCConstrain,maxTPCConstrain);
205 fTPCConstrain->SetBinEdges(2,binsPt);
206 fTPCConstrain->GetAxis(0)->SetTitle("(#phi-#phi_{ref})/#sigma");
207 fTPCConstrain->GetAxis(1)->SetTitle("phi (rad)");
208 fTPCConstrain->GetAxis(2)->SetTitle("p_{T}");
209 fTPCConstrain->GetAxis(3)->SetTitle("eta");
210 fTPCConstrain->Sumw2();
215 AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
217 AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object");
220 fAnalysisFolder = CreateFolder("folderMatch","Analysis Matching Folder");
222 // save merge status in object
223 fMergeTHnSparseObj = fgMergeTHnSparse;
229 //_____________________________________________________________________________
230 void AliPerformanceMatch::ProcessITSTPC(Int_t iTrack, AliESDEvent *const esdEvent, AliStack* /*const stack*/, AliESDtrack *const esdTrack)
233 // addition to standard analysis - check if ITS stand-alone tracks have a match in the TPC
234 // Origin: A. Kalwait
235 // Modified: J. Otwinowski
236 if(!esdEvent) return;
237 if(!esdTrack) return;
238 // if(!esdFriendTrack) return;
240 // ITS stand alone tracks with SPD layers
241 if(!(esdTrack->GetStatus() & AliESDtrack::kITSpureSA)) return;
242 if(!(esdTrack->GetStatus() & AliESDtrack::kITSrefit)) return;
243 if(esdTrack->GetNcls(0)<4) return;
244 if(!esdTrack->HasPointOnITSLayer(0) && !esdTrack->HasPointOnITSLayer(1)) return;
246 const AliESDVertex* vtxESD = esdEvent->GetPrimaryVertexTracks();
247 // const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTPC();
248 AliESDtrack* tpcTrack2 = NULL;
249 Bool_t hasMatch = kFALSE;
250 for (Int_t jTrack = 0; jTrack < esdEvent->GetNumberOfTracks(); jTrack++) {
251 // loop for all those tracks and check if the corresponding TPC track is found
252 if (jTrack==iTrack) continue;
253 AliESDtrack *trackTPC = esdEvent->GetTrack(jTrack);
254 if (!trackTPC) continue;
255 if (!trackTPC->GetTPCInnerParam()) continue;
256 if(!(trackTPC->GetStatus() & AliESDtrack::kTPCrefit)) continue;
258 // TPC nClust/track after first tracking pass
259 // if(trackTPC->GetTPCNclsIter1()<fCutsRC->GetMinNClustersTPC()) continue;
260 tpcTrack2 = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent, jTrack);
261 if(!tpcTrack2) continue;
262 if(!tpcTrack2->RelateToVertex(vtxESD,esdEvent->GetMagneticField(),100.)) { delete tpcTrack2; tpcTrack2=0; continue; }
264 if(!fCutsRC->AcceptTrack(tpcTrack2)) { delete tpcTrack2; tpcTrack2=0; continue; }
266 if (TMath::Abs(esdTrack->GetY() - tpcTrack2->GetY()) > 3) { delete tpcTrack2; tpcTrack2=0; continue; }
267 if (TMath::Abs(esdTrack->GetSnp() - tpcTrack2->GetSnp()) > 0.2) { delete tpcTrack2; tpcTrack2=0; continue; }
268 if (TMath::Abs(esdTrack->GetTgl() - tpcTrack2->GetTgl()) > 0.2) { delete tpcTrack2; tpcTrack2=0; continue; }
274 FillHistograms(tpcTrack2,esdTrack,hasMatch);
278 //_____________________________________________________________________________
279 void AliPerformanceMatch::ProcessTPCITS(AliStack* /*const stack*/, AliESDtrack *const esdTrack)
282 // Match TPC and ITS min-bias tracks
283 // at radius between detectors
285 if(!esdTrack) return;
286 // if(!esdFriendTrack) return;
288 Bool_t isTPC = kFALSE;
289 Bool_t isMatch = kFALSE;
292 if(esdTrack->Charge()==0) return;
293 if(!esdTrack->GetTPCInnerParam()) return;
294 if(!(esdTrack->GetStatus()&AliESDtrack::kTPCrefit)) return;
296 if(!fCutsRC->AcceptTrack(esdTrack)) return;
300 if( (esdTrack->GetStatus()&AliESDtrack::kITSrefit))
304 Double_t vecTrackingEff[5] = { static_cast<Double_t>(isMatch),esdTrack->Phi(), esdTrack->Pt(),esdTrack->Eta(),static_cast<Double_t>(esdTrack->GetITSclusters(0)) };
305 fTrackingEffHisto->Fill(vecTrackingEff);
309 //_____________________________________________________________________________
310 /*void AliPerformanceMatch::ProcessTPCTRD(AliStack* , AliESDtrack *const esdTrack, AliESDfriendTrack *const esdFriendTrack)
315 void AliPerformanceMatch::ProcessTPCConstrain(AliStack* /*const stack*/, AliESDEvent *const esdEvent, AliESDtrack *const esdTrack){
317 // Contrain TPC inner track to the vertex
318 // then compare to the global tracks
320 if(!esdTrack) return;
321 if(esdTrack->Charge()==0) return;
323 if(!esdTrack->GetTPCInnerParam()) return;
324 if(!(esdTrack->GetStatus()&AliESDtrack::kITSrefit)) return;
325 if(!(esdTrack->GetStatus()&AliESDtrack::kTPCrefit)) return;
326 if(!fCutsRC->AcceptTrack(esdTrack)) return;
328 Double_t x[3]; esdTrack->GetXYZ(x);
329 Double_t b[3]; AliTracker::GetBxByBz(x,b);
330 Bool_t isOK = kFALSE;
331 const AliESDVertex* vtxESD = esdEvent->GetPrimaryVertexTracks();
333 AliExternalTrackParam * TPCinner = (AliExternalTrackParam *)esdTrack->GetTPCInnerParam();
334 if(!TPCinner) return;
336 // isOK = TPCinner->Rotate(esdTrack->GetAlpha());
337 //isOK = TPCinner->PropagateTo(esdTrack->GetX(),esdEvent->GetMagneticField());
339 AliExternalTrackParam * TPCinnerC = new AliExternalTrackParam(*TPCinner);
341 isOK = TPCinnerC->ConstrainToVertex(vtxESD, b);
343 // transform to the track reference frame
344 isOK = TPCinnerC->Rotate(esdTrack->GetAlpha());
345 isOK = TPCinnerC->PropagateTo(esdTrack->GetX(),esdEvent->GetMagneticField());
349 Double_t sigmaPhi=0,deltaPhi=0,pullPhi=0;
350 deltaPhi = TPCinnerC->GetSnp() - esdTrack->GetSnp();
351 sigmaPhi = TMath::Sqrt(esdTrack->GetSigmaSnp2()+TPCinnerC->GetSigmaSnp2());
353 pullPhi = deltaPhi/sigmaPhi;
355 Double_t vTPCConstrain[4] = {pullPhi,esdTrack->Phi(),esdTrack->Pt(),esdTrack->Eta()};
356 fTPCConstrain->Fill(vTPCConstrain);
363 //_____________________________________________________________________________
364 void AliPerformanceMatch::FillHistograms(AliESDtrack *const refParam, AliESDtrack *const param, Bool_t isRec)
367 // fill performance histograms
368 // (TPC always as reference)
371 if(!refParam) return;
376 // Deltas (dy,dz,dphi,dtheta,dpt)
378 Float_t delta[5] = {0};
380 delta[0] = param->GetY()-refParam->GetY();
381 delta[1] = param->GetZ()-refParam->GetZ();
382 delta[2] = TMath::ATan2(param->Py(),param->Px())-TMath::ATan2(refParam->Py(),refParam->Px());
383 delta[3] = TMath::ATan2(param->Pz(),param->Pt())-TMath::ATan2(refParam->Pz(),refParam->Pt());
384 if(refParam->Pt()) delta[4] = (param->Pt()-refParam->Pt())/refParam->Pt();
387 // Pulls (y,z,snp,tanl,1/pt)
389 Float_t sigma[5] = {0};
390 Float_t pull[5] = {0};
392 sigma[0] = TMath::Sqrt(param->GetSigmaY2()+refParam->GetSigmaY2());
393 sigma[1] = TMath::Sqrt(param->GetSigmaZ2()+refParam->GetSigmaZ2());
394 sigma[2] = TMath::Sqrt(param->GetSigmaSnp2()+refParam->GetSigmaSnp2());
395 sigma[3] = TMath::Sqrt(param->GetSigmaTgl2()+refParam->GetSigmaTgl2());
396 sigma[4] = TMath::Sqrt(param->GetSigma1Pt2()+refParam->GetSigma1Pt2());
397 if(sigma[0]) pull[0] = delta[0] / sigma[0];
398 if(sigma[1]) pull[1] = delta[1] / sigma[1];
399 if(sigma[2]) pull[2] = (param->GetSnp()-refParam->GetSnp()) / sigma[2];
400 if(sigma[3]) pull[3] = (param->GetTgl()-refParam->GetTgl()) / sigma[3];
401 if(sigma[4]) pull[4] = (param->OneOverPt()-refParam->OneOverPt()) / sigma[4];
405 Double_t vResolHisto[9] = {delta[0],delta[1],delta[2],delta[3],delta[4],refParam->Phi(),refParam->Eta(),refParam->Pt(),static_cast<Double_t>(isRec)};
407 fResolHisto->Fill(vResolHisto);
409 Double_t vPullHisto[9] = {pull[0],pull[1],pull[2],pull[3],pull[4],refParam->Phi(),refParam->Eta(),refParam->OneOverPt(),static_cast<Double_t>(isRec)};
411 fPullHisto->Fill(vPullHisto);
414 //_____________________________________________________________________________
415 void AliPerformanceMatch::Exec(AliMCEvent* const mcEvent, AliESDEvent *const esdEvent, AliESDfriend */*const esdFriend*/, const Bool_t bUseMC, const Bool_t /*bUseESDfriend*/)
417 // Process comparison information
421 Error("Exec","esdEvent not available");
424 AliHeader* header = 0;
425 AliGenEventHeader* genHeader = 0;
432 Error("Exec","mcEvent not available");
435 // get MC event header
436 header = mcEvent->Header();
438 Error("Exec","Header not available");
442 stack = mcEvent->Stack();
444 Error("Exec","Stack not available");
448 genHeader = header->GenEventHeader();
450 Error("Exec","Could not retrieve genHeader from Header");
453 genHeader->PrimaryVertex(vtxMC);
457 /* if(bUseESDfriend) {
459 Error("Exec","esdFriend not available");
465 if(!bUseMC && GetTriggerClass()) {
466 Bool_t isEventTriggered = esdEvent->IsTriggerClassFired(GetTriggerClass());
467 if(!isEventTriggered) return;
470 // get TPC event vertex
471 const AliESDVertex *vtxESD = esdEvent->GetPrimaryVertexTPC();
472 if(vtxESD && (vtxESD->GetStatus()<=0)) return;
475 for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++)
477 AliESDtrack *track = esdEvent->GetTrack(iTrack);
480 /*AliESDfriendTrack *friendTrack=0;
482 friendTrack=esdFriend->GetTrack(iTrack);
483 if(!friendTrack) continue;
486 if(GetAnalysisMode() == 0){
487 if(!IsUseTOFBunchCrossing()){
488 ProcessTPCITS(stack,track);
491 if( track->GetTOFBunchCrossing(esdEvent->GetMagneticField())==0) {
492 ProcessTPCITS(stack,track);
495 /* else if(GetAnalysisMode() == 2) ProcessTPCTRD(stack,track,friendTrack);*/
496 else if(GetAnalysisMode() == 1) {ProcessITSTPC(iTrack,esdEvent,stack,track);}
497 else if(GetAnalysisMode() == 2){
498 if(!IsUseTOFBunchCrossing()){
499 ProcessTPCConstrain(stack,esdEvent,track);
502 if( track->GetTOFBunchCrossing(esdEvent->GetMagneticField())==0) {
503 ProcessTPCConstrain(stack,esdEvent,track);
507 printf("ERROR: AnalysisMode %d \n",fAnalysisMode);
514 //_____________________________________________________________________________
515 TH1F* AliPerformanceMatch::MakeResol(TH2F * his, Int_t integ, Bool_t type, Int_t cut){
516 // Create resolution histograms
519 if (!gPad) new TCanvas;
520 hisr = AliTreeDraw::CreateResHistoII(his,&hism,integ,kTRUE,cut);
521 if (type) return hism;
526 //_____________________________________________________________________________
527 void AliPerformanceMatch::Analyse() {
528 // Analyse comparison information and store output histograms
529 // in the folder "folderMatch"
533 TH1::AddDirectory(kFALSE);
538 TObjArray *aFolderObj = new TObjArray;
540 // write results in the folder
541 // TCanvas * c = new TCanvas("Phi resol Tan","Phi resol Tan");
547 if(GetAnalysisMode()==1) {
549 fResolHisto->GetAxis(8)->SetRangeUser(1.0,2.0); // only reconstructed
550 fPullHisto->GetAxis(8)->SetRangeUser(1.0,2.0); // only reconstructed
551 for(Int_t i=0; i<5; i++)
553 for(Int_t j=5; j<8; j++)
555 //if(j!=8) fResolHisto->GetAxis(8)->SetRangeUser(-0.9,0.89); // eta window
556 if(j!=6) fResolHisto->GetAxis(6)->SetRangeUser(0.0,1.49); // eta window
557 else fResolHisto->GetAxis(6)->SetRangeUser(-1.5,1.49);
558 fResolHisto->GetAxis(7)->SetRangeUser(0.01,10.); // pt threshold
561 AddProjection(aFolderObj, "match", fResolHisto, i, j, &selString);
564 if(j!=6) fPullHisto->GetAxis(6)->SetRangeUser(0.0,1.49); // eta window
565 else fPullHisto->GetAxis(6)->SetRangeUser(-1.5,1.49); // eta window
566 fPullHisto->GetAxis(7)->SetRangeUser(0.01,10.); // pt threshold
568 AddProjection(aFolderObj, "match", fPullHisto, i, j, &selString);
574 // TPC efficiency wrt ITS
576 if(GetAnalysisMode()==0) {
577 selString = "trackingeff";
578 AddProjection(aFolderObj, "match", fTrackingEffHisto, 0, &selString);
580 for(Int_t i=1; i<5; i++)
582 // all ITS standalone tracks
583 fTrackingEffHisto->GetAxis(0)->SetRange(1,fTrackingEffHisto->GetAxis(0)->GetNbins());
584 selString = "trackingeff_all";
585 AddProjection(aFolderObj, "match", fTrackingEffHisto, i, 3,&selString);
587 // TPC tracks which has matching with TPC
588 fTrackingEffHisto->GetAxis(0)->SetRange(2,2);
589 selString = "trackingeff_tpc";
590 AddProjection(aFolderObj, "match", fTrackingEffHisto, i, 3,&selString);
595 //TPC constrained to vertex
597 if(GetAnalysisMode()==2) {
599 // for(Int_t i=1; i<4; i++)
600 AddProjection(aFolderObj, "constrain", fTPCConstrain, 0, 1, 2, &selString);
601 AddProjection(aFolderObj, "constrain", fTPCConstrain, 0, 1, 3, &selString);
602 AddProjection(aFolderObj, "constrain", fTPCConstrain, 0, 2, 3, &selString);
605 printf("exportToFolder\n");
606 fAnalysisFolder = ExportToFolder(aFolderObj);
608 // delete only TObjArray
609 if(fFolderObj) delete fFolderObj;
610 fFolderObj = aFolderObj;
616 //_____________________________________________________________________________
617 TFolder* AliPerformanceMatch::ExportToFolder(TObjArray * array)
619 // recreate folder avery time and export objects to new one
621 AliPerformanceMatch * comp=this;
622 TFolder *folder = comp->GetAnalysisFolder();
625 TFolder *newFolder = 0;
627 Int_t size = array->GetSize();
630 // get name and title from old folder
631 name = folder->GetName();
632 title = folder->GetTitle();
638 newFolder = CreateFolder(name.Data(),title.Data());
639 newFolder->SetOwner();
641 // add objects to folder
643 newFolder->Add(array->At(i));
651 //_____________________________________________________________________________
652 TFolder* AliPerformanceMatch::CreateFolder(TString name,TString title) {
653 // create folder for analysed histograms
656 folder = new TFolder(name.Data(),title.Data());
661 //_____________________________________________________________________________
662 Long64_t AliPerformanceMatch::Merge(TCollection* const list)
664 // Merge list of objects (needed by PROOF)
672 Bool_t merge = ((fgUseMergeTHnSparse && fgMergeTHnSparse) || (!fgUseMergeTHnSparse && fMergeTHnSparseObj));
674 TIterator* iter = list->MakeIterator();
676 TObjArray* objArrayList = 0;
677 objArrayList = new TObjArray();
679 // collection of generated histograms
681 while((obj = iter->Next()) != 0)
683 AliPerformanceMatch* entry = dynamic_cast<AliPerformanceMatch*>(obj);
684 if (entry == 0) continue;
686 if ((fResolHisto) && (entry->fResolHisto)) { fResolHisto->Add(entry->fResolHisto); }
687 if ((fPullHisto) && (entry->fPullHisto)) { fPullHisto->Add(entry->fPullHisto); }
688 if ((fTrackingEffHisto) && (entry->fTrackingEffHisto)) { fTrackingEffHisto->Add(entry->fTrackingEffHisto); }
690 if ((fTPCConstrain) && (entry->fTPCConstrain)) { fTPCConstrain->Add(entry->fTPCConstrain); }
692 // the analysisfolder is only merged if present
693 if (entry->fFolderObj) { objArrayList->Add(entry->fFolderObj); }
697 if (fFolderObj) { fFolderObj->Merge(objArrayList); }
698 // to signal that track histos were not merged: reset
699 if (!merge) { fResolHisto->Reset(); fPullHisto->Reset(); fTrackingEffHisto->Reset(); fTPCConstrain->Reset(); }
701 if (objArrayList) delete objArrayList; objArrayList=0;