]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/SPECTRA/ChargedHadrons/dNdPt/AlidNdPtHelper.cxx
88e3a2c1700e1cc21557c6c8302577b783f134b9
[u/mrichter/AliRoot.git] / PWGLF / SPECTRA / ChargedHadrons / dNdPt / AlidNdPtHelper.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15 //
16 // static dNdPt helper functions
17 //
18 // basic functionality to select events and tracks 
19 // for dNdPt analysis
20 //
21 // Origin: Jan Fiete Grosse-Oetringhaus
22 // Modified and Extended: Jacek Otwinowski 19/11/2009
23 //
24 // 
25
26 #include <TROOT.h>
27 #include <TCanvas.h>
28 #include <TF1.h>
29 #include <TH1.h>
30 #include <TH2.h>
31 #include <TH3.h>
32
33 #include <AliHeader.h>
34 #include <AliStack.h>
35 #include <AliLog.h>
36 #include <AliESD.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"
47
48 //____________________________________________________________________
49 ClassImp(AlidNdPtHelper)
50
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)
53 {
54   // Get the vertex from the ESD and returns it if the vertex is valid
55   //
56   // Second argument decides which vertex is used (this selects
57   // also the quality criteria that are applied)
58
59   if(!aEsd) 
60   { 
61     ::Error("AlidNdPtHelper::GetVertex()","esd event is NULL");
62     return NULL;  
63   }
64  
65   if(!evtCuts || !accCuts || !trackCuts) 
66   { 
67     ::Error("AlidNdPtHelper::GetVertex()","cuts not available");
68     return NULL;  
69   }
70
71   const AliESDVertex* vertex = 0;
72   AliESDVertex *initVertex = 0;
73   if (analysisMode == kSPD || 
74       analysisMode == kTPCSPDvtx || analysisMode == kTPCSPDvtxUpdate || analysisMode == kTPCITSHybrid)
75   {
76     vertex = aEsd->GetPrimaryVertexSPD();
77     if (debug)
78       Printf("AlidNdPtHelper::GetVertex: Returning SPD vertex");
79   }  
80   else if (analysisMode == kTPCITS  || analysisMode == kTPCTrackSPDvtx || analysisMode == kTPCTrackSPDvtxUpdate || 
81            analysisMode == kTPCITSHybridTrackSPDvtx || analysisMode == kTPCITSHybridTrackSPDvtxDCArPt || 
82            analysisMode == kITSStandAloneTrackSPDvtx ||  analysisMode == kITSStandAloneTPCTrackSPDvtx)
83   {
84     vertex = aEsd->GetPrimaryVertexTracks();
85     if(!vertex) return NULL;
86     if(vertex->GetNContributors()<1) {
87       // SPD vertex
88       vertex = aEsd->GetPrimaryVertexSPD();
89     }
90   }
91   else if (analysisMode == kTPC)
92   {
93     if(bRedoTPC) {
94
95       Double_t kBz = aEsd->GetMagneticField();
96       AliVertexerTracks vertexer(kBz);
97
98       if(bUseMeanVertex) {
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();
104       }
105
106       Double_t maxDCAr = accCuts->GetMaxDCAr();
107       Double_t maxDCAz = accCuts->GetMaxDCAz();
108       Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
109
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);
112
113       // TPC track preselection
114       Int_t ntracks = aEsd->GetNumberOfTracks();
115       TObjArray array(ntracks);
116       UShort_t *id = new UShort_t[ntracks];
117
118
119       Int_t count=0;
120       for (Int_t i=0;i <ntracks; i++) {
121         AliESDtrack *t = aEsd->GetTrack(i);
122         if (!t) continue;
123         if (t->Charge() == 0) continue;
124         if (!t->GetTPCInnerParam()) continue;
125         if (t->GetTPCNcls()<vertexer.GetMinClusters()) continue;
126         AliExternalTrackParam  *tpcTrack  = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
127         if(tpcTrack) { 
128           array.AddLast(tpcTrack);
129           id[count] = (UShort_t)t->GetID();
130           count++;
131         }
132       } 
133       AliESDVertex *vTPC = vertexer.VertexForSelectedTracks(&array,id, kTRUE, kTRUE, bUseMeanVertex);
134       if(!vTPC) { 
135         delete [] id; id=NULL;
136         return 0;
137       }
138       
139       // set recreated TPC vertex
140       aEsd->SetPrimaryVertexTPC(vTPC);
141
142       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); i++) {
143         AliESDtrack *t = aEsd->GetTrack(i);
144         if(!t) continue;
145
146         Double_t x[3]; t->GetXYZ(x);
147         Double_t b[3]; AliTracker::GetBxByBz(x,b);
148         t->RelateToVertexTPCBxByBz(vTPC, b, kVeryBig);
149       }
150       
151       delete vTPC;
152       array.Delete();
153       delete [] id; id=NULL;
154
155     }
156     vertex = aEsd->GetPrimaryVertexTPC();
157     if (debug)
158      Printf("AlidNdPtHelper::GetVertex: Returning vertex from tracks");
159    }
160    else
161      Printf("AlidNdPtHelper::GetVertex: ERROR: Invalid second argument %d", analysisMode);
162
163     if (!vertex) {
164      if (debug)
165       Printf("AlidNdPtHelper::GetVertex: No vertex found in ESD");
166       return 0;
167     }
168
169   if (debug)
170   {
171     Printf("AlidNdPtHelper::GetVertex: Returning valid vertex: %s", vertex->GetTitle());
172     vertex->Print();
173   }
174   
175   if(initVertex) delete initVertex; initVertex=NULL;
176   return vertex;
177 }
178
179 //____________________________________________________________________
180 Bool_t AlidNdPtHelper::TestRecVertex(const AliESDVertex* vertex, const AliESDVertex* vertexSPD, AnalysisMode analysisMode, Bool_t debug)
181 {
182   // Checks if a vertex meets the needed quality criteria
183   if(!vertex) return kFALSE;
184   if(!vertex->GetStatus()) return kFALSE;
185
186   Float_t requiredZResolution = -1;
187   if (analysisMode == kSPD || analysisMode == kTPCITS || 
188       analysisMode == kTPCSPDvtx || analysisMode == kTPCSPDvtxUpdate || analysisMode == kTPCITSHybrid ||
189       analysisMode == kTPCTrackSPDvtx || analysisMode == kTPCTrackSPDvtxUpdate || analysisMode == kTPCITSHybridTrackSPDvtx || analysisMode == kTPCITSHybridTrackSPDvtxDCArPt
190            || analysisMode == kITSStandAloneTrackSPDvtx ||  analysisMode == kITSStandAloneTPCTrackSPDvtx)
191   {
192     requiredZResolution = 1000;
193   }
194   else if (analysisMode == kTPC)
195     requiredZResolution = 10.;
196
197   // check resolution
198   Double_t zRes = vertex->GetZRes();
199
200   if (zRes > requiredZResolution) {
201     if (debug)
202       Printf("AlidNdPtHelper::TestVertex: Resolution too poor %f (required: %f", zRes, requiredZResolution);
203     return kFALSE;
204   }
205
206   // always check for SPD vertex
207   if(!vertexSPD) return kFALSE;
208   if(!vertexSPD->GetStatus()) return kFALSE;
209   if (vertexSPD->IsFromVertexerZ())
210   {
211     if (vertexSPD->GetDispersion() > 0.02) 
212     {
213       if (debug)
214         Printf("AliPWG0Helper::TestVertex: Delta Phi too large in Vertexer Z: %f (required: %f", vertex->GetDispersion(), 0.02);
215       return kFALSE;
216     }
217   }
218
219
220   return kTRUE;
221 }
222
223 //____________________________________________________________________
224 Bool_t AlidNdPtHelper::IsGoodImpPar(const AliESDtrack *const track)
225 {
226 //
227 // check whether particle has good DCAr(Pt) impact
228 // parameter. Only for TPC+ITS tracks (7*sigma cut)
229 // Origin: Andrea Dainese
230 //
231
232 Float_t d0z0[2],covd0z0[3];
233 track->GetImpactParameters(d0z0,covd0z0);
234 Float_t sigma= 0.0050+0.0060/TMath::Power(track->Pt(),0.9);
235 Float_t d0max = 7.*sigma;
236 if(TMath::Abs(d0z0[0]) < d0max) return kTRUE;
237
238 return kFALSE;
239 }
240
241 //____________________________________________________________________
242 Bool_t AlidNdPtHelper::IsPrimaryParticle(AliStack* const stack, Int_t idx, ParticleMode particleMode)
243 {
244 //
245 // check primary particles 
246 // depending on the particle mode
247 //
248   if(!stack) return kFALSE;
249
250   TParticle* particle = stack->Particle(idx);
251   if (!particle) return  kFALSE;
252
253   // only charged particles
254   Double_t charge = particle->GetPDG()->Charge()/3.;
255   if (TMath::Abs(charge) < 0.001) return kFALSE;
256
257   Int_t pdg = TMath::Abs(particle->GetPdgCode());
258
259   // physical primary
260   Bool_t prim = stack->IsPhysicalPrimary(idx);
261
262   if(particleMode==kMCPion) {
263     if(prim && pdg==kPiPlus) return kTRUE;
264     else return kFALSE;
265   } 
266
267   if (particleMode==kMCKaon) {
268     if(prim && pdg==kKPlus) return kTRUE;
269     else return kFALSE;
270   }
271     
272   if (particleMode==kMCProton) {
273     if(prim && pdg==kProton) return kTRUE;
274     else return kFALSE;
275   }
276
277   if(particleMode==kMCRest) {
278     if(prim && pdg!=kPiPlus && pdg!=kKPlus && pdg!=kProton) return kTRUE;
279     else return kFALSE;
280   }
281
282 return prim;
283 }
284
285 //____________________________________________________________________
286 Bool_t AlidNdPtHelper::IsCosmicTrack(AliESDtrack *const track1, AliESDtrack *const track2)
287 {
288 //
289 // check cosmic tracks
290 //
291   if(!track1) return kFALSE;
292   if(!track2) return kFALSE;
293
294   //
295   // cosmic tracks in TPC
296   //
297   //if( TMath::Abs( track1->Theta() - track2->Theta() ) < 0.004  && 
298   //  ((TMath::Abs(track1->Phi()-track2->Phi()) - TMath::Pi() )<0.004) && 
299   //if( track1->Pt() > 4.0 && (TMath::Abs(track1->Phi()-track2->Phi())-TMath::Pi())<0.1 
300   //    && (TMath::Abs(track1->Eta()+track2->Eta())-0.1) < 0.0 && (track1->Charge()+track2->Charge()) == 0)
301
302   //Float_t scaleF= 6.0;
303   if ( track1->Pt() > 4 && track2->Pt() > 4 && 
304        //(TMath::Abs(track1->GetSnp()-track2->GetSnp())-2.) < scaleF * TMath::Sqrt(track1->GetSigmaSnp2()+track2->GetSigmaSnp2()) &&
305        //TMath::Abs(track1->GetTgl()-track2->GetTgl())   < scaleF * TMath::Sqrt(track1->GetSigmaTgl2()+track2->GetSigmaTgl2()) &&
306        //TMath::Abs(track1->OneOverPt()-track2->OneOverPt()) < scaleF * TMath::Sqrt(track1->GetSigma1Pt2()+track2->GetSigma1Pt2()) && 
307        (track1->Charge()+track2->Charge()) == 0 && 
308        track1->Eta()*track2->Eta()<0.0 && TMath::Abs(track1->Eta()+track2->Eta())<0.03 &&
309        TMath::Abs(TMath::Abs(track1->Phi()-track2->Phi())-TMath::Pi())<0.1
310      )  
311   {
312     printf("COSMIC  candidate \n");
313
314     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());
315     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());
316     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()); 
317     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())); 
318     return kTRUE;
319   }
320      
321 return kFALSE; 
322 }
323
324 //____________________________________________________________________
325 void AlidNdPtHelper::PrintConf(AnalysisMode analysisMode, AliTriggerAnalysis::Trigger trigger)
326 {
327   //
328   // Prints the given configuration
329   //
330
331   TString str(">>>> Running with ");
332
333   switch (analysisMode)
334   {
335     case kInvalid: str += "invalid setting"; break;
336     case kSPD : str += "SPD-only"; break;
337     case kTPC : str += "TPC-only"; break;
338     case kTPCITS : str += "Global tracking"; break;
339     case kTPCSPDvtx : str += "TPC tracking + SPD event vertex"; break;
340     case kTPCSPDvtxUpdate : str += "TPC tracks updated with SPD event vertex point"; break;
341     case kTPCTrackSPDvtx : str += "TPC tracking + Tracks event vertex or SPD event vertex"; break;
342     case kTPCTrackSPDvtxUpdate : str += "TPC tracks updated with Track or SPD event vertex point"; break;
343     case kTPCITSHybrid : str += "TPC tracking + ITS refit + >1 SPD cluster"; break;
344     case kTPCITSHybridTrackSPDvtx : str += "TPC tracking + ITS refit + >1 SPD cluster + Tracks event vertex or SPD event vertex"; break;
345     case kTPCITSHybridTrackSPDvtxDCArPt : str += "TPC tracking + ITS refit + >1 SPD cluster + Tracks event vertex or SPD event vertex + DCAr(pt)"; break;
346     case kITSStandAloneTrackSPDvtx : str += "ITS standalone + Tracks event vertex or SPD event vertex + DCAr(pt)"; break;
347     case kITSStandAloneTPCTrackSPDvtx : str += "kITSStandAloneTPCTrackSPDvtx + TPC stand alone track  + Tracks event vertex or SPD event vertex + DCAr(pt)"; break;
348     case kMCRec : str += "TPC tracking + Replace rec. with MC values"; break;
349   }
350   str += " and trigger ";
351
352   str += AliTriggerAnalysis::GetTriggerName(trigger);
353
354   str += " <<<<";
355
356   Printf("%s", str.Data());
357 }
358
359 //____________________________________________________________________
360 Int_t AlidNdPtHelper::ConvertPdgToPid(const TParticle *const particle) {
361 //
362 // Convert Abs(pdg) to pid 
363 // (0 - e, 1 - muons, 2 - pions, 3 - kaons, 4 - protons, 5 -all rest)
364 //
365 Int_t pid=-1;
366
367   if (TMath::Abs(particle->GetPdgCode()) == kElectron)         { pid = 0; }
368   else if (TMath::Abs(particle->GetPdgCode()) == kMuonMinus) { pid = 1; }
369   else if (TMath::Abs(particle->GetPdgCode()) == kPiPlus)    { pid = 2; }
370   else if (TMath::Abs(particle->GetPdgCode()) == kKPlus)     { pid = 3; }
371   else if (TMath::Abs(particle->GetPdgCode()) == kProton)    { pid = 4; }
372   else                                                       { pid = 5; }
373
374 return pid;
375 }
376
377 //_____________________________________________________________________________
378 TH1F* AlidNdPtHelper::CreateResHisto(TH2F* const hRes2, TH1F **phMean, Int_t integ,  Bool_t drawBinFits, Int_t minHistEntries)
379 {
380 //
381 // Create mean and resolution 
382 // histograms
383 //
384   TVirtualPad* currentPad = gPad;
385   TAxis* axis = hRes2->GetXaxis();
386   Int_t nBins = axis->GetNbins();
387   //Bool_t overflowBinFits = kFALSE;
388   TH1F* hRes, *hMean;
389   if (axis->GetXbins()->GetSize()){
390     hRes = new TH1F("hRes", "", nBins, axis->GetXbins()->GetArray());
391     hMean = new TH1F("hMean", "", nBins, axis->GetXbins()->GetArray());
392   }
393   else{
394     hRes = new TH1F("hRes", "", nBins, axis->GetXmin(), axis->GetXmax());
395     hMean = new TH1F("hMean", "", nBins, axis->GetXmin(), axis->GetXmax());
396
397   }
398   hRes->SetStats(false);
399   hRes->SetOption("E");
400   hRes->SetMinimum(0.);
401   //
402   hMean->SetStats(false);
403   hMean->SetOption("E");
404  
405   // create the fit function
406   TF1 * fitFunc = new TF1("G","[0]*exp(-(x-[1])*(x-[1])/(2.*[2]*[2]))",-3,3);
407   
408   fitFunc->SetLineWidth(2);
409   fitFunc->SetFillStyle(0);
410   // create canvas for fits
411   TCanvas* canBinFits = NULL;
412   //Int_t nPads = (overflowBinFits) ? nBins+2 : nBins;
413   Int_t nPads = nBins;
414   Int_t nx = Int_t(sqrt(nPads-1.));// + 1;
415   Int_t ny = (nPads-1) / nx + 1;
416   if (drawBinFits) {
417     canBinFits = (TCanvas*)gROOT->FindObject("canBinFits");
418     if (canBinFits) delete canBinFits;
419     canBinFits = new TCanvas("canBinFits", "fits of bins", 200, 100, 500, 700);
420     canBinFits->Divide(nx, ny);
421   }
422
423   // loop over x bins and fit projection
424   //Int_t dBin = ((overflowBinFits) ? 1 : 0);
425   Int_t dBin = 0;
426   for (Int_t bin = 1-dBin; bin <= nBins+dBin; bin++) {
427     if (drawBinFits) canBinFits->cd(bin + dBin);
428     Int_t bin0=TMath::Max(bin-integ,0);
429     Int_t bin1=TMath::Min(bin+integ,nBins);
430     TH1D* hBin = hRes2->ProjectionY("hBin", bin0, bin1);
431     //    
432     if (hBin->GetEntries() > minHistEntries) {
433       fitFunc->SetParameters(hBin->GetMaximum(),hBin->GetMean(),hBin->GetRMS());
434       hBin->Fit(fitFunc,"s");
435       Double_t sigma = TMath::Abs(fitFunc->GetParameter(2));
436
437       if (sigma > 0.){
438         hRes->SetBinContent(bin, TMath::Abs(fitFunc->GetParameter(2)));
439         hMean->SetBinContent(bin, fitFunc->GetParameter(1));    
440       }
441       else{
442         hRes->SetBinContent(bin, 0.);
443         hMean->SetBinContent(bin,0);
444       }
445       hRes->SetBinError(bin, fitFunc->GetParError(2));
446       hMean->SetBinError(bin, fitFunc->GetParError(1));
447       
448       //
449       //
450
451     } else {
452       hRes->SetBinContent(bin, 0.);
453       hRes->SetBinError(bin, 0.);
454       hMean->SetBinContent(bin, 0.);
455       hMean->SetBinError(bin, 0.);
456     }
457     
458
459     if (drawBinFits) {
460       char name[256];
461       if (bin == 0) {
462         snprintf(name,256, "%s < %.4g", axis->GetTitle(), axis->GetBinUpEdge(bin));
463       } else if (bin == nBins+1) {
464         snprintf(name,256, "%.4g < %s", axis->GetBinLowEdge(bin), axis->GetTitle());
465       } else {
466         snprintf(name,256, "%.4g < %s < %.4g", axis->GetBinLowEdge(bin),
467                 axis->GetTitle(), axis->GetBinUpEdge(bin));
468       }
469       canBinFits->cd(bin + dBin);
470       hBin->SetTitle(name);
471       hBin->SetStats(kTRUE);
472       hBin->DrawCopy("E");
473       canBinFits->Update();
474       canBinFits->Modified();
475       canBinFits->Update();
476     }
477     
478     delete hBin;
479   }
480
481   delete fitFunc;
482   currentPad->cd();
483   *phMean = hMean;
484   return hRes;
485 }
486
487 //_____________________________________________________________________________
488 TH1F* AlidNdPtHelper::MakeResol(TH2F * his, Int_t integ, Bool_t type, Bool_t drawBins, Int_t minHistEntries){
489 // Create resolution histograms
490   
491      TH1F *hisr=0, *hism=0;
492      if (!gPad) new TCanvas;
493          hisr = CreateResHisto(his,&hism,integ,drawBins,minHistEntries);
494          if (type) return hism;
495          else return hisr;
496
497 return hisr;     
498 }
499
500 //_____________________________________________________________________________
501 TObjArray* AlidNdPtHelper::GetAllChargedTracks(AliESDEvent *esdEvent, AnalysisMode analysisMode)
502 {
503   //
504   // all charged TPC particles 
505   //
506   TObjArray *allTracks = new TObjArray();
507   if(!allTracks) return allTracks;
508
509   AliESDtrack *track=0;
510   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
511   { 
512     if(analysisMode == AlidNdPtHelper::kTPC) { 
513       //
514       // track must be deleted by user 
515       // esd track parameters are replaced by TPCinner
516       //
517       track = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent,iTrack);
518       if(!track) continue;
519     } 
520     else if (analysisMode == AlidNdPtHelper::kTPCSPDvtx || analysisMode == AlidNdPtHelper::kTPCSPDvtxUpdate)
521     {
522       //
523       // track must be deleted by the user 
524       // esd track parameters are replaced by TPCinner
525       //
526       track = AlidNdPtHelper::GetTPCOnlyTrackSPDvtx(esdEvent,iTrack,kFALSE);
527       if(!track) continue;
528     }
529     else if (analysisMode == AlidNdPtHelper::kTPCTrackSPDvtx || analysisMode == AlidNdPtHelper::kTPCTrackSPDvtxUpdate)
530     {
531       //
532       // track must be deleted by the user 
533       // esd track parameters are replaced by TPCinner
534       //
535       track = AlidNdPtHelper::GetTPCOnlyTrackTrackSPDvtx(esdEvent,iTrack,kFALSE);
536       if(!track) continue;
537     }
538     else if(analysisMode == AlidNdPtHelper::kTPCITSHybrid )
539     {
540       track = AlidNdPtHelper::GetTrackSPDvtx(esdEvent,iTrack,kFALSE);
541     }
542     else if(analysisMode == AlidNdPtHelper::kTPCITSHybridTrackSPDvtx || analysisMode == AlidNdPtHelper::kTPCITSHybridTrackSPDvtxDCArPt || analysisMode == AlidNdPtHelper::kITSStandAloneTrackSPDvtx || analysisMode ==AlidNdPtHelper::kITSStandAloneTPCTrackSPDvtx)
543     {
544       track = AlidNdPtHelper::GetTrackTrackSPDvtx(esdEvent,iTrack,kFALSE);
545     }
546     else 
547     {
548       track = esdEvent->GetTrack(iTrack);
549     }
550
551     if(!track) continue;
552
553     if(track->Charge()==0) { 
554       if(analysisMode == AlidNdPtHelper::kTPC || 
555          analysisMode == AlidNdPtHelper::kTPCSPDvtx || analysisMode == AlidNdPtHelper::kTPCTrackSPDvtx  ||
556          analysisMode == AlidNdPtHelper::kTPCSPDvtxUpdate || analysisMode == AlidNdPtHelper::kTPCTrackSPDvtxUpdate) 
557       {
558         delete track; continue; 
559       } else {
560         continue;
561       } 
562     }
563
564     allTracks->Add(track);
565   }
566
567   if(analysisMode == AlidNdPtHelper::kTPC || 
568      analysisMode == AlidNdPtHelper::kTPCSPDvtx || analysisMode == AlidNdPtHelper::kTPCTrackSPDvtx || 
569      analysisMode == AlidNdPtHelper::kTPCSPDvtxUpdate || analysisMode == AlidNdPtHelper::kTPCTrackSPDvtxUpdate) {
570      
571      allTracks->SetOwner(kTRUE);
572   }
573
574 return allTracks;
575 }
576
577 //_____________________________________________________________________________
578 AliESDtrack *AlidNdPtHelper::GetTPCOnlyTrackSPDvtx(const AliESDEvent* esdEvent, Int_t iTrack, Bool_t bUpdate)
579 {
580 //
581 // Create ESD tracks from TPCinner parameters.
582 // Propagte to DCA to SPD vertex.
583 // Update using SPD vertex point (parameter)
584 //
585 // It is user responsibility to delete these tracks
586 //
587
588   if (!esdEvent) return NULL;
589   if (!esdEvent->GetPrimaryVertexSPD() ) { return NULL; }
590   if (!esdEvent->GetPrimaryVertexSPD()->GetStatus() ) { return  NULL; }
591    
592   // 
593   AliESDtrack* track = esdEvent->GetTrack(iTrack);
594   if (!track)
595     return NULL;
596
597   Bool_t isOK = kFALSE;
598   Double_t x[3]; track->GetXYZ(x);
599   Double_t b[3]; AliTracker::GetBxByBz(x,b);
600
601   // create new ESD track
602   AliESDtrack *tpcTrack = new AliESDtrack();
603  
604   // relate TPC-only tracks (TPCinner) to SPD vertex
605   AliExternalTrackParam cParam;
606   if(bUpdate) {  
607     isOK = track->RelateToVertexTPCBxByBz(esdEvent->GetPrimaryVertexSPD(),b,kVeryBig,&cParam);
608     track->Set(cParam.GetX(),cParam.GetAlpha(),cParam.GetParameter(),cParam.GetCovariance());
609
610     // reject fake tracks
611     if(track->Pt() > 10000.)  {
612       ::Error("Exclude no physical tracks","pt>10000. GeV");
613       delete tpcTrack; 
614       return NULL;
615     }
616   }
617   else {
618     isOK = track->RelateToVertexTPCBxByBz(esdEvent->GetPrimaryVertexSPD(), b, kVeryBig);
619   }
620
621   // only true if we have a tpc track
622   if (!track->FillTPCOnlyTrack(*tpcTrack))
623   {
624     delete tpcTrack;
625     return NULL;
626   }
627   
628   if(!isOK) return NULL;
629
630 return tpcTrack;
631
632
633 //_____________________________________________________________________________
634 AliESDtrack *AlidNdPtHelper::GetTPCOnlyTrackTrackSPDvtx(const AliESDEvent* esdEvent, Int_t iTrack, Bool_t bUpdate)
635 {
636 //
637 // Create ESD tracks from TPCinner parameters.
638 // Propagte to DCA to Track or SPD vertex.
639 // Update using SPD vertex point (parameter)
640 //
641 // It is user responsibility to delete these tracks
642 //
643   if (!esdEvent) return NULL;
644   const AliESDVertex *vertex = esdEvent->GetPrimaryVertexTracks();
645   if(vertex->GetNContributors()<1) {
646     // SPD vertex
647     vertex = esdEvent->GetPrimaryVertexSPD();
648   }
649   if(!vertex) return NULL;
650  
651   // 
652   AliESDtrack* track = esdEvent->GetTrack(iTrack);
653   if (!track)
654     return NULL;
655
656   Bool_t isOK = kFALSE;
657   Double_t x[3]; track->GetXYZ(x);
658   Double_t b[3]; AliTracker::GetBxByBz(x,b);
659
660   // create new ESD track
661   AliESDtrack *tpcTrack = new AliESDtrack();
662  
663   // relate TPC-only tracks (TPCinner) to SPD vertex
664   AliExternalTrackParam cParam;
665   if(bUpdate) {  
666     isOK = track->RelateToVertexTPCBxByBz(vertex,b,kVeryBig,&cParam);
667     track->Set(cParam.GetX(),cParam.GetAlpha(),cParam.GetParameter(),cParam.GetCovariance());
668
669     // reject fake tracks
670     if(track->Pt() > 10000.)  {
671       ::Error("Exclude no physical tracks","pt>10000. GeV");
672       delete tpcTrack; 
673       return NULL;
674     }
675   }
676   else {
677     isOK = track->RelateToVertexTPCBxByBz(vertex, b, kVeryBig);
678   }
679
680   // only true if we have a tpc track
681   if (!track->FillTPCOnlyTrack(*tpcTrack))
682   {
683     delete tpcTrack;
684     return NULL;
685   }
686   
687   if(!isOK) return NULL;
688
689 return tpcTrack;
690
691
692 //_____________________________________________________________________________
693 AliESDtrack *AlidNdPtHelper::GetTrackSPDvtx(const AliESDEvent* esdEvent, Int_t iTrack, Bool_t bUpdate)
694 {
695 //
696 // Propagte track to DCA to SPD vertex.
697 // Update using SPD vertex point (parameter)
698 //
699   if (!esdEvent) return NULL;
700   if (!esdEvent->GetPrimaryVertexSPD() ) { return NULL; }
701   if (!esdEvent->GetPrimaryVertexSPD()->GetStatus() ) { return  NULL; }
702    
703   // 
704   AliESDtrack* track = esdEvent->GetTrack(iTrack);
705   if (!track)
706     return NULL;
707
708   Bool_t isOK = kFALSE;
709   Double_t x[3]; track->GetXYZ(x);
710   Double_t b[3]; AliTracker::GetBxByBz(x,b);
711
712   // relate tracks to SPD vertex
713   AliExternalTrackParam cParam;
714   if(bUpdate) {  
715     isOK = track->RelateToVertexBxByBz(esdEvent->GetPrimaryVertexSPD(),b,kVeryBig,&cParam);
716     track->Set(cParam.GetX(),cParam.GetAlpha(),cParam.GetParameter(),cParam.GetCovariance());
717
718     // reject fake tracks
719     if(track->Pt() > 10000.)  {
720       ::Error("Exclude no physical tracks","pt>10000. GeV");
721       return NULL;
722     }
723   }
724   else {
725     isOK = track->RelateToVertexBxByBz(esdEvent->GetPrimaryVertexSPD(), b, kVeryBig);
726   }
727  
728   if(!isOK) return NULL;
729
730 return track;
731
732
733 //_____________________________________________________________________________
734 AliESDtrack *AlidNdPtHelper::GetTrackTrackSPDvtx(const AliESDEvent* esdEvent, Int_t iTrack, Bool_t bUpdate)
735 {
736 //
737 // Propagte track to DCA to Track or SPD vertex.
738 // Update using SPD vertex point (parameter)
739 //
740   if (!esdEvent) return NULL;
741
742   const AliESDVertex *vertex = esdEvent->GetPrimaryVertexTracks();
743   if(vertex->GetNContributors()<1) {
744     // SPD vertex
745     vertex = esdEvent->GetPrimaryVertexSPD();
746   }
747   if(!vertex) return NULL;
748
749   // 
750   AliESDtrack* track = esdEvent->GetTrack(iTrack);
751   if (!track)
752     return NULL;
753
754   Bool_t isOK = kFALSE;
755   Double_t x[3]; track->GetXYZ(x);
756   Double_t b[3]; AliTracker::GetBxByBz(x,b);
757
758   // relate tracks to SPD vertex
759   AliExternalTrackParam cParam;
760   if(bUpdate) {  
761     isOK = track->RelateToVertexBxByBz(vertex,b,kVeryBig,&cParam);
762     track->Set(cParam.GetX(),cParam.GetAlpha(),cParam.GetParameter(),cParam.GetCovariance());
763
764     // reject fake tracks
765     if(track->Pt() > 10000.)  {
766       ::Error("Exclude no physical tracks","pt>10000. GeV");
767       return NULL;
768     }
769   }
770   else {
771     isOK = track->RelateToVertexBxByBz(vertex, b, kVeryBig);
772   }
773  
774   if(!isOK) return NULL;
775
776 return track;
777
778
779 //_____________________________________________________________________________
780 Bool_t AlidNdPtHelper::SelectEvent(const AliESDEvent* const esdEvent, AliESDtrackCuts* const esdTrackCuts) {
781 // select events with at least
782 // one reconstructed primary track in acceptance
783 // pT>0.5 GeV/c, |eta|<0.8 for cross section studies
784
785 if(!esdEvent) return kFALSE;
786 if(!esdTrackCuts) return kFALSE;
787
788   AliESDtrack *track=0;
789   Int_t count = 0;
790   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
791   { 
792     track = esdEvent->GetTrack(iTrack);
793     if(!track) continue;
794     if(track->Charge()==0) continue;
795     if(!esdTrackCuts->AcceptTrack(track)) continue;
796     if(track->Pt() < 0.5) continue;
797     if(TMath::Abs(track->Eta()) > 0.8) continue;
798
799     count++;
800   }
801
802   if(count > 0) return kTRUE;
803   else return kFALSE;
804
805 return kFALSE;
806 }
807
808 //_____________________________________________________________________________
809 Bool_t AlidNdPtHelper::SelectMCEvent(AliMCEvent* const mcEvent) {
810 //
811 // select events with at least
812 // one prompt (MC primary) track in acceptance
813 // pT>0.5 GeV/c, |eta|<0.8 for cross section studies
814 //
815
816 if(!mcEvent) return kFALSE;
817 AliStack* stack = mcEvent->Stack(); 
818 if(!stack) return kFALSE;
819
820   Int_t count = 0;
821   for (Int_t iMc = 0; iMc < stack->GetNtrack(); ++iMc) 
822   {
823     TParticle* particle = stack->Particle(iMc);
824     if (!particle) continue;
825
826     // only charged particles
827     if(!particle->GetPDG()) continue;
828     Double_t charge = particle->GetPDG()->Charge()/3.;
829     if(charge == 0) continue;
830
831     // physical primary
832     Bool_t prim = stack->IsPhysicalPrimary(iMc);
833     if(!prim) continue;
834
835     if(particle->Pt() < 0.5) continue;
836     if(TMath::Abs(particle->Eta()) > 0.8) continue;
837
838     count++;
839   }
840
841   if(count > 0) return kTRUE;
842   else return kFALSE;
843
844 return kFALSE;
845 }
846
847 //_____________________________________________________________________________
848 Int_t AlidNdPtHelper::GetTPCMBTrackMult(const AliESDEvent *const esdEvent,const AlidNdPtEventCuts *const evtCuts, const AlidNdPtAcceptanceCuts *const accCuts,const  AliESDtrackCuts *const trackCuts)
849 {
850   //
851   // get MB event track multiplicity
852   //
853   if(!esdEvent) 
854   { 
855     ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","esd event is NULL");
856     return 0;  
857   }
858  
859   if(!evtCuts || !accCuts || !trackCuts) 
860   { 
861     ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","cuts not available");
862     return 0;  
863   }
864
865   //
866   Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
867   Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
868   AliESDVertex vtx0(pos,err);
869
870   //
871   Float_t maxDCAr = accCuts->GetMaxDCAr();
872   Float_t maxDCAz = accCuts->GetMaxDCAz();
873   Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
874   //
875   Int_t ntracks = esdEvent->GetNumberOfTracks();
876   Double_t dca[2],cov[3];
877   Int_t mult=0;
878   for (Int_t i=0;i <ntracks; i++){
879     AliESDtrack *t = esdEvent->GetTrack(i);
880     if (!t) continue;
881     if (t->Charge() == 0) continue;
882     if (!t->GetTPCInnerParam()) continue;
883     if (t->GetTPCNcls()<minTPCClust) continue;
884     //
885     Double_t x[3]; t->GetXYZ(x);
886     Double_t b[3]; AliTracker::GetBxByBz(x,b);
887     const Double_t kMaxStep = 1;   //max step over the material
888
889     AliExternalTrackParam  *tpcTrack  = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
890     if(!tpcTrack) return 0;
891
892     if (!tpcTrack->PropagateToDCABxByBz(&vtx0,b,kMaxStep,dca,cov)) 
893     {
894       if(tpcTrack) delete tpcTrack; 
895       continue;
896     }
897     //
898     if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
899       if(tpcTrack) delete tpcTrack; 
900       continue;
901     }
902
903     mult++;    
904
905     if(tpcTrack) delete tpcTrack; 
906   }
907
908 return mult;  
909 }
910
911 //_____________________________________________________________________________
912 Int_t AlidNdPtHelper::GetTPCMBPrimTrackMult(const AliESDEvent *const esdEvent, AliStack *const  stack, const AlidNdPtEventCuts *const evtCuts, const AlidNdPtAcceptanceCuts *const accCuts, const AliESDtrackCuts *const trackCuts)
913 {
914   //
915   // get MB primary event track multiplicity
916   //
917   if(!esdEvent) 
918   { 
919     ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
920     return 0;  
921   }
922
923   if(!stack) 
924   { 
925     ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
926     return 0;  
927   }
928  
929   if(!evtCuts || !accCuts || !trackCuts) 
930   { 
931     ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","cuts not available");
932     return 0;  
933   }
934
935   //
936   Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
937   Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
938   AliESDVertex vtx0(pos,err);
939
940   //
941   Float_t maxDCAr = accCuts->GetMaxDCAr();
942   Float_t maxDCAz = accCuts->GetMaxDCAz();
943   Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
944
945   //
946   Int_t ntracks = esdEvent->GetNumberOfTracks();
947   Double_t dca[2],cov[3];
948   Int_t mult=0;
949   for (Int_t i=0;i <ntracks; i++){
950     AliESDtrack *t = esdEvent->GetTrack(i);
951     if (!t) continue;
952     if (t->Charge() == 0) continue;
953     if (!t->GetTPCInnerParam()) continue;
954     if (t->GetTPCNcls()<minTPCClust) continue;
955     //
956     Double_t x[3]; t->GetXYZ(x);
957     Double_t b[3]; AliTracker::GetBxByBz(x,b);
958     const Double_t kMaxStep = 1;   //max step over the material
959
960     AliExternalTrackParam  *tpcTrack  = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
961     if(!tpcTrack) return 0;
962
963     if (!tpcTrack->PropagateToDCABxByBz(&vtx0,b,kMaxStep,dca,cov)) 
964     {
965       if(tpcTrack) delete tpcTrack; 
966       continue;
967     }
968     //
969     if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
970       if(tpcTrack) delete tpcTrack; 
971       continue;
972     }
973
974     Int_t label = TMath::Abs(t->GetLabel());
975     TParticle *part = stack->Particle(label);
976     if(!part) { 
977       if(tpcTrack) delete tpcTrack; 
978       continue;
979     }
980     if(!stack->IsPhysicalPrimary(label)) 
981     { 
982       if(tpcTrack) delete tpcTrack; 
983       continue;
984     }
985
986     mult++;    
987
988     if(tpcTrack) delete tpcTrack; 
989   }
990
991 return mult;  
992 }
993
994
995 //_____________________________________________________________________________
996 Double_t AlidNdPtHelper::GetStrangenessCorrFactor(const Double_t pt)
997 {
998 // data driven correction factor for secondaries
999 // underestimated secondaries with strangeness in Pythia (A. Dainese)
1000
1001 //
1002 // pt=0.17; fact=1
1003 // pt=0.4; fact=1.07
1004 // pt=0.6; fact=1.25
1005 // pt>=1.2; fact=1.5
1006 //
1007
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);
1012 return 1.5;
1013
1014 }
1015
1016 //_____________________________________________________________________________
1017 Double_t AlidNdPtHelper::GetStrangenessCorrFactorPbPb(const Double_t pt)
1018 {
1019 // data driven correction factor for secondaries (PbPb)
1020
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);
1026 return 1.67;
1027
1028 }
1029
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)
1032 {
1033 //
1034 // linear interpolation
1035 //
1036   return ((y2-y1)/(x2-x1))*pt+(y2-(((y2-y1)/(x2-x1))*x2)); 
1037 }
1038
1039 //_____________________________________________________________________________
1040 Int_t AlidNdPtHelper::GetMCTrueTrackMult(AliMCEvent *const mcEvent, AlidNdPtEventCuts *const evtCuts, AlidNdPtAcceptanceCuts *const accCuts)
1041 {
1042   //
1043   // calculate mc event true track multiplicity
1044   //
1045   if(!mcEvent) return 0;
1046
1047   AliStack* stack = 0;
1048   Int_t mult = 0;
1049
1050   // MC particle stack
1051   stack = mcEvent->Stack();
1052   if (!stack) return 0;
1053
1054   //
1055   //printf("minZv %f, maxZv %f \n", evtCuts->GetMinZv(), evtCuts->GetMaxZv());
1056   //
1057
1058   Bool_t isEventOK = evtCuts->AcceptMCEvent(mcEvent);
1059   if(!isEventOK) return 0; 
1060
1061   Int_t nPart  = stack->GetNtrack();
1062   for (Int_t iMc = 0; iMc < nPart; ++iMc) 
1063   {
1064      TParticle* particle = stack->Particle(iMc);
1065      if (!particle)
1066      continue;
1067
1068      // only charged particles
1069      if(!particle->GetPDG()) continue;
1070      Double_t charge = particle->GetPDG()->Charge()/3.;
1071      if (TMath::Abs(charge) < 0.001)
1072      continue;
1073       
1074      // physical primary
1075      Bool_t prim = stack->IsPhysicalPrimary(iMc);
1076      if(!prim) continue;
1077
1078      // checked accepted without pt cut
1079      //if(accCuts->AcceptTrack(particle)) 
1080      if( particle->Eta() > accCuts->GetMinEta() && particle->Eta() < accCuts->GetMaxEta() ) 
1081      {
1082        mult++;
1083      }
1084   }
1085
1086 return mult;  
1087 }
1088
1089 //_______________________________________________________________________
1090 void  AlidNdPtHelper::PrintMCInfo(AliStack *const pStack,Int_t label)
1091 {
1092 // print information about particles in the stack
1093
1094   if(!pStack)return;
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());
1099   part->Print();
1100   TParticle* mother = part;
1101   Int_t imo = part->GetFirstMother();
1102   Int_t nprim = pStack->GetNprimary();
1103
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());
1107       mother->Print();
1108       imo =  mother->GetFirstMother();
1109  }
1110
1111  Printf("########################");
1112 }
1113
1114
1115 //_____________________________________________________________________________
1116 TH1* AlidNdPtHelper::GetContCorrHisto(TH1 *const hist) 
1117 {
1118 //
1119 // get contamination histogram
1120 //
1121  if(!hist) return 0;
1122
1123  Int_t nbins = hist->GetNbinsX();
1124  TH1 *hCont = (TH1D *)hist->Clone();
1125
1126  for(Int_t i=0; i<=nbins+1; i++) {
1127    Double_t binContent = hist->GetBinContent(i);
1128    Double_t binError = hist->GetBinError(i);
1129
1130    hCont->SetBinContent(i,1.-binContent);
1131    hCont->SetBinError(i,binError);
1132  }
1133
1134 return hCont;
1135 }
1136
1137
1138 //_____________________________________________________________________________
1139 TH1* AlidNdPtHelper::ScaleByBinWidth(TH1 *const hist) 
1140 {
1141 //
1142 // scale by bin width
1143 //
1144  if(!hist) return 0;
1145
1146  TH1 *hScale = (TH1D *)hist->Clone();
1147  hScale->Scale(1.,"width");
1148
1149 return hScale;
1150 }
1151
1152 //_____________________________________________________________________________
1153 TH1* AlidNdPtHelper::CalcRelativeDifference(const TH1 *const hist1, const TH1 *const hist2) 
1154 {
1155 //
1156 // calculate rel. difference 
1157 //
1158
1159  if(!hist1) return 0;
1160  if(!hist2) return 0;
1161
1162  TH1 *h1Clone = (TH1D *)hist1->Clone();
1163  h1Clone->Sumw2();
1164
1165  // (rec-mc)/mc
1166  h1Clone->Add(hist2,-1);
1167  h1Clone->Divide(hist2);
1168
1169 return h1Clone;
1170 }
1171
1172 //_____________________________________________________________________________
1173 TH1* AlidNdPtHelper::CalcRelativeDifferenceFun(const TH1 *const hist1, TF1 *const fun) 
1174 {
1175 //
1176 // calculate rel. difference
1177 // between histogram and function
1178 //
1179  if(!hist1) return 0;
1180  if(!fun) return 0;
1181
1182  TH1 *h1Clone = (TH1D *)hist1->Clone();
1183  h1Clone->Sumw2();
1184
1185  // 
1186  h1Clone->Add(fun,-1);
1187  h1Clone->Divide(hist1);
1188
1189 return h1Clone;
1190 }
1191
1192 //_____________________________________________________________________________
1193 TH1* AlidNdPtHelper::NormalizeToEvent(const TH2 *const hist1, const TH1 *const hist2) 
1194 {
1195 // normalise to event for a given multiplicity bin
1196 // return pt histogram 
1197
1198  if(!hist1) return 0;
1199  if(!hist2) return 0;
1200  char name[256];
1201
1202  Int_t nbinsX = hist1->GetNbinsX();
1203  //Int_t nbinsY = hist1->GetNbinsY();
1204
1205  TH1D *histNorm = 0;
1206  for(Int_t i=0; i<=nbinsX+1; i++) {
1207    snprintf(name,256,"mom_%d",i);
1208    TH1D *hist = (TH1D*)hist1->ProjectionY(name,i+1,i+1);
1209
1210    snprintf(name,256,"mom_norm");
1211    if(i==0) { 
1212      histNorm = (TH1D *)hist->Clone(name);
1213      histNorm->Reset();
1214    }
1215
1216    Double_t nbEvents = hist2->GetBinContent(i);
1217    if(!nbEvents) { nbEvents = 1.; };
1218
1219    hist->Scale(1./nbEvents);
1220    histNorm->Add(hist);
1221  }
1222
1223 return histNorm;
1224 }
1225
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; 
1231
1232 THnSparse *h =(THnSparse*)hist1->Clone(name);;
1233 h->Divide(hist1,hist2,1,1,"B");
1234
1235 return h;
1236 }
1237
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; 
1243
1244 TH2D *h =(TH2D*)hist1->Clone(name);;
1245 h->Divide(hist1,hist2,1,1,"B");
1246
1247 return h;
1248 }
1249
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; 
1255
1256 TH1D *h =(TH1D*)hist1->Clone(name);;
1257 h->Divide(hist1,hist2,1,1,"B");
1258
1259 return h;
1260 }
1261
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; 
1267
1268 THnSparse *hist =  GenerateCorrMatrix(hist1, hist2, name);
1269 if(!hist) return 0;
1270
1271 // only for non ZERO bins!!!!
1272
1273 Int_t* coord = new Int_t[hist->GetNdimensions()];
1274 memset(coord, 0, sizeof(Int_t) * hist->GetNdimensions());
1275
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);
1282   }
1283
1284 delete [] coord;
1285
1286 return hist;
1287 }
1288
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; 
1294
1295 TH2 *hist = GenerateCorrMatrix(hist1, hist2, name);
1296 if(!hist) return 0;
1297
1298 Int_t nBinsX = hist->GetNbinsX();
1299 Int_t nBinsY = hist->GetNbinsY();
1300
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);
1307   }
1308   }
1309
1310 return hist;
1311 }
1312
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; 
1318
1319 TH1 *hist = GenerateCorrMatrix(hist1, hist2, name);
1320 if(!hist) return 0;
1321
1322 Int_t nBinsX = hist->GetNbinsX();
1323
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);
1329   }
1330
1331 return hist;
1332 }
1333
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){
1336   //
1337   // TPC Z vertexer
1338   //
1339   if(!esdEvent)
1340   { 
1341     ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
1342     return NULL;  
1343   }
1344
1345   if(!evtCuts || !accCuts || !trackCuts) 
1346   { 
1347     ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
1348     return NULL;  
1349   }
1350
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);
1354
1355   Double_t maxDCAr = accCuts->GetMaxDCAr();
1356   Double_t maxDCAz = accCuts->GetMaxDCAz();
1357   Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
1358
1359   //
1360   Int_t ntracks = esdEvent->GetNumberOfTracks();
1361   TVectorD ztrack(ntracks);
1362   Double_t dca[2],cov[3];
1363   Int_t counter=0;
1364   for (Int_t i=0;i <ntracks; i++){
1365     AliESDtrack *t = esdEvent->GetTrack(i);
1366     if (!t) continue;
1367     if (!t->GetTPCInnerParam()) continue;
1368     if (t->GetTPCNcls()<minTPCClust) continue;
1369     //
1370
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
1374
1375     AliExternalTrackParam  *tpcTrack  = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
1376     if(!tpcTrack) return 0;
1377     if (!tpcTrack->PropagateToDCABxByBz(&vtx0,b,kMaxStep,dca,cov)) continue;
1378
1379     //
1380     if (TMath::Abs(dca[0])>maxDCAr) continue;
1381     //if (TMath::Sqrt(cov[0])>sigmaXYcut) continue;    
1382     if (TMath::Abs(tpcTrack->GetZ())>maxDCAz) continue;
1383
1384     ztrack[counter]=tpcTrack->GetZ();
1385     counter++;    
1386
1387     if(tpcTrack) delete tpcTrack;
1388   }
1389
1390   //
1391   // Find LTM z position
1392   //
1393   Double_t mean=0, sigma=0;
1394   if (counter<ntracksMin) return 0;
1395   //
1396   Int_t nused = TMath::Nint(counter*fraction);
1397   if (nused==counter) nused=counter-1;  
1398   if (nused>1){
1399     AliMathBase::EvaluateUni(counter, ztrack.GetMatrixArray(), mean,sigma, TMath::Nint(counter*fraction));
1400     sigma/=TMath::Sqrt(nused);
1401   }else{
1402     mean  = TMath::Mean(counter, ztrack.GetMatrixArray());
1403     sigma = TMath::RMS(counter, ztrack.GetMatrixArray());
1404     sigma/=TMath::Sqrt(counter-1);
1405   }
1406   vtxpos[2]=mean;
1407   vtxsigma[2]=sigma;
1408   const AliESDVertex* vertex = new AliESDVertex(vtxpos, vtxsigma);
1409   return vertex;
1410 }
1411
1412 //_____________________________________________________________________________
1413 Int_t  AlidNdPtHelper::GetSPDMBTrackMult(const AliESDEvent* const esdEvent, Float_t deltaThetaCut, Float_t deltaPhiCut) 
1414 {
1415   //
1416   // SPD track multiplicity
1417   //
1418
1419   // get tracklets
1420   const AliMultiplicity* mult = esdEvent->GetMultiplicity();
1421   if (!mult)
1422      return 0;
1423
1424   // get multiplicity from SPD tracklets
1425   Int_t inputCount = 0; 
1426   for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
1427   {
1428     //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
1429
1430      Float_t phi = mult->GetPhi(i);
1431      if (phi < 0)
1432        phi += TMath::Pi() * 2;
1433      Float_t deltaPhi = mult->GetDeltaPhi(i);
1434      Float_t deltaTheta = mult->GetDeltaTheta(i);
1435
1436      if (TMath::Abs(deltaPhi) > 1)
1437        printf("WARNING: Very high Delta Phi: %d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), deltaPhi);
1438
1439      if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
1440         continue;
1441
1442      if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
1443         continue;
1444       
1445      ++inputCount;
1446   }
1447
1448 return inputCount;
1449 }
1450
1451 //_____________________________________________________________________________
1452 Int_t  AlidNdPtHelper::GetSPDMBPrimTrackMult(const AliESDEvent* const esdEvent, AliStack* const stack, Float_t deltaThetaCut, Float_t deltaPhiCut) 
1453 {
1454   //
1455   // SPD track multiplicity
1456   //
1457
1458   // get tracklets
1459   const AliMultiplicity* mult = esdEvent->GetMultiplicity();
1460   if (!mult)
1461      return 0;
1462
1463   // get multiplicity from SPD tracklets
1464   Int_t inputCount = 0; 
1465   for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
1466   {
1467     //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
1468
1469      Float_t phi = mult->GetPhi(i);
1470      if (phi < 0)
1471        phi += TMath::Pi() * 2;
1472      Float_t deltaPhi = mult->GetDeltaPhi(i);
1473      Float_t deltaTheta = mult->GetDeltaTheta(i);
1474
1475      if (TMath::Abs(deltaPhi) > 1)
1476        printf("WARNING: Very high Delta Phi: %d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), deltaPhi);
1477
1478      if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
1479         continue;
1480
1481      if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
1482         continue;
1483
1484
1485      if (mult->GetLabel(i, 0) < 0 || mult->GetLabel(i, 0) != mult->GetLabel(i, 1) || 
1486          !stack->IsPhysicalPrimary(mult->GetLabel(i, 0)))
1487         continue;
1488
1489       
1490      ++inputCount;
1491   }
1492
1493 return inputCount;
1494 }
1495 //_____________________________________________________________________________
1496
1497 THnSparse* AlidNdPtHelper::RebinTHnSparse(const THnSparse* hist1, THnSparse* hist2, const Char_t* newname, Option_t* option)
1498 {
1499     THnSparse* htemp = 0;
1500     const THnSparse* hist = 0;
1501     TString opt = option;
1502     opt.ToLower();
1503     Bool_t calcErrors = kFALSE;
1504     Bool_t useRange = kFALSE;
1505     Bool_t overwrite = kFALSE;
1506     if (opt.Contains("e")) { calcErrors = kTRUE; } // calcluate correct errors (not implemented)
1507     if (opt.Contains("r")) { useRange = kTRUE; }   // use the axis range given in hist1
1508     if (opt.Contains("o")) { overwrite = kTRUE; }  // overwrite hist2 instead of creating a new one
1509     Int_t ndim  = hist1->GetNdimensions();
1510     if (ndim != hist2->GetNdimensions()) {
1511         printf("AlidNdPtHelper::RebinTHnSparse: ERROR: Histograms have different dimensions \n");
1512         return 0;
1513     }    
1514     Int_t* dims = new Int_t[ndim];
1515     for (Int_t i = 0; i < ndim; i++) { dims[i] = i; }
1516     if (useRange) {         
1517         htemp = hist1->Projection(ndim,dims,"e");
1518         hist = htemp; 
1519     } else { hist = hist1; }
1520     //THnSparse* hnew = hist2->Projection(ndim,dims,"o");
1521     //hnew->SetName(newname);
1522     THnSparse* hnew = 0;
1523     if (overwrite) { 
1524         hnew = hist2;
1525     } else {
1526         hnew = (THnSparse*) hist2->Clone(newname);
1527     }
1528     for (Int_t i = 0; i < ndim; i++) { hnew->GetAxis(i)->SetRange(); }
1529     hnew->SetTitle(hist1->GetTitle());
1530     hnew->Reset();
1531     hnew->Sumw2();
1532     Double_t content;
1533     Double_t error;
1534     Int_t* c = new Int_t[ndim];
1535     Double_t* x = new Double_t[ndim];
1536     Long64_t n = hist->GetNbins();
1537     for (Long64_t j = 0; j < n; j++) {
1538         content = hist->GetBinContent(j,c);
1539         error = hist->GetBinError(j);
1540         for (Int_t i = 0; i < ndim; i++) {
1541             x[i] = hist->GetAxis(i)->GetBinCenter(c[i]);
1542         }
1543         /* function IsInRange is protected, shit!
1544         if (useRange) {
1545             if (! hist1->IsInRange(c)) continue;
1546         }
1547         */
1548         if (calcErrors) {
1549            // implementation to be done
1550         } else {
1551            hnew->Fill(x,content);
1552         }
1553     }
1554     delete[] c; c=0;
1555     delete[] x; x=0;
1556     delete[] dims; dims=0;
1557     if (htemp) { delete htemp; htemp = 0;}
1558     return hnew;
1559 }
1560
1561