]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/SPECTRA/AliAnalysisTaskESDDedx.cxx
First dE/dx Analysis
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliAnalysisTaskESDDedx.cxx
1 /**************************************************************************
2  * Author: Boris Hippolyte.                                               *
3  * Contributors are mentioned in the code where appropriate.              *
4  *                                                                        *
5  * Permission to use, copy, modify and distribute this software and its   *
6  * documentation strictly for non-commercial purposes is hereby granted   *
7  * without fee, provided that the above copyright notice appears in all   *
8  * copies and that both the copyright notice and this permission notice   *
9  * appear in the supporting documentation. The authors make no claims     *
10  * about the suitability of this software for any purpose. It is          *
11  * provided "as is" without express or implied warranty.                  *
12  **************************************************************************/
13
14 //-----------------------------------------------------------------
15 //                 AliAnalysisTaskESDCheckV0 class
16 //            This task is for QAing the dE/dx from the ESD
17 //              Origin: B.H. Nov2007, hippolyt@in2p3.fr
18 //-----------------------------------------------------------------
19
20 //-----------------------------------------------------------------
21 //    ---> Next upgrades
22 //      1) extent the use of fMidPseudoRapidityFlag to primaries
23 //      2) discuss the GetSigmaToVertex(AliESDtrack*)
24 //      3) decide the kind of refit to be used
25 //-----------------------------------------------------------------
26
27 #include "TChain.h"
28 #include "TTree.h"
29 #include "TList.h"
30 #include "TH1F.h"
31 #include "TH2F.h"
32 #include "TCanvas.h"
33 #include "TVector3.h"
34
35 #include "AliAnalysisTask.h"
36 #include "AliAnalysisManager.h"
37
38 #include "AliESDEvent.h"
39 #include "AliESDInputHandler.h"
40
41 #include "AliESDv0.h"
42
43 #include "AliAnalysisTaskESDDedx.h"
44
45 ClassImp(AliAnalysisTaskESDDedx)
46
47 //________________________________________________________________________
48   AliAnalysisTaskESDDedx::AliAnalysisTaskESDDedx(const char    *rName,
49                                                  const Bool_t   rAllConstrainedFlag,
50                                                  const Bool_t   rMidPseudoRapidityFlag,
51                                                  const Bool_t   rSelTrackRemoveKink,
52                                                  const Int_t    rSelTrackMinClustersTPC,
53                                                  const Int_t    rSelTrackMinClustersITS,
54                                                  const Float_t  rSelTrackMaxChi2PerClusterTPC,
55                                                  const Float_t  rSelTrackMaxChi2PerClusterITS,
56                                                  const Double_t rSelTrackMaxCov11,
57                                                  const Double_t rSelTrackMaxCov22,
58                                                  const Double_t rSelTrackMaxCov33,
59                                                  const Double_t rSelTrackMaxCov44,
60                                                  const Double_t rSelTrackMaxCov55) 
61     : AliAnalysisTask(rName, ""), fESD(0), fListHist(), fHistPtot(0),
62
63     fHistMultiplicity(0), fHistTPCDedxVsMomentum(0), fHistITSDedxVsMomentum(0),
64     fHistMassK0(0), fHistMassLambda(0), fHistMassAntiLambda(0),
65     fHistTPCDedxVsMomPosK0(0), fHistTPCDedxVsMomNegK0(0),
66     fHistDiffInOutMomentum(0), fHistDiffPrimOutMomentum(0),
67     fHistDiffPrimMeanMomentum(0), fHistPercPrimMeanMomentum(0), fHistPrimEta(0),
68     fHistPercPrimMeanMomentumVsEta(0), fHistPercPrimMeanMomentumVsPrim(0),
69
70     fHistMultiplicityCuts(0), fHistTPCDedxVsMomentumCuts(0), fHistITSDedxVsMomentumCuts(0),
71     fHistMassK0Cuts(0), fHistMassLambdaCuts(0), fHistMassAntiLambdaCuts(0), 
72     fHistTPCDedxVsMomPosK0Cuts(0), fHistTPCDedxVsMomNegK0Cuts(0), 
73     fHistDiffInOutMomentumCuts(0), fHistDiffPrimOutMomentumCuts(0),
74     fHistDiffPrimMeanMomentumCuts(0), fHistPercPrimMeanMomentumCuts(0), fHistPrimEtaCuts(0),
75     fHistPercPrimMeanMomentumVsEtaCuts(0), fHistPercPrimMeanMomentumVsPrimCuts(0),
76
77     fAllConstrainedFlag(rAllConstrainedFlag), fMidPseudoRapidityFlag(rMidPseudoRapidityFlag),
78
79     fSelTrackRemoveKink(rSelTrackRemoveKink),
80     fSelTrackMinClustersTPC(rSelTrackMinClustersTPC), fSelTrackMinClustersITS(rSelTrackMinClustersITS),
81     fSelTrackMaxChi2PerClusterTPC(rSelTrackMaxChi2PerClusterTPC), fSelTrackMaxChi2PerClusterITS(rSelTrackMaxChi2PerClusterITS),
82     fSelTrackMaxCov11(rSelTrackMaxCov11), fSelTrackMaxCov22(rSelTrackMaxCov22), fSelTrackMaxCov33(rSelTrackMaxCov33),
83     fSelTrackMaxCov44(rSelTrackMaxCov44), fSelTrackMaxCov55(rSelTrackMaxCov55)
84 {
85   // Constructor
86
87   // Define input and output slots here
88   // Input slot #0 works with a TChain
89   DefineInput(0, TChain::Class());
90   // Output slot #0 writes into a TH1 container
91   DefineOutput(0, TH1F::Class());
92   // Output slot #1 writes into a TList container
93   DefineOutput(1, TList::Class());
94 }
95
96 //________________________________________________________________________
97 void AliAnalysisTaskESDDedx::ConnectInputData(Option_t *) 
98 {
99   // Connect ESD or AOD here
100   // Called once
101
102   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
103   if (!tree) {
104     Printf("ERROR: Could not read chain from input slot 0");
105   } else {
106     // Disable all branches and enable only the needed ones
107     // The next two lines are different when data produced as AliESDEvent is read
108     tree->SetBranchStatus("*", kFALSE);
109     tree->SetBranchStatus("fTracks.*", kTRUE);
110     tree->SetBranchStatus("fV0s.*", kTRUE);
111
112     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
113
114     if (!esdH) {
115       Printf("ERROR: Could not get ESDInputHandler");
116     } else
117       fESD = esdH->GetEvent();
118   }
119 }
120
121 //________________________________________________________________________
122 void AliAnalysisTaskESDDedx::CreateOutputObjects()
123 {
124   // Create histograms
125   // Called once
126   fListHist = new TList();
127   if (!fHistPtot) {
128     fHistPtot = new TH1F("fHistPtot", "P_{tot} distribution;P_{tot} (GeV/c);dN/dP_{tot} (c/GeV)", 15, 0.1, 3.1);
129     fHistPtot->SetMarkerStyle(kFullCircle);
130     fListHist->Add(fHistPtot);
131   }
132   if (!fHistMultiplicity) {
133     fHistMultiplicity = new TH1F("fHistMultiplicity", "Multiplicity distribution;Number of tracks;Events", 250, 0, 250);
134     fHistMultiplicity->SetMarkerStyle(kFullCircle);
135     fListHist->Add(fHistMultiplicity);
136   }
137   if (!fHistTPCDedxVsMomentum) {
138     fHistTPCDedxVsMomentum = new TH2F("h2TPCDedxVsMomentum","Bethe-Bloch Distribution for TPC;Momentum (GeV/c);dE/dx (A.U.)",1500,0,15,100,0,100);
139     fListHist->Add(fHistTPCDedxVsMomentum);
140   }
141   if (!fHistITSDedxVsMomentum) {
142     fHistITSDedxVsMomentum = new TH2F("h2ITSDedxVsMomentum","Bethe-Bloch Distribution for ITS;Momentum (GeV/c);dE/dx (A.U.)",1500,0,15,100,0,100);
143     fListHist->Add(fHistITSDedxVsMomentum);
144   }
145   if (!fHistMassK0) {
146     fHistMassK0 = new TH1F("h1MassK0","K^{0} candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});Counts",100,0.4,0.6);
147     fListHist->Add(fHistMassK0);
148   }
149   if (!fHistMassLambda) {
150     fHistMassLambda = new TH1F("h1MassLambda","#Lambda^{0} candidates;M(p#pi^{-}) (GeV/c^{2});Counts",75,1.05,1.2);
151     fListHist->Add(fHistMassLambda);
152   }
153   if (!fHistMassAntiLambda) {
154     fHistMassAntiLambda = new TH1F("h1MassAntiLambda","#bar{#Lambda}^{0} candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts",75,1.05,1.2);
155     fListHist->Add(fHistMassAntiLambda);
156   }
157   if (!fHistTPCDedxVsMomPosK0) {
158     fHistTPCDedxVsMomPosK0 = new TH2F("h2TPCDedxVsMomPosK0","Bethe-Bloch Distribution for TPC;Momentum (GeV/c);dE/dx (A.U.)",1500,0,15,100,0,100);
159     fListHist->Add(fHistTPCDedxVsMomPosK0);
160   }
161   if (!fHistTPCDedxVsMomNegK0) {
162     fHistTPCDedxVsMomNegK0 = new TH2F("h2TPCDedxVsMomNegK0","Bethe-Bloch Distribution for TPC;Momentum (GeV/c);dE/dx (A.U.)",1500,0,15,100,0,100);
163     fListHist->Add(fHistTPCDedxVsMomNegK0);
164   }
165   if (!fHistDiffInOutMomentum) {
166     fHistDiffInOutMomentum = new TH1F("h1DiffInOutMomentum","Momentum Difference In-Out;Momentum (GeV/c);Counts",200,-0.1,0.1); 
167     fListHist->Add(fHistDiffInOutMomentum);
168   }
169   if (!fHistDiffPrimOutMomentum) {
170     fHistDiffPrimOutMomentum = new TH1F("h1DiffPrimOutMomentum","Momentum Difference Prim-Out;Momentum (GeV/c);Counts",200,-0.1,0.1); 
171     fListHist->Add(fHistDiffPrimOutMomentum);
172   }
173   if (!fHistDiffPrimMeanMomentum) {
174     fHistDiffPrimMeanMomentum = new TH1F("h1DiffPrimMeanMomentum","Momentum Difference Prim-Mean;Momentum (GeV/c);Counts",200,-0.1,0.1); 
175     fListHist->Add(fHistDiffPrimMeanMomentum);
176   }
177   if (!fHistPercPrimMeanMomentum) {
178     fHistPercPrimMeanMomentum = new TH1F("h1PercPrimMeanMomentum","Momentum Percentage (Prim-Mean)/Prim;(Primary-Mean)/Primary (%);Counts",200,-10,10); 
179     fListHist->Add(fHistPercPrimMeanMomentum);
180   }
181   if (!fHistPrimEta) {
182     fHistPrimEta = new TH1F("h1PrimEta","Pseudorapidity Distribution (Primaries);#eta;Counts",300,-1.5,1.5); 
183     fListHist->Add(fHistPrimEta);
184   }
185   if (!fHistPercPrimMeanMomentumVsEta) {
186     fHistPercPrimMeanMomentumVsEta = new TH2F("h2PercPrimMeanMomentumVsEta","Momentum Percentage (Prim-Mean)/Prim vs #eta;#eta;(Primary-Mean)/Primary (%)",100,-1.5,1.5,100,-10,10); 
187     fListHist->Add(fHistPercPrimMeanMomentumVsEta);
188   }
189   if (!fHistPercPrimMeanMomentumVsPrim) {
190     fHistPercPrimMeanMomentumVsPrim = new TH2F("h2PercPrimMeanMomentumVsPrim","Momentum Percentage (Prim-Mean)/Prim vs Prim;Momentum (GeV/c);(Primary-Mean)/Primary (%)",1500,0,15,100,-10,10); 
191     fListHist->Add(fHistPercPrimMeanMomentumVsPrim);
192   }
193   // Histograms after selections
194   if (!fHistMultiplicityCuts) {
195     fHistMultiplicityCuts = new TH1F("h1HistMultiplicityCuts","Number of Good tracks;Number of tracks;Event",201,-0.5,200.5);
196     fListHist->Add(fHistMultiplicityCuts);
197   }
198   if (!fHistTPCDedxVsMomentumCuts) {
199     fHistTPCDedxVsMomentumCuts = new TH2F("h2TPCDedxVsMomentumCuts","Bethe-Bloch Distribution for TPC w/cuts;Momentum (GeV/c);dE/dx (A.U.)",1500,0,15,100,0,100);
200     fListHist->Add(fHistTPCDedxVsMomentumCuts);
201   }
202   if (!fHistITSDedxVsMomentumCuts) {
203     fHistITSDedxVsMomentumCuts = new TH2F("h2ITSDedxVsMomentumCuts","Bethe-Bloch Distribution for ITS w/cuts;Momentum (GeV/c);dE/dx (A.U.)",1500,0,15,100,0,100);
204     fListHist->Add(fHistITSDedxVsMomentumCuts);
205   }
206   if (!fHistMassK0Cuts) {
207     fHistMassK0Cuts = new TH1F("h1MassK0Cuts","K^{0} candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});Counts",100,0.4,0.6);
208     fListHist->Add(fHistMassK0Cuts);
209   }
210   if (!fHistMassLambdaCuts) {
211     fHistMassLambdaCuts = new TH1F("h1MassLambdaCuts","#Lambda^{0} candidates;M(p#pi^{-}) (GeV/c^{2});Counts",75,1.05,1.2);
212     fListHist->Add(fHistMassLambdaCuts);
213   }
214   if (!fHistMassAntiLambdaCuts) {
215     fHistMassAntiLambdaCuts = new TH1F("h1MassAntiLambdaCuts","#bar{#Lambda}^{0} candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts",75,1.05,1.2);
216     fListHist->Add(fHistMassAntiLambdaCuts);
217   }
218   if (!fHistTPCDedxVsMomPosK0Cuts) {
219     fHistTPCDedxVsMomPosK0Cuts = new TH2F("h2TPCDedxVsMomPosK0Cuts","Bethe-Bloch Distribution for TPC;Momentum (GeV/c);dE/dx (A.U.)",1500,0,15,100,0,100);
220     fListHist->Add(fHistTPCDedxVsMomPosK0Cuts);
221   }
222   if (!fHistTPCDedxVsMomNegK0Cuts) {
223     fHistTPCDedxVsMomNegK0Cuts = new TH2F("h2TPCDedxVsMomNegK0Cuts","Bethe-Bloch Distribution for TPC;Momentum (GeV/c);dE/dx (A.U.)",1500,0,15,100,0,100);
224     fListHist->Add(fHistTPCDedxVsMomNegK0Cuts);
225   }
226   if (!fHistDiffInOutMomentumCuts) {
227     fHistDiffInOutMomentumCuts = new TH1F("h1DiffInOutMomentumCuts","Momentum Difference In-Out;Momentum (GeV/c);Counts",200,-0.1,0.1); 
228     fListHist->Add(fHistDiffInOutMomentumCuts);
229   }
230   if (!fHistDiffPrimOutMomentumCuts) {
231     fHistDiffPrimOutMomentumCuts = new TH1F("h1DiffPrimOutMomentumCuts","Momentum Difference Prim-Out;Momentum (GeV/c);Counts",200,-0.1,0.1); 
232     fListHist->Add(fHistDiffPrimOutMomentumCuts);
233   }
234   if (!fHistDiffPrimMeanMomentumCuts) {
235     fHistDiffPrimMeanMomentumCuts = new TH1F("h1DiffPrimMeanMomentumCuts","Momentum Difference Prim-Mean;Momentum (GeV/c);Counts",200,-0.1,0.1); 
236     fListHist->Add(fHistDiffPrimMeanMomentumCuts);
237   }
238   if (!fHistPercPrimMeanMomentumCuts) {
239     fHistPercPrimMeanMomentumCuts = new TH1F("h1PercPrimMeanMomentumCuts","Momentum Percentage (Prim-Mean)/Prim;(Primary-Mean)/Primary (%);Counts",200,-10,10); 
240     fListHist->Add(fHistPercPrimMeanMomentumCuts);
241   }
242   if (!fHistPrimEtaCuts) {
243     fHistPrimEtaCuts = new TH1F("h1PrimEtaCuts","Pseudorapidity Distribution (Primaries);#eta;Counts",300,-1.5,1.5); 
244     fListHist->Add(fHistPrimEtaCuts);
245   }
246   if (!fHistPercPrimMeanMomentumVsEtaCuts) {
247     fHistPercPrimMeanMomentumVsEtaCuts = new TH2F("h2PercPrimMeanMomentumVsEtaCuts","Momentum Percentage (Prim-Mean)/Prim vs #eta;#eta;(Primary-Mean)/Primary (%)",100,-1.5,1.5,100,-10,10); 
248     fListHist->Add(fHistPercPrimMeanMomentumVsEtaCuts);
249   }
250   if (!fHistPercPrimMeanMomentumVsPrimCuts) {
251     fHistPercPrimMeanMomentumVsPrimCuts = new TH2F("h2PercPrimMeanMomentumVsPrimCuts","Momentum Percentage (Prim-Mean)/Prim vs Prim;Momentum (GeV/c);(Primary-Mean)/Primary (%)",1500,0,15,100,-10,10); 
252     fListHist->Add(fHistPercPrimMeanMomentumVsPrimCuts);
253   }
254 }
255
256 //________________________________________________________________________
257 void AliAnalysisTaskESDDedx::Exec(Option_t *) 
258 {
259   // Main loop
260   // Called for each event
261
262   if (!fESD) {
263     Printf("ERROR: fESD not available");
264     return;
265   }
266
267   Int_t    nGetTracks        = fESD->GetNumberOfTracks();
268   //  Printf("There are %d tracks in this event", nGetTracks);
269
270   Int_t    nAllTracks        = 0, nGoodTracks       = 0;
271   Double_t lTPCDedx          = 0, lITSDedx          = 0;
272
273   // Track loop
274   for (Int_t iTracks = 0; iTracks < nGetTracks; iTracks++) {
275     AliESDtrack* track = fESD->GetTrack(iTracks);
276     if (!track) {
277       Printf("ERROR: Could not retrieve track %d", iTracks);
278       continue;
279     }
280     fHistPtot->Fill(track->P()); // Fill a ptot spectrum before any selection
281
282     Double_t lPrimvtxMomemtum[3];
283     Bool_t lIsTrackConstrained = track->GetConstrainedPxPyPz(lPrimvtxMomemtum);
284     if (fAllConstrainedFlag && (!lIsTrackConstrained)) continue; // Constrained or not to Prim Vertex.
285
286     const AliExternalTrackParam *lInnerETP = track->GetInnerParam();
287     const AliExternalTrackParam *lOuterETP = track->GetOuterParam();
288     if ((!lInnerETP)||(!lOuterETP)) continue;  // No inner nor outer params for this track.
289
290     Double_t lInnerMomentum = 0, lOuterMomentum = 0, lMeanMomentum = 0;
291     if (lInnerETP) lInnerMomentum = lInnerETP->GetP(); // The "inner" momentum.
292     if (lOuterETP) lOuterMomentum = lOuterETP->GetP(); // The "outer" momentum.
293     if (lInnerETP&&lOuterETP) lMeanMomentum = (lInnerMomentum+lOuterMomentum)/2.; // The "mean" momentum 
294
295     lTPCDedx = track->GetTPCsignal();
296     lITSDedx = track->GetITSsignal();
297     // This needs to be discussed and implemented
298     // nSigmaToVertex = GetSigmaToVertex(track);
299
300     Float_t lPrimPseudoRap = 0, lPrimMomentum = 0;
301     if (lIsTrackConstrained){
302       TVector3 pPrimvtx(lPrimvtxMomemtum);
303       lPrimPseudoRap = pPrimvtx.Eta();
304       lPrimMomentum = pPrimvtx.Mag();
305     }
306
307     if (
308         ((!fAllConstrainedFlag)||(TMath::Abs(lPrimPseudoRap) < 0.5)) &&
309         (lTPCDedx)
310         )
311       // fAllConstrained = 0 -> !fAllConstrained = 1 so ok anyway
312       // fAllConstrained = 1 -> !fAllConstrained = 0 so it depends on eta condition 
313       {
314         nAllTracks++;
315         fHistTPCDedxVsMomentum->Fill(lInnerMomentum,lTPCDedx);
316         fHistDiffInOutMomentum->Fill(lInnerMomentum-lOuterMomentum);
317         if (lIsTrackConstrained){
318           fHistITSDedxVsMomentum->Fill(lPrimMomentum,lITSDedx);
319           fHistDiffPrimOutMomentum->Fill(lPrimMomentum-lOuterMomentum);
320           fHistDiffPrimMeanMomentum->Fill(lPrimMomentum-lMeanMomentum);
321           if (lPrimMomentum){
322             fHistPercPrimMeanMomentum->Fill(100.*(lPrimMomentum-lMeanMomentum)/lPrimMomentum);
323             fHistPercPrimMeanMomentumVsEta->Fill(lPrimPseudoRap,100.*(lPrimMomentum-lMeanMomentum)/lPrimMomentum);
324             fHistPercPrimMeanMomentumVsPrim->Fill(lPrimMomentum,100.*(lPrimMomentum-lMeanMomentum)/lPrimMomentum);
325           }
326           fHistPrimEta->Fill(lPrimPseudoRap);
327         }
328         if(
329            // (!fAllConstrainedFlag || (nSigmaToVertex < 3)) &&
330            (IsAccepted(track))
331            ){
332           fHistTPCDedxVsMomentumCuts->Fill(lInnerMomentum,lTPCDedx);
333           fHistDiffInOutMomentumCuts->Fill(lInnerMomentum-lOuterMomentum);
334           if (lIsTrackConstrained){
335             fHistITSDedxVsMomentumCuts->Fill(lPrimMomentum,lITSDedx);
336             fHistDiffPrimOutMomentumCuts->Fill(lPrimMomentum-lOuterMomentum);
337             fHistDiffPrimMeanMomentumCuts->Fill(lPrimMomentum-lMeanMomentum);
338             if (lPrimMomentum){
339               fHistPercPrimMeanMomentumCuts->Fill(100.*(lPrimMomentum-lMeanMomentum)/lPrimMomentum);
340               fHistPercPrimMeanMomentumVsEtaCuts->Fill(lPrimPseudoRap,100.*(lPrimMomentum-lMeanMomentum)/lPrimMomentum);
341               fHistPercPrimMeanMomentumVsPrimCuts->Fill(lPrimMomentum,100.*(lPrimMomentum-lMeanMomentum)/lPrimMomentum);
342             }
343             fHistPrimEtaCuts->Fill(lPrimPseudoRap);
344           }
345           nGoodTracks++;
346         }
347       }
348   } //track loop 
349   fHistMultiplicity->Fill(nGetTracks);
350   fHistMultiplicityCuts->Fill(nGoodTracks);
351
352   Int_t nv0s = 0;
353   nv0s = fESD->GetNumberOfV0s();
354
355   Int_t    lIndexTrackPos       = 0,  lIndexTrackNeg       = 0;
356   Double_t lTPCDedxPos          = 0,  lTPCDedxNeg          = 0;
357   Float_t lInvMassK0 = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
358   for (Int_t iV0 = 0; iV0 < nv0s; iV0++) 
359     {// This is the V0 loop
360       AliESDv0 *v0 = fESD->GetV0(iV0);
361       if (!v0) continue;
362
363       // Getting invariant mass infos directly from ESD
364       v0->ChangeMassHypothesis(310);
365       lInvMassK0 = v0->GetEffMass();
366       v0->ChangeMassHypothesis(3122);
367       lInvMassLambda = v0->GetEffMass();
368       v0->ChangeMassHypothesis(-3122);
369       lInvMassAntiLambda = v0->GetEffMass();
370
371       // Filling invariant mass histos for all candidates
372       fHistMassK0->Fill(lInvMassK0);
373       fHistMassLambda->Fill(lInvMassLambda);
374       fHistMassAntiLambda->Fill(lInvMassAntiLambda);
375
376       // Accessing the daughter track infos
377       lIndexTrackPos = TMath::Abs(v0->GetPindex());
378       lIndexTrackNeg = TMath::Abs(v0->GetNindex());
379       AliESDtrack* trackPos = fESD->GetTrack(lIndexTrackPos);
380       AliESDtrack* trackNeg = fESD->GetTrack(lIndexTrackNeg);
381       
382       Double_t lPrimvtxMomPos[3];
383       Bool_t lIsPosConstrained = trackPos->GetConstrainedPxPyPz(lPrimvtxMomPos);
384       if (fAllConstrainedFlag && (!lIsPosConstrained)) continue; // Constrained or not to Prim Vertex.
385
386       Double_t lPrimvtxMomNeg[3];
387       Bool_t lIsNegConstrained = trackNeg->GetConstrainedPxPyPz(lPrimvtxMomNeg);
388       if (fAllConstrainedFlag && (!lIsNegConstrained)) continue; // Constrained or not to Prim Vertex.
389
390       const AliExternalTrackParam *lInnerPosETP = trackPos->GetInnerParam();
391       const AliExternalTrackParam *lInnerNegETP = trackNeg->GetInnerParam();
392       if ((!lInnerPosETP)||(!lInnerNegETP)) continue; // No inner params for at least one track.
393
394       Double_t lInnerMomPos = 0, lInnerMomNeg = 0;
395       if (lInnerPosETP) lInnerMomPos = lInnerPosETP->GetP();
396       if (lInnerNegETP) lInnerMomNeg = lInnerNegETP->GetP();
397
398       lTPCDedxPos = trackPos->GetTPCsignal();
399       lTPCDedxNeg = trackNeg->GetTPCsignal();
400
401       Float_t lPrimPseudoRapPos = 0, lPrimPseudoRapNeg = 0;
402       if (lIsPosConstrained && lIsNegConstrained){
403         TVector3 pPrimPos(lPrimvtxMomPos);
404         TVector3 pPrimNeg(lPrimvtxMomNeg);
405         lPrimPseudoRapPos = pPrimPos.Eta();
406         lPrimPseudoRapNeg = pPrimNeg.Eta();
407       }
408
409       if (TMath::Abs(lInvMassK0-0.497)<0.01) {
410         fHistMassK0Cuts->Fill(lInvMassK0);
411         fHistTPCDedxVsMomPosK0->Fill(lInnerMomPos,lTPCDedxPos);
412         fHistTPCDedxVsMomNegK0->Fill(lInnerMomNeg,lTPCDedxNeg);
413         if (
414             (!fMidPseudoRapidityFlag) || ( TMath::Abs(lPrimPseudoRapPos) < 0.5) &&
415             (IsAccepted(trackPos))
416             )
417           fHistTPCDedxVsMomPosK0Cuts->Fill(lInnerMomPos,lTPCDedxPos);
418         if (
419             (!fMidPseudoRapidityFlag) || ( TMath::Abs(lPrimPseudoRapNeg) < 0.5) &&
420             (IsAccepted(trackNeg))
421             )
422           fHistTPCDedxVsMomNegK0Cuts->Fill(lInnerMomNeg,lTPCDedxNeg);
423       }
424       if (TMath::Abs(lInvMassLambda-1.115)<0.01) {
425         fHistMassLambdaCuts->Fill(lInvMassLambda);
426       }
427       if (TMath::Abs(lInvMassAntiLambda-1.115)<0.01) {
428         fHistMassAntiLambdaCuts->Fill(lInvMassAntiLambda);
429       }
430     }//V0 loop 
431   
432   // Post output data.
433   PostData(0, fHistPtot);
434   PostData(1, fListHist);
435 }      
436
437 //________________________________________________________________________
438 void AliAnalysisTaskESDDedx::Terminate(Option_t *) 
439 {
440   // Draw result to the screen
441   // Called once at the end of the query
442
443   fHistPtot = dynamic_cast<TH1F*> (GetOutputData(0));
444   if (!fHistPtot) {
445     Printf("ERROR: fHistPtot not available");
446     return;
447   }
448    
449   TCanvas *cTaskDedx = new TCanvas("AliAnalysisTaskESDDedx","Ptot",10,10,510,510);
450   cTaskDedx->cd(1)->SetLogy();
451   fHistPtot->DrawCopy("E");
452 }
453 //____________________________________________________________________//
454 Bool_t AliAnalysisTaskESDDedx::IsAccepted(AliESDtrack* track) {
455   // Checks if the track is excluded from the cuts
456
457   Int_t   fIdxInt[200];
458   Int_t   nClustersITS = track->GetITSclusters(fIdxInt);
459   Int_t   nClustersTPC = track->GetTPCclusters(fIdxInt);
460
461   Float_t chi2PerClusterITS = -1;
462   Float_t chi2PerClusterTPC = -1;
463   if (nClustersTPC!=0)
464     chi2PerClusterTPC = track->GetTPCchi2()/Float_t(nClustersTPC);
465
466   Double_t extCov[15];
467   track->GetExternalCovariance(extCov);
468
469   if (fSelTrackRemoveKink)
470     if (track->GetKinkIndex(0)>0) return kFALSE;
471   if (fSelTrackMinClustersTPC)
472     if (nClustersTPC < fSelTrackMinClustersTPC) return kFALSE;
473   if (fSelTrackMinClustersITS)
474     if (nClustersITS < fSelTrackMinClustersITS) return kFALSE;
475   if (fSelTrackMaxChi2PerClusterTPC)
476     if (chi2PerClusterTPC > fSelTrackMaxChi2PerClusterTPC) return kFALSE; 
477   if (fSelTrackMaxChi2PerClusterITS)
478     if (chi2PerClusterITS > fSelTrackMaxChi2PerClusterITS) return kFALSE; 
479
480   if (fSelTrackMaxCov11)
481     if (extCov[0]  > fSelTrackMaxCov11) return kFALSE;
482   if (fSelTrackMaxCov22)
483     if (extCov[2]  > fSelTrackMaxCov22) return kFALSE;
484   if (fSelTrackMaxCov33)
485     if (extCov[5]  > fSelTrackMaxCov33) return kFALSE;
486   if (fSelTrackMaxCov44)
487     if (extCov[9]  > fSelTrackMaxCov44) return kFALSE;
488   if(fSelTrackMaxCov55)
489     if (extCov[14] > fSelTrackMaxCov55) return kFALSE;
490   /*
491   if(fMaxSigmaToVertexFlag)
492     if(GetSigmaToVertex(track) > fMaxSigmaToVertex) return kFALSE;
493   if(fITSRefitFlag)
494     if ((track->GetStatus() & AliESDtrack::kITSrefit) == 0) return kFALSE;
495   if(fTPCRefitFlag)
496     if ((track->GetStatus() & AliESDtrack::kTPCrefit) == 0) return kFALSE;
497   if((Pt < fMinPt) || (Pt > fMaxPt)) return kFALSE;
498   */
499   return kTRUE;
500 }