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 // static dNdPt helper functions
18 // basic functionality to select events and tracks
21 // Origin: Jan Fiete Grosse-Oetringhaus
22 // Modified and Extended: Jacek Otwinowski 19/11/2009
33 #include <AliHeader.h>
37 #include <AliESDEvent.h>
38 #include <AliMCEvent.h>
39 #include <AliESDVertex.h>
40 #include <AliVertexerTracks.h>
41 #include <AliMathBase.h>
42 #include <AliESDtrackCuts.h>
43 #include "dNdPt/AlidNdPtEventCuts.h"
44 #include "dNdPt/AlidNdPtAcceptanceCuts.h"
45 #include "dNdPt/AlidNdPtHelper.h"
47 //____________________________________________________________________
48 ClassImp(AlidNdPtHelper)
50 //____________________________________________________________________
51 const AliESDVertex* AlidNdPtHelper::GetVertex(AliESDEvent* const aEsd, AlidNdPtEventCuts *const evtCuts, AlidNdPtAcceptanceCuts *const accCuts, AliESDtrackCuts *const trackCuts, AnalysisMode analysisMode, Bool_t debug, Bool_t bRedoTPC, Bool_t bUseMeanVertex)
53 // Get the vertex from the ESD and returns it if the vertex is valid
55 // Second argument decides which vertex is used (this selects
56 // also the quality criteria that are applied)
60 ::Error("AlidNdPtHelper::GetVertex()","esd event is NULL");
64 if(!evtCuts || !accCuts || !trackCuts)
66 ::Error("AlidNdPtHelper::GetVertex()","cuts not available");
70 const AliESDVertex* vertex = 0;
71 AliESDVertex *initVertex = 0;
72 if (analysisMode == kSPD || analysisMode == kTPCITS || analysisMode == kTPCSPDvtx || analysisMode == kTPCSPDvtxUpdate)
74 vertex = aEsd->GetPrimaryVertexSPD();
76 Printf("AlidNdPtHelper::GetVertex: Returning SPD vertex");
78 else if (analysisMode == kTPC)
82 Double_t kBz = aEsd->GetMagneticField();
83 AliVertexerTracks vertexer(kBz);
86 Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
87 Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
88 initVertex = new AliESDVertex(pos,err);
89 vertexer.SetVtxStart(initVertex);
90 vertexer.SetConstraintOn();
93 Double_t maxDCAr = accCuts->GetMaxDCAr();
94 Double_t maxDCAz = accCuts->GetMaxDCAz();
95 Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
97 //vertexer.SetTPCMode(Double_t dcacut=0.1, Double_t dcacutIter0=1.0, Double_t maxd0z0=5.0, Int_t minCls=10, Int_t mintrks=1, Double_t nsigma=3., Double_t mindetfitter=0.1, Double_t maxtgl=1.5, Double_t fidR=3., Double_t fidZ=30., Int_t finderAlgo=1, Int_t finderAlgoIter0=4);
98 vertexer.SetTPCMode(0.1,1.0,5.0,minTPCClust,1,3.,0.1,2.0,maxDCAr,maxDCAz,1,4);
100 // TPC track preselection
101 Int_t ntracks = aEsd->GetNumberOfTracks();
102 TObjArray array(ntracks);
103 UShort_t *id = new UShort_t[ntracks];
106 for (Int_t i=0;i <ntracks; i++) {
107 AliESDtrack *t = aEsd->GetTrack(i);
109 if (t->Charge() == 0) continue;
110 if (!t->GetTPCInnerParam()) continue;
111 if (t->GetTPCNcls()<vertexer.GetMinClusters()) continue;
112 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
114 array.AddLast(tpcTrack);
115 id[count] = (UShort_t)t->GetID();
119 AliESDVertex *vTPC = vertexer.VertexForSelectedTracks(&array,id, kTRUE, kTRUE, bUseMeanVertex);
121 // set recreated TPC vertex
122 aEsd->SetPrimaryVertexTPC(vTPC);
124 for (Int_t i=0; i<aEsd->GetNumberOfTracks(); i++) {
125 AliESDtrack *t = aEsd->GetTrack(i);
126 t->RelateToVertexTPC(vTPC, kBz, kVeryBig);
131 delete [] id; id=NULL;
134 vertex = aEsd->GetPrimaryVertexTPC();
136 Printf("AlidNdPtHelper::GetVertex: Returning vertex from tracks");
139 Printf("AlidNdPtHelper::GetVertex: ERROR: Invalid second argument %d", analysisMode);
143 Printf("AlidNdPtHelper::GetVertex: No vertex found in ESD");
149 Printf("AlidNdPtHelper::GetVertex: Returning valid vertex: %s", vertex->GetTitle());
153 if(initVertex) delete initVertex; initVertex=NULL;
157 //____________________________________________________________________
158 Bool_t AlidNdPtHelper::TestRecVertex(const AliESDVertex* vertex, AnalysisMode analysisMode, Bool_t debug)
160 // Checks if a vertex meets the needed quality criteria
161 if(!vertex) return kFALSE;
163 Float_t requiredZResolution = -1;
164 if (analysisMode == kSPD || analysisMode == kTPCITS || analysisMode == kTPCSPDvtx || analysisMode == kTPCSPDvtxUpdate)
166 requiredZResolution = 0.1;
168 else if (analysisMode == kTPC)
169 requiredZResolution = 10.;
171 // check Ncontributors
172 if (vertex->GetNContributors() <= 0) {
174 Printf("AlidNdPtHelper::GetVertex: NContributors() <= 0: %d",vertex->GetNContributors());
175 Printf("AlidNdPtHelper::GetVertex: NIndices(): %d",vertex->GetNIndices());
182 Double_t zRes = vertex->GetZRes();
183 if ((zRes-0.000000001) < 0.0) {
184 Printf("AlidNdPtHelper::GetVertex: UNEXPECTED: resolution is 0.");
188 if (zRes > requiredZResolution) {
190 Printf("AlidNdPtHelper::TestVertex: Resolution too poor %f (required: %f", zRes, requiredZResolution);
197 //____________________________________________________________________
198 Bool_t AlidNdPtHelper::IsPrimaryParticle(AliStack* const stack, Int_t idx, ParticleMode particleMode)
200 // check primary particles
201 // depending on the particle mode
203 if(!stack) return kFALSE;
205 TParticle* particle = stack->Particle(idx);
206 if (!particle) return kFALSE;
208 // only charged particles
209 Double_t charge = particle->GetPDG()->Charge()/3.;
210 if (TMath::Abs(charge) < 0.001) return kFALSE;
212 Int_t pdg = TMath::Abs(particle->GetPdgCode());
215 Bool_t prim = stack->IsPhysicalPrimary(idx);
217 if(particleMode==kMCPion) {
218 if(prim && pdg==kPiPlus) return kTRUE;
222 if (particleMode==kMCKaon) {
223 if(prim && pdg==kKPlus) return kTRUE;
227 if (particleMode==kMCProton) {
228 if(prim && pdg==kProton) return kTRUE;
235 //____________________________________________________________________
237 Bool_t AlidNdPtHelper::IsCosmicTrack(TObjArray *const allChargedTracks, AliESDtrack *const track1, Int_t trackIdx, AlidNdPtAcceptanceCuts *const accCuts, AliESDtrackCuts *const trackCuts)
240 // check cosmic tracks
242 if(!allChargedTracks) return kFALSE;
243 if(!track1) return kFALSE;
244 if(!accCuts) return kFALSE;
245 if(!trackCuts) return kFALSE;
247 Int_t entries = allChargedTracks->GetEntries();
248 for(Int_t i=0; i<entries;++i)
251 // exclude the same tracks
253 if(i == trackIdx) continue;
255 AliESDtrack *track2 = (AliESDtrack*)allChargedTracks->At(i);
256 if(!track2) continue;
257 if(track2->Charge()==0) continue;
259 if(track1->Pt() > 6. && track2->Pt() > 6. && (track1->Charge() + track2->Charge()) == 0 )
261 printf("track1->Theta() %f, track1->Eta() %f, track1->Phi() %f, track1->Charge() %d \n", track1->Theta(), track1->Eta(), track1->Phi(), track1->Charge());
262 printf("track2->Theta() %f, track2->Eta() %f, track2->Phi() %f, track2->Charge() %d \n", track2->Theta(), track2->Eta(), track2->Phi(), track2->Charge());
264 printf("deta %f, dphi %f, dq %d \n", track1->Eta()-track2->Eta(), track1->Phi()-track2->Phi(), track1->Charge()+track2->Charge());
269 // cosmic tracks in TPC
271 //if( TMath::Abs( track1->Theta() - track2->Theta() ) < 0.004 &&
272 // ((TMath::Abs(track1->Phi()-track2->Phi()) - TMath::Pi() )<0.004) &&
273 if( (track1->Pt()-track2->Pt()) < 0.1 && track1->Pt() > 4.0 && (track1->Charge()+track2->Charge()) == 0 )
275 //printf("COSMIC candidate \n");
276 printf("track1->Theta() %f, track1->Eta() %f, track1->Phi() %f, track1->Charge() %d \n", track1->Theta(), track1->Eta(), track1->Phi(), track1->Charge());
277 printf("track2->Theta() %f, track2->Eta() %f, track2->Phi() %f, track2->Charge() %d \n", track2->Theta(), track2->Eta(), track2->Phi(), track2->Charge());
278 printf("dtheta %f, deta %f, dphi %f, dq %d \n", track1->Theta()-track2->Theta(), track1->Eta()-track2->Eta(), track1->Phi()-track2->Phi(), track1->Charge()+track2->Charge());
287 //____________________________________________________________________
288 Bool_t AlidNdPtHelper::IsCosmicTrack(AliESDtrack *const track1, AliESDtrack *const track2)
291 // check cosmic tracks
293 if(!track1) return kFALSE;
294 if(!track2) return kFALSE;
297 if(track1->Pt() > 6. && track2->Pt() > 6. && (track1->Charge() + track2->Charge()) == 0 )
299 printf("track1->Theta() %f, track1->Eta() %f, track1->Phi() %f, track1->Charge() %d \n", track1->Theta(), track1->Eta(), track1->Phi(), track1->Charge());
300 printf("track2->Theta() %f, track2->Eta() %f, track2->Phi() %f, track2->Charge() %d \n", track2->Theta(), track2->Eta(), track2->Phi(), track2->Charge());
302 printf("deta %f, dphi %f, dq %d \n", track1->Eta()-track2->Eta(), track1->Phi()-track2->Phi(), track1->Charge()+track2->Charge());
308 // cosmic tracks in TPC
310 //if( TMath::Abs( track1->Theta() - track2->Theta() ) < 0.004 &&
311 // ((TMath::Abs(track1->Phi()-track2->Phi()) - TMath::Pi() )<0.004) &&
313 if( (track1->Pt()-track2->Pt()) < 0.2 && track1->Pt() > 3.0 &&
314 (track1->Charge()+track2->Charge()) == 0 )
315 //if( track1->Pt() > 6. || track2->Pt() > 6. )
317 //printf("COSMIC candidate \n");
319 printf("track1->Theta() %f, track1->Eta() %f, track1->Phi() %f, track1->Charge() %d \n", track1->Theta(), track1->Eta(), track1->Phi(), track1->Charge());
320 printf("track2->Theta() %f, track2->Eta() %f, track2->Phi() %f, track2->Charge() %d \n", track2->Theta(), track2->Eta(), track2->Phi(), track2->Charge());
321 printf("dtheta %f, deta %f, dphi %f, dq %d \n", track1->Theta()-track2->Theta(), track1->Eta()-track2->Eta(), track1->Phi()-track2->Phi(), track1->Charge()+track2->Charge());
329 //____________________________________________________________________
330 void AlidNdPtHelper::PrintConf(AnalysisMode analysisMode, AliTriggerAnalysis::Trigger trigger)
333 // Prints the given configuration
336 TString str(">>>> Running with ");
338 switch (analysisMode)
340 case kInvalid: str += "invalid setting"; break;
341 case kSPD : str += "SPD-only"; break;
342 case kTPC : str += "TPC-only"; break;
343 case kTPCITS : str += "Global tracking"; break;
344 case kTPCSPDvtx : str += "TPC tracking + SPD event vertex"; break;
345 case kTPCSPDvtxUpdate : str += "TPC tracks updated with SPD event vertex point"; break;
346 case kMCRec : str += "TPC tracking + Replace rec. with MC values"; break;
348 str += " and trigger ";
350 str += AliTriggerAnalysis::GetTriggerName(trigger);
354 Printf("%s", str.Data());
357 //____________________________________________________________________
358 Int_t AlidNdPtHelper::ConvertPdgToPid(TParticle *particle) {
360 // Convert Abs(pdg) to pid
361 // (0 - e, 1 - muons, 2 - pions, 3 - kaons, 4 - protons, 5 -all rest)
365 if (TMath::Abs(particle->GetPdgCode()) == kElectron) { pid = 0; }
366 else if (TMath::Abs(particle->GetPdgCode()) == kMuonMinus) { pid = 1; }
367 else if (TMath::Abs(particle->GetPdgCode()) == kPiPlus) { pid = 2; }
368 else if (TMath::Abs(particle->GetPdgCode()) == kKPlus) { pid = 3; }
369 else if (TMath::Abs(particle->GetPdgCode()) == kProton) { pid = 4; }
375 //_____________________________________________________________________________
376 TH1F* AlidNdPtHelper::CreateResHisto(TH2F* const hRes2, TH1F **phMean, Int_t integ, Bool_t drawBinFits, Int_t minHistEntries)
379 // Create mean and resolution
382 TVirtualPad* currentPad = gPad;
383 TAxis* axis = hRes2->GetXaxis();
384 Int_t nBins = axis->GetNbins();
385 Bool_t overflowBinFits = kFALSE;
387 if (axis->GetXbins()->GetSize()){
388 hRes = new TH1F("hRes", "", nBins, axis->GetXbins()->GetArray());
389 hMean = new TH1F("hMean", "", nBins, axis->GetXbins()->GetArray());
392 hRes = new TH1F("hRes", "", nBins, axis->GetXmin(), axis->GetXmax());
393 hMean = new TH1F("hMean", "", nBins, axis->GetXmin(), axis->GetXmax());
396 hRes->SetStats(false);
397 hRes->SetOption("E");
398 hRes->SetMinimum(0.);
400 hMean->SetStats(false);
401 hMean->SetOption("E");
403 // create the fit function
404 TF1 * fitFunc = new TF1("G","[0]*exp(-(x-[1])*(x-[1])/(2.*[2]*[2]))",-3,3);
406 fitFunc->SetLineWidth(2);
407 fitFunc->SetFillStyle(0);
408 // create canvas for fits
409 TCanvas* canBinFits = NULL;
410 Int_t nPads = (overflowBinFits) ? nBins+2 : nBins;
411 Int_t nx = Int_t(sqrt(nPads-1.));// + 1;
412 Int_t ny = (nPads-1) / nx + 1;
414 canBinFits = (TCanvas*)gROOT->FindObject("canBinFits");
415 if (canBinFits) delete canBinFits;
416 canBinFits = new TCanvas("canBinFits", "fits of bins", 200, 100, 500, 700);
417 canBinFits->Divide(nx, ny);
420 // loop over x bins and fit projection
421 Int_t dBin = ((overflowBinFits) ? 1 : 0);
422 for (Int_t bin = 1-dBin; bin <= nBins+dBin; bin++) {
423 if (drawBinFits) canBinFits->cd(bin + dBin);
424 Int_t bin0=TMath::Max(bin-integ,0);
425 Int_t bin1=TMath::Min(bin+integ,nBins);
426 TH1D* hBin = hRes2->ProjectionY("hBin", bin0, bin1);
428 if (hBin->GetEntries() > minHistEntries) {
429 fitFunc->SetParameters(hBin->GetMaximum(),hBin->GetMean(),hBin->GetRMS());
430 hBin->Fit(fitFunc,"s");
431 Double_t sigma = TMath::Abs(fitFunc->GetParameter(2));
434 hRes->SetBinContent(bin, TMath::Abs(fitFunc->GetParameter(2)));
435 hMean->SetBinContent(bin, fitFunc->GetParameter(1));
438 hRes->SetBinContent(bin, 0.);
439 hMean->SetBinContent(bin,0);
441 hRes->SetBinError(bin, fitFunc->GetParError(2));
442 hMean->SetBinError(bin, fitFunc->GetParError(1));
448 hRes->SetBinContent(bin, 0.);
449 hRes->SetBinError(bin, 0.);
450 hMean->SetBinContent(bin, 0.);
451 hMean->SetBinError(bin, 0.);
458 sprintf(name, "%s < %.4g", axis->GetTitle(), axis->GetBinUpEdge(bin));
459 } else if (bin == nBins+1) {
460 sprintf(name, "%.4g < %s", axis->GetBinLowEdge(bin), axis->GetTitle());
462 sprintf(name, "%.4g < %s < %.4g", axis->GetBinLowEdge(bin),
463 axis->GetTitle(), axis->GetBinUpEdge(bin));
465 canBinFits->cd(bin + dBin);
466 hBin->SetTitle(name);
467 hBin->SetStats(kTRUE);
469 canBinFits->Update();
470 canBinFits->Modified();
471 canBinFits->Update();
483 //_____________________________________________________________________________
484 TH1F* AlidNdPtHelper::MakeResol(TH2F * his, Int_t integ, Bool_t type, Bool_t drawBins, Int_t minHistEntries){
485 // Create resolution histograms
487 TH1F *hisr=0, *hism=0;
488 if (!gPad) new TCanvas;
489 hisr = CreateResHisto(his,&hism,integ,drawBins,minHistEntries);
490 if (type) return hism;
496 //_____________________________________________________________________________
497 TObjArray* AlidNdPtHelper::GetAllChargedTracks(AliESDEvent *esdEvent, AnalysisMode analysisMode)
500 // all charged TPC particles
502 TObjArray *allTracks = new TObjArray();
503 if(!allTracks) return allTracks;
505 AliESDtrack *track=0;
506 for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++)
508 if(analysisMode == AlidNdPtHelper::kTPC) {
509 // track must be deleted by user
510 track = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent,iTrack);
513 else if (analysisMode == AlidNdPtHelper::kTPCSPDvtx || AlidNdPtHelper::kTPCSPDvtxUpdate)
515 // track must be deleted by the user
516 track = AlidNdPtHelper::GetTPCOnlyTrackSPDvtx(esdEvent,iTrack,kFALSE);
520 track=esdEvent->GetTrack(iTrack);
524 if(track->Charge()==0) {
525 if(analysisMode == AlidNdPtHelper::kTPC || analysisMode == AlidNdPtHelper::kTPCSPDvtx ||
526 analysisMode == AlidNdPtHelper::kTPCSPDvtxUpdate)
528 delete track; continue;
534 allTracks->Add(track);
537 if(analysisMode == AlidNdPtHelper::kTPC || analysisMode == AlidNdPtHelper::kTPCSPDvtx ||
538 analysisMode == AlidNdPtHelper::kTPCSPDvtxUpdate) {
540 allTracks->SetOwner(kTRUE);
546 //_____________________________________________________________________________
547 AliESDtrack *AlidNdPtHelper::GetTPCOnlyTrackSPDvtx(AliESDEvent* esdEvent, Int_t iTrack, Bool_t bUpdate)
550 // Create ESD tracks from TPCinner parameters.
551 // Propagte to DCA to SPD vertex.
552 // Update using SPD vertex point (parameter)
554 // It is user responsibility to delete these tracks
557 if (!esdEvent) return NULL;
558 if (!esdEvent->GetPrimaryVertexSPD() ) { return NULL; }
559 if (!esdEvent->GetPrimaryVertexSPD()->GetStatus() ) { return NULL; }
562 AliESDtrack* track = esdEvent->GetTrack(iTrack);
566 // create new ESD track
567 AliESDtrack *tpcTrack = new AliESDtrack();
569 // relate TPC-only tracks (TPCinner) to SPD vertex
570 AliExternalTrackParam cParam;
572 track->RelateToVertexTPC(esdEvent->GetPrimaryVertexSPD(),esdEvent->GetMagneticField(),kVeryBig,&cParam);
573 track->Set(cParam.GetX(),cParam.GetAlpha(),cParam.GetParameter(),cParam.GetCovariance());
575 // reject fake tracks
576 if(track->Pt() > 10000.) {
577 ::Error("Exclude no physical tracks","pt>10000. GeV");
583 track->RelateToVertexTPC(esdEvent->GetPrimaryVertexSPD(), esdEvent->GetMagneticField(), kVeryBig);
586 // only true if we have a tpc track
587 if (!track->FillTPCOnlyTrack(*tpcTrack))
596 //_____________________________________________________________________________
597 Int_t AlidNdPtHelper::GetTPCMBTrackMult(AliESDEvent *const esdEvent, AlidNdPtEventCuts *const evtCuts, AlidNdPtAcceptanceCuts *const accCuts, AliESDtrackCuts *const trackCuts)
600 // get MB event track multiplicity
604 ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","esd event is NULL");
608 if(!evtCuts || !accCuts || !trackCuts)
610 ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","cuts not available");
615 Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
616 Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
617 AliESDVertex vtx0(pos,err);
620 Float_t maxDCAr = accCuts->GetMaxDCAr();
621 Float_t maxDCAz = accCuts->GetMaxDCAz();
622 Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
624 Int_t ntracks = esdEvent->GetNumberOfTracks();
625 Double_t dca[2],cov[3];
627 for (Int_t i=0;i <ntracks; i++){
628 AliESDtrack *t = esdEvent->GetTrack(i);
630 if (t->Charge() == 0) continue;
631 if (!t->GetTPCInnerParam()) continue;
632 if (t->GetTPCNcls()<minTPCClust) continue;
634 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
635 if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov))
637 if(tpcTrack) delete tpcTrack;
641 if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
642 if(tpcTrack) delete tpcTrack;
648 if(tpcTrack) delete tpcTrack;
654 //_____________________________________________________________________________
655 Int_t AlidNdPtHelper::GetTPCMBPrimTrackMult(AliESDEvent *const esdEvent, AliStack *const stack, AlidNdPtEventCuts *const evtCuts, AlidNdPtAcceptanceCuts *const accCuts, AliESDtrackCuts *const trackCuts)
658 // get MB primary event track multiplicity
662 ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
668 ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
672 if(!evtCuts || !accCuts || !trackCuts)
674 ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","cuts not available");
679 Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
680 Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
681 AliESDVertex vtx0(pos,err);
684 Float_t maxDCAr = accCuts->GetMaxDCAr();
685 Float_t maxDCAz = accCuts->GetMaxDCAz();
686 Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
689 Int_t ntracks = esdEvent->GetNumberOfTracks();
690 Double_t dca[2],cov[3];
692 for (Int_t i=0;i <ntracks; i++){
693 AliESDtrack *t = esdEvent->GetTrack(i);
695 if (t->Charge() == 0) continue;
696 if (!t->GetTPCInnerParam()) continue;
697 if (t->GetTPCNcls()<minTPCClust) continue;
699 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
700 if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov))
702 if(tpcTrack) delete tpcTrack;
706 if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
707 if(tpcTrack) delete tpcTrack;
711 Int_t label = TMath::Abs(t->GetLabel());
712 TParticle *part = stack->Particle(label);
714 if(tpcTrack) delete tpcTrack;
717 if(!stack->IsPhysicalPrimary(label))
719 if(tpcTrack) delete tpcTrack;
725 if(tpcTrack) delete tpcTrack;
735 //_____________________________________________________________________________
736 Int_t AlidNdPtHelper::GetMCTrueTrackMult(AliMCEvent *const mcEvent, AlidNdPtEventCuts *const evtCuts, AlidNdPtAcceptanceCuts *const accCuts)
739 // calculate mc event true track multiplicity
741 if(!mcEvent) return 0;
747 stack = mcEvent->Stack();
748 if (!stack) return 0;
750 Bool_t isEventOK = evtCuts->AcceptMCEvent(mcEvent);
751 if(!isEventOK) return 0;
753 Int_t nPart = stack->GetNtrack();
754 for (Int_t iMc = 0; iMc < nPart; ++iMc)
756 TParticle* particle = stack->Particle(iMc);
760 // only charged particles
761 Double_t charge = particle->GetPDG()->Charge()/3.;
762 if (TMath::Abs(charge) < 0.001)
766 Bool_t prim = stack->IsPhysicalPrimary(iMc);
770 if(accCuts->AcceptTrack(particle))
779 //_______________________________________________________________________
780 void AlidNdPtHelper::PrintMCInfo(AliStack *const pStack,Int_t label)
782 // print information about particles in the stack
785 label = TMath::Abs(label);
786 TParticle *part = pStack->Particle(label);
787 Printf("########################");
788 Printf("%s:%d %d UniqueID %d PDG %d P %3.3f",(char*)__FILE__,__LINE__,label,part->GetUniqueID(),part->GetPdgCode(),part->P());
790 TParticle* mother = part;
791 Int_t imo = part->GetFirstMother();
792 Int_t nprim = pStack->GetNprimary();
794 while((imo >= nprim)) {
795 mother = pStack->Particle(imo);
796 Printf("Mother %s:%d Label %d UniqueID %d PDG %d P %3.3f",(char*)__FILE__,__LINE__,imo,mother->GetUniqueID(),mother->GetPdgCode(),mother->P());
798 imo = mother->GetFirstMother();
801 Printf("########################");
805 //_____________________________________________________________________________
806 TH1* AlidNdPtHelper::GetContCorrHisto(TH1 *const hist)
809 // get contamination histogram
813 Int_t nbins = hist->GetNbinsX();
814 TH1 *hCont = (TH1D *)hist->Clone();
816 for(Int_t i=0; i<=nbins+1; i++) {
817 Double_t binContent = hist->GetBinContent(i);
818 Double_t binError = hist->GetBinError(i);
820 hCont->SetBinContent(i,1.-binContent);
821 hCont->SetBinError(i,binError);
828 //_____________________________________________________________________________
829 TH1* AlidNdPtHelper::ScaleByBinWidth(TH1 *const hist)
832 // scale by bin width
836 TH1 *hScale = (TH1D *)hist->Clone();
837 hScale->Scale(1.,"width");
842 //_____________________________________________________________________________
843 TH1* AlidNdPtHelper::CalcRelativeDifference(TH1 *const hist1, TH1 *const hist2)
846 // calculate rel. difference
852 TH1 *h1Clone = (TH1D *)hist1->Clone();
856 h1Clone->Add(hist2,-1);
857 h1Clone->Divide(hist2);
862 //_____________________________________________________________________________
863 TH1* AlidNdPtHelper::CalcRelativeDifferenceFun(TH1 *const hist1, TF1 *const fun)
866 // calculate rel. difference
867 // between histogram and function
872 TH1 *h1Clone = (TH1D *)hist1->Clone();
876 h1Clone->Add(fun,-1);
877 h1Clone->Divide(hist1);
882 //_____________________________________________________________________________
883 TH1* AlidNdPtHelper::NormalizeToEvent(TH2 *const hist1, TH1 *const hist2)
885 // normalise to event for a given multiplicity bin
886 // return pt histogram
892 Int_t nbinsX = hist1->GetNbinsX();
893 //Int_t nbinsY = hist1->GetNbinsY();
896 for(Int_t i=0; i<=nbinsX+1; i++) {
897 sprintf(name,"mom_%d",i);
898 TH1D *hist = (TH1D*)hist1->ProjectionY(name,i+1,i+1);
900 sprintf(name,"mom_norm");
902 histNorm = (TH1D *)hist->Clone(name);
906 Double_t nbEvents = hist2->GetBinContent(i);
907 if(!nbEvents) { nbEvents = 1.; };
909 hist->Scale(1./nbEvents);
916 //_____________________________________________________________________________
917 THnSparse* AlidNdPtHelper::GenerateCorrMatrix(THnSparse *const hist1, THnSparse *const hist2, char *const name) {
918 // generate correction matrix
919 if(!hist1 || !hist2) return 0;
921 THnSparse *h =(THnSparse*)hist1->Clone(name);;
922 h->Divide(hist1,hist2,1,1,"B");
927 //_____________________________________________________________________________
928 TH2* AlidNdPtHelper::GenerateCorrMatrix(TH2 *const hist1, TH2 *const hist2, char *const name) {
929 // generate correction matrix
930 if(!hist1 || !hist2) return 0;
932 TH2D *h =(TH2D*)hist1->Clone(name);;
933 h->Divide(hist1,hist2,1,1,"B");
938 //_____________________________________________________________________________
939 TH1* AlidNdPtHelper::GenerateCorrMatrix(TH1 *const hist1, TH1 *const hist2, char *const name) {
940 // generate correction matrix
941 if(!hist1 || !hist2) return 0;
943 TH1D *h =(TH1D*)hist1->Clone(name);;
944 h->Divide(hist1,hist2,1,1,"B");
949 //_____________________________________________________________________________
950 THnSparse* AlidNdPtHelper::GenerateContCorrMatrix(THnSparse *const hist1, THnSparse *const hist2, char *const name) {
951 // generate contamination correction matrix
952 if(!hist1 || !hist2) return 0;
954 THnSparse *hist = GenerateCorrMatrix(hist1, hist2, name);
957 // only for non ZERO bins!!!!
959 Int_t* coord = new Int_t[hist->GetNdimensions()];
960 memset(coord, 0, sizeof(Int_t) * hist->GetNdimensions());
962 for (Long64_t i = 0; i < hist->GetNbins(); ++i) {
963 Double_t v = hist->GetBinContent(i, coord);
964 hist->SetBinContent(coord, 1.0-v);
965 //printf("v %f, hist->GetBinContent(i, coord) %f \n",v,hist->GetBinContent(i, coord));
966 Double_t err = hist->GetBinError(coord);
967 hist->SetBinError(coord, err);
975 //_____________________________________________________________________________
976 TH2* AlidNdPtHelper::GenerateContCorrMatrix(TH2 *const hist1, TH2 *const hist2, char *const name) {
977 // generate contamination correction matrix
978 if(!hist1 || !hist2) return 0;
980 TH2 *hist = GenerateCorrMatrix(hist1, hist2, name);
983 Int_t nBinsX = hist->GetNbinsX();
984 Int_t nBinsY = hist->GetNbinsY();
986 for (Int_t i = 0; i < nBinsX+1; i++) {
987 for (Int_t j = 0; j < nBinsY+1; j++) {
988 Double_t cont = hist->GetBinContent(i,j);
989 hist->SetBinContent(i,j,1.-cont);
990 Double_t err = hist->GetBinError(i,j);
991 hist->SetBinError(i,j,err);
998 //_____________________________________________________________________________
999 TH1* AlidNdPtHelper::GenerateContCorrMatrix(TH1 *const hist1, TH1 *const hist2, char *const name) {
1000 // generate contamination correction matrix
1001 if(!hist1 || !hist2) return 0;
1003 TH1 *hist = GenerateCorrMatrix(hist1, hist2, name);
1006 Int_t nBinsX = hist->GetNbinsX();
1008 for (Int_t i = 0; i < nBinsX+1; i++) {
1009 Double_t cont = hist->GetBinContent(i);
1010 hist->SetBinContent(i,1.-cont);
1011 Double_t err = hist->GetBinError(i);
1012 hist->SetBinError(i,err);
1018 //_____________________________________________________________________________
1019 const AliESDVertex* AlidNdPtHelper::GetTPCVertexZ(AliESDEvent* const esdEvent, AlidNdPtEventCuts *const evtCuts, AlidNdPtAcceptanceCuts *const accCuts, AliESDtrackCuts *const trackCuts, Float_t fraction, Int_t ntracksMin){
1025 ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
1029 if(!evtCuts || !accCuts || !trackCuts)
1031 ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
1035 Double_t vtxpos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
1036 Double_t vtxsigma[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
1037 AliESDVertex vtx0(vtxpos,vtxsigma);
1039 Double_t maxDCAr = accCuts->GetMaxDCAr();
1040 Double_t maxDCAz = accCuts->GetMaxDCAz();
1041 Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
1044 Int_t ntracks = esdEvent->GetNumberOfTracks();
1045 TVectorD ztrack(ntracks);
1046 Double_t dca[2],cov[3];
1048 for (Int_t i=0;i <ntracks; i++){
1049 AliESDtrack *t = esdEvent->GetTrack(i);
1051 if (!t->GetTPCInnerParam()) continue;
1052 if (t->GetTPCNcls()<minTPCClust) continue;
1054 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
1055 if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov)) continue;
1058 if (TMath::Abs(dca[0])>maxDCAr) continue;
1059 //if (TMath::Sqrt(cov[0])>sigmaXYcut) continue;
1060 if (TMath::Abs(tpcTrack->GetZ())>maxDCAz) continue;
1062 ztrack[counter]=tpcTrack->GetZ();
1065 if(tpcTrack) delete tpcTrack;
1069 // Find LTM z position
1071 Double_t mean=0, sigma=0;
1072 if (counter<ntracksMin) return 0;
1074 Int_t nused = TMath::Nint(counter*fraction);
1075 if (nused==counter) nused=counter-1;
1077 AliMathBase::EvaluateUni(counter, ztrack.GetMatrixArray(), mean,sigma, TMath::Nint(counter*fraction));
1078 sigma/=TMath::Sqrt(nused);
1080 mean = TMath::Mean(counter, ztrack.GetMatrixArray());
1081 sigma = TMath::RMS(counter, ztrack.GetMatrixArray());
1082 sigma/=TMath::Sqrt(counter-1);
1086 const AliESDVertex* vertex = new AliESDVertex(vtxpos, vtxsigma);
1090 //_____________________________________________________________________________
1091 Int_t AlidNdPtHelper::GetSPDMBTrackMult(AliESDEvent* const esdEvent, Float_t deltaThetaCut, Float_t deltaPhiCut)
1094 // SPD track multiplicity
1098 const AliMultiplicity* mult = esdEvent->GetMultiplicity();
1102 // get multiplicity from SPD tracklets
1103 Int_t inputCount = 0;
1104 for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
1106 //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
1108 Float_t phi = mult->GetPhi(i);
1110 phi += TMath::Pi() * 2;
1111 Float_t deltaPhi = mult->GetDeltaPhi(i);
1112 Float_t deltaTheta = mult->GetDeltaTheta(i);
1114 if (TMath::Abs(deltaPhi) > 1)
1115 printf("WARNING: Very high Delta Phi: %d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), deltaPhi);
1117 if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
1120 if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
1129 //_____________________________________________________________________________
1130 Int_t AlidNdPtHelper::GetSPDMBPrimTrackMult(AliESDEvent* const esdEvent, AliStack* const stack, Float_t deltaThetaCut, Float_t deltaPhiCut)
1133 // SPD track multiplicity
1137 const AliMultiplicity* mult = esdEvent->GetMultiplicity();
1141 // get multiplicity from SPD tracklets
1142 Int_t inputCount = 0;
1143 for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
1145 //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
1147 Float_t phi = mult->GetPhi(i);
1149 phi += TMath::Pi() * 2;
1150 Float_t deltaPhi = mult->GetDeltaPhi(i);
1151 Float_t deltaTheta = mult->GetDeltaTheta(i);
1153 if (TMath::Abs(deltaPhi) > 1)
1154 printf("WARNING: Very high Delta Phi: %d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), deltaPhi);
1156 if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
1159 if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
1163 if (mult->GetLabel(i, 0) < 0 || mult->GetLabel(i, 0) != mult->GetLabel(i, 1) ||
1164 !stack->IsPhysicalPrimary(mult->GetLabel(i, 0)))