factoring out AliTriggerAnalysis class
[u/mrichter/AliRoot.git] / PWG0 / 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 #include <TROOT.h>
17 #include <TParticle.h>
18 #include <TParticlePDG.h>
19 #include <TH1.h>
20 #include <TH2.h>
21 #include <TH3.h>
22 #include <TCanvas.h>
23 #include <TList.h>
24 #include <TTree.h>
25 #include <TBranch.h>
26 #include <TLeaf.h>
27 #include <TArrayI.h>
28 #include <TF1.h>
29
30 #include <AliHeader.h>
31 #include <AliStack.h>
32 #include <AliLog.h>
33
34 #include <AliLog.h>
35 #include <AliESD.h>
36 #include <AliESDEvent.h>
37 #include <AliMCEvent.h>
38 #include <AliESDVertex.h>
39 #include <AliVertexerTracks.h>
40
41 #include <AliGenEventHeader.h>
42 #include <AliGenPythiaEventHeader.h>
43 #include <AliGenCocktailEventHeader.h>
44 #include <AliGenDPMjetEventHeader.h>
45
46 #include <AliMathBase.h>
47 #include <AliESDtrackCuts.h>
48 #include "dNdPt/AlidNdPtEventCuts.h"
49 #include "dNdPt/AlidNdPtAcceptanceCuts.h"
50 #include "dNdPt/AlidNdPtHelper.h"
51
52 //____________________________________________________________________
53 ClassImp(AlidNdPtHelper)
54
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)
57 {
58   // Get the vertex from the ESD and returns it if the vertex is valid
59   //
60   // Second argument decides which vertex is used (this selects
61   // also the quality criteria that are applied)
62
63   if(!aEsd) 
64   { 
65     ::Error("AlidNdPtHelper::GetVertex()","esd event is NULL");
66     return NULL;  
67   }
68  
69   if(!evtCuts || !accCuts || !trackCuts) 
70   { 
71     ::Error("AlidNdPtHelper::GetVertex()","cuts not available");
72     return NULL;  
73   }
74
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 )
80   {
81     vertex = aEsd->GetPrimaryVertexSPD();
82     if (debug)
83       Printf("AlidNdPtHelper::GetVertex: Returning SPD vertex");
84   }
85   else if (analysisMode == kTPC || analysisMode == kMCRec || 
86            analysisMode == kMCPion || analysisMode == kMCKaon || analysisMode == kMCProton || 
87            analysisMode == kPlus || analysisMode == kMinus)
88   {
89     if(bRedoTPC) {
90
91       Double_t kBz = aEsd->GetMagneticField();
92       AliVertexerTracks vertexer(kBz);
93
94       if(bUseMeanVertex) {
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();
100       }
101
102       Double_t maxDCAr = accCuts->GetMaxDCAr();
103       Double_t maxDCAz = accCuts->GetMaxDCAz();
104       Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
105
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);
108
109       // TPC track preselection
110       Int_t ntracks = aEsd->GetNumberOfTracks();
111       TObjArray array(ntracks);
112       UShort_t *id = new UShort_t[ntracks];
113
114       Int_t count=0;
115       for (Int_t i=0;i <ntracks; i++) {
116         AliESDtrack *t = aEsd->GetTrack(i);
117         if (!t) continue;
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()));
122         if(tpcTrack) { 
123           array.AddLast(tpcTrack);
124           id[count] = (UShort_t)t->GetID();
125           count++;
126         }
127       } 
128       AliESDVertex *vTPC = vertexer.VertexForSelectedTracks(&array,id, kTRUE, kTRUE, bUseMeanVertex);
129       
130       // set recreated TPC vertex
131       aEsd->SetPrimaryVertexTPC(vTPC);
132
133       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); i++) {
134         AliESDtrack *t = aEsd->GetTrack(i);
135         t->RelateToVertexTPC(vTPC, kBz, kVeryBig);
136       }
137       
138       delete vTPC;
139       array.Delete();
140       delete [] id; id=NULL;
141
142     }
143     vertex = aEsd->GetPrimaryVertexTPC();
144     if (debug)
145       Printf("AlidNdPtHelper::GetVertex: Returning vertex from tracks");
146     }
147       else
148        Printf("AlidNdPtHelper::GetVertex: ERROR: Invalid second argument %d", analysisMode);
149
150     if (!vertex) {
151      if (debug)
152       Printf("AlidNdPtHelper::GetVertex: No vertex found in ESD");
153       return 0;
154     }
155
156   if (debug)
157   {
158     Printf("AlidNdPtHelper::GetVertex: Returning valid vertex: %s", vertex->GetTitle());
159     vertex->Print();
160   }
161   
162   if(initVertex) delete initVertex; initVertex=NULL;
163   return vertex;
164 }
165
166 //____________________________________________________________________
167 Bool_t AlidNdPtHelper::TestRecVertex(const AliESDVertex* vertex, AnalysisMode analysisMode, Bool_t debug)
168 {
169   // Checks if a vertex meets the needed quality criteria
170   if(!vertex) return kFALSE;
171
172   Float_t requiredZResolution = -1;
173   if (analysisMode == kSPD || analysisMode == kTPCITS || analysisMode == kTPCSPDvtx)
174   {
175     requiredZResolution = 0.1;
176   }
177   else if (analysisMode == kTPC || analysisMode == kMCRec || 
178            analysisMode == kMCPion || analysisMode == kMCKaon || 
179            analysisMode == kMCProton || analysisMode ==kPlus || analysisMode ==kMinus)
180     requiredZResolution = 10.;
181
182   // check Ncontributors
183   if (vertex->GetNContributors() <= 0) {
184     if (debug){
185       Printf("AlidNdPtHelper::GetVertex: NContributors() <= 0: %d",vertex->GetNContributors());
186       Printf("AlidNdPtHelper::GetVertex: NIndices(): %d",vertex->GetNIndices());
187       vertex->Print();
188     }
189     return kFALSE;
190   }
191
192   // check resolution
193   Double_t zRes = vertex->GetZRes();
194   if (zRes == 0) {
195     Printf("AlidNdPtHelper::GetVertex: UNEXPECTED: resolution is 0.");
196     return kFALSE;
197   }
198
199   if (zRes > requiredZResolution) {
200     if (debug)
201       Printf("AlidNdPtHelper::TestVertex: Resolution too poor %f (required: %f", zRes, requiredZResolution);
202     return kFALSE;
203   }
204
205   return kTRUE;
206 }
207
208 //____________________________________________________________________
209 Bool_t AlidNdPtHelper::IsPrimaryParticle(AliStack* stack, Int_t idx, AnalysisMode analysisMode)
210 {
211 // check primary particles 
212 // depending on the analysis mode
213 //
214   if(!stack) return kFALSE;
215
216   TParticle* particle = stack->Particle(idx);
217   if (!particle) return  kFALSE;
218
219   // only charged particles
220   Double_t charge = particle->GetPDG()->Charge()/3.;
221   if (charge == 0.0) return kFALSE;
222
223   Int_t pdg = TMath::Abs(particle->GetPdgCode());
224
225   // physical primary
226   Bool_t prim = stack->IsPhysicalPrimary(idx);
227
228   if(analysisMode==kMCPion) {
229     if(prim && pdg==kPiPlus) return kTRUE;
230     else return kFALSE;
231   } 
232
233   if (analysisMode==kMCKaon) {
234     if(prim && pdg==kKPlus) return kTRUE;
235     else return kFALSE;
236   }
237     
238   if (analysisMode==kMCProton) {
239     if(prim && pdg==kProton) return kTRUE;
240     else return kFALSE;
241   }
242
243 return prim;
244 }
245
246 //____________________________________________________________________
247 void AlidNdPtHelper::PrintConf(AnalysisMode analysisMode, AliTriggerAnalysis::Trigger trigger)
248 {
249   //
250   // Prints the given configuration
251   //
252
253   TString str(">>>> Running with ");
254
255   switch (analysisMode)
256   {
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;
268   }
269   str += " and trigger ";
270
271   str += AliTriggerAnalysis::GetTriggerName(trigger);
272
273   str += " <<<<";
274
275   Printf("%s", str.Data());
276 }
277
278 //____________________________________________________________________
279 Int_t AlidNdPtHelper::ConvertPdgToPid(TParticle *particle) {
280 //
281 // Convert Abs(pdg) to pid 
282 // (0 - e, 1 - muons, 2 - pions, 3 - kaons, 4 - protons, 5 -all rest)
283 //
284 Int_t pid=-1;
285
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; }
291   else                                                       { pid = 5; }
292
293 return pid;
294 }
295
296 //_____________________________________________________________________________
297 TH1F* AlidNdPtHelper::CreateResHisto(TH2F* hRes2, TH1F **phMean, Int_t integ,  Bool_t drawBinFits, Int_t minHistEntries)
298 {
299 //
300 // Create mean and resolution 
301 // histograms
302 //
303   TVirtualPad* currentPad = gPad;
304   TAxis* axis = hRes2->GetXaxis();
305   Int_t nBins = axis->GetNbins();
306   Bool_t overflowBinFits = kFALSE;
307   TH1F* hRes, *hMean;
308   if (axis->GetXbins()->GetSize()){
309     hRes = new TH1F("hRes", "", nBins, axis->GetXbins()->GetArray());
310     hMean = new TH1F("hMean", "", nBins, axis->GetXbins()->GetArray());
311   }
312   else{
313     hRes = new TH1F("hRes", "", nBins, axis->GetXmin(), axis->GetXmax());
314     hMean = new TH1F("hMean", "", nBins, axis->GetXmin(), axis->GetXmax());
315
316   }
317   hRes->SetStats(false);
318   hRes->SetOption("E");
319   hRes->SetMinimum(0.);
320   //
321   hMean->SetStats(false);
322   hMean->SetOption("E");
323  
324   // create the fit function
325   TF1 * fitFunc = new TF1("G","[0]*exp(-(x-[1])*(x-[1])/(2.*[2]*[2]))",-3,3);
326   
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;
334   if (drawBinFits) {
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);
339   }
340
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);
348     //    
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));
353
354       if (sigma > 0.){
355         hRes->SetBinContent(bin, TMath::Abs(fitFunc->GetParameter(2)));
356         hMean->SetBinContent(bin, fitFunc->GetParameter(1));    
357       }
358       else{
359         hRes->SetBinContent(bin, 0.);
360         hMean->SetBinContent(bin,0);
361       }
362       hRes->SetBinError(bin, fitFunc->GetParError(2));
363       hMean->SetBinError(bin, fitFunc->GetParError(1));
364       
365       //
366       //
367
368     } else {
369       hRes->SetBinContent(bin, 0.);
370       hRes->SetBinError(bin, 0.);
371       hMean->SetBinContent(bin, 0.);
372       hMean->SetBinError(bin, 0.);
373     }
374     
375
376     if (drawBinFits) {
377       char name[256];
378       if (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());
382       } else {
383         sprintf(name, "%.4g < %s < %.4g", axis->GetBinLowEdge(bin),
384                 axis->GetTitle(), axis->GetBinUpEdge(bin));
385       }
386       canBinFits->cd(bin + dBin);
387       hBin->SetTitle(name);
388       hBin->SetStats(kTRUE);
389       hBin->DrawCopy("E");
390       canBinFits->Update();
391       canBinFits->Modified();
392       canBinFits->Update();
393     }
394     
395     delete hBin;
396   }
397
398   delete fitFunc;
399   currentPad->cd();
400   *phMean = hMean;
401   return hRes;
402 }
403
404 //_____________________________________________________________________________
405 TH1F* AlidNdPtHelper::MakeResol(TH2F * his, Int_t integ, Bool_t type, Bool_t drawBins, Int_t minHistEntries){
406 // Create resolution histograms
407   
408      TH1F *hisr=0, *hism=0;
409      if (!gPad) new TCanvas;
410          hisr = CreateResHisto(his,&hism,integ,drawBins,minHistEntries);
411          if (type) return hism;
412          else return hisr;
413
414 return hisr;     
415 }
416
417 //_____________________________________________________________________________
418 TObjArray* AlidNdPtHelper::GetAllChargedTracks(AliESDEvent *esdEvent, AnalysisMode analysisMode)
419 {
420   //
421   // all charged TPC particles 
422   //
423   TObjArray *allTracks = new TObjArray();
424   if(!allTracks) return allTracks;
425
426   AliESDtrack *track=0;
427   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
428   { 
429     if(analysisMode == AlidNdPtHelper::kTPC || analysisMode == AlidNdPtHelper::kTPCSPDvtx || 
430        analysisMode == AlidNdPtHelper::kMCRec || analysisMode == kMCPion || analysisMode == kMCKaon || 
431        analysisMode == kMCProton || analysisMode ==kPlus || analysisMode ==kMinus) { 
432
433       // track must be deleted by the user 
434       track = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent,iTrack);
435     } else {
436       track=esdEvent->GetTrack(iTrack);
437     }
438     if(!track) continue;
439
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) {
444
445         delete track; continue; 
446       } else {
447         continue;
448       } 
449     }
450
451     allTracks->Add(track);
452   }
453
454   if(analysisMode == AlidNdPtHelper::kTPC || analysisMode == AlidNdPtHelper::kTPCSPDvtx || 
455      analysisMode == AlidNdPtHelper::kMCRec || analysisMode == kMCPion || analysisMode == kMCKaon || 
456      analysisMode == kMCProton || analysisMode ==kPlus || analysisMode ==kMinus) {
457      
458      allTracks->SetOwner(kTRUE);
459   }
460
461 return allTracks;
462 }
463
464 //_____________________________________________________________________________
465 Int_t AlidNdPtHelper::GetTPCMBTrackMult(AliESDEvent *esdEvent, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts)
466 {
467   //
468   // get MB event track multiplicity
469   //
470   if(!esdEvent) 
471   { 
472     ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","esd event is NULL");
473     return 0;  
474   }
475  
476   if(!evtCuts || !accCuts || !trackCuts) 
477   { 
478     ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","cuts not available");
479     return 0;  
480   }
481
482   //
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);
486
487   //
488   Float_t maxDCAr = accCuts->GetMaxDCAr();
489   Float_t maxDCAz = accCuts->GetMaxDCAz();
490   Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
491   //
492   Int_t ntracks = esdEvent->GetNumberOfTracks();
493   Double_t dca[2],cov[3];
494   Int_t mult=0;
495   for (Int_t i=0;i <ntracks; i++){
496     AliESDtrack *t = esdEvent->GetTrack(i);
497     if (!t) continue;
498     if (t->Charge() == 0) continue;
499     if (!t->GetTPCInnerParam()) continue;
500     if (t->GetTPCNcls()<minTPCClust) continue;
501     //
502     AliExternalTrackParam  *tpcTrack  = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
503     if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov)) 
504     {
505       if(tpcTrack) delete tpcTrack; 
506       continue;
507     }
508     //
509     if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
510       if(tpcTrack) delete tpcTrack; 
511       continue;
512     }
513
514     mult++;    
515
516     if(tpcTrack) delete tpcTrack; 
517   }
518
519 return mult;  
520 }
521
522 //_____________________________________________________________________________
523 Int_t AlidNdPtHelper::GetTPCMBPrimTrackMult(AliESDEvent *esdEvent, AliStack * stack, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts)
524 {
525   //
526   // get MB primary event track multiplicity
527   //
528   if(!esdEvent) 
529   { 
530     ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
531     return 0;  
532   }
533
534   if(!stack) 
535   { 
536     ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
537     return 0;  
538   }
539  
540   if(!evtCuts || !accCuts || !trackCuts) 
541   { 
542     ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","cuts not available");
543     return 0;  
544   }
545
546   //
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);
550
551   //
552   Float_t maxDCAr = accCuts->GetMaxDCAr();
553   Float_t maxDCAz = accCuts->GetMaxDCAz();
554   Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
555
556   //
557   Int_t ntracks = esdEvent->GetNumberOfTracks();
558   Double_t dca[2],cov[3];
559   Int_t mult=0;
560   for (Int_t i=0;i <ntracks; i++){
561     AliESDtrack *t = esdEvent->GetTrack(i);
562     if (!t) continue;
563     if (t->Charge() == 0) continue;
564     if (!t->GetTPCInnerParam()) continue;
565     if (t->GetTPCNcls()<minTPCClust) continue;
566     //
567     AliExternalTrackParam  *tpcTrack  = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
568     if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov)) 
569     {
570       if(tpcTrack) delete tpcTrack; 
571       continue;
572     }
573     //
574     if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
575       if(tpcTrack) delete tpcTrack; 
576       continue;
577     }
578
579     Int_t label = TMath::Abs(t->GetLabel());
580     TParticle *part = stack->Particle(label);
581     if(!part) { 
582       if(tpcTrack) delete tpcTrack; 
583       continue;
584     }
585     if(!stack->IsPhysicalPrimary(label)) 
586     { 
587       if(tpcTrack) delete tpcTrack; 
588       continue;
589     }
590
591     mult++;    
592
593     if(tpcTrack) delete tpcTrack; 
594   }
595
596 return mult;  
597 }
598
599
600
601
602
603 //_____________________________________________________________________________
604 Int_t AlidNdPtHelper::GetMCTrueTrackMult(AliMCEvent *mcEvent, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts)
605 {
606   //
607   // calculate mc event true track multiplicity
608   //
609   if(!mcEvent) return 0;
610
611   AliStack* stack = 0;
612   Int_t mult = 0;
613
614   // MC particle stack
615   stack = mcEvent->Stack();
616   if (!stack) return 0;
617
618   Bool_t isEventOK = evtCuts->AcceptMCEvent(mcEvent);
619   if(!isEventOK) return 0; 
620
621   Int_t nPart  = stack->GetNtrack();
622   for (Int_t iMc = 0; iMc < nPart; ++iMc) 
623   {
624      TParticle* particle = stack->Particle(iMc);
625      if (!particle)
626      continue;
627
628      // only charged particles
629      Double_t charge = particle->GetPDG()->Charge()/3.;
630      if (charge == 0.0)
631      continue;
632       
633      // physical primary
634      Bool_t prim = stack->IsPhysicalPrimary(iMc);
635      if(!prim) continue;
636
637      // checked accepted
638      if(accCuts->AcceptTrack(particle)) 
639      {
640        mult++;
641      }
642   }
643
644 return mult;  
645 }
646
647 //_______________________________________________________________________
648 void  AlidNdPtHelper::PrintMCInfo(AliStack *pStack,Int_t label)
649 {
650 // print information about particles in the stack
651
652   if(!pStack)return;
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());
657   part->Print();
658   TParticle* mother = part;
659   Int_t imo = part->GetFirstMother();
660   Int_t nprim = pStack->GetNprimary();
661
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());
665       mother->Print();
666       imo =  mother->GetFirstMother();
667  }
668
669  Printf("########################");
670 }
671
672
673 //_____________________________________________________________________________
674 TH1* AlidNdPtHelper::GetContCorrHisto(TH1 *hist) 
675 {
676 //
677 // get contamination histogram
678 //
679  if(!hist) return 0;
680
681  Int_t nbins = hist->GetNbinsX();
682  TH1 *h_cont = (TH1D *)hist->Clone();
683
684  for(Int_t i=0; i<=nbins+1; i++) {
685    Double_t binContent = hist->GetBinContent(i);
686    Double_t binError = hist->GetBinError(i);
687
688    h_cont->SetBinContent(i,1.-binContent);
689    h_cont->SetBinError(i,binError);
690  }
691
692 return h_cont;
693 }
694
695
696 //_____________________________________________________________________________
697 TH1* AlidNdPtHelper::ScaleByBinWidth(TH1 *hist) 
698 {
699 //
700 // scale by bin width
701 //
702  if(!hist) return 0;
703
704  TH1 *h_scale = (TH1D *)hist->Clone();
705  h_scale->Scale(1.,"width");
706
707 return h_scale;
708 }
709
710 //_____________________________________________________________________________
711 TH1* AlidNdPtHelper::CalcRelativeDifference(TH1 *hist1, TH1 *hist2) 
712 {
713 //
714 // calculate rel. difference 
715 //
716
717  if(!hist1) return 0;
718  if(!hist2) return 0;
719
720  TH1 *h1_clone = (TH1D *)hist1->Clone();
721  h1_clone->Sumw2();
722
723  // (rec-mc)/mc
724  h1_clone->Add(hist2,-1);
725  h1_clone->Divide(hist2);
726
727 return h1_clone;
728 }
729
730 //_____________________________________________________________________________
731 TH1* AlidNdPtHelper::CalcRelativeDifferenceFun(TH1 *hist1, TF1 *fun) 
732 {
733 //
734 // calculate rel. difference
735 // between histogram and function
736 //
737  if(!hist1) return 0;
738  if(!fun) return 0;
739
740  TH1 *h1_clone = (TH1D *)hist1->Clone();
741  h1_clone->Sumw2();
742
743  // 
744  h1_clone->Add(fun,-1);
745  h1_clone->Divide(hist1);
746
747 return h1_clone;
748 }
749
750 //_____________________________________________________________________________
751 TH1* AlidNdPtHelper::NormalizeToEvent(TH2 *hist1, TH1 *hist2) 
752 {
753 // normalise to event for a given multiplicity bin
754 // return pt histogram 
755
756  if(!hist1) return 0;
757  if(!hist2) return 0;
758  char name[256];
759
760  Int_t nbinsX = hist1->GetNbinsX();
761  //Int_t nbinsY = hist1->GetNbinsY();
762
763  TH1D *hist_norm = 0;
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);
767
768    sprintf(name,"mom_norm");
769    if(i==0) { 
770      hist_norm = (TH1D *)hist->Clone(name);
771      hist_norm->Reset();
772    }
773
774    Double_t nbEvents = hist2->GetBinContent(i);
775    if(!nbEvents) { nbEvents = 1.; };
776
777    hist->Scale(1./nbEvents);
778    hist_norm->Add(hist);
779  }
780
781 return hist_norm;
782 }
783
784 //_____________________________________________________________________________
785 THnSparse* AlidNdPtHelper::GenerateCorrMatrix(THnSparse *hist1, THnSparse *hist2, char *name) {
786 // generate correction matrix
787 if(!hist1 || !hist2) return 0; 
788
789 THnSparse *h =(THnSparse*)hist1->Clone(name);;
790 h->Divide(hist1,hist2,1,1,"B");
791
792 return h;
793 }
794
795 //_____________________________________________________________________________
796 TH2* AlidNdPtHelper::GenerateCorrMatrix(TH2 *hist1, TH2 *hist2, char *name) {
797 // generate correction matrix
798 if(!hist1 || !hist2) return 0; 
799
800 TH2D *h =(TH2D*)hist1->Clone(name);;
801 h->Divide(hist1,hist2,1,1,"B");
802
803 return h;
804 }
805
806 //_____________________________________________________________________________
807 TH1* AlidNdPtHelper::GenerateCorrMatrix(TH1 *hist1, TH1 *hist2, char *name) {
808 // generate correction matrix
809 if(!hist1 || !hist2) return 0; 
810
811 TH1D *h =(TH1D*)hist1->Clone(name);;
812 h->Divide(hist1,hist2,1,1,"B");
813
814 return h;
815 }
816
817 //_____________________________________________________________________________
818 THnSparse* AlidNdPtHelper::GenerateContCorrMatrix(THnSparse *hist1, THnSparse *hist2, char *name) {
819 // generate contamination correction matrix
820 if(!hist1 || !hist2) return 0; 
821
822 THnSparse *hist =  GenerateCorrMatrix(hist1, hist2, name);
823 if(!hist) return 0;
824
825 // only for non ZERO bins!!!!
826
827 Int_t* coord = new Int_t[hist->GetNdimensions()];
828 memset(coord, 0, sizeof(Int_t) * hist->GetNdimensions());
829
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);
836   }
837
838 delete [] coord;
839
840 return hist;
841 }
842
843 //_____________________________________________________________________________
844 TH2* AlidNdPtHelper::GenerateContCorrMatrix(TH2 *hist1, TH2 *hist2, char *name) {
845 // generate contamination correction matrix
846 if(!hist1 || !hist2) return 0; 
847
848 TH2 *hist = GenerateCorrMatrix(hist1, hist2, name);
849 if(!hist) return 0;
850
851 Int_t nBinsX = hist->GetNbinsX();
852 Int_t nBinsY = hist->GetNbinsY();
853
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);
860   }
861   }
862
863 return hist;
864 }
865
866 //_____________________________________________________________________________
867 TH1* AlidNdPtHelper::GenerateContCorrMatrix(TH1 *hist1, TH1 *hist2, char *name) {
868 // generate contamination correction matrix
869 if(!hist1 || !hist2) return 0; 
870
871 TH1 *hist = GenerateCorrMatrix(hist1, hist2, name);
872 if(!hist) return 0;
873
874 Int_t nBinsX = hist->GetNbinsX();
875
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);
881   }
882
883 return hist;
884 }
885
886 //_____________________________________________________________________________
887 const AliESDVertex* AlidNdPtHelper::GetTPCVertexZ(AliESDEvent* esdEvent, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts, Float_t fraction, Int_t ntracksMin){
888   //
889   // TPC Z vertexer
890   //
891   if(!esdEvent)
892   { 
893     ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
894     return NULL;  
895   }
896
897   if(!evtCuts || !accCuts || !trackCuts) 
898   { 
899     ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
900     return NULL;  
901   }
902
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);
906
907   Double_t maxDCAr = accCuts->GetMaxDCAr();
908   Double_t maxDCAz = accCuts->GetMaxDCAz();
909   Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
910
911   //
912   Int_t ntracks = esdEvent->GetNumberOfTracks();
913   TVectorD ztrack(ntracks);
914   Double_t dca[2],cov[3];
915   Int_t counter=0;
916   for (Int_t i=0;i <ntracks; i++){
917     AliESDtrack *t = esdEvent->GetTrack(i);
918     if (!t) continue;
919     if (!t->GetTPCInnerParam()) continue;
920     if (t->GetTPCNcls()<minTPCClust) continue;
921     //
922     AliExternalTrackParam  *tpcTrack  = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
923     if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov)) continue;
924
925     //
926     if (TMath::Abs(dca[0])>maxDCAr) continue;
927     //if (TMath::Sqrt(cov[0])>sigmaXYcut) continue;    
928     if (TMath::Abs(tpcTrack->GetZ())>maxDCAz) continue;
929
930     ztrack[counter]=tpcTrack->GetZ();
931     counter++;    
932
933     if(tpcTrack) delete tpcTrack;
934   }
935
936   //
937   // Find LTM z position
938   //
939   Double_t mean=0, sigma=0;
940   if (counter<ntracksMin) return 0;
941   //
942   Int_t nused = TMath::Nint(counter*fraction);
943   if (nused==counter) nused=counter-1;  
944   if (nused>1){
945     AliMathBase::EvaluateUni(counter, ztrack.GetMatrixArray(), mean,sigma, TMath::Nint(counter*fraction));
946     sigma/=TMath::Sqrt(nused);
947   }else{
948     mean  = TMath::Mean(counter, ztrack.GetMatrixArray());
949     sigma = TMath::RMS(counter, ztrack.GetMatrixArray());
950     sigma/=TMath::Sqrt(counter-1);
951   }
952   vtxpos[2]=mean;
953   vtxsigma[2]=sigma;
954   const AliESDVertex* vertex = new AliESDVertex(vtxpos, vtxsigma);
955   return vertex;
956 }
957
958 //_____________________________________________________________________________
959 Int_t  AlidNdPtHelper::GetSPDMBTrackMult(AliESDEvent* esdEvent, Float_t deltaThetaCut, Float_t deltaPhiCut) 
960 {
961   //
962   // SPD track multiplicity
963   //
964
965   // get tracklets
966   const AliMultiplicity* mult = esdEvent->GetMultiplicity();
967   if (!mult)
968      return 0;
969
970   // get multiplicity from SPD tracklets
971   Int_t inputCount = 0; 
972   for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
973   {
974     //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
975
976      Float_t phi = mult->GetPhi(i);
977      if (phi < 0)
978        phi += TMath::Pi() * 2;
979      Float_t deltaPhi = mult->GetDeltaPhi(i);
980      Float_t deltaTheta = mult->GetDeltaTheta(i);
981
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);
984
985      if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
986         continue;
987
988      if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
989         continue;
990       
991      ++inputCount;
992   }
993
994 return inputCount;
995 }
996
997 //_____________________________________________________________________________
998 Int_t  AlidNdPtHelper::GetSPDMBPrimTrackMult(AliESDEvent* esdEvent, AliStack* stack, Float_t deltaThetaCut, Float_t deltaPhiCut) 
999 {
1000   //
1001   // SPD track multiplicity
1002   //
1003
1004   // get tracklets
1005   const AliMultiplicity* mult = esdEvent->GetMultiplicity();
1006   if (!mult)
1007      return 0;
1008
1009   // get multiplicity from SPD tracklets
1010   Int_t inputCount = 0; 
1011   for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
1012   {
1013     //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
1014
1015      Float_t phi = mult->GetPhi(i);
1016      if (phi < 0)
1017        phi += TMath::Pi() * 2;
1018      Float_t deltaPhi = mult->GetDeltaPhi(i);
1019      Float_t deltaTheta = mult->GetDeltaTheta(i);
1020
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);
1023
1024      if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
1025         continue;
1026
1027      if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
1028         continue;
1029
1030
1031      if (mult->GetLabel(i, 0) < 0 || mult->GetLabel(i, 0) != mult->GetLabel(i, 1) || 
1032          !stack->IsPhysicalPrimary(mult->GetLabel(i, 0)))
1033         continue;
1034
1035       
1036      ++inputCount;
1037   }
1038
1039 return inputCount;
1040 }
1041
1042