]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGHF/hfe/AliAnalysisTaskElecHadronCorrel.cxx
Compatibility with the Root trunk
[u/mrichter/AliRoot.git] / PWGHF / hfe / AliAnalysisTaskElecHadronCorrel.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
17 ////////////////////////////////////////////////////////////////////////
18 //                                                                    //
19 //  Task for Heavy Flavour Electron-Hadron DeltaPhi Correlation       //
20 //  Non-Photonic Electron identified with Invariant mass              //
21 //  analysis methos in function  SelectPhotonicElectron               //
22 //  DeltaPhi calculated in function  ElectronHadCorrel                // 
23 //                                                                    //
24 //  Author: Deepa Thomas (Utrecht University)                         //
25 //                                                                    //
26 ////////////////////////////////////////////////////////////////////////
27
28 #include "TChain.h"
29 #include "TTree.h"
30 #include "TH2F.h"
31 #include "TMath.h"
32 #include "TCanvas.h"
33 #include "THnSparse.h"
34 #include "TLorentzVector.h"
35 #include "TString.h"
36 #include "TFile.h"
37
38 #include "AliAnalysisTask.h"
39 #include "AliAnalysisManager.h"
40
41 #include "AliESDEvent.h"
42 #include "AliESDHandler.h"
43 #include "AliAODEvent.h"
44 #include "AliAODHandler.h"
45
46 #include "AliAnalysisTaskElecHadronCorrel.h"
47 #include "TGeoGlobalMagField.h"
48 #include "AliLog.h"
49 #include "AliAnalysisTaskSE.h"
50 #include "TRefArray.h"
51 #include "TVector.h"
52 #include "AliESDInputHandler.h"
53 #include "AliESDpid.h"
54 #include "AliESDtrackCuts.h"
55 #include "AliPhysicsSelection.h"
56 #include "AliESDCaloCluster.h"
57 #include "AliAODCaloCluster.h"
58 #include "AliESDCaloTrigger.h"
59 #include "AliEMCALRecoUtils.h"
60 #include "AliEMCALGeometry.h"
61 #include "AliGeomManager.h"
62 #include "stdio.h"
63 #include "TGeoManager.h"
64 #include "iostream"
65 #include "fstream"
66
67 //#include "AliEventPoolManager.h"
68 #include "AliAnalysisTaskPhiCorrelations.h"
69
70 #include "AliCentrality.h"
71 #include "AliEMCALTrack.h"
72 //#include "AliEMCALTracker.h"
73 #include "AliMagF.h"
74
75 #include "AliKFParticle.h"
76 #include "AliKFVertex.h"
77
78 #include "AliPID.h"
79 #include "AliPIDResponse.h"
80 #include "AliHFEcontainer.h"
81 #include "AliHFEcuts.h"
82 #include "AliHFEpid.h"
83 #include "AliHFEpidBase.h"
84 #include "AliHFEpidQAmanager.h"
85 #include "AliHFEtools.h"
86 #include "AliCFContainer.h"
87 #include "AliCFManager.h"
88
89 using std::cout;
90 using std::endl;
91
92 ClassImp(AliAnalysisTaskElecHadronCorrel)
93 //________________________________________________________________________
94   AliAnalysisTaskElecHadronCorrel::AliAnalysisTaskElecHadronCorrel(const char *name) 
95   : AliAnalysisTaskSE(name)
96   ,fESD(0)
97   ,fGeom(0)
98   ,fOutputList(0)
99   ,fTrackCuts1(new AliESDtrackCuts)
100   ,fTrackCuts2(new AliESDtrackCuts)
101   ,fCuts(0)
102   ,fIdentifiedAsOutInz(kFALSE)
103   ,fPassTheEventCut(kFALSE)
104   ,fRejectKinkMother(kFALSE)
105   ,fVz(0.0)
106   ,fCFM(0)      
107   ,fPID(0)
108   ,fPIDqa(0)           
109   ,fOpeningAngleCut(0.1)
110   ,fInvmassCut(0.01)    
111 //  ,fPoolMgr(0x0)  
112   ,fNoEvents(0)
113   ,fTrkpt(0)
114   ,fTrkEovPBef(0)        
115   ,fTrkEovPBefHad(0)     
116   ,fTrkEovPAft(0)       
117   ,fTrkEovPAftOwn(0)    
118   ,fdEdxBef(0)   
119   ,fdEdxAft(0)   
120   ,fdEdxAftOwn(0)        
121   ,fInvmassLS(0)                
122   ,fInvmassULS(0)               
123   ,fOpeningAngleLS(0)   
124   ,fOpeningAngleULS(0)  
125   ,fSemiIncElecDphi(0)  
126   ,fPhotElecDphi(0)     
127   ,fInclusiveElecDphi(0)        
128   ,fDphiULSMassLow(0)   
129   ,fDphiLSMassLow(0)
130   ,fDphiULSMassLowNoPartner(0)   
131   ,fDphiLSMassLowNoPartner(0)
132   ,fPhotoElecPt(0)
133   ,fSemiInclElecPt(0)
134   ,fInclusiveElecPt(0)
135   ,fULSElecPt(0)
136   ,fLSElecPt(0)  
137   ,fTrackPtBefTrkCuts(0)         
138   ,fTrackPtAftTrkCuts(0)
139   ,fTPCnsigma(0)
140   ,fTPCnsigmaAft(0)
141   ,fTPCnsigmaAftOwn(0)
142   ,fNCellv1(0)
143   ,fClsEv1(0)
144   ,fNClusv1(0)
145   ,fKFParticleP(0)
146   ,fKFParticleE(0)
147   ,fInvmassLS1(0)       
148   ,fInvmassULS1(0)
149   ,fcentrality(0)     
150   ,fElecPhi(0)  
151   ,fElecPhiTPC(0)  
152   ,fElecPhiTPCEovP(0)  
153   ,fHadronPhi(0)  
154   ,fTrackHFEcuts(0)
155   ,fTrakPhiSPD1(0)
156   ,fTrakPhiSPD2(0)
157   ,fTrakPhiSPDOr(0)
158   ,fTrakPhiSPDAnd(0)
159   ,fTrackHFEcutsITS(0)  
160 /*  ,fNoMixedEvents(0)
161   ,fMixStat(0)       
162   ,fMixStat1(0)        
163   ,fMixedIncElecDphi(0)  
164   ,fMixedPhotElecDphi(0)
165   ,fMixedSemiIncElecDphi(0)  
166   ,fMixedDphiULSMassLow(0)  
167   ,fMixedDphiLSMassLow(0)  
168 */
169 //  ,fSparseElectron(0)  
170 //  ,fvalueElectron(0)   
171 {
172   //Named constructor
173
174   fPID = new AliHFEpid("hfePid");
175 //  fvalueElectron = new Double_t[8];
176
177   // Define input and output slots here
178   // Input slot #0 works with a TChain
179   DefineInput(0, TChain::Class());
180   // Output slot #0 id reserved by the base class for AOD
181   // Output slot #1 writes into a TH1 container
182   // DefineOutput(1, TH1I::Class());
183   DefineOutput(1, TList::Class());
184   //  DefineOutput(3, TTree::Class());
185 }
186
187 //________________________________________________________________________
188 AliAnalysisTaskElecHadronCorrel::AliAnalysisTaskElecHadronCorrel() 
189   : AliAnalysisTaskSE("DefaultAnalysis_AliAnalysisElecHadCorrel")
190   ,fESD(0)
191   ,fGeom(0)  
192   ,fOutputList(0)
193   ,fTrackCuts1(new AliESDtrackCuts)
194   ,fTrackCuts2(new AliESDtrackCuts)
195   ,fCuts(0)
196   ,fIdentifiedAsOutInz(kFALSE)
197   ,fPassTheEventCut(kFALSE)
198   ,fRejectKinkMother(kFALSE)
199   ,fVz(0.0)
200   ,fCFM(0)      
201   ,fPID(0)       
202   ,fPIDqa(0)           
203   ,fOpeningAngleCut(0.1)
204   ,fInvmassCut(0.01)    
205 //  ,fPoolMgr(0x0)    
206   ,fNoEvents(0)
207   ,fTrkpt(0)
208   ,fTrkEovPBef(0)        
209   ,fTrkEovPBefHad(0)     
210   ,fTrkEovPAft(0)        
211   ,fTrkEovPAftOwn(0)     
212   ,fdEdxBef(0)   
213   ,fdEdxAft(0)   
214   ,fdEdxAftOwn(0)        
215   ,fInvmassLS(0)                
216   ,fInvmassULS(0)               
217   ,fOpeningAngleLS(0)   
218   ,fOpeningAngleULS(0)  
219   ,fSemiIncElecDphi(0)  
220   ,fPhotElecDphi(0)     
221   ,fInclusiveElecDphi(0)        
222   ,fDphiULSMassLow(0)   
223   ,fDphiLSMassLow(0)
224   ,fDphiULSMassLowNoPartner(0)   
225   ,fDphiLSMassLowNoPartner(0)
226   ,fPhotoElecPt(0)
227   ,fSemiInclElecPt(0)
228   ,fInclusiveElecPt(0)
229   ,fULSElecPt(0)
230   ,fLSElecPt(0)  
231   ,fTrackPtBefTrkCuts(0)         
232   ,fTrackPtAftTrkCuts(0)                  
233   ,fTPCnsigma(0)        
234   ,fTPCnsigmaAft(0)     
235   ,fTPCnsigmaAftOwn(0)  
236   ,fNCellv1(0)  
237   ,fClsEv1(0)
238   ,fNClusv1(0)
239   ,fKFParticleP(0)
240   ,fKFParticleE(0)
241   ,fInvmassLS1(0)   
242   ,fInvmassULS1(0)  
243   ,fcentrality(0)     
244   ,fElecPhi(0)
245   ,fElecPhiTPC(0)
246   ,fElecPhiTPCEovP(0)  
247   ,fHadronPhi(0)
248   ,fTrackHFEcuts(0)
249   ,fTrakPhiSPD1(0)
250   ,fTrakPhiSPD2(0)
251   ,fTrakPhiSPDOr(0)
252   ,fTrakPhiSPDAnd(0)
253   ,fTrackHFEcutsITS(0)  
254 /*  ,fNoMixedEvents(0)
255   ,fMixStat(0)      
256   ,fMixStat1(0)     
257   ,fMixedIncElecDphi(0)  
258   ,fMixedPhotElecDphi(0)
259   ,fMixedSemiIncElecDphi(0)
260   ,fMixedDphiULSMassLow(0) 
261   ,fMixedDphiLSMassLow(0)      
262 */
263     //  ,fSparseElectron(0)  
264     //    ,fvalueElectron(0)  
265 {
266   //Default constructor
267   fPID = new AliHFEpid("hfePid");
268   //  fvalueElectron = new Double_t[8];
269
270   // Constructor
271   // Define input and output slots here
272   // Input slot #0 works with a TChain
273   DefineInput(0, TChain::Class());
274   // Output slot #0 id reserved by the base class for AOD
275   // Output slot #1 writes into a TH1 container
276   // DefineOutput(1, TH1I::Class());
277   DefineOutput(1, TList::Class());
278   //DefineOutput(3, TTree::Class());
279 }
280 //_________________________________________
281
282 AliAnalysisTaskElecHadronCorrel::~AliAnalysisTaskElecHadronCorrel()
283 {
284   //Destructor 
285
286   delete fOutputList;
287   delete fGeom;
288   delete fPID;
289   delete fCFM;
290   delete fPIDqa;
291   delete fTrackCuts1;
292   delete fTrackCuts2;
293   //   delete fSparseElectron;
294   //   delete []fvalueElectron;
295 }
296 //_________________________________________
297
298 void AliAnalysisTaskElecHadronCorrel::UserExec(Option_t*)
299 {
300   //Main loop
301   //Called for each event
302
303   // create pointer to event
304   fESD = dynamic_cast<AliESDEvent*>(InputEvent());
305   if (!fESD) {
306     printf("ERROR: fESD not available\n");
307     return;
308   }
309
310   if(!fCuts){
311     AliError("HFE cuts not available");
312     return;
313   }
314
315   if(!fPID->IsInitialized()){ 
316     // Initialize PID with the given run number
317     AliWarning("PID not initialised, get from Run no");
318     fPID->InitializePID(fESD->GetRunNumber());
319   }
320
321   //-------trigger selection
322   UInt_t res = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
323   if (res==0)
324     return;
325
326   //    if( (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & AliVEvent::kFastOnly) )
327   //            return;
328
329   if(!(((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & (AliVEvent::kCentral))) return;
330
331   AliCentrality *fCentrality = (AliCentrality*)fESD->GetCentrality();
332
333   Float_t centvalue = fCentrality->GetCentralityPercentile("V0M");
334   fcentrality->Fill(centvalue);    
335   cout << "cent val" << centvalue <<endl;
336   if(centvalue<0 || centvalue>10) return;
337
338   cout << "event no : " <<fESD->GetRunNumber() <<endl;
339   Int_t fNOtrks =  fESD->GetNumberOfTracks();
340   const AliESDVertex *pVtx = fESD->GetPrimaryVertex();
341
342   Double_t pVtxZ = -999;
343   pVtxZ = pVtx->GetZ();
344
345   // Event cut
346   //    if(!fCFM->CheckEventCuts(AliHFEcuts::kEventStepReconstructed, fESD)) return;
347
348   if(TMath::Abs(pVtxZ)>10) return;
349   fNoEvents->Fill(0);
350
351   if(fNOtrks<2) return;
352
353   AliPIDResponse *pidResponse = fInputHandler->GetPIDResponse();
354   if(!pidResponse){
355     AliDebug(1, "Using default PID Response");
356      pidResponse = AliHFEtools::GetDefaultPID(kFALSE, fInputEvent->IsA() == AliAODEvent::Class()); 
357    }
358
359    fPID->SetPIDResponse(pidResponse);
360
361    fCFM->SetRecEventInfo(fESD);
362 /*
363    //Event mixing
364    AliEventPool* pool = fPoolMgr->GetEventPool(centvalue, pVtxZ); // Get the buffer associated with the current centrality and z-vtx
365    if (!pool)
366      AliFatal(Form("No pool found for centrality = %f, zVtx = %f", centvalue, pVtxZ));
367
368    TObjArray* tracksClone = CloneAndReduceTrackList();
369    tracksClone->SetOwner();
370    pool->UpdatePool(tracksClone);
371 */
372    // Track loop 
373    for (Int_t iTracks = 0; iTracks < fESD->GetNumberOfTracks(); iTracks++) {
374      AliESDtrack* track = fESD->GetTrack(iTracks);
375      if (!track) {
376        printf("ERROR: Could not receive track %d\n", iTracks);
377        continue;
378      }
379
380      if(track->Pt()<1) continue;
381
382      fTrackPtBefTrkCuts->Fill(track->Pt());             
383
384      // RecKine: ITSTPC cuts  
385      if(!ProcessCutStep(AliHFEcuts::kStepRecKineITSTPC, track)) continue;
386
387      //RecKink
388      if(fRejectKinkMother) { // Quick and dirty fix to reject both kink mothers and daughters
389        if(track->GetKinkIndex(0) != 0) continue;
390      } 
391
392      // RecPrim
393      if(!ProcessCutStep(AliHFEcuts::kStepRecPrim, track)) continue;
394
395      // HFE cuts: TPC PID cleanup
396      if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsTPC, track)) continue;
397
398      fTrackHFEcuts->Fill(track->Phi());
399
400      //track phi distribution for diff ITS layer hit
401      if(track->HasPointOnITSLayer(0)) fTrakPhiSPD1->Fill(track->Phi());
402      if(track->HasPointOnITSLayer(1)) fTrakPhiSPD2->Fill(track->Phi());
403
404      if(track->HasPointOnITSLayer(0) || track->HasPointOnITSLayer(1)) fTrakPhiSPDOr->Fill(track->Phi());
405      if(track->HasPointOnITSLayer(0) && track->HasPointOnITSLayer(1)) fTrakPhiSPDAnd->Fill(track->Phi());
406
407      // HFEcuts: ITS layers cuts
408      if(!ProcessCutStep(AliHFEcuts::kStepHFEcutsITS, track)) continue;
409
410      fTrackHFEcutsITS->Fill(track->Phi());
411
412      fTrackPtAftTrkCuts->Fill(track->Pt());             
413
414      Double_t fClsE = -999, p = -999, fEovP=-999, pt = -999, dEdx=-999, fTPCnSigma=0;
415      pt = track->Pt();
416      p = track->P();
417      dEdx = track->GetTPCsignal();
418      fTPCnSigma = fPID->GetPIDResponse() ? fPID->GetPIDResponse()->NumberOfSigmasTPC(track, AliPID::kElectron) : 1000;
419
420      //TPC electron phi
421      if(fTPCnSigma >= -2 && fTPCnSigma <= 2){
422        fElecPhiTPC->Fill(track->Phi());
423      }
424
425      //eta cut (-0.7,0.7)
426      if(track->Eta() < -0.7 || track->Eta() > 0.7) continue;
427
428      // Track extrapolation to EMCAL
429      Int_t fClsId = track->GetEMCALcluster();
430      if(fClsId <0) continue;
431      AliESDCaloCluster *cluster = fESD->GetCaloCluster(fClsId);
432      if(TMath::Abs(cluster->GetTrackDx())>0.05 || TMath::Abs(cluster->GetTrackDz())>0.05) continue;    
433      fdEdxBef->Fill(p,dEdx);
434      fTPCnsigma->Fill(p,fTPCnSigma);
435
436      fTrkpt->Fill(pt);
437      fClsE = cluster->E();
438      cout << "cluster E = " << fClsE <<endl;
439      fEovP = fClsE/p;
440      /*
441         fvalueElectron[0] = pt;
442         fvalueElectron[1] = p;
443         fvalueElectron[2] = fTPCnSigma;
444         fvalueElectron[3] = dEdx;
445         fvalueElectron[4] = fEovP;
446         fvalueElectron[5] = cluster->GetM20();
447         fvalueElectron[6] = cluster->GetM02();
448         fvalueElectron[7] = cluster->GetDispersion();
449
450         fSparseElectron->Fill(fvalueElectron);
451       */
452      if(fTPCnSigma >= -2 && fTPCnSigma <= 2){
453        fTrkEovPBef->Fill(pt,fEovP);
454      }
455      if(fTPCnSigma < -3.5)fTrkEovPBefHad->Fill(pt,fEovP);
456      /*
457         Int_t pidpassed = 0;
458      //--- track accepted, do PID
459      AliHFEpidObject hfetrack;
460      hfetrack.SetAnalysisType(AliHFEpidObject::kESDanalysis);
461      hfetrack.SetRecTrack(track);
462      hfetrack.SetPbPb();
463      if(fPID->IsSelected(&hfetrack, NULL, "", fPIDqa)) pidpassed = 1;
464
465      if(pidpassed==1){
466      cout << "dedx, E/p :  "<< dEdx << ", " << fEovP <<endl;
467      fTrkEovPAft->Fill(pt,fEovP);
468      fdEdxAft->Fill(p,dEdx);
469      fTPCnsigmaAft->Fill(p,fTPCnSigma);
470      }
471       */
472
473      //Electron id with TPC and E/p
474      if(fTPCnSigma >= -2 && fTPCnSigma <= 2 && fEovP >= 0.8 && fEovP <=1.2) {
475        fElecPhiTPCEovP->Fill(track->Phi());
476
477        //Electron id with shower shape  
478        if(cluster->GetM20()<0.2 && cluster->GetM02()< 0.5 && cluster->GetDispersion()<1){
479          fElecPhi->Fill(track->Phi());
480          fTrkEovPAftOwn->Fill(pt,fEovP);
481          fdEdxAftOwn->Fill(p,dEdx);
482          fTPCnsigmaAftOwn->Fill(p,fTPCnSigma);
483
484          Bool_t fFlagPhotonicElec = kFALSE;
485          // select photonic electron
486          SelectPhotonicElectron(iTracks,track,fFlagPhotonicElec);
487          //Inclusive electron-hadron correlation
488          ElectronHadCorrel(iTracks, track, fInclusiveElecDphi);
489          fInclusiveElecPt->Fill(pt);
490         // MixedEvent(track,fMixedIncElecDphi);
491          
492          // photonic electron
493          if(fFlagPhotonicElec){
494          //Electron hadron correlation
495          ElectronHadCorrel(iTracks, track, fPhotElecDphi);
496          fPhotoElecPt->Fill(pt);
497         // MixedEvent(track,fMixedPhotElecDphi);
498          }
499
500          // Semi inclusive electron 
501          if(!fFlagPhotonicElec){
502            //Electron hadron correlation
503            ElectronHadCorrel(iTracks, track, fSemiIncElecDphi);
504            fSemiInclElecPt->Fill(pt);
505           // MixedEvent(track,fMixedSemiIncElecDphi);
506          }
507          
508        }
509      }
510    }
511
512    //EMC clusters  
513    Int_t clsNo = fESD->GetNumberOfCaloClusters();
514    fNClusv1->Fill(clsNo); 
515    for(Int_t iclus=0; iclus<clsNo ; iclus++){ 
516      AliESDCaloCluster* clus = fESD->GetCaloCluster(iclus);
517      if(!clus->IsEMCAL()) continue; 
518      fNCellv1->Fill(clus->GetNCells());
519      fClsEv1->Fill(clus->E());  
520    }
521
522
523    PostData(1, fOutputList);
524 }
525 //_________________________________________
526 void AliAnalysisTaskElecHadronCorrel::UserCreateOutputObjects()
527 {
528   //Create histograms
529   //  TGeoManager::Import("geometry.root");
530   //  fGeom =  AliEMCALGeometry::GetInstance("EMCAL_COMPLETEV1");
531
532   //--------Initialize PID
533   fPID->SetHasMCData(kFALSE);
534   if(!fPID->GetNumberOfPIDdetectors()) 
535   {
536     fPID->AddDetector("TPC", 0);
537     fPID->AddDetector("EMCAL", 1);
538   }
539
540   fPID->SortDetectors(); 
541   fPIDqa = new AliHFEpidQAmanager();
542   fPIDqa->Initialize(fPID);
543
544   //--------Initialize correction Framework and Cuts
545   fCFM = new AliCFManager;
546   const Int_t kNcutSteps = AliHFEcuts::kNcutStepsMCTrack + AliHFEcuts::kNcutStepsRecTrack + AliHFEcuts::kNcutStepsDETrack;
547   fCFM->SetNStepParticle(kNcutSteps);
548   for(Int_t istep = 0; istep < kNcutSteps; istep++)
549     fCFM->SetParticleCutsList(istep, NULL);
550
551   if(!fCuts){
552     AliWarning("Cuts not available. Default cuts will be used");
553     fCuts = new AliHFEcuts;
554     fCuts->CreateStandardCuts();
555   }
556   fCuts->Initialize(fCFM);
557 /*
558   //Mixed event initialising
559   Int_t trackDepth = 2000;
560   Int_t poolsize   = 1000;
561
562   Int_t nCentralityBins  = 5;
563   Double_t CentralityBins[] = {0,2,4,6,8,10};
564
565   Int_t nZvtxBins  = 2;
566   Double_t vertexBins[] = {-10,0,10};
567
568   fPoolMgr = new AliEventPoolManager(poolsize, trackDepth, nCentralityBins, (Double_t*) CentralityBins, nZvtxBins, (Double_t*) vertexBins);
569 */
570   //---------Output Tlist
571   fOutputList = new TList();
572   fOutputList->SetOwner();
573   fOutputList->Add(fPIDqa->MakeList("PIDQA"));
574
575   fNoEvents = new TH1F("fNoEvents","",1,0,1) ;
576   fOutputList->Add(fNoEvents);
577
578   fcentrality = new TH1F("fcentrality","centrality", 100,0,100);
579   fOutputList->Add(fcentrality);
580
581   fTrkpt = new TH1F("fTrkpt","track pt",1000,0,50);
582   fOutputList->Add(fTrkpt);
583
584   fTrackPtBefTrkCuts = new TH1F("fTrackPtBefTrkCuts","track pt before track cuts",1000,0,50);
585   fOutputList->Add(fTrackPtBefTrkCuts);
586
587   fTrackPtAftTrkCuts = new TH1F("fTrackPtAftTrkCuts","track pt after track cuts",1000,0,50);
588   fOutputList->Add(fTrackPtAftTrkCuts);
589
590   fTPCnsigma = new TH2F("fTPCnsigma", "TPC - n sigma",1000,0,50,200,-10,10);
591   fOutputList->Add(fTPCnsigma);
592
593   fTPCnsigmaAft = new TH2F("fTPCnsigmaAft", "TPC - n sigma after hfepid",1000,0,50,200,-10,10);
594   fOutputList->Add(fTPCnsigmaAft);
595
596   fTPCnsigmaAftOwn = new TH2F("fTPCnsigmaAftOwn", "TPC - n sigma after own pid",1000,0,50,200,-10,10);
597   fOutputList->Add(fTPCnsigmaAftOwn);
598
599   fTrkEovPBef = new TH2F("fTrkEovPBef","track E/p before HFE pid",1000,0,50,100,0,2);
600   fOutputList->Add(fTrkEovPBef);
601
602   fTrkEovPBefHad = new TH2F("fTrkEovPBefHad","track E/p for TPCnsig < 3.5",1000,0,50,100,0,2);
603   fOutputList->Add(fTrkEovPBefHad);
604
605   fTrkEovPAft = new TH2F("fTrkEovPAft","track E/p after HFE pid",1000,0,50,100,0,2);
606   fOutputList->Add(fTrkEovPAft);
607
608   fTrkEovPAftOwn = new TH2F("fTrkEovPAftOwn","track E/p after own pid",1000,0,50,100,0,2);
609   fOutputList->Add(fTrkEovPAftOwn);
610
611   fdEdxBef = new TH2F("fdEdxBef","track dEdx vs p before HFE pid",1000,0,50,150,0,150);
612   fOutputList->Add(fdEdxBef);
613
614   fdEdxAft = new TH2F("fdEdxAft","track dEdx vs p after HFE pid",1000,0,50,150,0,150);
615   fOutputList->Add(fdEdxAft);
616
617   fdEdxAftOwn = new TH2F("fdEdxAftOwn","track dEdx vs p own HFE pid",1000,0,50,150,0,150);
618   fOutputList->Add(fdEdxAftOwn);
619
620   fElecPhi = new TH1F("fElecPhi", "Electron phi",1000,0,6.28);
621   fOutputList->Add(fElecPhi);
622
623   fElecPhiTPC = new TH1F("fElecPhiTPC", "Electron phi after TPC cut",1000,0,6.28);
624   fOutputList->Add(fElecPhiTPC);
625
626   fElecPhiTPCEovP = new TH1F("fElecPhiTPCEovP", "Electron phi after TPC and E/p cut",1000,0,6.28);
627   fOutputList->Add(fElecPhiTPCEovP);
628
629   fHadronPhi = new TH1F("fHadronPhi", "Hadron phi",1000,0,6.28);
630   fOutputList->Add(fHadronPhi);
631
632   fTrackHFEcuts = new TH1F("fTrackHFEcuts","Track phi for HFE cuts",1000,0,6.28);
633   fOutputList->Add(fTrackHFEcuts);
634
635   fTrakPhiSPD1 = new TH1F("fTrakPhiSPD1","Track phi for hit in SPD layer 1",1000,0,6.28);
636   fOutputList->Add(fTrakPhiSPD1);
637
638   fTrakPhiSPD2 = new TH1F("fTrakPhiSPD2","Track phi for hit in SPD layer 2",1000,0,6.28);
639   fOutputList->Add(fTrakPhiSPD2);
640
641   fTrakPhiSPDOr = new TH1F("fTrakPhiSPDOr","Track phi for hit in any SPD layer",1000,0,6.28);
642   fOutputList->Add(fTrakPhiSPDOr);
643
644   fTrakPhiSPDAnd = new TH1F("fTrakPhiSPDAnd","Track phi for hit in both SPD layer",1000,0,6.28);
645   fOutputList->Add(fTrakPhiSPDAnd);
646
647   fTrackHFEcutsITS = new TH1F("fTrackHFEcutsITS","Track phi for HFE cuts + ITS HFE cuts",1000,0,6.28);
648   fOutputList->Add(fTrackHFEcutsITS);
649
650   //  fInvmassLS = new TH1F("fInvmassLS", "Inv mass of LS (e,e) if mass cal is correct; mass(GeV/c^2); counts;", 1000,0,1.0);
651   //  fOutputList->Add(fInvmassLS);
652
653   //  fInvmassULS = new TH1F("fInvmassULS", "Inv mass of ULS (e,e); mass(GeV/c^2) if mass cal is correct; counts;", 1000,0,1.0);
654   //  fOutputList->Add(fInvmassULS);
655
656   fOpeningAngleLS = new TH1F("fOpeningAngleLS","Opening angle for LS pairs",100,0,1);
657   fOutputList->Add(fOpeningAngleLS);
658
659   fOpeningAngleULS = new TH1F("fOpeningAngleULS","Opening angle for ULS pairs",100,0,1);
660   fOutputList->Add(fOpeningAngleULS);
661
662   fSemiIncElecDphi = new TH2F("fSemiIncElecDphi", "Semi Inclusive elec-had Dphi correlation",200,0,20,100,-1.6,4.75);
663   fOutputList->Add(fSemiIncElecDphi);
664
665   fPhotElecDphi = new TH2F("fPhotElecDphi", "Photon elec-had Dphi correlation",200,0,20,100,-1.6,4.75);
666   fOutputList->Add(fPhotElecDphi);
667
668   fInclusiveElecDphi = new TH2F("fInclusiveElecDphi", "Inclusive elec-had Dphi correlation",200,0,20,100,-1.6,4.75);
669   fOutputList->Add(fInclusiveElecDphi);
670
671   fDphiULSMassLow = new TH2F("fDphiULSMassLow", "e-h Dphi ULS, mass<0.01",200,0,20,100,-1.6,4.75);
672   fOutputList->Add(fDphiULSMassLow);
673
674   fDphiLSMassLow = new TH2F("fDphiLSMassLow", "e-h Dphi LS, mass<0.01",200,0,20,100,-1.6,4.75);
675   fOutputList->Add(fDphiLSMassLow);
676
677   fDphiULSMassLowNoPartner = new TH2F("fDphiULSMassLowNoPartner", "e-h Dphi ULS with no partner, mass<mass cut,",200,0,20,100,-1.6,4.75);
678   fOutputList->Add(fDphiULSMassLowNoPartner);
679
680   fDphiLSMassLowNoPartner = new TH2F("fDphiLSMassLowNoPartner", "e-h Dphi LS with no partner, mass<mass cut",200,0,20,100,-1.6,4.75);
681   fOutputList->Add(fDphiLSMassLowNoPartner);
682
683   fPhotoElecPt = new TH1F("fPhotoElecPt", "photonic electron pt",1000,0,100);
684   fOutputList->Add(fPhotoElecPt);
685
686   fSemiInclElecPt = new TH1F("fSemiInclElecPt", "Semi-inclusive electron pt",1000,0,100);
687   fOutputList->Add(fSemiInclElecPt);
688
689   fInclusiveElecPt = new TH1F("fInclElecPt", "Inclusive electron pt",1000,0,100);
690   fOutputList->Add(fInclusiveElecPt);
691
692   fULSElecPt = new TH1F("fULSElecPt", "ULS electron pt",1000,0,100);
693   fOutputList->Add(fULSElecPt);
694
695   fLSElecPt = new TH1F("fLSElecPt", "LS electron pt",1000,0,100);
696   fOutputList->Add(fLSElecPt);
697
698   fNCellv1 = new TH1F("fNCellv1","Ncell in clus (v1); NCell; count",100,0,100) ;
699   fOutputList->Add(fNCellv1);
700
701   fClsEv1 = new TH1F("fClsEv1", "Clus E(v1); Cls E; count",1000,0,100); 
702   fOutputList->Add(fClsEv1); 
703
704   fNClusv1 = new TH1F("fNClusv1","Nclus in event (v1); NClus; count",500,0,500) ; 
705   fOutputList->Add(fNClusv1);
706
707   fKFParticleP = new TH1F("fKFParticleP","KFparticle rec P; P(GeV/c)",1000,0,50);
708   fOutputList->Add(fKFParticleP);
709
710   fKFParticleE = new TH1F("fKFParticleE", "KfParticle rec E; E; count",1000,0,100); 
711   fOutputList->Add(fKFParticleE);
712
713   fInvmassLS1 = new TH1F("fInvmassLS1", "Inv mass of LS (e,e); mass(GeV/c^2); counts;", 1000,0,1.0);
714   fOutputList->Add(fInvmassLS1);
715
716   fInvmassULS1 = new TH1F("fInvmassULS1", "Inv mass of ULS (e,e); mass(GeV/c^2); counts;", 1000,0,1.0);
717   fOutputList->Add(fInvmassULS1);
718 /*
719   fNoMixedEvents = new TH1F("fNoMixedEvents","",1,0,1) ;
720   fOutputList->Add(fNoMixedEvents);
721
722   fMixStat = new TH2F("fMixStat","no of events in pool  vs multiplicity;Nevent in pool;N hadrons",200,0,200,500,0,500);
723   fOutputList->Add(fMixStat);                                                             
724
725   fMixStat1 = new TH2F("fMixStat1","no of events in pool  vs zvtx;Nevents in pool;zvtx",200,0,200,200,-11,11);
726   fOutputList->Add(fMixStat1);
727
728   fMixedIncElecDphi = new TH2F("fMixedIncElecDphi", "Mixed event - Inclusive elec-had Dphi correlation",200,0,20,100,-1.6,4.75);
729   fOutputList->Add(fMixedIncElecDphi);
730
731   fMixedSemiIncElecDphi = new TH2F("fMixedSemiIncElecDphi", "Mixed event - Semi Inclusive elec-had Dphi correlation",200,0,20,100,-1.6,4.75);
732   fOutputList->Add(fMixedSemiIncElecDphi);
733
734   fMixedPhotElecDphi = new TH2F("fMixedPhotElecDphi", "Mixed event - Photo elec-had Dphi correlation",200,0,20,100,-1.6,4.75);
735   fOutputList->Add(fMixedPhotElecDphi);
736
737   fMixedDphiULSMassLow = new TH2F("fMixedDphiULSMassLow", "Mixed event - ULS mass < cut elec-had Dphi correlation",200,0,20,100,-1.6,4.75);
738   fOutputList->Add(fMixedDphiULSMassLow);
739
740   fMixedDphiLSMassLow = new TH2F("fMixedDphiLSMassLow", "Mixed event - LS mass < cut elec-had Dphi correlation",200,0,20,100,-1.6,4.75);
741   fOutputList->Add(fMixedDphiLSMassLow);
742  */
743   /*
744      Int_t binsv1[8]={1000,1000,200,150,100,100,100,100}; //pt, p, TPCnsig, dEdx, E/p, M20, M02, dispersion 
745      Double_t xminv1[8]={0,0,-10,0,0,0,0,0};
746      Double_t xmaxv1[8]={50,50,10,150,2,2,2,2};
747      fSparseElectron = new THnSparseD ("Electron","Electron",8,binsv1,xminv1,xmaxv1);
748      fOutputList->Add(fSparseElectron);
749    */
750   PostData(1,fOutputList);
751 }
752
753 //________________________________________________________________________
754 void AliAnalysisTaskElecHadronCorrel::Terminate(Option_t *)
755 {
756   // Info("Terminate");
757   AliAnalysisTaskSE::Terminate();
758 }
759
760 //________________________________________________________________________
761 Bool_t AliAnalysisTaskElecHadronCorrel::ProcessCutStep(Int_t cutStep, AliVParticle *track)
762 {
763   // Check single track cuts for a given cut step
764   const Int_t kMCOffset = AliHFEcuts::kNcutStepsMCTrack;
765   if(!fCFM->CheckParticleCuts(cutStep + kMCOffset, track)) return kFALSE;
766   return kTRUE;
767 }
768 //_________________________________________
769 void AliAnalysisTaskElecHadronCorrel::SelectPhotonicElectron(Int_t itrack, AliESDtrack *track, Bool_t &fFlagPhotonicElec)
770 {
771   //Identify non-heavy flavour electrons using Invariant mass method
772
773   fTrackCuts1->SetAcceptKinkDaughters(kFALSE);
774   fTrackCuts1->SetRequireTPCRefit(kTRUE);
775   fTrackCuts1->SetEtaRange(-0.9,0.9);
776   fTrackCuts1->SetRequireSigmaToVertex(kTRUE);
777   fTrackCuts1->SetMaxChi2PerClusterTPC(3.5);
778   fTrackCuts1->SetMinNClustersTPC(80);
779
780   //  const AliESDVertex *pVtx = fESD->GetPrimaryVertex();
781
782   Bool_t flagPhotonicElec = kFALSE;
783
784   for(Int_t jTracks = itrack+1; jTracks<fESD->GetNumberOfTracks(); jTracks++){
785     AliESDtrack* trackAsso = fESD->GetTrack(jTracks);
786     if (!trackAsso) {
787       printf("ERROR: Could not receive track %d\n", jTracks);
788       continue;
789     }
790
791     Double_t dEdxAsso = -999., ptAsso=-999., openingAngle = -999.;
792     Double_t mass=-999., width = -999;
793     Bool_t fFlagLS=kFALSE, fFlagULS=kFALSE;
794
795     dEdxAsso = trackAsso->GetTPCsignal();
796     ptAsso = trackAsso->Pt();
797     Int_t chargeAsso = trackAsso->Charge();
798     Int_t charge = track->Charge();
799
800     if(ptAsso <0.3) continue;
801     if(!fTrackCuts1->AcceptTrack(trackAsso)) continue;
802     if(dEdxAsso <70 || dEdxAsso>100) continue; //11a pass1
803
804     Int_t fPDGe1 = 11; Int_t fPDGe2 = 11;
805     if(charge>0) fPDGe1 = -11;
806     if(chargeAsso>0) fPDGe2 = -11;
807
808     if(charge == chargeAsso) fFlagLS = kTRUE;
809     if(charge != chargeAsso) fFlagULS = kTRUE;
810
811     AliKFParticle ge1 = AliKFParticle(*track, fPDGe1);
812     AliKFParticle ge2 = AliKFParticle(*trackAsso, fPDGe2);
813     AliKFParticle recg(ge1, ge2);
814     /*
815        Double_t recP2=-999.0, recP=-999.0, recE=-999.0, m2=-999.0, m=-999.0;
816        recP2= (recg.GetPx()*recg.GetPx() + recg.GetPy()*recg.GetPy() + recg.GetPz()*recg.GetPz());
817        recP = TMath::Sqrt(recP2);
818        recE = recg.GetE();
819        fKFParticleP->Fill(recP);
820        fKFParticleE->Fill(recE);
821        m2 = (recg.GetE()*recg.GetE() - recg.GetPx()*recg.GetPx() - recg.GetPy()*recg.GetPy() - recg.GetPz()*recg.GetPz());
822        m = TMath::Sqrt(m2);
823      */
824     if(recg.GetNDF()<1) continue;
825     Double_t chi2recg = recg.GetChi2()/recg.GetNDF();
826     if(TMath::Sqrt(TMath::Abs(chi2recg))>3.) continue;
827
828     openingAngle = ge1.GetAngle(ge2);
829     if(fFlagLS) fOpeningAngleLS->Fill(openingAngle);
830     if(fFlagULS) fOpeningAngleULS->Fill(openingAngle);
831
832  //   if(openingAngle > fOpeningAngleCut) continue;
833
834     Int_t MassCorrect;
835     MassCorrect = recg.GetMass(mass,width);
836     cout << "mass = " << mass <<endl;
837
838     if(fFlagLS) fInvmassLS1->Fill(mass);
839     if(fFlagULS) fInvmassULS1->Fill(mass);
840
841     //    if(MassCorrect==1){
842 //    if(fFlagLS) fInvmassLS->Fill(m);
843 //    if(fFlagULS) fInvmassULS->Fill(m);
844
845     if(mass<fInvmassCut){
846       if(fFlagULS)
847       {
848         ElectronHadCorrel(itrack,track,fDphiULSMassLow);
849         fULSElecPt->Fill(track->Pt());
850        // MixedEvent(track,fMixedDphiULSMassLow);
851       }
852       if(fFlagLS)
853       {
854         ElectronHadCorrel(itrack,track,fDphiLSMassLow);
855         fLSElecPt->Fill(track->Pt());
856        // MixedEvent(track,fMixedDphiLSMassLow);
857       }
858       if(fFlagLS) ElectronHadCorrelNoPartner(itrack,jTracks,track,fDphiLSMassLowNoPartner);
859       if(fFlagULS) ElectronHadCorrelNoPartner(itrack,jTracks,track,fDphiULSMassLowNoPartner);
860     }
861
862     if(mass<fInvmassCut && fFlagULS && !flagPhotonicElec){
863       flagPhotonicElec = kTRUE;
864     }
865     //   }
866
867   }
868   fFlagPhotonicElec = flagPhotonicElec;
869
870 }
871 //_________________________________________
872 void AliAnalysisTaskElecHadronCorrel::ElectronHadCorrel(Int_t itrack, AliESDtrack *track, TH2F *DphiPt)
873 {
874   //Construct Delta Phi between electrons and hadrons
875
876   fTrackCuts2->SetAcceptKinkDaughters(kFALSE);
877   fTrackCuts2->SetRequireTPCRefit(kTRUE);
878   fTrackCuts2->SetRequireITSRefit(kTRUE);
879   fTrackCuts2->SetEtaRange(-0.9,0.9);
880   fTrackCuts2->SetRequireSigmaToVertex(kTRUE);
881   fTrackCuts2->SetMaxChi2PerClusterTPC(3.5);
882   fTrackCuts2->SetMinNClustersTPC(80);
883
884   for(Int_t ktracks = 0; ktracks<fESD->GetNumberOfTracks(); ktracks++){
885     AliESDtrack* trackHad = fESD->GetTrack(ktracks);
886     if (!trackHad) {
887       printf("ERROR: Could not receive track %d\n", ktracks);
888       continue;
889     }
890     if(ktracks == itrack) continue; //do not select the same electron
891
892     Double_t ptHad= -999, pHad=-999., dEdxHad = -999;
893     Double_t ptEle = -999;
894     Double_t phiEle = -999, phiHad = -999, Dphi = -999;
895     Double_t pi = 3.14;
896
897     dEdxHad = trackHad->GetTPCsignal();
898     ptHad = trackHad->Pt();
899     pHad = trackHad->P();
900     ptEle = track->Pt();
901
902     if(ptHad <2) continue;
903     if(ptHad > ptEle) continue;
904     if(!fTrackCuts2->AcceptTrack(trackHad)) continue;
905    
906     fHadronPhi->Fill(trackHad->Phi());
907     
908     phiEle = track->Phi();
909     phiHad = trackHad->Phi();
910     Dphi = phiEle - phiHad;
911     if (Dphi > 3*pi/2)
912       Dphi = Dphi - 2*pi;
913     if (Dphi < -pi/2)
914       Dphi = Dphi + 2*pi;
915
916     DphiPt->Fill(ptEle,Dphi);
917
918   }
919 }
920 //_________________________________________
921 void AliAnalysisTaskElecHadronCorrel::ElectronHadCorrelNoPartner(Int_t itrack,Int_t jtrack, AliESDtrack *track, TH2F *DphiPtNew)
922 {
923   //Construct Delta Phi between electrons and hadrons for electrons from invariant mass calculation excluding associated track
924
925   fTrackCuts2->SetAcceptKinkDaughters(kFALSE);
926   fTrackCuts2->SetRequireTPCRefit(kTRUE);
927   fTrackCuts2->SetRequireITSRefit(kTRUE);
928   fTrackCuts2->SetEtaRange(-0.9,0.9);
929   fTrackCuts2->SetRequireSigmaToVertex(kTRUE);
930   fTrackCuts2->SetMaxChi2PerClusterTPC(3.5);
931   fTrackCuts2->SetMinNClustersTPC(80);
932
933   for(Int_t ktracks = 0; ktracks<fESD->GetNumberOfTracks(); ktracks++){
934     AliESDtrack* trackHad = fESD->GetTrack(ktracks);
935     if (!trackHad) {
936       printf("ERROR: Could not receive track %d\n", ktracks);
937       continue;
938     }
939     if(ktracks == itrack || ktracks == jtrack) continue; //do not select the same electron and associated track from inv mass cal
940
941
942     Double_t ptHad= -999, pHad=-999., dEdxHad = -999;
943     Double_t ptEle = -999;
944     Double_t phiEle = -999, phiHad = -999, Dphi = -999;
945     Double_t pi = 3.14;
946
947     dEdxHad = trackHad->GetTPCsignal();
948     ptHad = trackHad->Pt();
949     pHad = trackHad->P();
950     ptEle = track->Pt();
951
952     if(ptHad <2) continue;
953     if(ptHad > ptEle) continue;
954     if(!fTrackCuts2->AcceptTrack(trackHad)) continue;
955
956     phiEle = track->Phi();
957     phiHad = trackHad->Phi();
958     Dphi = phiEle - phiHad;
959     if (Dphi > 3*pi/2)
960       Dphi = Dphi - 2*pi;
961     if (Dphi < -pi/2)
962       Dphi = Dphi + 2*pi;
963
964     DphiPtNew->Fill(ptEle,Dphi);
965   }
966 }
967 /*
968 //_________________________________________
969 void AliAnalysisTaskElecHadronCorrel::MixedEvent(AliESDtrack *track, TH2F *DphiPt)
970 {
971
972   const AliESDVertex *pVtx = fESD->GetPrimaryVertex();
973   Double_t zVtx;
974   zVtx = pVtx->GetZ();
975
976
977   AliCentrality *fCentrality = (AliCentrality*)fESD->GetCentrality();
978   Double_t centvalue = fCentrality->GetCentralityPercentile("V0M");
979
980   AliEventPool* pool = fPoolMgr->GetEventPool(centvalue, zVtx); // Get the buffer associated with the current centrality and z-vtx
981   if (!pool)
982     AliFatal(Form("No pool found for centrality = %f, zVtx = %f", centvalue, zVtx));
983
984   //  pool->PrintInfo();
985   if (pool->GetCurrentNEvents() >= 5) // start mixing when 5 events are in the buffer
986   {
987     Int_t nMix = pool->GetCurrentNEvents();
988     fNoMixedEvents->Fill(0);
989     fMixStat->Fill(pool->GetCurrentNEvents(),centvalue);
990     fMixStat1->Fill(pool->GetCurrentNEvents(),zVtx);
991
992     cout << "nMix = " << nMix << " tracks in pool = " << pool->NTracksInPool() << endl;
993     for (Int_t jMix=0; jMix<pool->GetCurrentNEvents(); jMix++)  // mix with each event in the buffer
994     {
995       TObjArray* bgTracks = pool->GetEvent(jMix);
996       for (Int_t i=0;i<bgTracks->GetEntriesFast(); i++)
997       {
998         AliVParticle* mixtrk = (AliVParticle*) bgTracks->At(i);
999
1000         Double_t mixtrkPhi = -999;
1001         Double_t ptEle = -999;
1002         Double_t phiEle = -999, Dphi = -999;
1003         Double_t pi = 3.14;
1004         Double_t ptmixtrk = -999;
1005
1006         ptEle = track->Pt();
1007         if(ptmixtrk > ptEle) continue;
1008
1009         mixtrkPhi = mixtrk->Phi();
1010         phiEle = track->Phi();
1011         Dphi = phiEle - mixtrkPhi;
1012
1013         if (Dphi > 3*pi/2)
1014           Dphi = Dphi - 2*pi;
1015         if (Dphi < -pi/2)
1016           Dphi = Dphi + 2*pi;
1017         cout << "dphi ; " << Dphi <<endl;
1018         DphiPt->Fill(ptEle,Dphi);
1019       }
1020     }
1021
1022   }
1023
1024 }
1025 //___________________________________________
1026 TObjArray*  AliAnalysisTaskElecHadronCorrel::CloneAndReduceTrackList()
1027 {
1028   // clones a track list by using AliDPhiBasicParticle which uses much less memory (used for event mixing)
1029
1030   fTrackCuts2->SetAcceptKinkDaughters(kFALSE);
1031   fTrackCuts2->SetRequireTPCRefit(kTRUE);
1032   fTrackCuts2->SetRequireITSRefit(kTRUE);
1033   fTrackCuts2->SetEtaRange(-0.9,0.9);
1034   fTrackCuts2->SetRequireSigmaToVertex(kTRUE);
1035   fTrackCuts2->SetMaxChi2PerClusterTPC(3.5);
1036   fTrackCuts2->SetMinNClustersTPC(80);
1037
1038   TObjArray* tracksClone = new TObjArray;
1039   tracksClone->SetOwner(kTRUE);
1040
1041   for(Int_t ktracks = 0; ktracks<fESD->GetNumberOfTracks(); ktracks++){
1042     AliESDtrack* track = fESD->GetTrack(ktracks);
1043     if (!track) {
1044       printf("ERROR: Could not receive track %d\n", ktracks);
1045       continue;
1046     }
1047
1048     //   if(ktracks == iTrack) continue;
1049     Double_t ptHad= -999, pHad=-999.;
1050     ptHad = track->Pt();
1051     pHad = track->P();
1052     
1053     if(ptHad <2) continue;
1054     if(!fTrackCuts2->AcceptTrack(track)) continue;
1055
1056     cout << "hadrntrk"<<endl;
1057     AliVParticle* particle = (AliVParticle*) fESD->GetTrack(ktracks);
1058     tracksClone->Add(new AliDPhiBasicParticle(particle->Eta(), particle->Phi(), particle->Pt(), particle->Charge()));
1059   }
1060
1061   return tracksClone;
1062 }
1063 */