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 **************************************************************************/
17 #include <TParticle.h>
18 #include <TParticlePDG.h>
29 #include <TLorentzVector.h>
31 #include <AliHeader.h>
36 #include <AliESDEvent.h>
37 #include <AliMCEvent.h>
38 #include <AliESDVertex.h>
39 #include <AliVertexerTracks.h>
41 #include <AliGenEventHeader.h>
42 #include <AliGenPythiaEventHeader.h>
43 #include <AliGenCocktailEventHeader.h>
44 #include <AliGenDPMjetEventHeader.h>
46 #include <AliMathBase.h>
47 #include <AliESDtrackCuts.h>
48 #include "dNdPt/AlidNdPtEventCuts.h"
49 #include "dNdPt/AlidNdPtAcceptanceCuts.h"
50 #include "dNdPt/AlidNdPtHelper.h"
52 //____________________________________________________________________
53 ClassImp(AlidNdPtHelper)
55 //____________________________________________________________________
56 const AliESDVertex* AlidNdPtHelper::GetVertex(AliESDEvent* aEsd, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts, AnalysisMode analysisMode, Bool_t debug, Bool_t bRedoTPC, Bool_t bUseMeanVertex)
58 // Get the vertex from the ESD and returns it if the vertex is valid
60 // Second argument decides which vertex is used (this selects
61 // also the quality criteria that are applied)
65 ::Error("AlidNdPtHelper::GetVertex()","esd event is NULL");
69 if(!evtCuts || !accCuts || !trackCuts)
71 ::Error("AlidNdPtHelper::GetVertex()","cuts not available");
75 const AliESDVertex* vertex = 0;
76 AliESDVertex *initVertex = 0;
77 if (analysisMode == kSPD || analysisMode == kTPCITS || analysisMode == kTPCSPDvtx || analysisMode == kTPCSPDvtxUpdate)
79 vertex = aEsd->GetPrimaryVertexSPD();
81 Printf("AlidNdPtHelper::GetVertex: Returning SPD vertex");
83 else if (analysisMode == kTPC)
87 Double_t kBz = aEsd->GetMagneticField();
88 AliVertexerTracks vertexer(kBz);
91 Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
92 Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
93 initVertex = new AliESDVertex(pos,err);
94 vertexer.SetVtxStart(initVertex);
95 vertexer.SetConstraintOn();
98 Double_t maxDCAr = accCuts->GetMaxDCAr();
99 Double_t maxDCAz = accCuts->GetMaxDCAz();
100 Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
102 //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);
103 vertexer.SetTPCMode(0.1,1.0,5.0,minTPCClust,1,3.,0.1,2.0,maxDCAr,maxDCAz,1,4);
105 // TPC track preselection
106 Int_t ntracks = aEsd->GetNumberOfTracks();
107 TObjArray array(ntracks);
108 UShort_t *id = new UShort_t[ntracks];
111 for (Int_t i=0;i <ntracks; i++) {
112 AliESDtrack *t = aEsd->GetTrack(i);
114 if (t->Charge() == 0) continue;
115 if (!t->GetTPCInnerParam()) continue;
116 if (t->GetTPCNcls()<vertexer.GetMinClusters()) continue;
117 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
119 array.AddLast(tpcTrack);
120 id[count] = (UShort_t)t->GetID();
124 AliESDVertex *vTPC = vertexer.VertexForSelectedTracks(&array,id, kTRUE, kTRUE, bUseMeanVertex);
126 // set recreated TPC vertex
127 aEsd->SetPrimaryVertexTPC(vTPC);
129 for (Int_t i=0; i<aEsd->GetNumberOfTracks(); i++) {
130 AliESDtrack *t = aEsd->GetTrack(i);
131 t->RelateToVertexTPC(vTPC, kBz, kVeryBig);
136 delete [] id; id=NULL;
139 vertex = aEsd->GetPrimaryVertexTPC();
141 Printf("AlidNdPtHelper::GetVertex: Returning vertex from tracks");
144 Printf("AlidNdPtHelper::GetVertex: ERROR: Invalid second argument %d", analysisMode);
148 Printf("AlidNdPtHelper::GetVertex: No vertex found in ESD");
154 Printf("AlidNdPtHelper::GetVertex: Returning valid vertex: %s", vertex->GetTitle());
158 if(initVertex) delete initVertex; initVertex=NULL;
162 //____________________________________________________________________
163 Bool_t AlidNdPtHelper::TestRecVertex(const AliESDVertex* vertex, AnalysisMode analysisMode, Bool_t debug)
165 // Checks if a vertex meets the needed quality criteria
166 if(!vertex) return kFALSE;
168 Float_t requiredZResolution = -1;
169 if (analysisMode == kSPD || analysisMode == kTPCITS || analysisMode == kTPCSPDvtx || analysisMode == kTPCSPDvtxUpdate)
171 requiredZResolution = 0.1;
173 else if (analysisMode == kTPC)
174 requiredZResolution = 10.;
176 // check Ncontributors
177 if (vertex->GetNContributors() <= 0) {
179 Printf("AlidNdPtHelper::GetVertex: NContributors() <= 0: %d",vertex->GetNContributors());
180 Printf("AlidNdPtHelper::GetVertex: NIndices(): %d",vertex->GetNIndices());
187 Double_t zRes = vertex->GetZRes();
189 Printf("AlidNdPtHelper::GetVertex: UNEXPECTED: resolution is 0.");
193 if (zRes > requiredZResolution) {
195 Printf("AlidNdPtHelper::TestVertex: Resolution too poor %f (required: %f", zRes, requiredZResolution);
202 //____________________________________________________________________
203 Bool_t AlidNdPtHelper::IsPrimaryParticle(AliStack* stack, Int_t idx, ParticleMode particleMode)
205 // check primary particles
206 // depending on the particle mode
208 if(!stack) return kFALSE;
210 TParticle* particle = stack->Particle(idx);
211 if (!particle) return kFALSE;
213 // only charged particles
214 Double_t charge = particle->GetPDG()->Charge()/3.;
215 if (charge == 0.0) return kFALSE;
217 Int_t pdg = TMath::Abs(particle->GetPdgCode());
220 Bool_t prim = stack->IsPhysicalPrimary(idx);
222 if(particleMode==kMCPion) {
223 if(prim && pdg==kPiPlus) return kTRUE;
227 if (particleMode==kMCKaon) {
228 if(prim && pdg==kKPlus) return kTRUE;
232 if (particleMode==kMCProton) {
233 if(prim && pdg==kProton) return kTRUE;
240 //____________________________________________________________________
241 Bool_t AlidNdPtHelper::IsCosmicTrack(TObjArray *allChargedTracks, AliESDtrack *track1, Int_t trackIdx, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts)
244 // check cosmic tracks
246 if(!allChargedTracks) return kFALSE;
247 if(!track1) return kFALSE;
248 if(!accCuts) return kFALSE;
249 if(!trackCuts) return kFALSE;
251 Int_t entries = allChargedTracks->GetEntries();
252 for(Int_t i=0; i<entries;++i)
255 // exclude the same tracks
257 if(i == trackIdx) continue;
259 AliESDtrack *track2 = (AliESDtrack*)allChargedTracks->At(i);
260 if(!track2) continue;
261 if(track2->Charge()==0) continue;
264 if(track1->Pt() > 6. && track2->Pt() > 6. && (track1->Charge() + track2->Charge()) == 0 )
266 printf("track1->Theta() %f, track1->Eta() %f, track1->Phi() %f, track1->Charge() %d \n", track1->Theta(), track1->Eta(), track1->Phi(), track1->Charge());
267 printf("track2->Theta() %f, track2->Eta() %f, track2->Phi() %f, track2->Charge() %d \n", track2->Theta(), track2->Eta(), track2->Phi(), track2->Charge());
269 printf("deta %f, dphi %f, dq %d \n", track1->Eta()-track2->Eta(), track1->Phi()-track2->Phi(), track1->Charge()+track2->Charge());
275 // cosmic tracks in TPC
277 //if( TMath::Abs( track1->Theta() - track2->Theta() ) < 0.004 &&
278 // ((TMath::Abs(track1->Phi()-track2->Phi()) - TMath::Pi() )<0.004) &&
279 if( (track1->Pt()-track2->Pt()) < 0.1 && track1->Pt() > 6.0 &&
280 (track1->Charge()+track2->Charge()) == 0 )
282 //printf("COSMIC candidate \n");
283 printf("track1->Theta() %f, track1->Eta() %f, track1->Phi() %f, track1->Charge() %d \n", track1->Theta(), track1->Eta(), track1->Phi(), track1->Charge());
284 printf("track2->Theta() %f, track2->Eta() %f, track2->Phi() %f, track2->Charge() %d \n", track2->Theta(), track2->Eta(), track2->Phi(), track2->Charge());
285 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());
293 //____________________________________________________________________
294 void AlidNdPtHelper::PrintConf(AnalysisMode analysisMode, AliTriggerAnalysis::Trigger trigger)
297 // Prints the given configuration
300 TString str(">>>> Running with ");
302 switch (analysisMode)
304 case kInvalid: str += "invalid setting"; break;
305 case kSPD : str += "SPD-only"; break;
306 case kTPC : str += "TPC-only"; break;
307 case kTPCITS : str += "Global tracking"; break;
308 case kTPCSPDvtx : str += "TPC tracking + SPD event vertex"; break;
309 case kTPCSPDvtxUpdate : str += "TPC tracks updated with SPD event vertex point"; break;
310 case kMCRec : str += "TPC tracking + Replace rec. with MC values"; break;
312 str += " and trigger ";
314 str += AliTriggerAnalysis::GetTriggerName(trigger);
318 Printf("%s", str.Data());
321 //____________________________________________________________________
322 Int_t AlidNdPtHelper::ConvertPdgToPid(TParticle *particle) {
324 // Convert Abs(pdg) to pid
325 // (0 - e, 1 - muons, 2 - pions, 3 - kaons, 4 - protons, 5 -all rest)
329 if (TMath::Abs(particle->GetPdgCode()) == kElectron) { pid = 0; }
330 else if (TMath::Abs(particle->GetPdgCode()) == kMuonMinus) { pid = 1; }
331 else if (TMath::Abs(particle->GetPdgCode()) == kPiPlus) { pid = 2; }
332 else if (TMath::Abs(particle->GetPdgCode()) == kKPlus) { pid = 3; }
333 else if (TMath::Abs(particle->GetPdgCode()) == kProton) { pid = 4; }
339 //_____________________________________________________________________________
340 TH1F* AlidNdPtHelper::CreateResHisto(TH2F* hRes2, TH1F **phMean, Int_t integ, Bool_t drawBinFits, Int_t minHistEntries)
343 // Create mean and resolution
346 TVirtualPad* currentPad = gPad;
347 TAxis* axis = hRes2->GetXaxis();
348 Int_t nBins = axis->GetNbins();
349 Bool_t overflowBinFits = kFALSE;
351 if (axis->GetXbins()->GetSize()){
352 hRes = new TH1F("hRes", "", nBins, axis->GetXbins()->GetArray());
353 hMean = new TH1F("hMean", "", nBins, axis->GetXbins()->GetArray());
356 hRes = new TH1F("hRes", "", nBins, axis->GetXmin(), axis->GetXmax());
357 hMean = new TH1F("hMean", "", nBins, axis->GetXmin(), axis->GetXmax());
360 hRes->SetStats(false);
361 hRes->SetOption("E");
362 hRes->SetMinimum(0.);
364 hMean->SetStats(false);
365 hMean->SetOption("E");
367 // create the fit function
368 TF1 * fitFunc = new TF1("G","[0]*exp(-(x-[1])*(x-[1])/(2.*[2]*[2]))",-3,3);
370 fitFunc->SetLineWidth(2);
371 fitFunc->SetFillStyle(0);
372 // create canvas for fits
373 TCanvas* canBinFits = NULL;
374 Int_t nPads = (overflowBinFits) ? nBins+2 : nBins;
375 Int_t nx = Int_t(sqrt(nPads-1.));// + 1;
376 Int_t ny = (nPads-1) / nx + 1;
378 canBinFits = (TCanvas*)gROOT->FindObject("canBinFits");
379 if (canBinFits) delete canBinFits;
380 canBinFits = new TCanvas("canBinFits", "fits of bins", 200, 100, 500, 700);
381 canBinFits->Divide(nx, ny);
384 // loop over x bins and fit projection
385 Int_t dBin = ((overflowBinFits) ? 1 : 0);
386 for (Int_t bin = 1-dBin; bin <= nBins+dBin; bin++) {
387 if (drawBinFits) canBinFits->cd(bin + dBin);
388 Int_t bin0=TMath::Max(bin-integ,0);
389 Int_t bin1=TMath::Min(bin+integ,nBins);
390 TH1D* hBin = hRes2->ProjectionY("hBin", bin0, bin1);
392 if (hBin->GetEntries() > minHistEntries) {
393 fitFunc->SetParameters(hBin->GetMaximum(),hBin->GetMean(),hBin->GetRMS());
394 hBin->Fit(fitFunc,"s");
395 Double_t sigma = TMath::Abs(fitFunc->GetParameter(2));
398 hRes->SetBinContent(bin, TMath::Abs(fitFunc->GetParameter(2)));
399 hMean->SetBinContent(bin, fitFunc->GetParameter(1));
402 hRes->SetBinContent(bin, 0.);
403 hMean->SetBinContent(bin,0);
405 hRes->SetBinError(bin, fitFunc->GetParError(2));
406 hMean->SetBinError(bin, fitFunc->GetParError(1));
412 hRes->SetBinContent(bin, 0.);
413 hRes->SetBinError(bin, 0.);
414 hMean->SetBinContent(bin, 0.);
415 hMean->SetBinError(bin, 0.);
422 sprintf(name, "%s < %.4g", axis->GetTitle(), axis->GetBinUpEdge(bin));
423 } else if (bin == nBins+1) {
424 sprintf(name, "%.4g < %s", axis->GetBinLowEdge(bin), axis->GetTitle());
426 sprintf(name, "%.4g < %s < %.4g", axis->GetBinLowEdge(bin),
427 axis->GetTitle(), axis->GetBinUpEdge(bin));
429 canBinFits->cd(bin + dBin);
430 hBin->SetTitle(name);
431 hBin->SetStats(kTRUE);
433 canBinFits->Update();
434 canBinFits->Modified();
435 canBinFits->Update();
447 //_____________________________________________________________________________
448 TH1F* AlidNdPtHelper::MakeResol(TH2F * his, Int_t integ, Bool_t type, Bool_t drawBins, Int_t minHistEntries){
449 // Create resolution histograms
451 TH1F *hisr=0, *hism=0;
452 if (!gPad) new TCanvas;
453 hisr = CreateResHisto(his,&hism,integ,drawBins,minHistEntries);
454 if (type) return hism;
460 //_____________________________________________________________________________
461 TObjArray* AlidNdPtHelper::GetAllChargedTracks(AliESDEvent *esdEvent, AnalysisMode analysisMode)
464 // all charged TPC particles
466 TObjArray *allTracks = new TObjArray();
467 if(!allTracks) return allTracks;
469 AliESDtrack *track=0;
470 for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++)
472 if(analysisMode == AlidNdPtHelper::kTPC) {
473 // track must be deleted by user
474 track = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent,iTrack);
477 else if (analysisMode == AlidNdPtHelper::kTPCSPDvtx || AlidNdPtHelper::kTPCSPDvtxUpdate)
479 // track must be deleted by the user
480 track = AlidNdPtHelper::GetTPCOnlyTrackSPDvtx(esdEvent,iTrack,kFALSE);
484 track=esdEvent->GetTrack(iTrack);
488 if(track->Charge()==0) {
489 if(analysisMode == AlidNdPtHelper::kTPC || analysisMode == AlidNdPtHelper::kTPCSPDvtx ||
490 analysisMode == AlidNdPtHelper::kTPCSPDvtxUpdate)
492 delete track; continue;
498 allTracks->Add(track);
501 if(analysisMode == AlidNdPtHelper::kTPC || analysisMode == AlidNdPtHelper::kTPCSPDvtx ||
502 analysisMode == AlidNdPtHelper::kTPCSPDvtxUpdate) {
504 allTracks->SetOwner(kTRUE);
510 //_____________________________________________________________________________
511 AliESDtrack *AlidNdPtHelper::GetTPCOnlyTrackSPDvtx(AliESDEvent* esdEvent, Int_t iTrack, Bool_t bUpdate)
514 // Create ESD tracks from TPCinner parameters.
515 // Propagte to DCA to SPD vertex.
516 // Update using SPD vertex point (parameter)
518 // It is user responsibility to delete these tracks
521 if (!esdEvent) return NULL;
522 if (!esdEvent->GetPrimaryVertexSPD() ) { return NULL; }
523 if (!esdEvent->GetPrimaryVertexSPD()->GetStatus() ) { return NULL; }
526 AliESDtrack* track = esdEvent->GetTrack(iTrack);
530 // create new ESD track
531 AliESDtrack *tpcTrack = new AliESDtrack();
533 // relate TPC-only tracks (TPCinner) to SPD vertex
534 AliExternalTrackParam cParam;
536 track->RelateToVertexTPC(esdEvent->GetPrimaryVertexSPD(),esdEvent->GetMagneticField(),kVeryBig,&cParam);
537 track->Set(cParam.GetX(),cParam.GetAlpha(),cParam.GetParameter(),cParam.GetCovariance());
539 // reject fake tracks
540 if(track->Pt() > 10000.) {
541 ::Error("Exclude no physical tracks","pt>10000. GeV");
547 track->RelateToVertexTPC(esdEvent->GetPrimaryVertexSPD(), esdEvent->GetMagneticField(), kVeryBig);
550 // only true if we have a tpc track
551 if (!track->FillTPCOnlyTrack(*tpcTrack))
560 //_____________________________________________________________________________
561 Int_t AlidNdPtHelper::GetTPCMBTrackMult(AliESDEvent *esdEvent, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts)
564 // get MB event track multiplicity
568 ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","esd event is NULL");
572 if(!evtCuts || !accCuts || !trackCuts)
574 ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","cuts not available");
579 Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
580 Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
581 AliESDVertex vtx0(pos,err);
584 Float_t maxDCAr = accCuts->GetMaxDCAr();
585 Float_t maxDCAz = accCuts->GetMaxDCAz();
586 Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
588 Int_t ntracks = esdEvent->GetNumberOfTracks();
589 Double_t dca[2],cov[3];
591 for (Int_t i=0;i <ntracks; i++){
592 AliESDtrack *t = esdEvent->GetTrack(i);
594 if (t->Charge() == 0) continue;
595 if (!t->GetTPCInnerParam()) continue;
596 if (t->GetTPCNcls()<minTPCClust) continue;
598 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
599 if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov))
601 if(tpcTrack) delete tpcTrack;
605 if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
606 if(tpcTrack) delete tpcTrack;
612 if(tpcTrack) delete tpcTrack;
618 //_____________________________________________________________________________
619 Int_t AlidNdPtHelper::GetTPCMBPrimTrackMult(AliESDEvent *esdEvent, AliStack * stack, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts)
622 // get MB primary event track multiplicity
626 ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
632 ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
636 if(!evtCuts || !accCuts || !trackCuts)
638 ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","cuts not available");
643 Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
644 Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
645 AliESDVertex vtx0(pos,err);
648 Float_t maxDCAr = accCuts->GetMaxDCAr();
649 Float_t maxDCAz = accCuts->GetMaxDCAz();
650 Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
653 Int_t ntracks = esdEvent->GetNumberOfTracks();
654 Double_t dca[2],cov[3];
656 for (Int_t i=0;i <ntracks; i++){
657 AliESDtrack *t = esdEvent->GetTrack(i);
659 if (t->Charge() == 0) continue;
660 if (!t->GetTPCInnerParam()) continue;
661 if (t->GetTPCNcls()<minTPCClust) continue;
663 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
664 if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov))
666 if(tpcTrack) delete tpcTrack;
670 if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
671 if(tpcTrack) delete tpcTrack;
675 Int_t label = TMath::Abs(t->GetLabel());
676 TParticle *part = stack->Particle(label);
678 if(tpcTrack) delete tpcTrack;
681 if(!stack->IsPhysicalPrimary(label))
683 if(tpcTrack) delete tpcTrack;
689 if(tpcTrack) delete tpcTrack;
699 //_____________________________________________________________________________
700 Int_t AlidNdPtHelper::GetMCTrueTrackMult(AliMCEvent *mcEvent, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts)
703 // calculate mc event true track multiplicity
705 if(!mcEvent) return 0;
711 stack = mcEvent->Stack();
712 if (!stack) return 0;
714 Bool_t isEventOK = evtCuts->AcceptMCEvent(mcEvent);
715 if(!isEventOK) return 0;
717 Int_t nPart = stack->GetNtrack();
718 for (Int_t iMc = 0; iMc < nPart; ++iMc)
720 TParticle* particle = stack->Particle(iMc);
724 // only charged particles
725 Double_t charge = particle->GetPDG()->Charge()/3.;
730 Bool_t prim = stack->IsPhysicalPrimary(iMc);
734 if(accCuts->AcceptTrack(particle))
743 //_______________________________________________________________________
744 void AlidNdPtHelper::PrintMCInfo(AliStack *pStack,Int_t label)
746 // print information about particles in the stack
749 label = TMath::Abs(label);
750 TParticle *part = pStack->Particle(label);
751 Printf("########################");
752 Printf("%s:%d %d UniqueID %d PDG %d P %3.3f",(char*)__FILE__,__LINE__,label,part->GetUniqueID(),part->GetPdgCode(),part->P());
754 TParticle* mother = part;
755 Int_t imo = part->GetFirstMother();
756 Int_t nprim = pStack->GetNprimary();
758 while((imo >= nprim)) {
759 mother = pStack->Particle(imo);
760 Printf("Mother %s:%d Label %d UniqueID %d PDG %d P %3.3f",(char*)__FILE__,__LINE__,imo,mother->GetUniqueID(),mother->GetPdgCode(),mother->P());
762 imo = mother->GetFirstMother();
765 Printf("########################");
769 //_____________________________________________________________________________
770 TH1* AlidNdPtHelper::GetContCorrHisto(TH1 *hist)
773 // get contamination histogram
777 Int_t nbins = hist->GetNbinsX();
778 TH1 *h_cont = (TH1D *)hist->Clone();
780 for(Int_t i=0; i<=nbins+1; i++) {
781 Double_t binContent = hist->GetBinContent(i);
782 Double_t binError = hist->GetBinError(i);
784 h_cont->SetBinContent(i,1.-binContent);
785 h_cont->SetBinError(i,binError);
792 //_____________________________________________________________________________
793 TH1* AlidNdPtHelper::ScaleByBinWidth(TH1 *hist)
796 // scale by bin width
800 TH1 *h_scale = (TH1D *)hist->Clone();
801 h_scale->Scale(1.,"width");
806 //_____________________________________________________________________________
807 TH1* AlidNdPtHelper::CalcRelativeDifference(TH1 *hist1, TH1 *hist2)
810 // calculate rel. difference
816 TH1 *h1_clone = (TH1D *)hist1->Clone();
820 h1_clone->Add(hist2,-1);
821 h1_clone->Divide(hist2);
826 //_____________________________________________________________________________
827 TH1* AlidNdPtHelper::CalcRelativeDifferenceFun(TH1 *hist1, TF1 *fun)
830 // calculate rel. difference
831 // between histogram and function
836 TH1 *h1_clone = (TH1D *)hist1->Clone();
840 h1_clone->Add(fun,-1);
841 h1_clone->Divide(hist1);
846 //_____________________________________________________________________________
847 TH1* AlidNdPtHelper::NormalizeToEvent(TH2 *hist1, TH1 *hist2)
849 // normalise to event for a given multiplicity bin
850 // return pt histogram
856 Int_t nbinsX = hist1->GetNbinsX();
857 //Int_t nbinsY = hist1->GetNbinsY();
860 for(Int_t i=0; i<=nbinsX+1; i++) {
861 sprintf(name,"mom_%d",i);
862 TH1D *hist = (TH1D*)hist1->ProjectionY(name,i+1,i+1);
864 sprintf(name,"mom_norm");
866 hist_norm = (TH1D *)hist->Clone(name);
870 Double_t nbEvents = hist2->GetBinContent(i);
871 if(!nbEvents) { nbEvents = 1.; };
873 hist->Scale(1./nbEvents);
874 hist_norm->Add(hist);
880 //_____________________________________________________________________________
881 THnSparse* AlidNdPtHelper::GenerateCorrMatrix(THnSparse *hist1, THnSparse *hist2, char *name) {
882 // generate correction matrix
883 if(!hist1 || !hist2) return 0;
885 THnSparse *h =(THnSparse*)hist1->Clone(name);;
886 h->Divide(hist1,hist2,1,1,"B");
891 //_____________________________________________________________________________
892 TH2* AlidNdPtHelper::GenerateCorrMatrix(TH2 *hist1, TH2 *hist2, char *name) {
893 // generate correction matrix
894 if(!hist1 || !hist2) return 0;
896 TH2D *h =(TH2D*)hist1->Clone(name);;
897 h->Divide(hist1,hist2,1,1,"B");
902 //_____________________________________________________________________________
903 TH1* AlidNdPtHelper::GenerateCorrMatrix(TH1 *hist1, TH1 *hist2, char *name) {
904 // generate correction matrix
905 if(!hist1 || !hist2) return 0;
907 TH1D *h =(TH1D*)hist1->Clone(name);;
908 h->Divide(hist1,hist2,1,1,"B");
913 //_____________________________________________________________________________
914 THnSparse* AlidNdPtHelper::GenerateContCorrMatrix(THnSparse *hist1, THnSparse *hist2, char *name) {
915 // generate contamination correction matrix
916 if(!hist1 || !hist2) return 0;
918 THnSparse *hist = GenerateCorrMatrix(hist1, hist2, name);
921 // only for non ZERO bins!!!!
923 Int_t* coord = new Int_t[hist->GetNdimensions()];
924 memset(coord, 0, sizeof(Int_t) * hist->GetNdimensions());
926 for (Long64_t i = 0; i < hist->GetNbins(); ++i) {
927 Double_t v = hist->GetBinContent(i, coord);
928 hist->SetBinContent(coord, 1.0-v);
929 //printf("v %f, hist->GetBinContent(i, coord) %f \n",v,hist->GetBinContent(i, coord));
930 Double_t err = hist->GetBinError(coord);
931 hist->SetBinError(coord, err);
939 //_____________________________________________________________________________
940 TH2* AlidNdPtHelper::GenerateContCorrMatrix(TH2 *hist1, TH2 *hist2, char *name) {
941 // generate contamination correction matrix
942 if(!hist1 || !hist2) return 0;
944 TH2 *hist = GenerateCorrMatrix(hist1, hist2, name);
947 Int_t nBinsX = hist->GetNbinsX();
948 Int_t nBinsY = hist->GetNbinsY();
950 for (Int_t i = 0; i < nBinsX+1; i++) {
951 for (Int_t j = 0; j < nBinsY+1; j++) {
952 Double_t cont = hist->GetBinContent(i,j);
953 hist->SetBinContent(i,j,1.-cont);
954 Double_t err = hist->GetBinError(i,j);
955 hist->SetBinError(i,j,err);
962 //_____________________________________________________________________________
963 TH1* AlidNdPtHelper::GenerateContCorrMatrix(TH1 *hist1, TH1 *hist2, char *name) {
964 // generate contamination correction matrix
965 if(!hist1 || !hist2) return 0;
967 TH1 *hist = GenerateCorrMatrix(hist1, hist2, name);
970 Int_t nBinsX = hist->GetNbinsX();
972 for (Int_t i = 0; i < nBinsX+1; i++) {
973 Double_t cont = hist->GetBinContent(i);
974 hist->SetBinContent(i,1.-cont);
975 Double_t err = hist->GetBinError(i);
976 hist->SetBinError(i,err);
982 //_____________________________________________________________________________
983 const AliESDVertex* AlidNdPtHelper::GetTPCVertexZ(AliESDEvent* esdEvent, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts, Float_t fraction, Int_t ntracksMin){
989 ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
993 if(!evtCuts || !accCuts || !trackCuts)
995 ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
999 Double_t vtxpos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
1000 Double_t vtxsigma[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
1001 AliESDVertex vtx0(vtxpos,vtxsigma);
1003 Double_t maxDCAr = accCuts->GetMaxDCAr();
1004 Double_t maxDCAz = accCuts->GetMaxDCAz();
1005 Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
1008 Int_t ntracks = esdEvent->GetNumberOfTracks();
1009 TVectorD ztrack(ntracks);
1010 Double_t dca[2],cov[3];
1012 for (Int_t i=0;i <ntracks; i++){
1013 AliESDtrack *t = esdEvent->GetTrack(i);
1015 if (!t->GetTPCInnerParam()) continue;
1016 if (t->GetTPCNcls()<minTPCClust) continue;
1018 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
1019 if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov)) continue;
1022 if (TMath::Abs(dca[0])>maxDCAr) continue;
1023 //if (TMath::Sqrt(cov[0])>sigmaXYcut) continue;
1024 if (TMath::Abs(tpcTrack->GetZ())>maxDCAz) continue;
1026 ztrack[counter]=tpcTrack->GetZ();
1029 if(tpcTrack) delete tpcTrack;
1033 // Find LTM z position
1035 Double_t mean=0, sigma=0;
1036 if (counter<ntracksMin) return 0;
1038 Int_t nused = TMath::Nint(counter*fraction);
1039 if (nused==counter) nused=counter-1;
1041 AliMathBase::EvaluateUni(counter, ztrack.GetMatrixArray(), mean,sigma, TMath::Nint(counter*fraction));
1042 sigma/=TMath::Sqrt(nused);
1044 mean = TMath::Mean(counter, ztrack.GetMatrixArray());
1045 sigma = TMath::RMS(counter, ztrack.GetMatrixArray());
1046 sigma/=TMath::Sqrt(counter-1);
1050 const AliESDVertex* vertex = new AliESDVertex(vtxpos, vtxsigma);
1054 //_____________________________________________________________________________
1055 Int_t AlidNdPtHelper::GetSPDMBTrackMult(AliESDEvent* esdEvent, Float_t deltaThetaCut, Float_t deltaPhiCut)
1058 // SPD track multiplicity
1062 const AliMultiplicity* mult = esdEvent->GetMultiplicity();
1066 // get multiplicity from SPD tracklets
1067 Int_t inputCount = 0;
1068 for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
1070 //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
1072 Float_t phi = mult->GetPhi(i);
1074 phi += TMath::Pi() * 2;
1075 Float_t deltaPhi = mult->GetDeltaPhi(i);
1076 Float_t deltaTheta = mult->GetDeltaTheta(i);
1078 if (TMath::Abs(deltaPhi) > 1)
1079 printf("WARNING: Very high Delta Phi: %d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), deltaPhi);
1081 if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
1084 if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
1093 //_____________________________________________________________________________
1094 Int_t AlidNdPtHelper::GetSPDMBPrimTrackMult(AliESDEvent* esdEvent, AliStack* stack, Float_t deltaThetaCut, Float_t deltaPhiCut)
1097 // SPD track multiplicity
1101 const AliMultiplicity* mult = esdEvent->GetMultiplicity();
1105 // get multiplicity from SPD tracklets
1106 Int_t inputCount = 0;
1107 for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
1109 //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
1111 Float_t phi = mult->GetPhi(i);
1113 phi += TMath::Pi() * 2;
1114 Float_t deltaPhi = mult->GetDeltaPhi(i);
1115 Float_t deltaTheta = mult->GetDeltaTheta(i);
1117 if (TMath::Abs(deltaPhi) > 1)
1118 printf("WARNING: Very high Delta Phi: %d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), deltaPhi);
1120 if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
1123 if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
1127 if (mult->GetLabel(i, 0) < 0 || mult->GetLabel(i, 0) != mult->GetLabel(i, 1) ||
1128 !stack->IsPhysicalPrimary(mult->GetLabel(i, 0)))