]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGGA/GammaConv/AliAnalysisTaskMaterial.cxx
- changes to the Material budget, resolution and QA code
[u/mrichter/AliRoot.git] / PWGGA / GammaConv / AliAnalysisTaskMaterial.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                                                                                                                                              *
4  * Authors: Friederike Bock                                                                                                                     *
5  * Version 1.0                                                                                                                                                          *
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 // QA Task for V0 Reader V1
19 //---------------------------------------------
20 ////////////////////////////////////////////////
21
22 #include "AliAnalysisTaskMaterial.h"
23 #include "TChain.h"
24 #include "AliAnalysisManager.h"
25 #include "TParticle.h"
26 #include "TVectorF.h"
27 #include "AliPIDResponse.h"
28 #include "AliESDtrackCuts.h"
29 #include "TFile.h"
30
31 class iostream;
32
33 using namespace std;
34
35 ClassImp(AliAnalysisTaskMaterial)
36
37 AliAnalysisTaskMaterial::AliAnalysisTaskMaterial() : AliAnalysisTaskSE(),
38    fV0Reader(NULL),
39    fConversionGammas(NULL),
40    fConversionCuts(NULL),
41    fOutputList(NULL),
42    fEventList(NULL),
43    fRecGammaList(NULL),
44    fAllMCGammaList(NULL),
45    fAllMCConvGammaList(NULL),
46    fTreeEvent(NULL),
47    fTreeMaterialRec(NULL),
48    fTreeMaterialAllGamma(NULL),
49    fTreeMaterialConvGamma(NULL),
50    fPrimVtxZ(0.),
51    fNContrVtx(0),
52    fNESDtracksEta09(0),
53    fNESDtracksEta0914(0),
54    fNESDtracksEta14(0),
55    fGammaMCPt(0.),
56    fGammaMCTheta(0.),
57    fGammaMCConvPt(0.),
58    fGammaMCConvTheta(0.),
59    fMCConvCords(5),
60    fMCConvDaughterProp(4),
61    fGammaPt(0.),
62    fGammaTheta(0.),
63    fGammaChi2NDF(0.),
64    fRecCords(5),
65    fDaughterProp(4),
66    fKind(0),      
67    fIsHeavyIon(0),
68    fIsMC(kFALSE),
69    fESDEvent(NULL),
70    fMCEvent(NULL)
71 {
72
73    
74 }
75
76
77 //________________________________________________________________________
78 AliAnalysisTaskMaterial::AliAnalysisTaskMaterial(const char *name) : AliAnalysisTaskSE(name),
79    fV0Reader(NULL),
80    fConversionGammas(NULL),
81    fConversionCuts(NULL),
82    fOutputList(NULL),
83    fEventList(NULL),
84    fRecGammaList(NULL),
85    fAllMCGammaList(NULL),
86    fAllMCConvGammaList(NULL),
87    fTreeEvent(NULL),
88    fTreeMaterialRec(NULL),
89    fTreeMaterialAllGamma(NULL),
90    fTreeMaterialConvGamma(NULL),
91    fPrimVtxZ(0.),
92    fNContrVtx(0),
93    fNESDtracksEta09(0),
94    fNESDtracksEta0914(0),
95    fNESDtracksEta14(0),
96    fGammaMCPt(0.),
97    fGammaMCTheta(0.),
98    fGammaMCConvPt(0.),
99    fGammaMCConvTheta(0.),
100    fMCConvCords(5),
101    fMCConvDaughterProp(4),
102    fGammaPt(0.),
103    fGammaTheta(0.),
104    fGammaChi2NDF(0.),
105    fRecCords(5),
106    fDaughterProp(4),
107    fKind(0),      
108    fIsHeavyIon(0),
109    fIsMC(kFALSE),
110    fESDEvent(NULL),
111    fMCEvent(NULL)
112 {
113         // Default constructor
114
115         
116         DefineInput(0, TChain::Class());
117         DefineOutput(1, TList::Class());
118 }
119
120 //________________________________________________________________________
121 AliAnalysisTaskMaterial::~AliAnalysisTaskMaterial()
122 {
123         // default deconstructor
124 }
125 //________________________________________________________________________
126 void AliAnalysisTaskMaterial::UserCreateOutputObjects()
127 {
128         // Create User Output Objects
129
130         if(fOutputList != NULL){
131                 delete fOutputList;
132                 fOutputList = NULL;
133         }
134         if(fOutputList == NULL){
135                 fOutputList = new TList();
136                 fOutputList->SetOwner(kTRUE);
137         }
138         
139         fEventList = new TList();
140         fEventList->SetName("EventList");
141         fEventList->SetOwner(kTRUE);
142         fOutputList->Add(fEventList);
143         
144         fTreeEvent = new TTree("Event","Event");   
145         fTreeEvent->Branch("primVtxZ",&fPrimVtxZ,"fPrimVtxZ/F");
146         fTreeEvent->Branch("nContrVtx",&fNContrVtx,"fNContrVtx/I");
147         fTreeEvent->Branch("nGoodTracksEta09",&fNESDtracksEta09,"fNESDtracksEta09/I");
148         fTreeEvent->Branch("nGoodTracksEta14",&fNESDtracksEta14,"fNESDtracksEta14/I");
149         fEventList->Add(fTreeEvent);
150         
151         fRecGammaList= new TList();
152         fRecGammaList->SetName("RecGammaList");
153         fRecGammaList->SetOwner(kTRUE);
154         fOutputList->Add(fRecGammaList);
155                 
156         fTreeMaterialRec = new TTree("ConvPointRec","ConvPointRec");   
157         fTreeMaterialRec->Branch("nGoodTracksEta09",&fNESDtracksEta09,"fNESDtracksEta09/I");
158         fTreeMaterialRec->Branch("recCords",&fRecCords);
159         fTreeMaterialRec->Branch("daughterProp",&fDaughterProp);
160         fTreeMaterialRec->Branch("pt",&fGammaPt,"fGammaPt/F");
161         fTreeMaterialRec->Branch("theta",&fGammaTheta,"fGammaTheta/F");
162         fTreeMaterialRec->Branch("chi2ndf",&fGammaChi2NDF,"fGammaChi2NDF/F");
163         if (fIsMC) {
164                 fTreeMaterialRec->Branch("kind",&fKind,"fKind/b");
165         }   
166         fRecGammaList->Add(fTreeMaterialRec);
167         
168 //      if (fIsMC) {
169 //              fAllMCGammaList = new TList();
170 //              fAllMCGammaList->SetName("AllMCGammaList");
171 //              fAllMCGammaList->SetOwner(kTRUE);
172 //              fOutputList->Add(fAllMCGammaList);
173 //              
174 //              fTreeMaterialAllGamma = new TTree("AllGamma","AllGamma");   
175 //              fTreeMaterialAllGamma->Branch("pt",&fGammaMCPt,"fGammaMCPt/F");
176 //              fTreeMaterialAllGamma->Branch("theta",&fGammaMCTheta,"fGammaMCTheta/F");
177 //              fAllMCGammaList->Add(fTreeMaterialAllGamma);
178 //              
179 //              fAllMCConvGammaList = new TList();
180 //              fAllMCConvGammaList->SetName("AllMCGammaConvList");
181 //              fAllMCConvGammaList->SetOwner(kTRUE);
182 //              fOutputList->Add(fAllMCConvGammaList);
183 // 
184 //      //       fMCConvCords = new Float_t[5];
185 //      //       fMCConvDaughterProp = new Float_t[4];
186 // 
187 //              
188 //              fTreeMaterialConvGamma = new TTree("ConvGammaMC","ConvGammaMC");   
189 //              fTreeMaterialConvGamma->Branch("Cords",&fMCConvCords);
190 //              fTreeMaterialConvGamma->Branch("daughterProp",&fMCConvDaughterProp);
191 //              fTreeMaterialConvGamma->Branch("Pt",&fGammaMCConvPt,"fGammaMCConvPt/F");
192 //              fTreeMaterialConvGamma->Branch("Theta",&fGammaMCConvTheta,"fGammaMCConvTheta/F");   
193 //              fAllMCConvGammaList->Add(fTreeMaterialConvGamma);
194 //      }
195         
196         PostData(1, fOutputList);
197    
198 }
199
200 //________________________________________________________________________
201 void AliAnalysisTaskMaterial::UserExec(Option_t *){
202
203         fV0Reader=(AliV0ReaderV1*)AliAnalysisManager::GetAnalysisManager()->GetTask("V0ReaderV1");
204
205         Int_t eventQuality = ((AliConversionCuts*)fV0Reader->GetConversionCuts())->GetEventQuality();
206                 if(eventQuality != 0){// Event Not Accepted
207                 return;
208         }
209         fESDEvent = (AliESDEvent*) InputEvent();
210         if (fESDEvent==NULL) return;
211         if(MCEvent()){
212                 fMCEvent = MCEvent();
213         }
214  
215         if(MCEvent()){
216                 // Process MC Particle
217                 if(fConversionCuts->GetSignalRejection() != 0){
218 //              if(fESDEvent->IsA()==AliESDEvent::Class()){
219                         fConversionCuts->GetNotRejectedParticles(       fConversionCuts->GetSignalRejection(),
220                                                                                                                 fConversionCuts->GetAcceptedHeader(),
221                                                                                                                 fMCEvent);
222 //              }
223 //              else if(fInputEvent->IsA()==AliAODEvent::Class()){
224 //                      fConversionCuts->GetNotRejectedParticles(       fConversionCuts->GetSignalRejection(),
225 //                                                                                                              fConversionCuts->GetAcceptedHeader(),
226 //                                                                                                              fInputEvent);
227 //              }
228                 }
229         }
230    
231         if(fIsHeavyIon > 0 && !fConversionCuts->IsCentralitySelected(fESDEvent)) return;
232         fNESDtracksEta09 = CountTracks09(); // Estimate Event Multiplicity
233         fNESDtracksEta0914 = CountTracks0914(); // Estimate Event Multiplicity
234         fNESDtracksEta14 = fNESDtracksEta09 + fNESDtracksEta0914;
235         if(fESDEvent){
236                 if(fESDEvent->GetPrimaryVertexTracks()->GetNContributors()>0) {
237                         fNContrVtx = fESDEvent->GetPrimaryVertexTracks()->GetNContributors();
238                 } else {
239                         fNContrVtx = 0;
240                 }       
241 //              else if(fESDEvent->GetPrimaryVertexTracks()->GetNContributors()<1) {
242 //                      if(fESDEvent->GetPrimaryVertexSPD()->GetNContributors()>0) {
243 //                              fNContrVtx = fESDEvent->GetPrimaryVertexSPD()->GetNContributors();
244 //                      } else if(fESDEvent->GetPrimaryVertexSPD()->GetNContributors()<1) {
245 //                              fNContrVtx = 0;
246 //                      }
247 //              }
248         }
249         fPrimVtxZ = fESDEvent->GetPrimaryVertex()->GetZ();
250         
251 //      if (fIsHeavyIon == 2){
252 //              if (!(fNESDtracksEta09 > 20 && fNESDtracksEta09 < 80)) return;
253 //      }       
254         
255         
256         if (fTreeEvent){
257                 fTreeEvent->Fill();
258         }
259                 
260         fConversionGammas=fV0Reader->GetReconstructedGammas();
261         ProcessPhotons();
262 //      if(MCEvent()){
263 //              ProcessMCPhotons();
264 //      }
265         PostData(1, fOutputList);
266 }
267
268 // ///________________________________________________________________________
269 // void AliAnalysisTaskMaterial::FillMCTree(Int_t stackPos){
270 //      AliStack *MCStack = fMCEvent->Stack();
271 //      TParticle* candidate = (TParticle *)MCStack->Particle(stackPos);
272 //      
273 //      if(fConversionCuts->PhotonIsSelectedMC(candidate,MCStack,kFALSE)){
274 //              fGammaMCPt = candidate->Pt();
275 //              fGammaMCTheta = candidate->Theta();
276 //              if (fTreeMaterialAllGamma){
277 //                      fTreeMaterialAllGamma->Fill();
278 //              }       
279 //      }
280 //      if(fConversionCuts->PhotonIsSelectedMC(candidate,MCStack,kTRUE)){
281 //              fGammaMCConvPt = candidate->Pt();
282 //              fGammaMCConvTheta = candidate->Theta();
283 //              TParticle* daughter1 = (TParticle *)MCStack->Particle(candidate->GetFirstDaughter()); 
284 //              TParticle* daughter2 = (TParticle *)MCStack->Particle(candidate->GetLastDaughter()); 
285 //              fMCConvCords(0) = (Float_t)daughter1->Vx();
286 //              fMCConvCords(1) = (Float_t)daughter1->Vy();
287 //              fMCConvCords(2) = (Float_t)daughter1->Vz();
288 //              fMCConvCords(3) = (Float_t)daughter1->R();
289 //              fMCConvCords(4) = (Float_t)daughter1->Phi();
290 //              
291 //              fMCConvDaughterProp(0) = (Float_t)daughter1->Pt();
292 //              fMCConvDaughterProp(1) = (Float_t)daughter1->Theta();
293 //              fMCConvDaughterProp(2) = (Float_t)daughter2->Pt();
294 //              fMCConvDaughterProp(3) = (Float_t)daughter2->Theta();      
295 //                      
296 //              if (fTreeMaterialConvGamma){
297 //                      fTreeMaterialConvGamma->Fill();
298 //              }
299 //      } // Converted MC Gamma         
300 // }
301 // 
302 // ///________________________________________________________________________
303 // void AliAnalysisTaskMaterial::ProcessMCPhotons(){
304 //      // Loop over all primary MC particle
305 //      AliStack *ffMCStack = fMCEvent->Stack();
306 //      for(Int_t i = 0; i < ffMCStack->GetNprimary(); i++) {
307 //              TParticle* particle = (TParticle *)ffMCStack->Particle(i);
308 //              if (!particle) continue;
309 //              Int_t isMCFromMBHeader = -1;
310 //              if(fConversionCuts->GetSignalRejection() != 0){
311 //                      isMCFromMBHeader
312 //                              = fConversionCuts->IsParticleFromBGEvent(i, ffMCStack, fESDEvent);
313 //                      if(isMCFromMBHeader == 0) continue;
314 //              }               
315 //              if (particle->GetPdgCode() == 111 && particle->GetFirstDaughter() >= ffMCStack->GetNprimary()){
316 // //                   cout << "Undecayed pi0 found with mother: " << particle->GetMother(0) << endl;
317 //                      for (Int_t j = 0; j < 2 ; j++){
318 //                              FillMCTree(particle->GetDaughter(j));
319 //                      }
320 //              } else {
321 //                      FillMCTree(i);
322 //              }       
323 //      }       
324 // }
325
326 ///________________________________________________________________________
327 void AliAnalysisTaskMaterial::ProcessPhotons(){
328
329                 // Fill Histograms for QA and MC
330         for(Int_t firstGammaIndex=0;firstGammaIndex<fConversionGammas->GetEntriesFast();firstGammaIndex++){
331                 AliAODConversionPhoton *gamma=dynamic_cast<AliAODConversionPhoton*>(fConversionGammas->At(firstGammaIndex));
332                 if (gamma ==NULL) continue;
333                 if(!fConversionCuts->PhotonIsSelected(gamma,fESDEvent)) continue;
334
335                 fNESDtracksEta09 = CountTracks09();
336                 fGammaPt = gamma->GetPhotonPt();
337                 fGammaTheta = gamma->GetPhotonTheta();
338                 fGammaChi2NDF = gamma->GetChi2perNDF();
339                 fRecCords(0) = (Float_t)gamma->GetConversionX();
340                 fRecCords(1) = (Float_t)gamma->GetConversionY();
341                 fRecCords(2) = (Float_t)gamma->GetConversionZ();
342                 fRecCords(3) = (Float_t)gamma->GetConversionRadius();
343                 fRecCords(4) = (Float_t)gamma->GetPhotonPhi();
344                 
345                 AliESDtrack * negTrack = fConversionCuts->GetESDTrack(fESDEvent, gamma->GetTrackLabelNegative());
346                 AliESDtrack * posTrack = fConversionCuts->GetESDTrack(fESDEvent, gamma->GetTrackLabelPositive());
347                 fDaughterProp(0) = (Float_t)posTrack->Pt();
348                 fDaughterProp(1) = (Float_t)posTrack->Theta();
349                 fDaughterProp(2) = (Float_t)negTrack->Pt();
350                 fDaughterProp(3) = (Float_t)negTrack->Theta();
351
352                 fKind = 9;      
353                 
354                 if(fMCEvent){
355 //                      cout << "generating MC stack"<< endl;
356                         AliStack *fMCStack = fMCEvent->Stack();
357                         if (!fMCStack) continue;
358                         TParticle *posDaughter = gamma->GetPositiveMCDaughter(fMCStack);
359                         TParticle *negDaughter = gamma->GetNegativeMCDaughter(fMCStack);
360 //                      cout << "generate Daughters: "<<posDaughter << "\t" << negDaughter << endl;
361                 
362                         if(fMCStack && fConversionCuts->GetSignalRejection() != 0){
363                                 Int_t isPosFromMBHeader
364                                 = fConversionCuts->IsParticleFromBGEvent(gamma->GetMCLabelPositive(), fMCStack, fESDEvent);
365                                 Int_t isNegFromMBHeader
366                                 = fConversionCuts->IsParticleFromBGEvent(gamma->GetMCLabelNegative(), fMCStack, fESDEvent);
367                                 if( (isNegFromMBHeader < 1) || (isPosFromMBHeader < 1)) continue;
368                         }
369                                 
370                         if(posDaughter == NULL || negDaughter == NULL){ 
371                                 fKind = 9; // garbage
372 //                              cout << "one of the daughters not available" << endl;
373                         } else if(posDaughter->GetMother(0) != negDaughter->GetMother(0) || (posDaughter->GetMother(0) == negDaughter->GetMother(0) && posDaughter->GetMother(0) ==-1)){ 
374                                 // Not Same Mother == Combinatorial Bck
375                                 fKind = 1;
376 //                              cout << "not the same mother" << endl;
377                                 Int_t pdgCodePos; 
378                                 if (posDaughter->GetPdgCode()) pdgCodePos = posDaughter->GetPdgCode(); else continue;
379                                 Int_t pdgCodeNeg; 
380                                 if (negDaughter->GetPdgCode()) pdgCodeNeg = negDaughter->GetPdgCode(); else continue;
381 //                              cout << "PDG codes daughters: " << pdgCodePos << "\t" << pdgCodeNeg << endl;
382                                 if(TMath::Abs(pdgCodePos)==11 && TMath::Abs(pdgCodeNeg)==11)
383                                         fKind = 10; //Electron Combinatorial
384                                 if(TMath::Abs(pdgCodePos)==11 && TMath::Abs(pdgCodeNeg)==11 && (posDaughter->GetMother(0) == negDaughter->GetMother(0) && posDaughter->GetMother(0) ==-1))
385                                         fKind = 15; //direct Electron Combinatorial
386                                 
387                                 if(TMath::Abs(pdgCodePos)==211 && TMath::Abs(pdgCodeNeg)==211)
388                                         fKind = 11; //Pion Combinatorial
389                                 if((TMath::Abs(pdgCodePos)==211 && TMath::Abs(pdgCodeNeg)==2212) ||
390                                         (TMath::Abs(pdgCodePos)==2212 && TMath::Abs(pdgCodeNeg)==211))
391                                         fKind = 12; //Pion, Proton Combinatorics
392                                 if((TMath::Abs(pdgCodePos)==211 && TMath::Abs(pdgCodeNeg)==11) ||
393                                         (TMath::Abs(pdgCodePos)==11 && TMath::Abs(pdgCodeNeg)==211))
394                                         fKind = 13; //Pion, Electron Combinatorics
395                                 if (TMath::Abs(pdgCodePos)==321 || TMath::Abs(pdgCodeNeg)==321) 
396                                         fKind = 14; //Kaon combinatorics
397
398                         } else {
399 //                              cout << "same mother" << endl;
400                                 Int_t pdgCodePos; 
401                                 if (posDaughter->GetPdgCode()) pdgCodePos = posDaughter->GetPdgCode(); else continue;
402                                 Int_t pdgCodeNeg; 
403                                 if (negDaughter->GetPdgCode()) pdgCodeNeg = negDaughter->GetPdgCode(); else continue;
404 //                              cout << "PDG codes daughters: " << pdgCodePos << "\t" << pdgCodeNeg << endl;
405                                 Int_t pdgCode; 
406                                 if (gamma->GetMCParticle(fMCStack)->GetPdgCode()) pdgCode = gamma->GetMCParticle(fMCStack)->GetPdgCode(); else continue;
407 //                              cout << "PDG code: " << pdgCode << endl;
408                                 if(TMath::Abs(pdgCodePos)!=11 || TMath::Abs(pdgCodeNeg)!=11)
409                                         fKind = 2; // combinatorics from hadronic decays
410                                 else if ( !(pdgCodeNeg==pdgCodePos)){
411                                         TParticle *truePhotonCanditate = gamma->GetMCParticle(fMCStack);
412                                         Int_t motherLabelPhoton = truePhotonCanditate->GetMother(0);
413                                         if(pdgCode == 111) 
414                                                 fKind = 3; // pi0 Dalitz
415                                         else if (pdgCode == 221) 
416                                                 fKind = 4; // eta Dalitz
417                                         else if (!(negDaughter->GetUniqueID() != 5 || posDaughter->GetUniqueID() !=5)){
418                                                 if(pdgCode == 22 && motherLabelPhoton < fMCStack->GetNprimary()){
419                                                         fKind = 0; // primary photons
420                                                 } else if (pdgCode == 22){
421                                                         fKind = 5; //secondary photons
422                                                 }               
423                                         } else  fKind = 9; //garbage
424                                 } else fKind = 9; //garbage
425                         }                                       
426                 }
427                                 
428                 if (fTreeMaterialRec){
429                         fTreeMaterialRec->Fill();
430                 }
431         }
432 }
433
434 //________________________________________________________________________
435 Int_t AliAnalysisTaskMaterial::CountTracks09(){
436         Int_t fNumberOfESDTracks = 0;
437         if(fInputEvent->IsA()==AliESDEvent::Class()){
438         // Using standard function for setting Cuts
439                 
440                 AliStack *fMCStack = NULL;
441                 if (fMCEvent){
442                         fMCStack= fMCEvent->Stack();
443                         if (!fMCStack) return 0;
444                 }       
445                                 
446                 Bool_t selectPrimaries=kTRUE;
447                 AliESDtrackCuts *EsdTrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(selectPrimaries);
448                 EsdTrackCuts->SetMaxDCAToVertexZ(2);
449                 EsdTrackCuts->SetEtaRange(-0.9, 0.9);
450                 EsdTrackCuts->SetPtRange(0.15);
451                 
452                 for(Int_t iTracks = 0; iTracks < fInputEvent->GetNumberOfTracks(); iTracks++){
453                         AliESDtrack* curTrack = (AliESDtrack*) fInputEvent->GetTrack(iTracks);
454                         if(!curTrack) continue;
455                         if(EsdTrackCuts->AcceptTrack(curTrack) ){
456                                 if (fMCEvent){
457                                         if(fMCStack && fConversionCuts->GetSignalRejection() != 0){
458                                                 Int_t isFromMBHeader
459                                                 = fConversionCuts->IsParticleFromBGEvent(abs(curTrack->GetLabel()), fMCStack, fESDEvent);
460                                                 if( (isFromMBHeader < 1) ) continue;
461                                         }                                       
462                                 }       
463                                 fNumberOfESDTracks++;
464                         }       
465                 }
466                 delete EsdTrackCuts;
467                 EsdTrackCuts=0x0;
468         } 
469         return fNumberOfESDTracks;
470 }
471
472 //________________________________________________________________________
473 Int_t AliAnalysisTaskMaterial::CountTracks0914(){
474         Int_t fNumberOfESDTracks = 0;
475         if(fInputEvent->IsA()==AliESDEvent::Class()){
476                 // Using standard function for setting Cuts
477                 
478                 AliStack *fMCStack = NULL;
479                 if (fMCEvent){
480                         fMCStack= fMCEvent->Stack();
481                         if (!fMCStack) return 0;
482                 }       
483
484                 AliESDtrackCuts *EsdTrackCuts = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
485                 EsdTrackCuts->SetMaxDCAToVertexZ(5);
486                 EsdTrackCuts->SetEtaRange(0.9, 1.4);
487                 EsdTrackCuts->SetPtRange(0.15);
488                 
489                 for(Int_t iTracks = 0; iTracks < fInputEvent->GetNumberOfTracks(); iTracks++){
490                         AliESDtrack* curTrack = (AliESDtrack*) fInputEvent->GetTrack(iTracks);
491                         if(!curTrack) continue;
492                         if(EsdTrackCuts->AcceptTrack(curTrack) ){
493                                 if (fMCEvent){
494                                         if(fMCStack && fConversionCuts->GetSignalRejection() != 0){
495                                                 Int_t isFromMBHeader
496                                                 = fConversionCuts->IsParticleFromBGEvent(abs(curTrack->GetLabel()), fMCStack, fESDEvent);
497                                                 if( (isFromMBHeader < 1) ) continue;
498                                         }                                       
499                                 }       
500                                 fNumberOfESDTracks++;
501                         }
502                 }
503                 EsdTrackCuts->SetEtaRange(-1.4, -0.9);
504                 for(Int_t iTracks = 0; iTracks < fESDEvent->GetNumberOfTracks(); iTracks++){
505                         AliESDtrack* curTrack = fESDEvent->GetTrack(iTracks);
506                         if(!curTrack) continue;
507                         if(EsdTrackCuts->AcceptTrack(curTrack) ){
508                                 if (fMCEvent){
509                                         if(fMCStack && fConversionCuts->GetSignalRejection() != 0){
510                                                 Int_t isFromMBHeader
511                                                 = fConversionCuts->IsParticleFromBGEvent(abs(curTrack->GetLabel()), fMCStack, fESDEvent);
512                                                 if( (isFromMBHeader < 1) ) continue;
513                                         }                                       
514                                 }       
515                                 fNumberOfESDTracks++;
516                         }       
517                 }
518                 delete EsdTrackCuts;
519                 EsdTrackCuts=0x0;
520         } 
521         
522         return fNumberOfESDTracks;
523 }
524
525 //________________________________________________________________________
526 void AliAnalysisTaskMaterial::Terminate(Option_t *)
527 {
528 //    if (fStreamMaterial){
529 //       fStreamMaterial->GetFile()->Write();
530 //    }
531 }