fix coding rules violation
[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 #include <TLorentzVector.h>
30
31 #include <AliHeader.h>
32 #include <AliStack.h>
33 #include <AliLog.h>
34
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 || analysisMode == kTPCSPDvtxUpdate)
78   {
79     vertex = aEsd->GetPrimaryVertexSPD();
80     if (debug)
81       Printf("AlidNdPtHelper::GetVertex: Returning SPD vertex");
82   }
83   else if (analysisMode == kTPC)
84   {
85     if(bRedoTPC) {
86
87       Double_t kBz = aEsd->GetMagneticField();
88       AliVertexerTracks vertexer(kBz);
89
90       if(bUseMeanVertex) {
91          Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
92          Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
93          initVertex = new AliESDVertex(pos,err);
94          vertexer.SetVtxStart(initVertex);
95          vertexer.SetConstraintOn();
96       }
97
98       Double_t maxDCAr = accCuts->GetMaxDCAr();
99       Double_t maxDCAz = accCuts->GetMaxDCAz();
100       Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
101
102       //vertexer.SetTPCMode(Double_t dcacut=0.1, Double_t dcacutIter0=1.0, Double_t maxd0z0=5.0, Int_t minCls=10, Int_t mintrks=1, Double_t nsigma=3., Double_t mindetfitter=0.1, Double_t maxtgl=1.5, Double_t fidR=3., Double_t fidZ=30., Int_t finderAlgo=1, Int_t finderAlgoIter0=4);
103       vertexer.SetTPCMode(0.1,1.0,5.0,minTPCClust,1,3.,0.1,2.0,maxDCAr,maxDCAz,1,4);
104
105       // TPC track preselection
106       Int_t ntracks = aEsd->GetNumberOfTracks();
107       TObjArray array(ntracks);
108       UShort_t *id = new UShort_t[ntracks];
109
110       Int_t count=0;
111       for (Int_t i=0;i <ntracks; i++) {
112         AliESDtrack *t = aEsd->GetTrack(i);
113         if (!t) continue;
114         if (t->Charge() == 0) continue;
115         if (!t->GetTPCInnerParam()) continue;
116         if (t->GetTPCNcls()<vertexer.GetMinClusters()) continue;
117         AliExternalTrackParam  *tpcTrack  = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
118         if(tpcTrack) { 
119           array.AddLast(tpcTrack);
120           id[count] = (UShort_t)t->GetID();
121           count++;
122         }
123       } 
124       AliESDVertex *vTPC = vertexer.VertexForSelectedTracks(&array,id, kTRUE, kTRUE, bUseMeanVertex);
125       
126       // set recreated TPC vertex
127       aEsd->SetPrimaryVertexTPC(vTPC);
128
129       for (Int_t i=0; i<aEsd->GetNumberOfTracks(); i++) {
130         AliESDtrack *t = aEsd->GetTrack(i);
131         t->RelateToVertexTPC(vTPC, kBz, kVeryBig);
132       }
133       
134       delete vTPC;
135       array.Delete();
136       delete [] id; id=NULL;
137
138     }
139     vertex = aEsd->GetPrimaryVertexTPC();
140     if (debug)
141       Printf("AlidNdPtHelper::GetVertex: Returning vertex from tracks");
142     }
143       else
144        Printf("AlidNdPtHelper::GetVertex: ERROR: Invalid second argument %d", analysisMode);
145
146     if (!vertex) {
147      if (debug)
148       Printf("AlidNdPtHelper::GetVertex: No vertex found in ESD");
149       return 0;
150     }
151
152   if (debug)
153   {
154     Printf("AlidNdPtHelper::GetVertex: Returning valid vertex: %s", vertex->GetTitle());
155     vertex->Print();
156   }
157   
158   if(initVertex) delete initVertex; initVertex=NULL;
159   return vertex;
160 }
161
162 //____________________________________________________________________
163 Bool_t AlidNdPtHelper::TestRecVertex(const AliESDVertex* vertex, AnalysisMode analysisMode, Bool_t debug)
164 {
165   // Checks if a vertex meets the needed quality criteria
166   if(!vertex) return kFALSE;
167
168   Float_t requiredZResolution = -1;
169   if (analysisMode == kSPD || analysisMode == kTPCITS || analysisMode == kTPCSPDvtx || analysisMode == kTPCSPDvtxUpdate)
170   {
171     requiredZResolution = 0.1;
172   }
173   else if (analysisMode == kTPC)
174     requiredZResolution = 10.;
175
176   // check Ncontributors
177   if (vertex->GetNContributors() <= 0) {
178     if (debug){
179       Printf("AlidNdPtHelper::GetVertex: NContributors() <= 0: %d",vertex->GetNContributors());
180       Printf("AlidNdPtHelper::GetVertex: NIndices(): %d",vertex->GetNIndices());
181       vertex->Print();
182     }
183     return kFALSE;
184   }
185
186   // check resolution
187   Double_t zRes = vertex->GetZRes();
188   if (zRes == 0) {
189     Printf("AlidNdPtHelper::GetVertex: UNEXPECTED: resolution is 0.");
190     return kFALSE;
191   }
192
193   if (zRes > requiredZResolution) {
194     if (debug)
195       Printf("AlidNdPtHelper::TestVertex: Resolution too poor %f (required: %f", zRes, requiredZResolution);
196     return kFALSE;
197   }
198
199   return kTRUE;
200 }
201
202 //____________________________________________________________________
203 Bool_t AlidNdPtHelper::IsPrimaryParticle(AliStack* stack, Int_t idx, ParticleMode particleMode)
204 {
205 // check primary particles 
206 // depending on the particle mode
207 //
208   if(!stack) return kFALSE;
209
210   TParticle* particle = stack->Particle(idx);
211   if (!particle) return  kFALSE;
212
213   // only charged particles
214   Double_t charge = particle->GetPDG()->Charge()/3.;
215   if (charge == 0.0) return kFALSE;
216
217   Int_t pdg = TMath::Abs(particle->GetPdgCode());
218
219   // physical primary
220   Bool_t prim = stack->IsPhysicalPrimary(idx);
221
222   if(particleMode==kMCPion) {
223     if(prim && pdg==kPiPlus) return kTRUE;
224     else return kFALSE;
225   } 
226
227   if (particleMode==kMCKaon) {
228     if(prim && pdg==kKPlus) return kTRUE;
229     else return kFALSE;
230   }
231     
232   if (particleMode==kMCProton) {
233     if(prim && pdg==kProton) return kTRUE;
234     else return kFALSE;
235   }
236
237 return prim;
238 }
239
240 //____________________________________________________________________
241 /*
242 Bool_t AlidNdPtHelper::IsCosmicTrack(TObjArray *allChargedTracks, AliESDtrack *track1, Int_t trackIdx, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts)
243 {
244 //
245 // check cosmic tracks
246 //
247   if(!allChargedTracks) return kFALSE;
248   if(!track1) return kFALSE;
249   if(!accCuts) return kFALSE;
250   if(!trackCuts) return kFALSE;
251
252   Int_t entries = allChargedTracks->GetEntries();
253   for(Int_t i=0; i<entries;++i) 
254   {
255       //
256       // exclude the same tracks
257       //
258       if(i == trackIdx) continue;
259
260       AliESDtrack *track2 = (AliESDtrack*)allChargedTracks->At(i);
261       if(!track2) continue;
262       if(track2->Charge()==0) continue;
263
264       if(track1->Pt() > 6. && track2->Pt() > 6. && (track1->Charge() + track2->Charge()) == 0 ) 
265       {
266         printf("track1->Theta() %f, track1->Eta() %f, track1->Phi() %f, track1->Charge() %d  \n", track1->Theta(), track1->Eta(), track1->Phi(), track1->Charge());
267         printf("track2->Theta() %f, track2->Eta() %f, track2->Phi() %f, track2->Charge() %d  \n", track2->Theta(), track2->Eta(), track2->Phi(), track2->Charge());
268
269         printf("deta %f, dphi %f, dq %d  \n", track1->Eta()-track2->Eta(), track1->Phi()-track2->Phi(), track1->Charge()+track2->Charge()); 
270
271       }
272
273       //
274       // cosmic tracks in TPC
275       //
276       //if( TMath::Abs( track1->Theta() - track2->Theta() ) < 0.004  && 
277       //  ((TMath::Abs(track1->Phi()-track2->Phi()) - TMath::Pi() )<0.004) && 
278       if( (track1->Pt()-track2->Pt()) < 0.1 && track1->Pt() > 4.0 && (track1->Charge()+track2->Charge()) == 0 )
279       {
280         //printf("COSMIC  candidate \n");
281         printf("track1->Theta() %f, track1->Eta() %f, track1->Phi() %f, track1->Charge() %d  \n", track1->Theta(), track1->Eta(), track1->Phi(), track1->Charge());
282         printf("track2->Theta() %f, track2->Eta() %f, track2->Phi() %f, track2->Charge() %d  \n", track2->Theta(), track2->Eta(), track2->Phi(), track2->Charge());
283         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()); 
284         return kTRUE;
285       }
286    }
287      
288 return kFALSE; 
289 }
290 */
291
292 //____________________________________________________________________
293 Bool_t AlidNdPtHelper::IsCosmicTrack(AliESDtrack *track1, AliESDtrack *track2, Int_t /*trackIdx*/, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts)
294 {
295 //
296 // check cosmic tracks
297 //
298   if(!track1) return kFALSE;
299   if(!track2) return kFALSE;
300   if(!accCuts) return kFALSE;
301   if(!trackCuts) return kFALSE;
302
303   /*
304   if(track1->Pt() > 6. && track2->Pt() > 6. && (track1->Charge() + track2->Charge()) == 0 ) 
305   {
306         printf("track1->Theta() %f, track1->Eta() %f, track1->Phi() %f, track1->Charge() %d  \n", track1->Theta(), track1->Eta(), track1->Phi(), track1->Charge());
307         printf("track2->Theta() %f, track2->Eta() %f, track2->Phi() %f, track2->Charge() %d  \n", track2->Theta(), track2->Eta(), track2->Phi(), track2->Charge());
308
309         printf("deta %f, dphi %f, dq %d  \n", track1->Eta()-track2->Eta(), track1->Phi()-track2->Phi(), track1->Charge()+track2->Charge()); 
310
311       }
312       */
313
314       //
315       // cosmic tracks in TPC
316       //
317       //if( TMath::Abs( track1->Theta() - track2->Theta() ) < 0.004  && 
318       //  ((TMath::Abs(track1->Phi()-track2->Phi()) - TMath::Pi() )<0.004) && 
319
320       if( (track1->Pt()-track2->Pt()) < 0.2 && track1->Pt() > 3.0 && 
321              (track1->Charge()+track2->Charge()) == 0 )
322       {
323         //printf("COSMIC  candidate \n");
324         /*
325         printf("track1->Theta() %f, track1->Eta() %f, track1->Phi() %f, track1->Charge() %d  \n", track1->Theta(), track1->Eta(), track1->Phi(), track1->Charge());
326         printf("track2->Theta() %f, track2->Eta() %f, track2->Phi() %f, track2->Charge() %d  \n", track2->Theta(), track2->Eta(), track2->Phi(), track2->Charge());
327         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()); 
328         */
329         return kTRUE;
330       }
331      
332 return kFALSE; 
333 }
334
335 //____________________________________________________________________
336 void AlidNdPtHelper::PrintConf(AnalysisMode analysisMode, AliTriggerAnalysis::Trigger trigger)
337 {
338   //
339   // Prints the given configuration
340   //
341
342   TString str(">>>> Running with ");
343
344   switch (analysisMode)
345   {
346     case kInvalid: str += "invalid setting"; break;
347     case kSPD : str += "SPD-only"; break;
348     case kTPC : str += "TPC-only"; break;
349     case kTPCITS : str += "Global tracking"; break;
350     case kTPCSPDvtx : str += "TPC tracking + SPD event vertex"; break;
351     case kTPCSPDvtxUpdate : str += "TPC tracks updated with SPD event vertex point"; break;
352     case kMCRec : str += "TPC tracking + Replace rec. with MC values"; break;
353   }
354   str += " and trigger ";
355
356   str += AliTriggerAnalysis::GetTriggerName(trigger);
357
358   str += " <<<<";
359
360   Printf("%s", str.Data());
361 }
362
363 //____________________________________________________________________
364 Int_t AlidNdPtHelper::ConvertPdgToPid(TParticle *particle) {
365 //
366 // Convert Abs(pdg) to pid 
367 // (0 - e, 1 - muons, 2 - pions, 3 - kaons, 4 - protons, 5 -all rest)
368 //
369 Int_t pid=-1;
370
371   if (TMath::Abs(particle->GetPdgCode()) == kElectron)         { pid = 0; }
372   else if (TMath::Abs(particle->GetPdgCode()) == kMuonMinus) { pid = 1; }
373   else if (TMath::Abs(particle->GetPdgCode()) == kPiPlus)    { pid = 2; }
374   else if (TMath::Abs(particle->GetPdgCode()) == kKPlus)     { pid = 3; }
375   else if (TMath::Abs(particle->GetPdgCode()) == kProton)    { pid = 4; }
376   else                                                       { pid = 5; }
377
378 return pid;
379 }
380
381 //_____________________________________________________________________________
382 TH1F* AlidNdPtHelper::CreateResHisto(TH2F* hRes2, TH1F **phMean, Int_t integ,  Bool_t drawBinFits, Int_t minHistEntries)
383 {
384 //
385 // Create mean and resolution 
386 // histograms
387 //
388   TVirtualPad* currentPad = gPad;
389   TAxis* axis = hRes2->GetXaxis();
390   Int_t nBins = axis->GetNbins();
391   Bool_t overflowBinFits = kFALSE;
392   TH1F* hRes, *hMean;
393   if (axis->GetXbins()->GetSize()){
394     hRes = new TH1F("hRes", "", nBins, axis->GetXbins()->GetArray());
395     hMean = new TH1F("hMean", "", nBins, axis->GetXbins()->GetArray());
396   }
397   else{
398     hRes = new TH1F("hRes", "", nBins, axis->GetXmin(), axis->GetXmax());
399     hMean = new TH1F("hMean", "", nBins, axis->GetXmin(), axis->GetXmax());
400
401   }
402   hRes->SetStats(false);
403   hRes->SetOption("E");
404   hRes->SetMinimum(0.);
405   //
406   hMean->SetStats(false);
407   hMean->SetOption("E");
408  
409   // create the fit function
410   TF1 * fitFunc = new TF1("G","[0]*exp(-(x-[1])*(x-[1])/(2.*[2]*[2]))",-3,3);
411   
412   fitFunc->SetLineWidth(2);
413   fitFunc->SetFillStyle(0);
414   // create canvas for fits
415   TCanvas* canBinFits = NULL;
416   Int_t nPads = (overflowBinFits) ? nBins+2 : nBins;
417   Int_t nx = Int_t(sqrt(nPads-1.));// + 1;
418   Int_t ny = (nPads-1) / nx + 1;
419   if (drawBinFits) {
420     canBinFits = (TCanvas*)gROOT->FindObject("canBinFits");
421     if (canBinFits) delete canBinFits;
422     canBinFits = new TCanvas("canBinFits", "fits of bins", 200, 100, 500, 700);
423     canBinFits->Divide(nx, ny);
424   }
425
426   // loop over x bins and fit projection
427   Int_t dBin = ((overflowBinFits) ? 1 : 0);
428   for (Int_t bin = 1-dBin; bin <= nBins+dBin; bin++) {
429     if (drawBinFits) canBinFits->cd(bin + dBin);
430     Int_t bin0=TMath::Max(bin-integ,0);
431     Int_t bin1=TMath::Min(bin+integ,nBins);
432     TH1D* hBin = hRes2->ProjectionY("hBin", bin0, bin1);
433     //    
434     if (hBin->GetEntries() > minHistEntries) {
435       fitFunc->SetParameters(hBin->GetMaximum(),hBin->GetMean(),hBin->GetRMS());
436       hBin->Fit(fitFunc,"s");
437       Double_t sigma = TMath::Abs(fitFunc->GetParameter(2));
438
439       if (sigma > 0.){
440         hRes->SetBinContent(bin, TMath::Abs(fitFunc->GetParameter(2)));
441         hMean->SetBinContent(bin, fitFunc->GetParameter(1));    
442       }
443       else{
444         hRes->SetBinContent(bin, 0.);
445         hMean->SetBinContent(bin,0);
446       }
447       hRes->SetBinError(bin, fitFunc->GetParError(2));
448       hMean->SetBinError(bin, fitFunc->GetParError(1));
449       
450       //
451       //
452
453     } else {
454       hRes->SetBinContent(bin, 0.);
455       hRes->SetBinError(bin, 0.);
456       hMean->SetBinContent(bin, 0.);
457       hMean->SetBinError(bin, 0.);
458     }
459     
460
461     if (drawBinFits) {
462       char name[256];
463       if (bin == 0) {
464         sprintf(name, "%s < %.4g", axis->GetTitle(), axis->GetBinUpEdge(bin));
465       } else if (bin == nBins+1) {
466         sprintf(name, "%.4g < %s", axis->GetBinLowEdge(bin), axis->GetTitle());
467       } else {
468         sprintf(name, "%.4g < %s < %.4g", axis->GetBinLowEdge(bin),
469                 axis->GetTitle(), axis->GetBinUpEdge(bin));
470       }
471       canBinFits->cd(bin + dBin);
472       hBin->SetTitle(name);
473       hBin->SetStats(kTRUE);
474       hBin->DrawCopy("E");
475       canBinFits->Update();
476       canBinFits->Modified();
477       canBinFits->Update();
478     }
479     
480     delete hBin;
481   }
482
483   delete fitFunc;
484   currentPad->cd();
485   *phMean = hMean;
486   return hRes;
487 }
488
489 //_____________________________________________________________________________
490 TH1F* AlidNdPtHelper::MakeResol(TH2F * his, Int_t integ, Bool_t type, Bool_t drawBins, Int_t minHistEntries){
491 // Create resolution histograms
492   
493      TH1F *hisr=0, *hism=0;
494      if (!gPad) new TCanvas;
495          hisr = CreateResHisto(his,&hism,integ,drawBins,minHistEntries);
496          if (type) return hism;
497          else return hisr;
498
499 return hisr;     
500 }
501
502 //_____________________________________________________________________________
503 TObjArray* AlidNdPtHelper::GetAllChargedTracks(AliESDEvent *esdEvent, AnalysisMode analysisMode)
504 {
505   //
506   // all charged TPC particles 
507   //
508   TObjArray *allTracks = new TObjArray();
509   if(!allTracks) return allTracks;
510
511   AliESDtrack *track=0;
512   for (Int_t iTrack = 0; iTrack < esdEvent->GetNumberOfTracks(); iTrack++) 
513   { 
514     if(analysisMode == AlidNdPtHelper::kTPC) { 
515       // track must be deleted by user 
516       track = AliESDtrackCuts::GetTPCOnlyTrack(esdEvent,iTrack);
517       if(!track) continue;
518     } 
519     else if (analysisMode == AlidNdPtHelper::kTPCSPDvtx || AlidNdPtHelper::kTPCSPDvtxUpdate)
520     {
521       // track must be deleted by the user 
522       track = AlidNdPtHelper::GetTPCOnlyTrackSPDvtx(esdEvent,iTrack,kFALSE);
523       if(!track) continue;
524     }
525     else {
526       track=esdEvent->GetTrack(iTrack);
527     }
528     if(!track) continue;
529
530     if(track->Charge()==0) { 
531       if(analysisMode == AlidNdPtHelper::kTPC || analysisMode == AlidNdPtHelper::kTPCSPDvtx || 
532          analysisMode == AlidNdPtHelper::kTPCSPDvtxUpdate) 
533       {
534         delete track; continue; 
535       } else {
536         continue;
537       } 
538     }
539
540     allTracks->Add(track);
541   }
542
543   if(analysisMode == AlidNdPtHelper::kTPC || analysisMode == AlidNdPtHelper::kTPCSPDvtx || 
544      analysisMode == AlidNdPtHelper::kTPCSPDvtxUpdate) {
545      
546      allTracks->SetOwner(kTRUE);
547   }
548
549 return allTracks;
550 }
551
552 //_____________________________________________________________________________
553 AliESDtrack *AlidNdPtHelper::GetTPCOnlyTrackSPDvtx(AliESDEvent* esdEvent, Int_t iTrack, Bool_t bUpdate)
554 {
555 //
556 // Create ESD tracks from TPCinner parameters.
557 // Propagte to DCA to SPD vertex.
558 // Update using SPD vertex point (parameter)
559 //
560 // It is user responsibility to delete these tracks
561 //
562
563   if (!esdEvent) return NULL;
564   if (!esdEvent->GetPrimaryVertexSPD() ) { return NULL; }
565   if (!esdEvent->GetPrimaryVertexSPD()->GetStatus() ) { return  NULL; }
566    
567   // 
568   AliESDtrack* track = esdEvent->GetTrack(iTrack);
569   if (!track)
570     return NULL;
571
572   // create new ESD track
573   AliESDtrack *tpcTrack = new AliESDtrack();
574  
575   // relate TPC-only tracks (TPCinner) to SPD vertex
576   AliExternalTrackParam cParam;
577   if(bUpdate) {  
578     track->RelateToVertexTPC(esdEvent->GetPrimaryVertexSPD(),esdEvent->GetMagneticField(),kVeryBig,&cParam);
579     track->Set(cParam.GetX(),cParam.GetAlpha(),cParam.GetParameter(),cParam.GetCovariance());
580
581     // reject fake tracks
582     if(track->Pt() > 10000.)  {
583       ::Error("Exclude no physical tracks","pt>10000. GeV");
584       delete tpcTrack; 
585       return NULL;
586     }
587   }
588   else {
589     track->RelateToVertexTPC(esdEvent->GetPrimaryVertexSPD(), esdEvent->GetMagneticField(), kVeryBig);
590   }
591
592   // only true if we have a tpc track
593   if (!track->FillTPCOnlyTrack(*tpcTrack))
594   {
595     delete tpcTrack;
596     return NULL;
597   }
598
599 return tpcTrack;
600
601
602 //_____________________________________________________________________________
603 Int_t AlidNdPtHelper::GetTPCMBTrackMult(AliESDEvent *esdEvent, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts)
604 {
605   //
606   // get MB event track multiplicity
607   //
608   if(!esdEvent) 
609   { 
610     ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","esd event is NULL");
611     return 0;  
612   }
613  
614   if(!evtCuts || !accCuts || !trackCuts) 
615   { 
616     ::Error("AlidNdPtHelper::GetTPCMBTrackMult()","cuts not available");
617     return 0;  
618   }
619
620   //
621   Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
622   Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
623   AliESDVertex vtx0(pos,err);
624
625   //
626   Float_t maxDCAr = accCuts->GetMaxDCAr();
627   Float_t maxDCAz = accCuts->GetMaxDCAz();
628   Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
629   //
630   Int_t ntracks = esdEvent->GetNumberOfTracks();
631   Double_t dca[2],cov[3];
632   Int_t mult=0;
633   for (Int_t i=0;i <ntracks; i++){
634     AliESDtrack *t = esdEvent->GetTrack(i);
635     if (!t) continue;
636     if (t->Charge() == 0) continue;
637     if (!t->GetTPCInnerParam()) continue;
638     if (t->GetTPCNcls()<minTPCClust) continue;
639     //
640     AliExternalTrackParam  *tpcTrack  = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
641     if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov)) 
642     {
643       if(tpcTrack) delete tpcTrack; 
644       continue;
645     }
646     //
647     if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
648       if(tpcTrack) delete tpcTrack; 
649       continue;
650     }
651
652     mult++;    
653
654     if(tpcTrack) delete tpcTrack; 
655   }
656
657 return mult;  
658 }
659
660 //_____________________________________________________________________________
661 Int_t AlidNdPtHelper::GetTPCMBPrimTrackMult(AliESDEvent *esdEvent, AliStack * stack, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts)
662 {
663   //
664   // get MB primary event track multiplicity
665   //
666   if(!esdEvent) 
667   { 
668     ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
669     return 0;  
670   }
671
672   if(!stack) 
673   { 
674     ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","esd event is NULL");
675     return 0;  
676   }
677  
678   if(!evtCuts || !accCuts || !trackCuts) 
679   { 
680     ::Error("AlidNdPtHelper::GetTPCMBPrimTrackMult()","cuts not available");
681     return 0;  
682   }
683
684   //
685   Double_t pos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
686   Double_t err[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
687   AliESDVertex vtx0(pos,err);
688
689   //
690   Float_t maxDCAr = accCuts->GetMaxDCAr();
691   Float_t maxDCAz = accCuts->GetMaxDCAz();
692   Float_t minTPCClust = trackCuts->GetMinNClusterTPC();
693
694   //
695   Int_t ntracks = esdEvent->GetNumberOfTracks();
696   Double_t dca[2],cov[3];
697   Int_t mult=0;
698   for (Int_t i=0;i <ntracks; i++){
699     AliESDtrack *t = esdEvent->GetTrack(i);
700     if (!t) continue;
701     if (t->Charge() == 0) continue;
702     if (!t->GetTPCInnerParam()) continue;
703     if (t->GetTPCNcls()<minTPCClust) continue;
704     //
705     AliExternalTrackParam  *tpcTrack  = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
706     if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov)) 
707     {
708       if(tpcTrack) delete tpcTrack; 
709       continue;
710     }
711     //
712     if (TMath::Abs(dca[0])>maxDCAr || TMath::Abs(dca[1])>maxDCAz) {
713       if(tpcTrack) delete tpcTrack; 
714       continue;
715     }
716
717     Int_t label = TMath::Abs(t->GetLabel());
718     TParticle *part = stack->Particle(label);
719     if(!part) { 
720       if(tpcTrack) delete tpcTrack; 
721       continue;
722     }
723     if(!stack->IsPhysicalPrimary(label)) 
724     { 
725       if(tpcTrack) delete tpcTrack; 
726       continue;
727     }
728
729     mult++;    
730
731     if(tpcTrack) delete tpcTrack; 
732   }
733
734 return mult;  
735 }
736
737
738
739
740
741 //_____________________________________________________________________________
742 Int_t AlidNdPtHelper::GetMCTrueTrackMult(AliMCEvent *mcEvent, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts)
743 {
744   //
745   // calculate mc event true track multiplicity
746   //
747   if(!mcEvent) return 0;
748
749   AliStack* stack = 0;
750   Int_t mult = 0;
751
752   // MC particle stack
753   stack = mcEvent->Stack();
754   if (!stack) return 0;
755
756   Bool_t isEventOK = evtCuts->AcceptMCEvent(mcEvent);
757   if(!isEventOK) return 0; 
758
759   Int_t nPart  = stack->GetNtrack();
760   for (Int_t iMc = 0; iMc < nPart; ++iMc) 
761   {
762      TParticle* particle = stack->Particle(iMc);
763      if (!particle)
764      continue;
765
766      // only charged particles
767      Double_t charge = particle->GetPDG()->Charge()/3.;
768      if (charge == 0.0)
769      continue;
770       
771      // physical primary
772      Bool_t prim = stack->IsPhysicalPrimary(iMc);
773      if(!prim) continue;
774
775      // checked accepted
776      if(accCuts->AcceptTrack(particle)) 
777      {
778        mult++;
779      }
780   }
781
782 return mult;  
783 }
784
785 //_______________________________________________________________________
786 void  AlidNdPtHelper::PrintMCInfo(AliStack *pStack,Int_t label)
787 {
788 // print information about particles in the stack
789
790   if(!pStack)return;
791   label = TMath::Abs(label);
792   TParticle *part = pStack->Particle(label);
793   Printf("########################");
794   Printf("%s:%d %d UniqueID %d PDG %d P %3.3f",(char*)__FILE__,__LINE__,label,part->GetUniqueID(),part->GetPdgCode(),part->P());
795   part->Print();
796   TParticle* mother = part;
797   Int_t imo = part->GetFirstMother();
798   Int_t nprim = pStack->GetNprimary();
799
800   while((imo >= nprim)) {
801       mother =  pStack->Particle(imo);
802       Printf("Mother %s:%d Label %d UniqueID %d PDG %d P %3.3f",(char*)__FILE__,__LINE__,imo,mother->GetUniqueID(),mother->GetPdgCode(),mother->P());
803       mother->Print();
804       imo =  mother->GetFirstMother();
805  }
806
807  Printf("########################");
808 }
809
810
811 //_____________________________________________________________________________
812 TH1* AlidNdPtHelper::GetContCorrHisto(TH1 *hist) 
813 {
814 //
815 // get contamination histogram
816 //
817  if(!hist) return 0;
818
819  Int_t nbins = hist->GetNbinsX();
820  TH1 *h_cont = (TH1D *)hist->Clone();
821
822  for(Int_t i=0; i<=nbins+1; i++) {
823    Double_t binContent = hist->GetBinContent(i);
824    Double_t binError = hist->GetBinError(i);
825
826    h_cont->SetBinContent(i,1.-binContent);
827    h_cont->SetBinError(i,binError);
828  }
829
830 return h_cont;
831 }
832
833
834 //_____________________________________________________________________________
835 TH1* AlidNdPtHelper::ScaleByBinWidth(TH1 *hist) 
836 {
837 //
838 // scale by bin width
839 //
840  if(!hist) return 0;
841
842  TH1 *h_scale = (TH1D *)hist->Clone();
843  h_scale->Scale(1.,"width");
844
845 return h_scale;
846 }
847
848 //_____________________________________________________________________________
849 TH1* AlidNdPtHelper::CalcRelativeDifference(TH1 *hist1, TH1 *hist2) 
850 {
851 //
852 // calculate rel. difference 
853 //
854
855  if(!hist1) return 0;
856  if(!hist2) return 0;
857
858  TH1 *h1_clone = (TH1D *)hist1->Clone();
859  h1_clone->Sumw2();
860
861  // (rec-mc)/mc
862  h1_clone->Add(hist2,-1);
863  h1_clone->Divide(hist2);
864
865 return h1_clone;
866 }
867
868 //_____________________________________________________________________________
869 TH1* AlidNdPtHelper::CalcRelativeDifferenceFun(TH1 *hist1, TF1 *fun) 
870 {
871 //
872 // calculate rel. difference
873 // between histogram and function
874 //
875  if(!hist1) return 0;
876  if(!fun) return 0;
877
878  TH1 *h1_clone = (TH1D *)hist1->Clone();
879  h1_clone->Sumw2();
880
881  // 
882  h1_clone->Add(fun,-1);
883  h1_clone->Divide(hist1);
884
885 return h1_clone;
886 }
887
888 //_____________________________________________________________________________
889 TH1* AlidNdPtHelper::NormalizeToEvent(TH2 *hist1, TH1 *hist2) 
890 {
891 // normalise to event for a given multiplicity bin
892 // return pt histogram 
893
894  if(!hist1) return 0;
895  if(!hist2) return 0;
896  char name[256];
897
898  Int_t nbinsX = hist1->GetNbinsX();
899  //Int_t nbinsY = hist1->GetNbinsY();
900
901  TH1D *hist_norm = 0;
902  for(Int_t i=0; i<=nbinsX+1; i++) {
903    sprintf(name,"mom_%d",i);
904    TH1D *hist = (TH1D*)hist1->ProjectionY(name,i+1,i+1);
905
906    sprintf(name,"mom_norm");
907    if(i==0) { 
908      hist_norm = (TH1D *)hist->Clone(name);
909      hist_norm->Reset();
910    }
911
912    Double_t nbEvents = hist2->GetBinContent(i);
913    if(!nbEvents) { nbEvents = 1.; };
914
915    hist->Scale(1./nbEvents);
916    hist_norm->Add(hist);
917  }
918
919 return hist_norm;
920 }
921
922 //_____________________________________________________________________________
923 THnSparse* AlidNdPtHelper::GenerateCorrMatrix(THnSparse *hist1, THnSparse *hist2, char *name) {
924 // generate correction matrix
925 if(!hist1 || !hist2) return 0; 
926
927 THnSparse *h =(THnSparse*)hist1->Clone(name);;
928 h->Divide(hist1,hist2,1,1,"B");
929
930 return h;
931 }
932
933 //_____________________________________________________________________________
934 TH2* AlidNdPtHelper::GenerateCorrMatrix(TH2 *hist1, TH2 *hist2, char *name) {
935 // generate correction matrix
936 if(!hist1 || !hist2) return 0; 
937
938 TH2D *h =(TH2D*)hist1->Clone(name);;
939 h->Divide(hist1,hist2,1,1,"B");
940
941 return h;
942 }
943
944 //_____________________________________________________________________________
945 TH1* AlidNdPtHelper::GenerateCorrMatrix(TH1 *hist1, TH1 *hist2, char *name) {
946 // generate correction matrix
947 if(!hist1 || !hist2) return 0; 
948
949 TH1D *h =(TH1D*)hist1->Clone(name);;
950 h->Divide(hist1,hist2,1,1,"B");
951
952 return h;
953 }
954
955 //_____________________________________________________________________________
956 THnSparse* AlidNdPtHelper::GenerateContCorrMatrix(THnSparse *hist1, THnSparse *hist2, char *name) {
957 // generate contamination correction matrix
958 if(!hist1 || !hist2) return 0; 
959
960 THnSparse *hist =  GenerateCorrMatrix(hist1, hist2, name);
961 if(!hist) return 0;
962
963 // only for non ZERO bins!!!!
964
965 Int_t* coord = new Int_t[hist->GetNdimensions()];
966 memset(coord, 0, sizeof(Int_t) * hist->GetNdimensions());
967
968   for (Long64_t i = 0; i < hist->GetNbins(); ++i) {
969     Double_t v = hist->GetBinContent(i, coord);
970     hist->SetBinContent(coord, 1.0-v);
971     //printf("v %f, hist->GetBinContent(i, coord) %f \n",v,hist->GetBinContent(i, coord));
972     Double_t err = hist->GetBinError(coord);
973     hist->SetBinError(coord, err);
974   }
975
976 delete [] coord;
977
978 return hist;
979 }
980
981 //_____________________________________________________________________________
982 TH2* AlidNdPtHelper::GenerateContCorrMatrix(TH2 *hist1, TH2 *hist2, char *name) {
983 // generate contamination correction matrix
984 if(!hist1 || !hist2) return 0; 
985
986 TH2 *hist = GenerateCorrMatrix(hist1, hist2, name);
987 if(!hist) return 0;
988
989 Int_t nBinsX = hist->GetNbinsX();
990 Int_t nBinsY = hist->GetNbinsY();
991
992   for (Int_t i = 0; i < nBinsX+1; i++) {
993   for (Int_t j = 0; j < nBinsY+1; j++) {
994      Double_t cont = hist->GetBinContent(i,j);
995      hist->SetBinContent(i,j,1.-cont);
996      Double_t err = hist->GetBinError(i,j);
997      hist->SetBinError(i,j,err);
998   }
999   }
1000
1001 return hist;
1002 }
1003
1004 //_____________________________________________________________________________
1005 TH1* AlidNdPtHelper::GenerateContCorrMatrix(TH1 *hist1, TH1 *hist2, char *name) {
1006 // generate contamination correction matrix
1007 if(!hist1 || !hist2) return 0; 
1008
1009 TH1 *hist = GenerateCorrMatrix(hist1, hist2, name);
1010 if(!hist) return 0;
1011
1012 Int_t nBinsX = hist->GetNbinsX();
1013
1014   for (Int_t i = 0; i < nBinsX+1; i++) {
1015      Double_t cont = hist->GetBinContent(i);
1016      hist->SetBinContent(i,1.-cont);
1017      Double_t err = hist->GetBinError(i);
1018      hist->SetBinError(i,err);
1019   }
1020
1021 return hist;
1022 }
1023
1024 //_____________________________________________________________________________
1025 const AliESDVertex* AlidNdPtHelper::GetTPCVertexZ(AliESDEvent* esdEvent, AlidNdPtEventCuts *evtCuts, AlidNdPtAcceptanceCuts *accCuts, AliESDtrackCuts *trackCuts, Float_t fraction, Int_t ntracksMin){
1026   //
1027   // TPC Z vertexer
1028   //
1029   if(!esdEvent)
1030   { 
1031     ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
1032     return NULL;  
1033   }
1034
1035   if(!evtCuts || !accCuts || !trackCuts) 
1036   { 
1037     ::Error("AlidNdPtHelper::GetTPCVertexZ()","cuts not available");
1038     return NULL;  
1039   }
1040
1041   Double_t vtxpos[3]={evtCuts->GetMeanXv(),evtCuts->GetMeanYv(),evtCuts->GetMeanZv()};
1042   Double_t vtxsigma[3]={evtCuts->GetSigmaMeanXv(),evtCuts->GetSigmaMeanYv(),evtCuts->GetSigmaMeanZv()};
1043   AliESDVertex vtx0(vtxpos,vtxsigma);
1044
1045   Double_t maxDCAr = accCuts->GetMaxDCAr();
1046   Double_t maxDCAz = accCuts->GetMaxDCAz();
1047   Int_t minTPCClust = trackCuts->GetMinNClusterTPC();
1048
1049   //
1050   Int_t ntracks = esdEvent->GetNumberOfTracks();
1051   TVectorD ztrack(ntracks);
1052   Double_t dca[2],cov[3];
1053   Int_t counter=0;
1054   for (Int_t i=0;i <ntracks; i++){
1055     AliESDtrack *t = esdEvent->GetTrack(i);
1056     if (!t) continue;
1057     if (!t->GetTPCInnerParam()) continue;
1058     if (t->GetTPCNcls()<minTPCClust) continue;
1059     //
1060     AliExternalTrackParam  *tpcTrack  = new AliExternalTrackParam(*(t->GetTPCInnerParam()));
1061     if (!tpcTrack->PropagateToDCA(&vtx0,esdEvent->GetMagneticField(),100.,dca,cov)) continue;
1062
1063     //
1064     if (TMath::Abs(dca[0])>maxDCAr) continue;
1065     //if (TMath::Sqrt(cov[0])>sigmaXYcut) continue;    
1066     if (TMath::Abs(tpcTrack->GetZ())>maxDCAz) continue;
1067
1068     ztrack[counter]=tpcTrack->GetZ();
1069     counter++;    
1070
1071     if(tpcTrack) delete tpcTrack;
1072   }
1073
1074   //
1075   // Find LTM z position
1076   //
1077   Double_t mean=0, sigma=0;
1078   if (counter<ntracksMin) return 0;
1079   //
1080   Int_t nused = TMath::Nint(counter*fraction);
1081   if (nused==counter) nused=counter-1;  
1082   if (nused>1){
1083     AliMathBase::EvaluateUni(counter, ztrack.GetMatrixArray(), mean,sigma, TMath::Nint(counter*fraction));
1084     sigma/=TMath::Sqrt(nused);
1085   }else{
1086     mean  = TMath::Mean(counter, ztrack.GetMatrixArray());
1087     sigma = TMath::RMS(counter, ztrack.GetMatrixArray());
1088     sigma/=TMath::Sqrt(counter-1);
1089   }
1090   vtxpos[2]=mean;
1091   vtxsigma[2]=sigma;
1092   const AliESDVertex* vertex = new AliESDVertex(vtxpos, vtxsigma);
1093   return vertex;
1094 }
1095
1096 //_____________________________________________________________________________
1097 Int_t  AlidNdPtHelper::GetSPDMBTrackMult(AliESDEvent* esdEvent, Float_t deltaThetaCut, Float_t deltaPhiCut) 
1098 {
1099   //
1100   // SPD track multiplicity
1101   //
1102
1103   // get tracklets
1104   const AliMultiplicity* mult = esdEvent->GetMultiplicity();
1105   if (!mult)
1106      return 0;
1107
1108   // get multiplicity from SPD tracklets
1109   Int_t inputCount = 0; 
1110   for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
1111   {
1112     //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
1113
1114      Float_t phi = mult->GetPhi(i);
1115      if (phi < 0)
1116        phi += TMath::Pi() * 2;
1117      Float_t deltaPhi = mult->GetDeltaPhi(i);
1118      Float_t deltaTheta = mult->GetDeltaTheta(i);
1119
1120      if (TMath::Abs(deltaPhi) > 1)
1121        printf("WARNING: Very high Delta Phi: %d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), deltaPhi);
1122
1123      if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
1124         continue;
1125
1126      if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
1127         continue;
1128       
1129      ++inputCount;
1130   }
1131
1132 return inputCount;
1133 }
1134
1135 //_____________________________________________________________________________
1136 Int_t  AlidNdPtHelper::GetSPDMBPrimTrackMult(AliESDEvent* esdEvent, AliStack* stack, Float_t deltaThetaCut, Float_t deltaPhiCut) 
1137 {
1138   //
1139   // SPD track multiplicity
1140   //
1141
1142   // get tracklets
1143   const AliMultiplicity* mult = esdEvent->GetMultiplicity();
1144   if (!mult)
1145      return 0;
1146
1147   // get multiplicity from SPD tracklets
1148   Int_t inputCount = 0; 
1149   for (Int_t i=0; i<mult->GetNumberOfTracklets(); ++i)
1150   {
1151     //printf("%d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), mult->GetDeltaPhi(i));
1152
1153      Float_t phi = mult->GetPhi(i);
1154      if (phi < 0)
1155        phi += TMath::Pi() * 2;
1156      Float_t deltaPhi = mult->GetDeltaPhi(i);
1157      Float_t deltaTheta = mult->GetDeltaTheta(i);
1158
1159      if (TMath::Abs(deltaPhi) > 1)
1160        printf("WARNING: Very high Delta Phi: %d %f %f %f\n", i, mult->GetTheta(i), mult->GetPhi(i), deltaPhi);
1161
1162      if (deltaThetaCut > 0. && TMath::Abs(deltaTheta) > deltaThetaCut)
1163         continue;
1164
1165      if (deltaPhiCut > 0. && TMath::Abs(deltaPhi) > deltaPhiCut)
1166         continue;
1167
1168
1169      if (mult->GetLabel(i, 0) < 0 || mult->GetLabel(i, 0) != mult->GetLabel(i, 1) || 
1170          !stack->IsPhysicalPrimary(mult->GetLabel(i, 0)))
1171         continue;
1172
1173       
1174      ++inputCount;
1175   }
1176
1177 return inputCount;
1178 }
1179
1180