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 <AliTracker.h>
44 #include "dNdPt/AlidNdPtEventCuts.h"
45 #include "dNdPt/AlidNdPtAcceptanceCuts.h"
46 #include "dNdPt/AlidNdPtHelper.h"
48 //____________________________________________________________________
49 ClassImp(AlidNdPtHelper)
51 //____________________________________________________________________
52 const AliESDVertex* AlidNdPtHelper::GetVertex(AliESDEvent* const aEsd, const AlidNdPtEventCuts *const evtCuts, const AlidNdPtAcceptanceCuts *const accCuts, const AliESDtrackCuts *const trackCuts, AnalysisMode analysisMode, Bool_t debug, Bool_t bRedoTPC, Bool_t bUseMeanVertex)
54 // Get the vertex from the ESD and returns it if the vertex is valid
56 // Second argument decides which vertex is used (this selects
57 // also the quality criteria that are applied)
61 ::Error("AlidNdPtHelper::GetVertex()","esd event is NULL");
65 if(!evtCuts || !accCuts || !trackCuts)
67 ::Error("AlidNdPtHelper::GetVertex()","cuts not available");
71 const AliESDVertex* vertex = 0;
72 AliESDVertex *initVertex = 0;
73 if (analysisMode == kSPD || analysisMode == kTPCITS ||
74 analysisMode == kTPCSPDvtx || analysisMode == kTPCSPDvtxUpdate || analysisMode == kTPCITSHybrid)
76 vertex = aEsd->GetPrimaryVertexSPD();
78 Printf("AlidNdPtHelper::GetVertex: Returning SPD vertex");
80 else if (analysisMode == kTPCTrackSPDvtx || analysisMode == kTPCTrackSPDvtxUpdate ||
81 analysisMode == kTPCITSHybridTrackSPDvtx || analysisMode == kTPCITSHybridTrackSPDvtxDCArPt ||
82 analysisMode == kITSStandAloneTrackSPDvtx || analysisMode == kITSStandAloneTPCTrackSPDvtx)
84 vertex = aEsd->GetPrimaryVertexTracks();
85 if(!vertex) return NULL;
86 if(vertex->GetNContributors()<1) {
88 vertex = aEsd->GetPrimaryVertexSPD();
91 else if (analysisMode == kTPC)
95 Double_t kBz = aEsd->GetMagneticField();
96 AliVertexerTracks vertexer(kBz);
99 Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
100 Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
101 initVertex = new AliESDVertex(pos,err);
102 vertexer.SetVtxStart(initVertex);
103 vertexer.SetConstraintOn();
106 Double_t maxDCAr = accCuts->GetMaxDCAr();
107 Double_t maxDCAz = accCuts->GetMaxDCAz();
108 Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
110 //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);
111 vertexer.SetTPCMode(0.1,1.0,5.0,minTPCClust,1,3.,0.1,2.0,maxDCAr,maxDCAz,1,4);
113 // TPC track preselection
114 Int_t ntracks = aEsd->GetNumberOfTracks();
115 TObjArray array(ntracks);
116 UShort_t *id = new UShort_t[ntracks];
119 for (Int_t i=0;i <ntracks; i++) {
120 AliESDtrack *t = aEsd->GetTrack(i);
122 if (t->Charge() == 0) continue;
123 if (!t->GetTPCInnerParam()) continue;
124 if (t->GetTPCNcls()<vertexer.GetMinClusters()) continue;
125 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
127 array.AddLast(tpcTrack);
128 id[count] = (UShort_t)t->GetID();
132 AliESDVertex *vTPC = vertexer.VertexForSelectedTracks(&array,id, kTRUE, kTRUE, bUseMeanVertex);
134 // set recreated TPC vertex
135 aEsd->SetPrimaryVertexTPC(vTPC);
137 for (Int_t i=0; i<aEsd->GetNumberOfTracks(); i++) {
138 AliESDtrack *t = aEsd->GetTrack(i);
141 Double_t x[3]; t->GetXYZ(x);
142 Double_t b[3]; AliTracker::GetBxByBz(x,b);
143 t->RelateToVertexTPCBxByBz(vTPC, b, kVeryBig);
148 delete [] id; id=NULL;
151 vertex = aEsd->GetPrimaryVertexTPC();
153 Printf("AlidNdPtHelper::GetVertex: Returning vertex from tracks");
156 Printf("AlidNdPtHelper::GetVertex: ERROR: Invalid second argument %d", analysisMode);
160 Printf("AlidNdPtHelper::GetVertex: No vertex found in ESD");
166 Printf("AlidNdPtHelper::GetVertex: Returning valid vertex: %s", vertex->GetTitle());
170 if(initVertex) delete initVertex; initVertex=NULL;
174 //____________________________________________________________________
175 Bool_t AlidNdPtHelper::TestRecVertex(const AliESDVertex* vertex, const AliESDVertex* vertexSPD, AnalysisMode analysisMode, Bool_t debug)
177 // Checks if a vertex meets the needed quality criteria
178 if(!vertex) return kFALSE;
179 if(!vertex->GetStatus()) return kFALSE;
181 Float_t requiredZResolution = -1;
182 if (analysisMode == kSPD || analysisMode == kTPCITS ||
183 analysisMode == kTPCSPDvtx || analysisMode == kTPCSPDvtxUpdate || analysisMode == kTPCITSHybrid ||
184 analysisMode == kTPCTrackSPDvtx || analysisMode == kTPCTrackSPDvtxUpdate || analysisMode == kTPCITSHybridTrackSPDvtx || analysisMode == kTPCITSHybridTrackSPDvtxDCArPt
185 || analysisMode == kITSStandAloneTrackSPDvtx || analysisMode == kITSStandAloneTPCTrackSPDvtx)
187 requiredZResolution = 1000;
189 else if (analysisMode == kTPC)
190 requiredZResolution = 10.;
193 Double_t zRes = vertex->GetZRes();
195 if (zRes > requiredZResolution) {
197 Printf("AlidNdPtHelper::TestVertex: Resolution too poor %f (required: %f", zRes, requiredZResolution);
201 // always check for SPD vertex
202 if(!vertexSPD) return kFALSE;
203 if(!vertexSPD->GetStatus()) return kFALSE;
204 if (vertexSPD->IsFromVertexerZ())
206 if (vertexSPD->GetDispersion() > 0.02)
209 Printf("AliPWG0Helper::TestVertex: Delta Phi too large in Vertexer Z: %f (required: %f", vertex->GetDispersion(), 0.02);
215 // check Ncontributors
216 if (vertex->GetNContributors() <= 0) {
218 Printf("AlidNdPtHelper::GetVertex: NContributors() <= 0: %d",vertex->GetNContributors());
219 Printf("AlidNdPtHelper::GetVertex: NIndices(): %d",vertex->GetNIndices());
229 //____________________________________________________________________
230 Bool_t AlidNdPtHelper::IsGoodImpPar(const AliESDtrack *const track)
233 // check whether particle has good DCAr(Pt) impact
234 // parameter. Only for TPC+ITS tracks (7*sigma cut)
235 // Origin: Andrea Dainese
238 Float_t d0z0[2],covd0z0[3];
239 track->GetImpactParameters(d0z0,covd0z0);
240 Float_t sigma= 0.0050+0.0060/TMath::Power(track->Pt(),0.9);
241 Float_t d0max = 7.*sigma;
242 if(TMath::Abs(d0z0[0]) < d0max) return kTRUE;
247 //____________________________________________________________________
248 Bool_t AlidNdPtHelper::IsPrimaryParticle(AliStack* const stack, Int_t idx, ParticleMode particleMode)
251 // check primary particles
252 // depending on the particle mode
254 if(!stack) return kFALSE;
256 TParticle* particle = stack->Particle(idx);
257 if (!particle) return kFALSE;
259 // only charged particles
260 Double_t charge = particle->GetPDG()->Charge()/3.;
261 if (TMath::Abs(charge) < 0.001) return kFALSE;
263 Int_t pdg = TMath::Abs(particle->GetPdgCode());
266 Bool_t prim = stack->IsPhysicalPrimary(idx);
268 if(particleMode==kMCPion) {
269 if(prim && pdg==kPiPlus) return kTRUE;
273 if (particleMode==kMCKaon) {
274 if(prim && pdg==kKPlus) return kTRUE;
278 if (particleMode==kMCProton) {
279 if(prim && pdg==kProton) return kTRUE;
283 if(particleMode==kMCRest) {
284 if(prim && pdg!=kPiPlus && pdg!=kKPlus && pdg!=kProton) return kTRUE;
291 //____________________________________________________________________
292 Bool_t AlidNdPtHelper::IsCosmicTrack(AliESDtrack *const track1, AliESDtrack *const track2)
295 // check cosmic tracks
297 if(!track1) return kFALSE;
298 if(!track2) return kFALSE;
301 // cosmic tracks in TPC
303 //if( TMath::Abs( track1->Theta() - track2->Theta() ) < 0.004 &&
304 // ((TMath::Abs(track1->Phi()-track2->Phi()) - TMath::Pi() )<0.004) &&
305 //if( track1->Pt() > 4.0 && (TMath::Abs(track1->Phi()-track2->Phi())-TMath::Pi())<0.1
306 // && (TMath::Abs(track1->Eta()+track2->Eta())-0.1) < 0.0 && (track1->Charge()+track2->Charge()) == 0)
308 //Float_t scaleF= 6.0;
309 if ( track1->Pt() > 4 && track2->Pt() > 4 &&
310 //(TMath::Abs(track1->GetSnp()-track2->GetSnp())-2.) < scaleF * TMath::Sqrt(track1->GetSigmaSnp2()+track2->GetSigmaSnp2()) &&
311 //TMath::Abs(track1->GetTgl()-track2->GetTgl()) < scaleF * TMath::Sqrt(track1->GetSigmaTgl2()+track2->GetSigmaTgl2()) &&
312 //TMath::Abs(track1->OneOverPt()-track2->OneOverPt()) < scaleF * TMath::Sqrt(track1->GetSigma1Pt2()+track2->GetSigma1Pt2()) &&
313 (track1->Charge()+track2->Charge()) == 0 &&
314 track1->Eta()*track2->Eta()<0.0 && TMath::Abs(track1->Eta()+track2->Eta())<0.03 &&
315 TMath::Abs(TMath::Abs(track1->Phi()-track2->Phi())-TMath::Pi())<0.1
318 printf("COSMIC candidate \n");
320 printf("track1->Pt() %f, track1->Theta() %f, track1->Eta() %f, track1->Phi() %f, track1->Charge() %d \n", track1->Pt(), track1->Theta(), track1->Eta(), track1->Phi(), track1->Charge());
321 printf("track2->Pt() %f, track2->Theta() %f, track2->Eta() %f, track2->Phi() %f, track2->Charge() %d \n", track2->Pt(), track2->Theta(), track2->Eta(), track2->Phi(), track2->Charge());
322 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());
323 printf("dsphi %f, errsphi %f, dtanl %f, errtanl %f \n", TMath::Abs(track1->GetSnp()-track2->GetSnp()), TMath::Sqrt(track1->GetSigmaSnp2()+track2->GetSigmaSnp2()), TMath::Abs(track1->GetTgl()-track2->GetTgl()), TMath::Sqrt(track1->GetSigmaTgl2()+track2->GetSigmaTgl2()));
330 //____________________________________________________________________
331 void AlidNdPtHelper::PrintConf(AnalysisMode analysisMode, AliTriggerAnalysis::Trigger trigger)
334 // Prints the given configuration
337 TString str(">>>> Running with ");
339 switch (analysisMode)
341 case kInvalid: str += "invalid setting"; break;
342 case kSPD : str += "SPD-only"; break;
343 case kTPC : str += "TPC-only"; break;
344 case kTPCITS : str += "Global tracking"; break;
345 case kTPCSPDvtx : str += "TPC tracking + SPD event vertex"; break;
346 case kTPCSPDvtxUpdate : str += "TPC tracks updated with SPD event vertex point"; break;
347 case kTPCTrackSPDvtx : str += "TPC tracking + Tracks event vertex or SPD event vertex"; break;
348 case kTPCTrackSPDvtxUpdate : str += "TPC tracks updated with Track or SPD event vertex point"; break;
349 case kTPCITSHybrid : str += "TPC tracking + ITS refit + >1 SPD cluster"; break;
350 case kTPCITSHybridTrackSPDvtx : str += "TPC tracking + ITS refit + >1 SPD cluster + Tracks event vertex or SPD event vertex"; break;
351 case kTPCITSHybridTrackSPDvtxDCArPt : str += "TPC tracking + ITS refit + >1 SPD cluster + Tracks event vertex or SPD event vertex + DCAr(pt)"; break;
352 case kITSStandAloneTrackSPDvtx : str += "ITS standalone + Tracks event vertex or SPD event vertex + DCAr(pt)"; break;
353 case kITSStandAloneTPCTrackSPDvtx : str += "kITSStandAloneTPCTrackSPDvtx + TPC stand alone track + Tracks event vertex or SPD event vertex + DCAr(pt)"; break;
354 case kMCRec : str += "TPC tracking + Replace rec. with MC values"; break;
356 str += " and trigger ";
358 str += AliTriggerAnalysis::GetTriggerName(trigger);
362 Printf("%s", str.Data());
365 //____________________________________________________________________
366 Int_t AlidNdPtHelper::ConvertPdgToPid(const TParticle *const particle) {
368 // Convert Abs(pdg) to pid
369 // (0 - e, 1 - muons, 2 - pions, 3 - kaons, 4 - protons, 5 -all rest)
373 if (TMath::Abs(particle->GetPdgCode()) == kElectron) { pid = 0; }
374 else if (TMath::Abs(particle->GetPdgCode()) == kMuonMinus) { pid = 1; }
375 else if (TMath::Abs(particle->GetPdgCode()) == kPiPlus) { pid = 2; }
376 else if (TMath::Abs(particle->GetPdgCode()) == kKPlus) { pid = 3; }
377 else if (TMath::Abs(particle->GetPdgCode()) == kProton) { pid = 4; }
383 //_____________________________________________________________________________
384 TH1F* AlidNdPtHelper::CreateResHisto(TH2F* const hRes2, TH1F **phMean, Int_t integ, Bool_t drawBinFits, Int_t minHistEntries)
387 // Create mean and resolution
390 TVirtualPad* currentPad = gPad;
391 TAxis* axis = hRes2->GetXaxis();
392 Int_t nBins = axis->GetNbins();
393 Bool_t overflowBinFits = kFALSE;
395 if (axis->GetXbins()->GetSize()){
396 hRes = new TH1F("hRes", "", nBins, axis->GetXbins()->GetArray());
397 hMean = new TH1F("hMean", "", nBins, axis->GetXbins()->GetArray());
400 hRes = new TH1F("hRes", "", nBins, axis->GetXmin(), axis->GetXmax());
401 hMean = new TH1F("hMean", "", nBins, axis->GetXmin(), axis->GetXmax());
404 hRes->SetStats(false);
405 hRes->SetOption("E");
406 hRes->SetMinimum(0.);
408 hMean->SetStats(false);
409 hMean->SetOption("E");
411 // create the fit function
412 TF1 * fitFunc = new TF1("G","[0]*exp(-(x-[1])*(x-[1])/(2.*[2]*[2]))",-3,3);
414 fitFunc->SetLineWidth(2);
415 fitFunc->SetFillStyle(0);
416 // create canvas for fits
417 TCanvas* canBinFits = NULL;
418 Int_t nPads = (overflowBinFits) ? nBins+2 : nBins;
419 Int_t nx = Int_t(sqrt(nPads-1.));// + 1;
420 Int_t ny = (nPads-1) / nx + 1;
422 canBinFits = (TCanvas*)gROOT->FindObject("canBinFits");
423 if (canBinFits) delete canBinFits;
424 canBinFits = new TCanvas("canBinFits", "fits of bins", 200, 100, 500, 700);
425 canBinFits->Divide(nx, ny);
428 // loop over x bins and fit projection
429 Int_t dBin = ((overflowBinFits) ? 1 : 0);
430 for (Int_t bin = 1-dBin; bin <= nBins+dBin; bin++) {
431 if (drawBinFits) canBinFits->cd(bin + dBin);
432 Int_t bin0=TMath::Max(bin-integ,0);
433 Int_t bin1=TMath::Min(bin+integ,nBins);
434 TH1D* hBin = hRes2->ProjectionY("hBin", bin0, bin1);
436 if (hBin->GetEntries() > minHistEntries) {
437 fitFunc->SetParameters(hBin->GetMaximum(),hBin->GetMean(),hBin->GetRMS());
438 hBin->Fit(fitFunc,"s");
439 Double_t sigma = TMath::Abs(fitFunc->GetParameter(2));
442 hRes->SetBinContent(bin, TMath::Abs(fitFunc->GetParameter(2)));
443 hMean->SetBinContent(bin, fitFunc->GetParameter(1));
446 hRes->SetBinContent(bin, 0.);
447 hMean->SetBinContent(bin,0);
449 hRes->SetBinError(bin, fitFunc->GetParError(2));
450 hMean->SetBinError(bin, fitFunc->GetParError(1));
456 hRes->SetBinContent(bin, 0.);
457 hRes->SetBinError(bin, 0.);
458 hMean->SetBinContent(bin, 0.);
459 hMean->SetBinError(bin, 0.);
466 sprintf(name, "%s < %.4g", axis->GetTitle(), axis->GetBinUpEdge(bin));
467 } else if (bin == nBins+1) {
468 sprintf(name, "%.4g < %s", axis->GetBinLowEdge(bin), axis->GetTitle());
470 sprintf(name, "%.4g < %s < %.4g", axis->GetBinLowEdge(bin),
471 axis->GetTitle(), axis->GetBinUpEdge(bin));
473 canBinFits->cd(bin + dBin);
474 hBin->SetTitle(name);
475 hBin->SetStats(kTRUE);
477 canBinFits->Update();
478 canBinFits->Modified();
479 canBinFits->Update();
491 //_____________________________________________________________________________
492 TH1F* AlidNdPtHelper::MakeResol(TH2F * his, Int_t integ, Bool_t type, Bool_t drawBins, Int_t minHistEntries){
493 // Create resolution histograms
495 TH1F *hisr=0, *hism=0;
496 if (!gPad) new TCanvas;
497 hisr = CreateResHisto(his,&hism,integ,drawBins,minHistEntries);
498 if (type) return hism;
504 //_____________________________________________________________________________
505 TObjArray* AlidNdPtHelper::GetAllChargedTracks(AliESDEvent *esdEvent, AnalysisMode analysisMode)
508 // all charged TPC particles
510 TObjArray *allTracks = new TObjArray();
511 if(!allTracks) return allTracks;
513 AliESDtrack *track=0;
514 for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++)
516 if(analysisMode == AlidNdPtHelper::kTPC) {
518 // track must be deleted by user
519 // esd track parameters are replaced by TPCinner
521 track = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent,iTrack);
524 else if (analysisMode == AlidNdPtHelper::kTPCSPDvtx || analysisMode == AlidNdPtHelper::kTPCSPDvtxUpdate)
527 // track must be deleted by the user
528 // esd track parameters are replaced by TPCinner
530 track = AlidNdPtHelper::GetTPCOnlyTrackSPDvtx(esdEvent,iTrack,kFALSE);
533 else if (analysisMode == AlidNdPtHelper::kTPCTrackSPDvtx || analysisMode == AlidNdPtHelper::kTPCTrackSPDvtxUpdate)
536 // track must be deleted by the user
537 // esd track parameters are replaced by TPCinner
539 track = AlidNdPtHelper::GetTPCOnlyTrackTrackSPDvtx(esdEvent,iTrack,kFALSE);
542 else if(analysisMode == AlidNdPtHelper::kTPCITSHybrid )
544 track = AlidNdPtHelper::GetTrackSPDvtx(esdEvent,iTrack,kFALSE);
546 else if(analysisMode == AlidNdPtHelper::kTPCITSHybridTrackSPDvtx || analysisMode == AlidNdPtHelper::kTPCITSHybridTrackSPDvtxDCArPt || analysisMode == AlidNdPtHelper::kITSStandAloneTrackSPDvtx || analysisMode ==AlidNdPtHelper::kITSStandAloneTPCTrackSPDvtx)
548 track = AlidNdPtHelper::GetTrackTrackSPDvtx(esdEvent,iTrack,kFALSE);
552 track = esdEvent->GetTrack(iTrack);
557 if(track->Charge()==0) {
558 if(analysisMode == AlidNdPtHelper::kTPC ||
559 analysisMode == AlidNdPtHelper::kTPCSPDvtx || analysisMode == AlidNdPtHelper::kTPCTrackSPDvtx ||
560 analysisMode == AlidNdPtHelper::kTPCSPDvtxUpdate || analysisMode == AlidNdPtHelper::kTPCTrackSPDvtxUpdate)
562 delete track; continue;
568 allTracks->Add(track);
571 if(analysisMode == AlidNdPtHelper::kTPC ||
572 analysisMode == AlidNdPtHelper::kTPCSPDvtx || analysisMode == AlidNdPtHelper::kTPCTrackSPDvtx ||
573 analysisMode == AlidNdPtHelper::kTPCSPDvtxUpdate || analysisMode == AlidNdPtHelper::kTPCTrackSPDvtxUpdate) {
575 allTracks->SetOwner(kTRUE);
581 //_____________________________________________________________________________
582 AliESDtrack *AlidNdPtHelper::GetTPCOnlyTrackSPDvtx(const AliESDEvent* esdEvent, Int_t iTrack, Bool_t bUpdate)
585 // Create ESD tracks from TPCinner parameters.
586 // Propagte to DCA to SPD vertex.
587 // Update using SPD vertex point (parameter)
589 // It is user responsibility to delete these tracks
592 if (!esdEvent) return NULL;
593 if (!esdEvent->GetPrimaryVertexSPD() ) { return NULL; }
594 if (!esdEvent->GetPrimaryVertexSPD()->GetStatus() ) { return NULL; }
597 AliESDtrack* track = esdEvent->GetTrack(iTrack);
601 Bool_t isOK = kFALSE;
602 Double_t x[3]; track->GetXYZ(x);
603 Double_t b[3]; AliTracker::GetBxByBz(x,b);
605 // create new ESD track
606 AliESDtrack *tpcTrack = new AliESDtrack();
608 // relate TPC-only tracks (TPCinner) to SPD vertex
609 AliExternalTrackParam cParam;
611 isOK = track->RelateToVertexTPCBxByBz(esdEvent->GetPrimaryVertexSPD(),b,kVeryBig,&cParam);
612 track->Set(cParam.GetX(),cParam.GetAlpha(),cParam.GetParameter(),cParam.GetCovariance());
614 // reject fake tracks
615 if(track->Pt() > 10000.) {
616 ::Error("Exclude no physical tracks","pt>10000. GeV");
622 isOK = track->RelateToVertexTPCBxByBz(esdEvent->GetPrimaryVertexSPD(), b, kVeryBig);
625 // only true if we have a tpc track
626 if (!track->FillTPCOnlyTrack(*tpcTrack))
632 if(!isOK) return NULL;
637 //_____________________________________________________________________________
638 AliESDtrack *AlidNdPtHelper::GetTPCOnlyTrackTrackSPDvtx(const AliESDEvent* esdEvent, Int_t iTrack, Bool_t bUpdate)
641 // Create ESD tracks from TPCinner parameters.
642 // Propagte to DCA to Track or SPD vertex.
643 // Update using SPD vertex point (parameter)
645 // It is user responsibility to delete these tracks
647 if (!esdEvent) return NULL;
648 const AliESDVertex *vertex = esdEvent->GetPrimaryVertexTracks();
649 if(vertex->GetNContributors()<1) {
651 vertex = esdEvent->GetPrimaryVertexSPD();
653 if(!vertex) return NULL;
656 AliESDtrack* track = esdEvent->GetTrack(iTrack);
660 Bool_t isOK = kFALSE;
661 Double_t x[3]; track->GetXYZ(x);
662 Double_t b[3]; AliTracker::GetBxByBz(x,b);
664 // create new ESD track
665 AliESDtrack *tpcTrack = new AliESDtrack();
667 // relate TPC-only tracks (TPCinner) to SPD vertex
668 AliExternalTrackParam cParam;
670 isOK = track->RelateToVertexTPCBxByBz(vertex,b,kVeryBig,&cParam);
671 track->Set(cParam.GetX(),cParam.GetAlpha(),cParam.GetParameter(),cParam.GetCovariance());
673 // reject fake tracks
674 if(track->Pt() > 10000.) {
675 ::Error("Exclude no physical tracks","pt>10000. GeV");
681 isOK = track->RelateToVertexTPCBxByBz(vertex, b, kVeryBig);
684 // only true if we have a tpc track
685 if (!track->FillTPCOnlyTrack(*tpcTrack))
691 if(!isOK) return NULL;
696 //_____________________________________________________________________________
697 AliESDtrack *AlidNdPtHelper::GetTrackSPDvtx(const AliESDEvent* esdEvent, Int_t iTrack, Bool_t bUpdate)
700 // Propagte track to DCA to SPD vertex.
701 // Update using SPD vertex point (parameter)
703 if (!esdEvent) return NULL;
704 if (!esdEvent->GetPrimaryVertexSPD() ) { return NULL; }
705 if (!esdEvent->GetPrimaryVertexSPD()->GetStatus() ) { return NULL; }
708 AliESDtrack* track = esdEvent->GetTrack(iTrack);
712 Bool_t isOK = kFALSE;
713 Double_t x[3]; track->GetXYZ(x);
714 Double_t b[3]; AliTracker::GetBxByBz(x,b);
716 // relate tracks to SPD vertex
717 AliExternalTrackParam cParam;
719 isOK = track->RelateToVertexBxByBz(esdEvent->GetPrimaryVertexSPD(),b,kVeryBig,&cParam);
720 track->Set(cParam.GetX(),cParam.GetAlpha(),cParam.GetParameter(),cParam.GetCovariance());
722 // reject fake tracks
723 if(track->Pt() > 10000.) {
724 ::Error("Exclude no physical tracks","pt>10000. GeV");
729 isOK = track->RelateToVertexBxByBz(esdEvent->GetPrimaryVertexSPD(), b, kVeryBig);
732 if(!isOK) return NULL;
737 //_____________________________________________________________________________
738 AliESDtrack *AlidNdPtHelper::GetTrackTrackSPDvtx(const AliESDEvent* esdEvent, Int_t iTrack, Bool_t bUpdate)
741 // Propagte track to DCA to Track or SPD vertex.
742 // Update using SPD vertex point (parameter)
744 if (!esdEvent) return NULL;
746 const AliESDVertex *vertex = esdEvent->GetPrimaryVertexTracks();
747 if(vertex->GetNContributors()<1) {
749 vertex = esdEvent->GetPrimaryVertexSPD();
751 if(!vertex) return NULL;
754 AliESDtrack* track = esdEvent->GetTrack(iTrack);
758 Bool_t isOK = kFALSE;
759 Double_t x[3]; track->GetXYZ(x);
760 Double_t b[3]; AliTracker::GetBxByBz(x,b);
762 // relate tracks to SPD vertex
763 AliExternalTrackParam cParam;
765 isOK = track->RelateToVertexBxByBz(vertex,b,kVeryBig,&cParam);
766 track->Set(cParam.GetX(),cParam.GetAlpha(),cParam.GetParameter(),cParam.GetCovariance());
768 // reject fake tracks
769 if(track->Pt() > 10000.) {
770 ::Error("Exclude no physical tracks","pt>10000. GeV");
775 isOK = track->RelateToVertexBxByBz(vertex, b, kVeryBig);
778 if(!isOK) return NULL;
783 //_____________________________________________________________________________
784 Bool_t AlidNdPtHelper::SelectEvent(const AliESDEvent* const esdEvent, AliESDtrackCuts* const esdTrackCuts) {
785 // select events with at least
786 // one reconstructed primary track in acceptance
787 // pT>0.5 GeV/c, |eta|<0.8 for cross section studies
789 if(!esdEvent) return kFALSE;
790 if(!esdTrackCuts) return kFALSE;
792 AliESDtrack *track=0;
794 for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++)
796 track = esdEvent->GetTrack(iTrack);
798 if(track->Charge()==0) continue;
799 if(!esdTrackCuts->AcceptTrack(track)) continue;
800 if(track->Pt() < 0.5) continue;
801 if(TMath::Abs(track->Eta()) > 0.8) continue;
806 if(count > 0) return kTRUE;
812 //_____________________________________________________________________________
813 Bool_t AlidNdPtHelper::SelectMCEvent(AliMCEvent* const mcEvent) {
815 // select events with at least
816 // one prompt (MC primary) track in acceptance
817 // pT>0.5 GeV/c, |eta|<0.8 for cross section studies
820 if(!mcEvent) return kFALSE;
821 AliStack* stack = mcEvent->Stack();
822 if(!stack) return kFALSE;
825 for (Int_t iMc = 0; iMc < stack->GetNtrack(); ++iMc)
827 TParticle* particle = stack->Particle(iMc);
828 if (!particle) continue;
830 // only charged particles
831 if(!particle->GetPDG()) continue;
832 Double_t charge = particle->GetPDG()->Charge()/3.;
833 if(charge == 0) continue;
836 Bool_t prim = stack->IsPhysicalPrimary(iMc);
839 if(particle->Pt() < 0.5) continue;
840 if(TMath::Abs(particle->Eta()) > 0.8) continue;
845 if(count > 0) return kTRUE;
851 //_____________________________________________________________________________
852 Int_t AlidNdPtHelper::GetTPCMBTrackMult(const AliESDEvent *const esdEvent,const AlidNdPtEventCuts *const evtCuts, const AlidNdPtAcceptanceCuts *const accCuts,const AliESDtrackCuts *const trackCuts)
855 // get MB event track multiplicity
859 ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","esd event is NULL");
863 if(!evtCuts || !accCuts || !trackCuts)
865 ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","cuts not available");
870 Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
871 Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
872 AliESDVertex vtx0(pos,err);
875 Float_t maxDCAr = accCuts->GetMaxDCAr();
876 Float_t maxDCAz = accCuts->GetMaxDCAz();
877 Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
879 Int_t ntracks = esdEvent->GetNumberOfTracks();
880 Double_t dca[2],cov[3];
882 for (Int_t i=0;i <ntracks; i++){
883 AliESDtrack *t = esdEvent->GetTrack(i);
885 if (t->Charge() == 0) continue;
886 if (!t->GetTPCInnerParam()) continue;
887 if (t->GetTPCNcls()<minTPCClust) continue;
889 Double_t x[3]; t->GetXYZ(x);
890 Double_t b[3]; AliTracker::GetBxByBz(x,b);
891 const Double_t kMaxStep = 1; //max step over the material
893 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
894 if (!tpcTrack->PropagateToDCABxByBz(&vtx0,b,kMaxStep,dca,cov))
896 if(tpcTrack) delete tpcTrack;
900 if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
901 if(tpcTrack) delete tpcTrack;
907 if(tpcTrack) delete tpcTrack;
913 //_____________________________________________________________________________
914 Int_t AlidNdPtHelper::GetTPCMBPrimTrackMult(const AliESDEvent *const esdEvent, AliStack *const stack, const AlidNdPtEventCuts *const evtCuts, const AlidNdPtAcceptanceCuts *const accCuts, const AliESDtrackCuts *const trackCuts)
917 // get MB primary event track multiplicity
921 ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
927 ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
931 if(!evtCuts || !accCuts || !trackCuts)
933 ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","cuts not available");
938 Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
939 Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
940 AliESDVertex vtx0(pos,err);
943 Float_t maxDCAr = accCuts->GetMaxDCAr();
944 Float_t maxDCAz = accCuts->GetMaxDCAz();
945 Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
948 Int_t ntracks = esdEvent->GetNumberOfTracks();
949 Double_t dca[2],cov[3];
951 for (Int_t i=0;i <ntracks; i++){
952 AliESDtrack *t = esdEvent->GetTrack(i);
954 if (t->Charge() == 0) continue;
955 if (!t->GetTPCInnerParam()) continue;
956 if (t->GetTPCNcls()<minTPCClust) continue;
958 Double_t x[3]; t->GetXYZ(x);
959 Double_t b[3]; AliTracker::GetBxByBz(x,b);
960 const Double_t kMaxStep = 1; //max step over the material
962 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
963 if (!tpcTrack->PropagateToDCABxByBz(&vtx0,b,kMaxStep,dca,cov))
965 if(tpcTrack) delete tpcTrack;
969 if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
970 if(tpcTrack) delete tpcTrack;
974 Int_t label = TMath::Abs(t->GetLabel());
975 TParticle *part = stack->Particle(label);
977 if(tpcTrack) delete tpcTrack;
980 if(!stack->IsPhysicalPrimary(label))
982 if(tpcTrack) delete tpcTrack;
988 if(tpcTrack) delete tpcTrack;
995 //_____________________________________________________________________________
996 Double_t AlidNdPtHelper::GetStrangenessCorrFactor(const Double_t pt)
998 // data driven correction factor for secondaries
999 // underestimated secondaries with strangeness in Pythia (A. Dainese)
1003 // pt=0.4; fact=1.07
1004 // pt=0.6; fact=1.25
1005 // pt>=1.2; fact=1.5
1008 if (pt <= 0.17) return 1.0;
1009 if (pt <= 0.4) return GetLinearInterpolationValue(0.17,1.0,0.4,1.07, pt);
1010 if (pt <= 0.6) return GetLinearInterpolationValue(0.4,1.07,0.6,1.25, pt);
1011 if (pt <= 1.2) return GetLinearInterpolationValue(0.6,1.25,1.2,1.5, pt);
1016 //_____________________________________________________________________________
1017 Double_t AlidNdPtHelper::GetStrangenessCorrFactorPbPb(const Double_t pt)
1019 // data driven correction factor for secondaries (PbPb)
1021 if (pt <= 0.25) return 1.0;
1022 if (pt <= 0.5) return GetLinearInterpolationValue(0.25,1.0,0.5,1.4, pt);
1023 if (pt <= 1.0) return GetLinearInterpolationValue(0.5,1.4,1.0,1.47, pt);
1024 if (pt <= 2.0) return GetLinearInterpolationValue(1.0,1.47,2.0,1.56, pt);
1025 if (pt <= 5.0) return GetLinearInterpolationValue(2.0,1.56,5.0,1.67, pt);
1030 //___________________________________________________________________________
1031 Double_t AlidNdPtHelper::GetLinearInterpolationValue(const Double_t x1,const Double_t y1,const Double_t x2,const Double_t y2, const Double_t pt)
1034 // linear interpolation
1036 return ((y2-y1)/(x2-x1))*pt+(y2-(((y2-y1)/(x2-x1))*x2));
1039 //_____________________________________________________________________________
1040 Int_t AlidNdPtHelper::GetMCTrueTrackMult(AliMCEvent *const mcEvent, AlidNdPtEventCuts *const evtCuts, AlidNdPtAcceptanceCuts *const accCuts)
1043 // calculate mc event true track multiplicity
1045 if(!mcEvent) return 0;
1047 AliStack* stack = 0;
1050 // MC particle stack
1051 stack = mcEvent->Stack();
1052 if (!stack) return 0;
1055 //printf("minZv %f, maxZv %f \n", evtCuts->GetMinZv(), evtCuts->GetMaxZv());
1058 Bool_t isEventOK = evtCuts->AcceptMCEvent(mcEvent);
1059 if(!isEventOK) return 0;
1061 Int_t nPart = stack->GetNtrack();
1062 for (Int_t iMc = 0; iMc < nPart; ++iMc)
1064 TParticle* particle = stack->Particle(iMc);
1068 // only charged particles
1069 if(!particle->GetPDG()) continue;
1070 Double_t charge = particle->GetPDG()->Charge()/3.;
1071 if (TMath::Abs(charge) < 0.001)
1075 Bool_t prim = stack->IsPhysicalPrimary(iMc);
1078 // checked accepted without pt cut
1079 //if(accCuts->AcceptTrack(particle))
1080 if( particle->Eta() > accCuts->GetMinEta() && particle->Eta() < accCuts->GetMaxEta() )
1089 //_______________________________________________________________________
1090 void AlidNdPtHelper::PrintMCInfo(AliStack *const pStack,Int_t label)
1092 // print information about particles in the stack
1095 label = TMath::Abs(label);
1096 TParticle *part = pStack->Particle(label);
1097 Printf("########################");
1098 Printf("%s:%d %d UniqueID %d PDG %d P %3.3f",(char*)__FILE__,__LINE__,label,part->GetUniqueID(),part->GetPdgCode(),part->P());
1100 TParticle* mother = part;
1101 Int_t imo = part->GetFirstMother();
1102 Int_t nprim = pStack->GetNprimary();
1104 while((imo >= nprim)) {
1105 mother = pStack->Particle(imo);
1106 Printf("Mother %s:%d Label %d UniqueID %d PDG %d P %3.3f",(char*)__FILE__,__LINE__,imo,mother->GetUniqueID(),mother->GetPdgCode(),mother->P());
1108 imo = mother->GetFirstMother();
1111 Printf("########################");
1115 //_____________________________________________________________________________
1116 TH1* AlidNdPtHelper::GetContCorrHisto(TH1 *const hist)
1119 // get contamination histogram
1123 Int_t nbins = hist->GetNbinsX();
1124 TH1 *hCont = (TH1D *)hist->Clone();
1126 for(Int_t i=0; i<=nbins+1; i++) {
1127 Double_t binContent = hist->GetBinContent(i);
1128 Double_t binError = hist->GetBinError(i);
1130 hCont->SetBinContent(i,1.-binContent);
1131 hCont->SetBinError(i,binError);
1138 //_____________________________________________________________________________
1139 TH1* AlidNdPtHelper::ScaleByBinWidth(TH1 *const hist)
1142 // scale by bin width
1146 TH1 *hScale = (TH1D *)hist->Clone();
1147 hScale->Scale(1.,"width");
1152 //_____________________________________________________________________________
1153 TH1* AlidNdPtHelper::CalcRelativeDifference(const TH1 *const hist1, const TH1 *const hist2)
1156 // calculate rel. difference
1159 if(!hist1) return 0;
1160 if(!hist2) return 0;
1162 TH1 *h1Clone = (TH1D *)hist1->Clone();
1166 h1Clone->Add(hist2,-1);
1167 h1Clone->Divide(hist2);
1172 //_____________________________________________________________________________
1173 TH1* AlidNdPtHelper::CalcRelativeDifferenceFun(const TH1 *const hist1, TF1 *const fun)
1176 // calculate rel. difference
1177 // between histogram and function
1179 if(!hist1) return 0;
1182 TH1 *h1Clone = (TH1D *)hist1->Clone();
1186 h1Clone->Add(fun,-1);
1187 h1Clone->Divide(hist1);
1192 //_____________________________________________________________________________
1193 TH1* AlidNdPtHelper::NormalizeToEvent(const TH2 *const hist1, const TH1 *const hist2)
1195 // normalise to event for a given multiplicity bin
1196 // return pt histogram
1198 if(!hist1) return 0;
1199 if(!hist2) return 0;
1202 Int_t nbinsX = hist1->GetNbinsX();
1203 //Int_t nbinsY = hist1->GetNbinsY();
1206 for(Int_t i=0; i<=nbinsX+1; i++) {
1207 sprintf(name,"mom_%d",i);
1208 TH1D *hist = (TH1D*)hist1->ProjectionY(name,i+1,i+1);
1210 sprintf(name,"mom_norm");
1212 histNorm = (TH1D *)hist->Clone(name);
1216 Double_t nbEvents = hist2->GetBinContent(i);
1217 if(!nbEvents) { nbEvents = 1.; };
1219 hist->Scale(1./nbEvents);
1220 histNorm->Add(hist);
1226 //_____________________________________________________________________________
1227 //THnSparse* AlidNdPtHelper::GenerateCorrMatrix(THnSparse *const hist1, const THnSparse *const hist2, char *const name) {
1228 THnSparse* AlidNdPtHelper::GenerateCorrMatrix(THnSparse *const hist1, const THnSparse *const hist2, const char *name) {
1229 // generate correction matrix
1230 if(!hist1 || !hist2) return 0;
1232 THnSparse *h =(THnSparse*)hist1->Clone(name);;
1233 h->Divide(hist1,hist2,1,1,"B");
1238 //_____________________________________________________________________________
1239 //TH2* AlidNdPtHelper::GenerateCorrMatrix(TH2 *const hist1, TH2 *const hist2, char *const name) {
1240 TH2* AlidNdPtHelper::GenerateCorrMatrix(TH2 *const hist1, TH2 *const hist2, const char *name) {
1241 // generate correction matrix
1242 if(!hist1 || !hist2) return 0;
1244 TH2D *h =(TH2D*)hist1->Clone(name);;
1245 h->Divide(hist1,hist2,1,1,"B");
1250 //_____________________________________________________________________________
1251 //TH1* AlidNdPtHelper::GenerateCorrMatrix(TH1 *const hist1, TH1 *const hist2, char *const name) {
1252 TH1* AlidNdPtHelper::GenerateCorrMatrix(TH1 *const hist1, TH1 *const hist2, const char* name) {
1253 // generate correction matrix
1254 if(!hist1 || !hist2) return 0;
1256 TH1D *h =(TH1D*)hist1->Clone(name);;
1257 h->Divide(hist1,hist2,1,1,"B");
1262 //_____________________________________________________________________________
1263 //THnSparse* AlidNdPtHelper::GenerateContCorrMatrix(THnSparse *const hist1, THnSparse *const hist2, char *const name) {
1264 THnSparse* AlidNdPtHelper::GenerateContCorrMatrix(THnSparse *const hist1, const THnSparse *const hist2, const char* name) {
1265 // generate contamination correction matrix
1266 if(!hist1 || !hist2) return 0;
1268 THnSparse *hist = GenerateCorrMatrix(hist1, hist2, name);
1271 // only for non ZERO bins!!!!
1273 Int_t* coord = new Int_t[hist->GetNdimensions()];
1274 memset(coord, 0, sizeof(Int_t) * hist->GetNdimensions());
1276 for (Long64_t i = 0; i < hist->GetNbins(); ++i) {
1277 Double_t v = hist->GetBinContent(i, coord);
1278 hist->SetBinContent(coord, 1.0-v);
1279 //printf("v %f, hist->GetBinContent(i, coord) %f \n",v,hist->GetBinContent(i, coord));
1280 Double_t err = hist->GetBinError(coord);
1281 hist->SetBinError(coord, err);
1289 //_____________________________________________________________________________
1290 //TH2* AlidNdPtHelper::GenerateContCorrMatrix(TH2 *const hist1, TH2 *const hist2, char *const name) {
1291 TH2* AlidNdPtHelper::GenerateContCorrMatrix(TH2 *const hist1, TH2 *const hist2, const char* name) {
1292 // generate contamination correction matrix
1293 if(!hist1 || !hist2) return 0;
1295 TH2 *hist = GenerateCorrMatrix(hist1, hist2, name);
1298 Int_t nBinsX = hist->GetNbinsX();
1299 Int_t nBinsY = hist->GetNbinsY();
1301 for (Int_t i = 0; i < nBinsX+1; i++) {
1302 for (Int_t j = 0; j < nBinsY+1; j++) {
1303 Double_t cont = hist->GetBinContent(i,j);
1304 hist->SetBinContent(i,j,1.-cont);
1305 Double_t err = hist->GetBinError(i,j);
1306 hist->SetBinError(i,j,err);
1313 //_____________________________________________________________________________
1314 //TH1* AlidNdPtHelper::GenerateContCorrMatrix(TH1 *const hist1, TH1 *const hist2, char *const name) {
1315 TH1* AlidNdPtHelper::GenerateContCorrMatrix(TH1 *const hist1, TH1 *const hist2, const char* name) {
1316 // generate contamination correction matrix
1317 if(!hist1 || !hist2) return 0;
1319 TH1 *hist = GenerateCorrMatrix(hist1, hist2, name);
1322 Int_t nBinsX = hist->GetNbinsX();
1324 for (Int_t i = 0; i < nBinsX+1; i++) {
1325 Double_t cont = hist->GetBinContent(i);
1326 hist->SetBinContent(i,1.-cont);
1327 Double_t err = hist->GetBinError(i);
1328 hist->SetBinError(i,err);
1334 //_____________________________________________________________________________
1335 const AliESDVertex* AlidNdPtHelper::GetTPCVertexZ(const AliESDEvent* const esdEvent, const AlidNdPtEventCuts *const evtCuts, const AlidNdPtAcceptanceCuts *const accCuts, const AliESDtrackCuts *const trackCuts, Float_t fraction, Int_t ntracksMin){
1341 ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
1345 if(!evtCuts || !accCuts || !trackCuts)
1347 ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
1351 Double_t vtxpos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
1352 Double_t vtxsigma[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
1353 AliESDVertex vtx0(vtxpos,vtxsigma);
1355 Double_t maxDCAr = accCuts->GetMaxDCAr();
1356 Double_t maxDCAz = accCuts->GetMaxDCAz();
1357 Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
1360 Int_t ntracks = esdEvent->GetNumberOfTracks();
1361 TVectorD ztrack(ntracks);
1362 Double_t dca[2],cov[3];
1364 for (Int_t i=0;i <ntracks; i++){
1365 AliESDtrack *t = esdEvent->GetTrack(i);
1367 if (!t->GetTPCInnerParam()) continue;
1368 if (t->GetTPCNcls()<minTPCClust) continue;
1371 Double_t x[3]; t->GetXYZ(x);
1372 Double_t b[3]; AliTracker::GetBxByBz(x,b);
1373 const Double_t kMaxStep = 1; //max step over the material
1375 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
1376 if (!tpcTrack->PropagateToDCABxByBz(&vtx0,b,kMaxStep,dca,cov)) continue;
1379 if (TMath::Abs(dca[0])>maxDCAr) continue;
1380 //if (TMath::Sqrt(cov[0])>sigmaXYcut) continue;
1381 if (TMath::Abs(tpcTrack->GetZ())>maxDCAz) continue;
1383 ztrack[counter]=tpcTrack->GetZ();
1386 if(tpcTrack) delete tpcTrack;
1390 // Find LTM z position
1392 Double_t mean=0, sigma=0;
1393 if (counter<ntracksMin) return 0;
1395 Int_t nused = TMath::Nint(counter*fraction);
1396 if (nused==counter) nused=counter-1;
1398 AliMathBase::EvaluateUni(counter, ztrack.GetMatrixArray(), mean,sigma, TMath::Nint(counter*fraction));
1399 sigma/=TMath::Sqrt(nused);
1401 mean = TMath::Mean(counter, ztrack.GetMatrixArray());
1402 sigma = TMath::RMS(counter, ztrack.GetMatrixArray());
1403 sigma/=TMath::Sqrt(counter-1);
1407 const AliESDVertex* vertex = new AliESDVertex(vtxpos, vtxsigma);
1411 //_____________________________________________________________________________
1412 Int_t AlidNdPtHelper::GetSPDMBTrackMult(const AliESDEvent* const esdEvent, Float_t deltaThetaCut, Float_t deltaPhiCut)
1415 // SPD track multiplicity
1419 const AliMultiplicity* mult = esdEvent->GetMultiplicity();
1423 // get multiplicity from SPD tracklets
1424 Int_t inputCount = 0;
1425 for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
1427 //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
1429 Float_t phi = mult->GetPhi(i);
1431 phi += TMath::Pi() * 2;
1432 Float_t deltaPhi = mult->GetDeltaPhi(i);
1433 Float_t deltaTheta = mult->GetDeltaTheta(i);
1435 if (TMath::Abs(deltaPhi) > 1)
1436 printf("WARNING: Very high Delta Phi: %d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), deltaPhi);
1438 if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
1441 if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
1450 //_____________________________________________________________________________
1451 Int_t AlidNdPtHelper::GetSPDMBPrimTrackMult(const AliESDEvent* const esdEvent, AliStack* const stack, Float_t deltaThetaCut, Float_t deltaPhiCut)
1454 // SPD track multiplicity
1458 const AliMultiplicity* mult = esdEvent->GetMultiplicity();
1462 // get multiplicity from SPD tracklets
1463 Int_t inputCount = 0;
1464 for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
1466 //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
1468 Float_t phi = mult->GetPhi(i);
1470 phi += TMath::Pi() * 2;
1471 Float_t deltaPhi = mult->GetDeltaPhi(i);
1472 Float_t deltaTheta = mult->GetDeltaTheta(i);
1474 if (TMath::Abs(deltaPhi) > 1)
1475 printf("WARNING: Very high Delta Phi: %d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), deltaPhi);
1477 if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
1480 if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
1484 if (mult->GetLabel(i, 0) < 0 || mult->GetLabel(i, 0) != mult->GetLabel(i, 1) ||
1485 !stack->IsPhysicalPrimary(mult->GetLabel(i, 0)))