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>
30 #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 ||
78 analysisMode == kMCRec || analysisMode == kMCPion || analysisMode == kMCKaon ||
79 analysisMode == kMCProton || analysisMode ==kPlus || analysisMode ==kMinus )
81 vertex = aEsd->GetPrimaryVertexSPD();
83 Printf("AlidNdPtHelper::GetVertex: Returning SPD vertex");
85 else if (analysisMode == kTPC || analysisMode == kMCRec ||
86 analysisMode == kMCPion || analysisMode == kMCKaon || analysisMode == kMCProton ||
87 analysisMode == kPlus || analysisMode == kMinus)
91 Double_t kBz = aEsd->GetMagneticField();
92 AliVertexerTracks vertexer(kBz);
95 Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
96 Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
97 initVertex = new AliESDVertex(pos,err);
98 vertexer.SetVtxStart(initVertex);
99 vertexer.SetConstraintOn();
102 Double_t maxDCAr = accCuts->GetMaxDCAr();
103 Double_t maxDCAz = accCuts->GetMaxDCAz();
104 Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
106 //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);
107 vertexer.SetTPCMode(0.1,1.0,5.0,minTPCClust,1,3.,0.1,2.0,maxDCAr,maxDCAz,1,4);
109 // TPC track preselection
110 Int_t ntracks = aEsd->GetNumberOfTracks();
111 TObjArray array(ntracks);
112 UShort_t *id = new UShort_t[ntracks];
115 for (Int_t i=0;i <ntracks; i++) {
116 AliESDtrack *t = aEsd->GetTrack(i);
118 if (t->Charge() == 0) continue;
119 if (!t->GetTPCInnerParam()) continue;
120 if (t->GetTPCNcls()<vertexer.GetMinClusters()) continue;
121 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
123 array.AddLast(tpcTrack);
124 id[count] = (UShort_t)t->GetID();
128 AliESDVertex *vTPC = vertexer.VertexForSelectedTracks(&array,id, kTRUE, kTRUE, bUseMeanVertex);
130 // set recreated TPC vertex
131 aEsd->SetPrimaryVertexTPC(vTPC);
133 for (Int_t i=0; i<aEsd->GetNumberOfTracks(); i++) {
134 AliESDtrack *t = aEsd->GetTrack(i);
135 t->RelateToVertexTPC(vTPC, kBz, kVeryBig);
140 delete [] id; id=NULL;
143 vertex = aEsd->GetPrimaryVertexTPC();
145 Printf("AlidNdPtHelper::GetVertex: Returning vertex from tracks");
148 Printf("AlidNdPtHelper::GetVertex: ERROR: Invalid second argument %d", analysisMode);
152 Printf("AlidNdPtHelper::GetVertex: No vertex found in ESD");
158 Printf("AlidNdPtHelper::GetVertex: Returning valid vertex: %s", vertex->GetTitle());
162 if(initVertex) delete initVertex; initVertex=NULL;
166 //____________________________________________________________________
167 Bool_t AlidNdPtHelper::TestRecVertex(const AliESDVertex* vertex, AnalysisMode analysisMode, Bool_t debug)
169 // Checks if a vertex meets the needed quality criteria
170 if(!vertex) return kFALSE;
172 Float_t requiredZResolution = -1;
173 if (analysisMode == kSPD || analysisMode == kTPCITS || analysisMode == kTPCSPDvtx)
175 requiredZResolution = 0.1;
177 else if (analysisMode == kTPC || analysisMode == kMCRec ||
178 analysisMode == kMCPion || analysisMode == kMCKaon ||
179 analysisMode == kMCProton || analysisMode ==kPlus || analysisMode ==kMinus)
180 requiredZResolution = 10.;
182 // check Ncontributors
183 if (vertex->GetNContributors() <= 0) {
185 Printf("AlidNdPtHelper::GetVertex: NContributors() <= 0: %d",vertex->GetNContributors());
186 Printf("AlidNdPtHelper::GetVertex: NIndices(): %d",vertex->GetNIndices());
193 Double_t zRes = vertex->GetZRes();
195 Printf("AlidNdPtHelper::GetVertex: UNEXPECTED: resolution is 0.");
199 if (zRes > requiredZResolution) {
201 Printf("AlidNdPtHelper::TestVertex: Resolution too poor %f (required: %f", zRes, requiredZResolution);
208 //____________________________________________________________________
209 Bool_t AlidNdPtHelper::IsPrimaryParticle(AliStack* stack, Int_t idx, AnalysisMode analysisMode)
211 // check primary particles
212 // depending on the analysis mode
214 if(!stack) return kFALSE;
216 TParticle* particle = stack->Particle(idx);
217 if (!particle) return kFALSE;
219 // only charged particles
220 Double_t charge = particle->GetPDG()->Charge()/3.;
221 if (charge == 0.0) return kFALSE;
223 Int_t pdg = TMath::Abs(particle->GetPdgCode());
226 Bool_t prim = stack->IsPhysicalPrimary(idx);
228 if(analysisMode==kMCPion) {
229 if(prim && pdg==kPiPlus) return kTRUE;
233 if (analysisMode==kMCKaon) {
234 if(prim && pdg==kKPlus) return kTRUE;
238 if (analysisMode==kMCProton) {
239 if(prim && pdg==kProton) return kTRUE;
246 //____________________________________________________________________
247 void AlidNdPtHelper::PrintConf(AnalysisMode analysisMode, AliTriggerAnalysis::Trigger trigger)
250 // Prints the given configuration
253 TString str(">>>> Running with ");
255 switch (analysisMode)
257 case kInvalid: str += "invalid setting"; break;
258 case kSPD : str += "SPD-only"; break;
259 case kTPC : str += "TPC-only"; break;
260 case kTPCITS : str += "Global tracking"; break;
261 case kTPCSPDvtx : str += "TPC tracking + SPD event vertex"; break;
262 case kMCRec : str += "TPC tracking + Replace rec. with MC values"; break;
263 case kMCPion : str += "TPC tracking + only pion MC tracks"; break;
264 case kMCKaon : str += "TPC tracking + only kaon MC tracks"; break;
265 case kMCProton : str += "TPC tracking + only proton MC tracks"; break;
266 case kPlus: str += "TPC tracking + only positive charged tracks"; break;
267 case kMinus : str += "TPC tracking + only negative charge tracks"; break;
269 str += " and trigger ";
271 str += AliTriggerAnalysis::GetTriggerName(trigger);
275 Printf("%s", str.Data());
278 //____________________________________________________________________
279 Int_t AlidNdPtHelper::ConvertPdgToPid(TParticle *particle) {
281 // Convert Abs(pdg) to pid
282 // (0 - e, 1 - muons, 2 - pions, 3 - kaons, 4 - protons, 5 -all rest)
286 if (TMath::Abs(particle->GetPdgCode()) == kElectron) { pid = 0; }
287 else if (TMath::Abs(particle->GetPdgCode()) == kMuonMinus) { pid = 1; }
288 else if (TMath::Abs(particle->GetPdgCode()) == kPiPlus) { pid = 2; }
289 else if (TMath::Abs(particle->GetPdgCode()) == kKPlus) { pid = 3; }
290 else if (TMath::Abs(particle->GetPdgCode()) == kProton) { pid = 4; }
296 //_____________________________________________________________________________
297 TH1F* AlidNdPtHelper::CreateResHisto(TH2F* hRes2, TH1F **phMean, Int_t integ, Bool_t drawBinFits, Int_t minHistEntries)
300 // Create mean and resolution
303 TVirtualPad* currentPad = gPad;
304 TAxis* axis = hRes2->GetXaxis();
305 Int_t nBins = axis->GetNbins();
306 Bool_t overflowBinFits = kFALSE;
308 if (axis->GetXbins()->GetSize()){
309 hRes = new TH1F("hRes", "", nBins, axis->GetXbins()->GetArray());
310 hMean = new TH1F("hMean", "", nBins, axis->GetXbins()->GetArray());
313 hRes = new TH1F("hRes", "", nBins, axis->GetXmin(), axis->GetXmax());
314 hMean = new TH1F("hMean", "", nBins, axis->GetXmin(), axis->GetXmax());
317 hRes->SetStats(false);
318 hRes->SetOption("E");
319 hRes->SetMinimum(0.);
321 hMean->SetStats(false);
322 hMean->SetOption("E");
324 // create the fit function
325 TF1 * fitFunc = new TF1("G","[0]*exp(-(x-[1])*(x-[1])/(2.*[2]*[2]))",-3,3);
327 fitFunc->SetLineWidth(2);
328 fitFunc->SetFillStyle(0);
329 // create canvas for fits
330 TCanvas* canBinFits = NULL;
331 Int_t nPads = (overflowBinFits) ? nBins+2 : nBins;
332 Int_t nx = Int_t(sqrt(nPads-1.));// + 1;
333 Int_t ny = (nPads-1) / nx + 1;
335 canBinFits = (TCanvas*)gROOT->FindObject("canBinFits");
336 if (canBinFits) delete canBinFits;
337 canBinFits = new TCanvas("canBinFits", "fits of bins", 200, 100, 500, 700);
338 canBinFits->Divide(nx, ny);
341 // loop over x bins and fit projection
342 Int_t dBin = ((overflowBinFits) ? 1 : 0);
343 for (Int_t bin = 1-dBin; bin <= nBins+dBin; bin++) {
344 if (drawBinFits) canBinFits->cd(bin + dBin);
345 Int_t bin0=TMath::Max(bin-integ,0);
346 Int_t bin1=TMath::Min(bin+integ,nBins);
347 TH1D* hBin = hRes2->ProjectionY("hBin", bin0, bin1);
349 if (hBin->GetEntries() > minHistEntries) {
350 fitFunc->SetParameters(hBin->GetMaximum(),hBin->GetMean(),hBin->GetRMS());
351 hBin->Fit(fitFunc,"s");
352 Double_t sigma = TMath::Abs(fitFunc->GetParameter(2));
355 hRes->SetBinContent(bin, TMath::Abs(fitFunc->GetParameter(2)));
356 hMean->SetBinContent(bin, fitFunc->GetParameter(1));
359 hRes->SetBinContent(bin, 0.);
360 hMean->SetBinContent(bin,0);
362 hRes->SetBinError(bin, fitFunc->GetParError(2));
363 hMean->SetBinError(bin, fitFunc->GetParError(1));
369 hRes->SetBinContent(bin, 0.);
370 hRes->SetBinError(bin, 0.);
371 hMean->SetBinContent(bin, 0.);
372 hMean->SetBinError(bin, 0.);
379 sprintf(name, "%s < %.4g", axis->GetTitle(), axis->GetBinUpEdge(bin));
380 } else if (bin == nBins+1) {
381 sprintf(name, "%.4g < %s", axis->GetBinLowEdge(bin), axis->GetTitle());
383 sprintf(name, "%.4g < %s < %.4g", axis->GetBinLowEdge(bin),
384 axis->GetTitle(), axis->GetBinUpEdge(bin));
386 canBinFits->cd(bin + dBin);
387 hBin->SetTitle(name);
388 hBin->SetStats(kTRUE);
390 canBinFits->Update();
391 canBinFits->Modified();
392 canBinFits->Update();
404 //_____________________________________________________________________________
405 TH1F* AlidNdPtHelper::MakeResol(TH2F * his, Int_t integ, Bool_t type, Bool_t drawBins, Int_t minHistEntries){
406 // Create resolution histograms
408 TH1F *hisr=0, *hism=0;
409 if (!gPad) new TCanvas;
410 hisr = CreateResHisto(his,&hism,integ,drawBins,minHistEntries);
411 if (type) return hism;
417 //_____________________________________________________________________________
418 TObjArray* AlidNdPtHelper::GetAllChargedTracks(AliESDEvent *esdEvent, AnalysisMode analysisMode)
421 // all charged TPC particles
423 TObjArray *allTracks = new TObjArray();
424 if(!allTracks) return allTracks;
426 AliESDtrack *track=0;
427 for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++)
429 if(analysisMode == AlidNdPtHelper::kTPC || analysisMode == AlidNdPtHelper::kTPCSPDvtx ||
430 analysisMode == AlidNdPtHelper::kMCRec || analysisMode == kMCPion || analysisMode == kMCKaon ||
431 analysisMode == kMCProton || analysisMode ==kPlus || analysisMode ==kMinus) {
433 // track must be deleted by the user
434 track = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent,iTrack);
436 track=esdEvent->GetTrack(iTrack);
440 if(track->Charge()==0) {
441 if(analysisMode == AlidNdPtHelper::kTPC || analysisMode == AlidNdPtHelper::kTPCSPDvtx ||
442 analysisMode == AlidNdPtHelper::kMCRec || analysisMode == kMCPion || analysisMode == kMCKaon ||
443 analysisMode == kMCProton || analysisMode ==kPlus || analysisMode ==kMinus) {
445 delete track; continue;
451 allTracks->Add(track);
454 if(analysisMode == AlidNdPtHelper::kTPC || analysisMode == AlidNdPtHelper::kTPCSPDvtx ||
455 analysisMode == AlidNdPtHelper::kMCRec || analysisMode == kMCPion || analysisMode == kMCKaon ||
456 analysisMode == kMCProton || analysisMode ==kPlus || analysisMode ==kMinus) {
458 allTracks->SetOwner(kTRUE);
464 //_____________________________________________________________________________
465 Int_t AlidNdPtHelper::GetTPCMBTrackMult(AliESDEvent *esdEvent, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts)
468 // get MB event track multiplicity
472 ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","esd event is NULL");
476 if(!evtCuts || !accCuts || !trackCuts)
478 ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","cuts not available");
483 Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
484 Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
485 AliESDVertex vtx0(pos,err);
488 Float_t maxDCAr = accCuts->GetMaxDCAr();
489 Float_t maxDCAz = accCuts->GetMaxDCAz();
490 Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
492 Int_t ntracks = esdEvent->GetNumberOfTracks();
493 Double_t dca[2],cov[3];
495 for (Int_t i=0;i <ntracks; i++){
496 AliESDtrack *t = esdEvent->GetTrack(i);
498 if (t->Charge() == 0) continue;
499 if (!t->GetTPCInnerParam()) continue;
500 if (t->GetTPCNcls()<minTPCClust) continue;
502 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
503 if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov))
505 if(tpcTrack) delete tpcTrack;
509 if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
510 if(tpcTrack) delete tpcTrack;
516 if(tpcTrack) delete tpcTrack;
522 //_____________________________________________________________________________
523 Int_t AlidNdPtHelper::GetTPCMBPrimTrackMult(AliESDEvent *esdEvent, AliStack * stack, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts)
526 // get MB primary event track multiplicity
530 ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
536 ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
540 if(!evtCuts || !accCuts || !trackCuts)
542 ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","cuts not available");
547 Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
548 Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
549 AliESDVertex vtx0(pos,err);
552 Float_t maxDCAr = accCuts->GetMaxDCAr();
553 Float_t maxDCAz = accCuts->GetMaxDCAz();
554 Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
557 Int_t ntracks = esdEvent->GetNumberOfTracks();
558 Double_t dca[2],cov[3];
560 for (Int_t i=0;i <ntracks; i++){
561 AliESDtrack *t = esdEvent->GetTrack(i);
563 if (t->Charge() == 0) continue;
564 if (!t->GetTPCInnerParam()) continue;
565 if (t->GetTPCNcls()<minTPCClust) continue;
567 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
568 if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov))
570 if(tpcTrack) delete tpcTrack;
574 if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
575 if(tpcTrack) delete tpcTrack;
579 Int_t label = TMath::Abs(t->GetLabel());
580 TParticle *part = stack->Particle(label);
582 if(tpcTrack) delete tpcTrack;
585 if(!stack->IsPhysicalPrimary(label))
587 if(tpcTrack) delete tpcTrack;
593 if(tpcTrack) delete tpcTrack;
603 //_____________________________________________________________________________
604 Int_t AlidNdPtHelper::GetMCTrueTrackMult(AliMCEvent *mcEvent, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts)
607 // calculate mc event true track multiplicity
609 if(!mcEvent) return 0;
615 stack = mcEvent->Stack();
616 if (!stack) return 0;
618 Bool_t isEventOK = evtCuts->AcceptMCEvent(mcEvent);
619 if(!isEventOK) return 0;
621 Int_t nPart = stack->GetNtrack();
622 for (Int_t iMc = 0; iMc < nPart; ++iMc)
624 TParticle* particle = stack->Particle(iMc);
628 // only charged particles
629 Double_t charge = particle->GetPDG()->Charge()/3.;
634 Bool_t prim = stack->IsPhysicalPrimary(iMc);
638 if(accCuts->AcceptTrack(particle))
647 //_______________________________________________________________________
648 void AlidNdPtHelper::PrintMCInfo(AliStack *pStack,Int_t label)
650 // print information about particles in the stack
653 label = TMath::Abs(label);
654 TParticle *part = pStack->Particle(label);
655 Printf("########################");
656 Printf("%s:%d %d UniqueID %d PDG %d P %3.3f",(char*)__FILE__,__LINE__,label,part->GetUniqueID(),part->GetPdgCode(),part->P());
658 TParticle* mother = part;
659 Int_t imo = part->GetFirstMother();
660 Int_t nprim = pStack->GetNprimary();
662 while((imo >= nprim)) {
663 mother = pStack->Particle(imo);
664 Printf("Mother %s:%d Label %d UniqueID %d PDG %d P %3.3f",(char*)__FILE__,__LINE__,imo,mother->GetUniqueID(),mother->GetPdgCode(),mother->P());
666 imo = mother->GetFirstMother();
669 Printf("########################");
673 //_____________________________________________________________________________
674 TH1* AlidNdPtHelper::GetContCorrHisto(TH1 *hist)
677 // get contamination histogram
681 Int_t nbins = hist->GetNbinsX();
682 TH1 *h_cont = (TH1D *)hist->Clone();
684 for(Int_t i=0; i<=nbins+1; i++) {
685 Double_t binContent = hist->GetBinContent(i);
686 Double_t binError = hist->GetBinError(i);
688 h_cont->SetBinContent(i,1.-binContent);
689 h_cont->SetBinError(i,binError);
696 //_____________________________________________________________________________
697 TH1* AlidNdPtHelper::ScaleByBinWidth(TH1 *hist)
700 // scale by bin width
704 TH1 *h_scale = (TH1D *)hist->Clone();
705 h_scale->Scale(1.,"width");
710 //_____________________________________________________________________________
711 TH1* AlidNdPtHelper::CalcRelativeDifference(TH1 *hist1, TH1 *hist2)
714 // calculate rel. difference
720 TH1 *h1_clone = (TH1D *)hist1->Clone();
724 h1_clone->Add(hist2,-1);
725 h1_clone->Divide(hist2);
730 //_____________________________________________________________________________
731 TH1* AlidNdPtHelper::CalcRelativeDifferenceFun(TH1 *hist1, TF1 *fun)
734 // calculate rel. difference
735 // between histogram and function
740 TH1 *h1_clone = (TH1D *)hist1->Clone();
744 h1_clone->Add(fun,-1);
745 h1_clone->Divide(hist1);
750 //_____________________________________________________________________________
751 TH1* AlidNdPtHelper::NormalizeToEvent(TH2 *hist1, TH1 *hist2)
753 // normalise to event for a given multiplicity bin
754 // return pt histogram
760 Int_t nbinsX = hist1->GetNbinsX();
761 //Int_t nbinsY = hist1->GetNbinsY();
764 for(Int_t i=0; i<=nbinsX+1; i++) {
765 sprintf(name,"mom_%d",i);
766 TH1D *hist = (TH1D*)hist1->ProjectionY(name,i+1,i+1);
768 sprintf(name,"mom_norm");
770 hist_norm = (TH1D *)hist->Clone(name);
774 Double_t nbEvents = hist2->GetBinContent(i);
775 if(!nbEvents) { nbEvents = 1.; };
777 hist->Scale(1./nbEvents);
778 hist_norm->Add(hist);
784 //_____________________________________________________________________________
785 THnSparse* AlidNdPtHelper::GenerateCorrMatrix(THnSparse *hist1, THnSparse *hist2, char *name) {
786 // generate correction matrix
787 if(!hist1 || !hist2) return 0;
789 THnSparse *h =(THnSparse*)hist1->Clone(name);;
790 h->Divide(hist1,hist2,1,1,"B");
795 //_____________________________________________________________________________
796 TH2* AlidNdPtHelper::GenerateCorrMatrix(TH2 *hist1, TH2 *hist2, char *name) {
797 // generate correction matrix
798 if(!hist1 || !hist2) return 0;
800 TH2D *h =(TH2D*)hist1->Clone(name);;
801 h->Divide(hist1,hist2,1,1,"B");
806 //_____________________________________________________________________________
807 TH1* AlidNdPtHelper::GenerateCorrMatrix(TH1 *hist1, TH1 *hist2, char *name) {
808 // generate correction matrix
809 if(!hist1 || !hist2) return 0;
811 TH1D *h =(TH1D*)hist1->Clone(name);;
812 h->Divide(hist1,hist2,1,1,"B");
817 //_____________________________________________________________________________
818 THnSparse* AlidNdPtHelper::GenerateContCorrMatrix(THnSparse *hist1, THnSparse *hist2, char *name) {
819 // generate contamination correction matrix
820 if(!hist1 || !hist2) return 0;
822 THnSparse *hist = GenerateCorrMatrix(hist1, hist2, name);
825 // only for non ZERO bins!!!!
827 Int_t* coord = new Int_t[hist->GetNdimensions()];
828 memset(coord, 0, sizeof(Int_t) * hist->GetNdimensions());
830 for (Long64_t i = 0; i < hist->GetNbins(); ++i) {
831 Double_t v = hist->GetBinContent(i, coord);
832 hist->SetBinContent(coord, 1.0-v);
833 //printf("v %f, hist->GetBinContent(i, coord) %f \n",v,hist->GetBinContent(i, coord));
834 Double_t err = hist->GetBinError(coord);
835 hist->SetBinError(coord, err);
843 //_____________________________________________________________________________
844 TH2* AlidNdPtHelper::GenerateContCorrMatrix(TH2 *hist1, TH2 *hist2, char *name) {
845 // generate contamination correction matrix
846 if(!hist1 || !hist2) return 0;
848 TH2 *hist = GenerateCorrMatrix(hist1, hist2, name);
851 Int_t nBinsX = hist->GetNbinsX();
852 Int_t nBinsY = hist->GetNbinsY();
854 for (Int_t i = 0; i < nBinsX+1; i++) {
855 for (Int_t j = 0; j < nBinsY+1; j++) {
856 Double_t cont = hist->GetBinContent(i,j);
857 hist->SetBinContent(i,j,1.-cont);
858 Double_t err = hist->GetBinError(i,j);
859 hist->SetBinError(i,j,err);
866 //_____________________________________________________________________________
867 TH1* AlidNdPtHelper::GenerateContCorrMatrix(TH1 *hist1, TH1 *hist2, char *name) {
868 // generate contamination correction matrix
869 if(!hist1 || !hist2) return 0;
871 TH1 *hist = GenerateCorrMatrix(hist1, hist2, name);
874 Int_t nBinsX = hist->GetNbinsX();
876 for (Int_t i = 0; i < nBinsX+1; i++) {
877 Double_t cont = hist->GetBinContent(i);
878 hist->SetBinContent(i,1.-cont);
879 Double_t err = hist->GetBinError(i);
880 hist->SetBinError(i,err);
886 //_____________________________________________________________________________
887 const AliESDVertex* AlidNdPtHelper::GetTPCVertexZ(AliESDEvent* esdEvent, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts, Float_t fraction, Int_t ntracksMin){
893 ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
897 if(!evtCuts || !accCuts || !trackCuts)
899 ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
903 Double_t vtxpos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
904 Double_t vtxsigma[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
905 AliESDVertex vtx0(vtxpos,vtxsigma);
907 Double_t maxDCAr = accCuts->GetMaxDCAr();
908 Double_t maxDCAz = accCuts->GetMaxDCAz();
909 Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
912 Int_t ntracks = esdEvent->GetNumberOfTracks();
913 TVectorD ztrack(ntracks);
914 Double_t dca[2],cov[3];
916 for (Int_t i=0;i <ntracks; i++){
917 AliESDtrack *t = esdEvent->GetTrack(i);
919 if (!t->GetTPCInnerParam()) continue;
920 if (t->GetTPCNcls()<minTPCClust) continue;
922 AliExternalTrackParam *tpcTrack = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
923 if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov)) continue;
926 if (TMath::Abs(dca[0])>maxDCAr) continue;
927 //if (TMath::Sqrt(cov[0])>sigmaXYcut) continue;
928 if (TMath::Abs(tpcTrack->GetZ())>maxDCAz) continue;
930 ztrack[counter]=tpcTrack->GetZ();
933 if(tpcTrack) delete tpcTrack;
937 // Find LTM z position
939 Double_t mean=0, sigma=0;
940 if (counter<ntracksMin) return 0;
942 Int_t nused = TMath::Nint(counter*fraction);
943 if (nused==counter) nused=counter-1;
945 AliMathBase::EvaluateUni(counter, ztrack.GetMatrixArray(), mean,sigma, TMath::Nint(counter*fraction));
946 sigma/=TMath::Sqrt(nused);
948 mean = TMath::Mean(counter, ztrack.GetMatrixArray());
949 sigma = TMath::RMS(counter, ztrack.GetMatrixArray());
950 sigma/=TMath::Sqrt(counter-1);
954 const AliESDVertex* vertex = new AliESDVertex(vtxpos, vtxsigma);
958 //_____________________________________________________________________________
959 Int_t AlidNdPtHelper::GetSPDMBTrackMult(AliESDEvent* esdEvent, Float_t deltaThetaCut, Float_t deltaPhiCut)
962 // SPD track multiplicity
966 const AliMultiplicity* mult = esdEvent->GetMultiplicity();
970 // get multiplicity from SPD tracklets
971 Int_t inputCount = 0;
972 for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
974 //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
976 Float_t phi = mult->GetPhi(i);
978 phi += TMath::Pi() * 2;
979 Float_t deltaPhi = mult->GetDeltaPhi(i);
980 Float_t deltaTheta = mult->GetDeltaTheta(i);
982 if (TMath::Abs(deltaPhi) > 1)
983 printf("WARNING: Very high Delta Phi: %d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), deltaPhi);
985 if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
988 if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
997 //_____________________________________________________________________________
998 Int_t AlidNdPtHelper::GetSPDMBPrimTrackMult(AliESDEvent* esdEvent, AliStack* stack, Float_t deltaThetaCut, Float_t deltaPhiCut)
1001 // SPD track multiplicity
1005 const AliMultiplicity* mult = esdEvent->GetMultiplicity();
1009 // get multiplicity from SPD tracklets
1010 Int_t inputCount = 0;
1011 for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
1013 //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
1015 Float_t phi = mult->GetPhi(i);
1017 phi += TMath::Pi() * 2;
1018 Float_t deltaPhi = mult->GetDeltaPhi(i);
1019 Float_t deltaTheta = mult->GetDeltaTheta(i);
1021 if (TMath::Abs(deltaPhi) > 1)
1022 printf("WARNING: Very high Delta Phi: %d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), deltaPhi);
1024 if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
1027 if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
1031 if (mult->GetLabel(i, 0) < 0 || mult->GetLabel(i, 0) != mult->GetLabel(i, 1) ||
1032 !stack->IsPhysicalPrimary(mult->GetLabel(i, 0)))