]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/SPECTRA/AliAnalysisTaskStrange.cxx
Update AliResonanceKinkLikeSign::ConnectInputData (M.Gheata)
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliAnalysisTaskStrange.cxx
1
2 /**************************************************************************
3  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4  *                                                                        *
5  * Author: The ALICE Off-line Project.                                    *
6  * Contributors are mentioned in the code where appropriate.              *
7  *                                                                        *
8  * Permission to use, copy, modify and distribute this software and its   *
9  * documentation strictly for non-commercial purposes is hereby granted   *
10  * without fee, provided that the above copyright notice appears in all   *
11  * copies and that both the copyright notice and this permission notice   *
12  * appear in the supporting documentation. The authors make no claims     *
13  * about the suitability of this software for any purpose. It is          *
14  * provided "as is" without express or implied warranty.                  *
15  **************************************************************************/
16
17 //-----------------------------------------------------------------
18 //                 AliAnalysisTaskStrange class
19 //       This task is for single strange study from ESD/AOD
20 //          Origin: H.Ricaud, Helene.Ricaud@IReS.in2p3.fr
21 //-----------------------------------------------------------------
22 #include "TList.h"
23 #include "TH1F.h"
24 #include "TH2F.h"
25 #include "TCanvas.h"
26 #include "TMath.h"
27
28 #include "AliAnalysisTaskSE.h"
29
30 #include "AliESDVertex.h"
31 #include "AliESDEvent.h"
32 #include "AliAODEvent.h"
33
34 #include "AliESDv0.h"
35 #include "AliESDtrack.h"
36 #include "AliAODv0.h"
37 #include "AliAODTrack.h"
38
39 #include "AliLog.h"
40
41 #include "AliAnalysisTaskStrange.h"
42
43 ClassImp(AliAnalysisTaskStrange)
44
45 //________________________________________________________________________
46 AliAnalysisTaskStrange::AliAnalysisTaskStrange() 
47   : AliAnalysisTaskSE(), fAnalysisType("ESD"), fCollidingSystems(0), fUseCut("infoCut"), fListHist(),
48     fHistPrimaryVertexPosX(0), fHistPrimaryVertexPosY(0), fHistPrimaryVertexPosZ(0),
49     fHistTrackMultiplicity(0), fHistV0Multiplicity(0),
50     fHistDcaPosToPrimVertex(0), fHistDcaNegToPrimVertex(0),
51     fHistDcaPosToPrimVertexZoom(0), fHistDcaNegToPrimVertexZoom(0),
52     fHistRadiusV0(0), fHistDecayLengthV0(0), fHistDcaV0Daughters(0), fHistChi2(0),
53     fHistCosPointAngle(0), fHistCosPointAngleZoom(0),
54     fHistV0MultiplicityOff(0),
55     fHistPtVsYK0sOff(0), fHistPtVsYLambdaOff(0), fHistPtVsYAntiLambdaOff(0),
56     fHistMassK0sOff(0), fHistMassLambdaOff(0), fHistMassAntiLambdaOff(0),
57     fHistMassVsRadiusK0sOff(0), fHistMassVsRadiusLambdaOff(0), fHistMassVsRadiusAntiLambdaOff(0),
58     fHistPtVsMassK0sOff(0), fHistPtVsMassLambdaOff(0), fHistPtVsMassAntiLambdaOff(0),
59     fHistArmenterosPodolanskiOff(0),
60     fHistV0MultiplicityOn(0),
61     fHistPtVsYK0sOn(0), fHistPtVsYLambdaOn(0), fHistPtVsYAntiLambdaOn(0),
62     fHistMassK0sOn(0), fHistMassLambdaOn(0), fHistMassAntiLambdaOn(0),
63     fHistMassVsRadiusK0sOn(0), fHistMassVsRadiusLambdaOn(0), fHistMassVsRadiusAntiLambdaOn(0),
64     fHistPtVsMassK0sOn(0), fHistPtVsMassLambdaOn(0), fHistPtVsMassAntiLambdaOn(0),
65     fHistArmenterosPodolanskiOn(0)
66 {
67   // Dummy constructor
68 }
69 //________________________________________________________________________
70 AliAnalysisTaskStrange::AliAnalysisTaskStrange(const char *name) 
71   : AliAnalysisTaskSE(name), fAnalysisType("ESD"), fCollidingSystems(0), fUseCut("infocut"), fListHist(),
72     fHistPrimaryVertexPosX(0), fHistPrimaryVertexPosY(0), fHistPrimaryVertexPosZ(0),
73     fHistTrackMultiplicity(0), fHistV0Multiplicity(0),
74     fHistDcaPosToPrimVertex(0), fHistDcaNegToPrimVertex(0),
75     fHistDcaPosToPrimVertexZoom(0), fHistDcaNegToPrimVertexZoom(0),
76     fHistRadiusV0(0), fHistDecayLengthV0(0), fHistDcaV0Daughters(0), fHistChi2(0),
77     fHistCosPointAngle(0), fHistCosPointAngleZoom(0),
78     fHistV0MultiplicityOff(0),
79     fHistPtVsYK0sOff(0), fHistPtVsYLambdaOff(0), fHistPtVsYAntiLambdaOff(0),
80     fHistMassK0sOff(0), fHistMassLambdaOff(0), fHistMassAntiLambdaOff(0),
81     fHistMassVsRadiusK0sOff(0), fHistMassVsRadiusLambdaOff(0), fHistMassVsRadiusAntiLambdaOff(0),
82     fHistPtVsMassK0sOff(0), fHistPtVsMassLambdaOff(0), fHistPtVsMassAntiLambdaOff(0),
83     fHistArmenterosPodolanskiOff(0),
84     fHistV0MultiplicityOn(0),
85     fHistPtVsYK0sOn(0), fHistPtVsYLambdaOn(0), fHistPtVsYAntiLambdaOn(0),
86     fHistMassK0sOn(0), fHistMassLambdaOn(0), fHistMassAntiLambdaOn(0),
87     fHistMassVsRadiusK0sOn(0), fHistMassVsRadiusLambdaOn(0), fHistMassVsRadiusAntiLambdaOn(0),
88     fHistPtVsMassK0sOn(0), fHistPtVsMassLambdaOn(0), fHistPtVsMassAntiLambdaOn(0),
89     fHistArmenterosPodolanskiOn(0)
90 {
91   // Constructor
92   // Define output slots only here
93   // Output slot #1 writes into a TList container
94   DefineOutput(1, TList::Class());
95 }
96
97 //________________________________________________________________________
98 void AliAnalysisTaskStrange::UserCreateOutputObjects() 
99 {
100   // Create histograms
101   // Called once
102
103   fListHist = new TList();
104
105   // Primary Vertex:
106   fHistPrimaryVertexPosX       = new TH1F("h1PrimaryVertexPosX", "Primary Vertex Position X;Primary Vertex Position X (cm);Events",100,-0.5,0.5);
107   fListHist->Add(fHistPrimaryVertexPosX);
108   fHistPrimaryVertexPosY       = new TH1F("h1PrimaryVertexPosY", "Primary Vertex Position Y;Primary Vertex Position Y (cm);Events",100,-0.5,0.5);
109   fListHist->Add(fHistPrimaryVertexPosY);
110   fHistPrimaryVertexPosZ       = new TH1F("h1PrimaryVertexPosZ", "Primary Vertex Position Z;Primary Vertex Position Z (cm);Events",200,-2.0,2.0);
111   fListHist->Add(fHistPrimaryVertexPosZ);
112
113   // Multiplicity:
114   if (!fHistTrackMultiplicity) {
115     if (fCollidingSystems)
116       fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 200, 0, 40000);
117     else
118       fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 250, 0, 250);
119     fListHist->Add(fHistTrackMultiplicity);
120   }
121   if (!fHistV0Multiplicity) {
122     if (fCollidingSystems)
123       fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
124     else
125       fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
126     fListHist->Add(fHistV0Multiplicity);
127   }
128
129   // Selection checks:
130   fHistDcaPosToPrimVertex      = new TH2F("h2DcaPosToPrimVertex", "Positive V0 daughter;dca(cm);Status",500,0,5,2,-0.5,1.5);
131   fListHist->Add(fHistDcaPosToPrimVertex);
132   fHistDcaNegToPrimVertex      = new TH2F("h2DcaNegToPrimVertex", "Negative V0 daughter;dca(cm);Status",500,0,5,2,-0.5,1.5);
133   fListHist->Add(fHistDcaNegToPrimVertex);
134   fHistDcaPosToPrimVertexZoom  = new TH2F("h2DcaPosToPrimVertexZoom", "Positive V0 daughter;dca(cm);Status",100,0,0.1,2,-0.5,1.5);
135   fListHist->Add(fHistDcaPosToPrimVertexZoom);
136   fHistDcaNegToPrimVertexZoom  = new TH2F("h2DcaNegToPrimVertexZoom", "Negative V0 daughter;dca(cm);Status",100,0,0.1,2,-0.5,1.5);
137   fListHist->Add(fHistDcaNegToPrimVertexZoom);
138   fHistRadiusV0                = new TH2F("h2RadiusV0", "Radius;Radius(cm);Status",1000,0,100,2,-0.5,1.5);
139   fListHist->Add(fHistRadiusV0);
140   fHistDecayLengthV0           = new TH2F("h2DecayLengthV0", "V0s decay Length;decay length(cm);Status", 200, 0, 100,2,-0.5,1.5);
141   fListHist->Add(fHistDecayLengthV0);
142   fHistDcaV0Daughters          = new TH2F("h2DcaV0Daughters", "DCA between daughters;dca(cm);Status", 160, 0, 4,2,-0.5,1.5);
143   fListHist->Add(fHistDcaV0Daughters);
144   fHistChi2                    = new TH2F("h2Chi2", "V0s chi2;chi2;Status", 33, 0, 33,2,-0.5,1.5);
145   fListHist->Add(fHistChi2);
146   fHistCosPointAngle           = new TH2F("h2CosPointAngle", "Cosine of V0's pointing angle", 100,0,1,2,-0.5,1.5);
147   fListHist->Add(fHistCosPointAngle);
148   fHistCosPointAngleZoom       = new TH2F("h2CosPointAngleZoom", "Cosine of V0's pointing angle", 100,0.9,1,2,-0.5,1.5);
149   fListHist->Add(fHistCosPointAngleZoom);
150
151   // bounds of histograms:
152   // Radius
153   const Double_t radius[10] = {0.0,2.5,2.9,3.9,7.6,15.0,23.9,37.8,42.8,100.0};
154   Int_t nBinRadius        = 9;
155
156   // V0 offline distributions
157   if (!fHistV0MultiplicityOff) {
158     if (fCollidingSystems)
159       fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
160     else
161       fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50); 
162     fListHist->Add(fHistV0MultiplicityOff);
163   }
164   // Pt vs rapidity:
165   fHistPtVsYK0sOff             = new TH2F("h2PtVsYK0sOff", "K^{0} Offline candidates;p_{t} (GeV/c);rapidity",30,0,15,30,-1.5,1.5);
166   fListHist->Add(fHistPtVsYK0sOff);
167   fHistPtVsYLambdaOff          = new TH2F("h2PtVsYLambdaOff", "#Lambda^{0} Offline candidates;p_{t} (GeV/c);rapidity",30,0,15,30,-1.5,1.5);
168   fListHist->Add(fHistPtVsYLambdaOff);
169   fHistPtVsYAntiLambdaOff      = new TH2F("h2PtVsYAntiLambdaOff", "#bar{#Lambda}^{0} Offline candidates;p_{t} (GeV/c);rapidity",30,0,15,30,-1.5,1.5);
170   fListHist->Add(fHistPtVsYAntiLambdaOff);
171   // Mass:
172   fHistMassK0sOff               = new TH1F("h1MassK0sOff", "K^{0} Offline candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});Counts", 100, 0.4, 0.6);
173   fListHist->Add(fHistMassK0sOff);
174   fHistMassLambdaOff            = new TH1F("h1MassLambdaOff", "#Lambda^{0} Offline candidates;M(p#pi^{-}) (GeV/c^{2});Counts", 75, 1.05, 1.2);
175   fListHist->Add(fHistMassLambdaOff);
176   fHistMassAntiLambdaOff          = new TH1F("h1MassAntiLambdaOff", "#bar{#Lambda}^{0} Offline candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts", 75, 1.05, 1.2);
177   fListHist->Add(fHistMassAntiLambdaOff);
178   // Mass vs radius:
179   fHistMassVsRadiusK0sOff           = new TH2F("h2MassVsRadiusK0sOff", "K^{0} Offline candidates;radius (cm);M(#pi^{+}#pi^{-}) (GeV/c^{2})",nBinRadius,radius, 200, 0.4, 0.6);
180   fListHist->Add(fHistMassVsRadiusK0sOff);
181   fHistMassVsRadiusLambdaOff       = new TH2F("h2MassVsRadiusLambdaOff", "#Lambda Offline candidates;radius (cm);M(p#pi^{-}) (GeV/c^{2})",nBinRadius,radius, 140, 1.06, 1.2);
182   fListHist->Add(fHistMassVsRadiusLambdaOff);
183   fHistMassVsRadiusAntiLambdaOff = new TH2F("h2MassVsRadiusAntiLambdaOff", "#bar{#Lambda} Offline candidates;radius (cm);M(#bar{p}#pi^{+}) (GeV/c^{2})",nBinRadius,radius, 140, 1.06, 1.2);
184   fListHist->Add(fHistMassVsRadiusAntiLambdaOff);
185   // Pt Vs Mass:
186   fHistPtVsMassK0sOff             = new TH2F("h2PtVsMassK0sOff","K^{0} Offline candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});p_{t} (GeV/c)",200, 0.4, 0.6,100,0,10);
187   fListHist->Add(fHistPtVsMassK0sOff);
188   fHistPtVsMassLambdaOff         = new TH2F("h2PtVsMassLambdaOff","#Lambda^{0} Offline candidates;M(p#pi^{-}) (GeV/c^{2});p_{t} (GeV/c)",140, 1.06, 1.2,100,0,10);
189   fListHist->Add(fHistPtVsMassLambdaOff);
190   fHistPtVsMassAntiLambdaOff     = new TH2F("h2PtVsMassAntiLambdaOff","#bar{#Lambda}^{0} Offline candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});p_{t} (GeV/c)",140, 1.06, 1.2,100,0,10);
191   fListHist->Add(fHistPtVsMassAntiLambdaOff);
192   //ArmenterosPodolanski:
193   fHistArmenterosPodolanskiOff   = new TH2F("h2ArmenterosPodolanskiOff","Armenteros-Podolanski Offline phase space;#alpha;p_{t} arm",100,-1.0,1.0,50,0,0.5);
194   fListHist->Add(fHistArmenterosPodolanskiOff);
195
196   // V0 on-the-fly distributions
197   if (!fHistV0MultiplicityOn) {
198     if (fCollidingSystems)
199       fHistV0MultiplicityOn = new TH1F("fHistV0MultiplicityOn", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
200     else
201       fHistV0MultiplicityOn = new TH1F("fHistV0MultiplicityOn", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
202     fListHist->Add(fHistV0MultiplicityOn);
203   }
204   // Pt vs rapidity:
205   fHistPtVsYK0sOn              = new TH2F("h2PtVsYK0sOn", "K^{0} Onthefly candidates;p_{t} (GeV/c);rapidity",30,0,15,30,-1.5,1.5);
206   fListHist->Add(fHistPtVsYK0sOn);
207   fHistPtVsYLambdaOn           = new TH2F("h2PtVsYLambdaOn", "#Lambda^{0} Onthefly candidates;p_{t} (GeV/c);rapidity",30,0,15,30,-1.5,1.5);
208   fListHist->Add(fHistPtVsYLambdaOn);
209   fHistPtVsYAntiLambdaOn       = new TH2F("h2PtVsYAntiLambdaOn", "#bar{#Lambda}^{0} Onthefly candidates;p_{t} (GeV/c);rapidity",30,0,15,30,-1.5,1.5);
210   fListHist->Add(fHistPtVsYAntiLambdaOn);
211   // Mass:
212   fHistMassK0sOn                = new TH1F("h1MassK0sOn", "K^{0} Onthefly candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});Counts", 100, 0.4, 0.6);
213   fListHist->Add(fHistMassK0sOn);
214   fHistMassLambdaOn            = new TH1F("h1MassLambdaOn", "#Lambda^{0} Onthefly candidates;M(p#pi^{-}) (GeV/c^{2});Counts", 75, 1.05, 1.2);
215   fListHist->Add(fHistMassLambdaOn);
216   fHistMassAntiLambdaOn        = new TH1F("h1MassAntiLambdaOn", "#bar{#Lambda}^{0} Onthefly candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts", 75, 1.05, 1.2);
217   fListHist->Add(fHistMassAntiLambdaOn);
218   // Mass vs radius:
219   fHistMassVsRadiusK0sOn         = new TH2F("h2MassVsRadiusK0sOn", "K^{0} Onthefly candidates;radius (cm);M(#pi^{+}#pi^{-}) (GeV/c^{2})",nBinRadius,radius, 200, 0.4, 0.6);
220   fListHist->Add(fHistMassVsRadiusK0sOn);
221   fHistMassVsRadiusLambdaOn     = new TH2F("h2MassVsRadiusLambdaOn", "#Lambda Onthefly candidates;radius (cm);M(p#pi^{-}) (GeV/c^{2})",nBinRadius,radius, 140, 1.06, 1.2);
222   fListHist->Add(fHistMassVsRadiusLambdaOn);
223   fHistMassVsRadiusAntiLambdaOn  = new TH2F("h2MassVsRadiusAntiLambdaOn", "#bar{#Lambda} Onthefly candidates;radius (cm);M(#bar{p}#pi^{+}) (GeV/c^{2})",nBinRadius,radius, 140, 1.06, 1.2);
224   fListHist->Add(fHistMassVsRadiusAntiLambdaOn);
225   // Pt Vs Mass:
226   fHistPtVsMassK0sOn              = new TH2F("h2PtVsMassK0sOn","K^{0} Onthefly candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});p_{t} (GeV/c)",200, 0.4, 0.6,100,0,10);
227   fListHist->Add(fHistPtVsMassK0sOn);
228   fHistPtVsMassLambdaOn          = new TH2F("h2PtVsMassLambdaOn","#Lambda^{0} Onthefly candidates;M(p#pi^{-}) (GeV/c^{2});p_{t} (GeV/c)",140, 1.06, 1.2,100,0,10);
229   fListHist->Add(fHistPtVsMassLambdaOn);
230   fHistPtVsMassAntiLambdaOn      = new TH2F("h2PtVsMassAntiLambdaOn","#bar{#Lambda}^{0} Onthefly candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});p_{t} (GeV/c)",140, 1.06, 1.2,100,0,10);
231   fListHist->Add(fHistPtVsMassAntiLambdaOn);
232   //ArmenterosPodolanski:
233   fHistArmenterosPodolanskiOn    = new TH2F("h2ArmenterosPodolanskiOn","Armenteros-Podolanski Onthefly phase space;#alpha;p_{t} arm",100,-1.0,1.0,50,0,0.5);
234   fListHist->Add(fHistArmenterosPodolanskiOn);
235 }
236
237 //________________________________________________________________________
238 void AliAnalysisTaskStrange::UserExec(Option_t *) 
239 {
240   // Main loop
241   // Called for each event
242   AliVEvent* lEvent = InputEvent();
243   if (!lEvent) {
244     Printf("ERROR: Event not available");
245     return;
246   }
247
248   if (!(lEvent->GetNumberOfTracks())) {
249     //Printf("Strange analysis task: There is no track in this event");
250     return;
251   }
252   fHistTrackMultiplicity->Fill(lEvent->GetNumberOfTracks());
253
254   Double_t tPrimaryVtxPosition[3];
255
256   Int_t nv0s = 0;
257   nv0s = lEvent->GetNumberOfV0s();
258   //Printf("Strange analysis task: There are %d v0s in this event",nv0s);
259
260   Int_t    lOnFlyStatus = 0, nv0sOn = 0, nv0sOff = 0;
261   Double_t lChi2V0 = 0;
262   Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0;
263   Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0;
264   Double_t lV0CosineOfPointingAngle = 0;
265   Double_t lV0Radius = 0;
266   Double_t lV0DecayLength = 0;
267   Double_t lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
268   Double_t lPt       = 0, lRapK0s = 0, lRapLambda = 0;
269   Double_t lAlphaV0  = 0, lPtArmV0 = 0;
270
271   Double_t  tV0Position[3];
272
273   Double_t lMagneticField      = 999;
274
275
276   //***********************
277   // ESD loop
278   //***********************
279
280   if(fAnalysisType == "ESD") {
281
282     const AliESDVertex *primaryVtx = ((AliESDEvent*)lEvent)->GetPrimaryVertex();
283     tPrimaryVtxPosition[0] = primaryVtx->GetXv();
284     tPrimaryVtxPosition[1] = primaryVtx->GetYv();
285     tPrimaryVtxPosition[2] = primaryVtx->GetZv();
286
287     fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
288     fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
289     fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
290
291     lMagneticField = ((AliESDEvent*)lEvent)->GetMagneticField();
292   
293
294     for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
295       {// This is the begining of the V0 loop  
296         AliESDv0 *v0 = ((AliESDEvent*)lEvent)->GetV0(iV0);
297         if (!v0) continue;
298
299         // AliAODtrack (V0 Daughters)
300         UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
301         UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
302
303         AliESDtrack *pTrack = ((AliESDEvent*)lEvent)->GetTrack(lKeyPos);
304         AliESDtrack *nTrack = ((AliESDEvent*)lEvent)->GetTrack(lKeyNeg);
305         if (!pTrack || !nTrack) {
306           Printf("ERROR: Could not retreive one of the daughter track");
307           continue;
308         }
309
310         // Remove like-sign
311         if ( pTrack->GetSign() == nTrack->GetSign()){
312           //cout<< "like sign, continue"<< endl;
313           continue;
314         } 
315
316         // Tracks quality cuts 
317         if ( ( (pTrack->GetTPCNcls()) < 80 ) || ( (nTrack->GetTPCNcls()) < 80 ) ) continue;
318         
319         // TPC refit condition (done during reconstruction for Offline but not for On-the-fly)
320         if( !(pTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;      
321         if( !(nTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
322
323         // DCA between daughter and Primary Vertex:
324         if (pTrack) lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],lMagneticField) );
325
326         if (nTrack) lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],lMagneticField) );
327
328         // VO's main characteristics:
329         lOnFlyStatus             = v0->GetOnFlyStatus();
330         lChi2V0                  = v0->GetChi2V0();
331         lDcaV0Daughters          = v0->GetDcaV0Daughters();
332         lDcaV0ToPrimVertex       = v0->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
333         lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1], tPrimaryVtxPosition[2]);
334         v0->GetXYZ(tV0Position[0], tV0Position[1], tV0Position[2]);
335         lV0Radius      = TMath::Sqrt(tV0Position[0]*tV0Position[0]+tV0Position[1]*tV0Position[1]);
336         lV0DecayLength = TMath::Sqrt(TMath::Power(tV0Position[0] - tPrimaryVtxPosition[0],2) +
337                                      TMath::Power(tV0Position[1] - tPrimaryVtxPosition[1],2) +
338                                      TMath::Power(tV0Position[2] - tPrimaryVtxPosition[2],2 ));
339
340         // Invariant mass
341         v0->ChangeMassHypothesis(310);
342         lInvMassK0s = v0->GetEffMass();
343         v0->ChangeMassHypothesis(3122);
344         lInvMassLambda = v0->GetEffMass();
345         v0->ChangeMassHypothesis(-3122);
346         lInvMassAntiLambda = v0->GetEffMass();
347
348         // Rapidity:
349         lRapK0s    = v0->Y(310);
350         lRapLambda = v0->Y(3122);
351         
352         // Pt:
353         lPt = v0->Pt();
354         
355         // Armenteros variables: !!
356         lAlphaV0      = v0->AlphaV0();
357         lPtArmV0      = v0->PtArmV0();
358         
359         // Selections:
360         if (fUseCut.Contains("yes")) {
361           if ( (lDcaPosToPrimVertex      < 0.036 )||
362                (lDcaNegToPrimVertex      < 0.036 )||
363                (lDcaV0Daughters          > 0.5 )  ||
364                (lV0CosineOfPointingAngle < 0.99)     
365                ) continue;
366         }
367     
368         // Filling histograms
369         fHistDcaPosToPrimVertex->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
370         fHistDcaNegToPrimVertex->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
371         fHistDcaPosToPrimVertexZoom->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
372         fHistDcaNegToPrimVertexZoom->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
373         fHistRadiusV0->Fill(lV0Radius,lOnFlyStatus);
374         fHistDecayLengthV0->Fill(lV0DecayLength,lOnFlyStatus);
375         fHistDcaV0Daughters->Fill(lDcaV0Daughters,lOnFlyStatus);
376         fHistChi2->Fill(lChi2V0,lOnFlyStatus);
377         fHistCosPointAngle->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
378         if (lV0CosineOfPointingAngle >= 0.9) fHistCosPointAngleZoom->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
379         if(!lOnFlyStatus){
380           nv0sOff++;
381           fHistPtVsYK0sOff->Fill(lPt,lRapK0s);
382           fHistPtVsYLambdaOff->Fill(lPt,lRapLambda);
383           fHistPtVsYAntiLambdaOff->Fill(lPt,lRapLambda);
384           fHistArmenterosPodolanskiOff->Fill(lAlphaV0,lPtArmV0);
385         }
386         else {
387           nv0sOn++;
388           fHistPtVsYK0sOn->Fill(lPt,lRapK0s);
389           fHistPtVsYLambdaOn->Fill(lPt,lRapLambda);
390           fHistPtVsYAntiLambdaOn->Fill(lPt,lRapLambda);
391           fHistArmenterosPodolanskiOn->Fill(lAlphaV0,lPtArmV0);
392         }
393         // K0s invariant mass histograms:
394         if (TMath::Abs(lRapK0s) < 1) {  
395           if(!lOnFlyStatus){
396             fHistMassK0sOff->Fill(lInvMassK0s);
397             fHistMassVsRadiusK0sOff->Fill(lV0Radius,lInvMassK0s);
398             fHistPtVsMassK0sOff->Fill(lInvMassK0s,lPt);
399           }
400           else {
401             fHistMassK0sOn->Fill(lInvMassK0s);
402             fHistMassVsRadiusK0sOn->Fill(lV0Radius,lInvMassK0s);
403             fHistPtVsMassK0sOn->Fill(lInvMassK0s,lPt);
404           }
405         }
406         // Lambda and AntiLambda invariant mass histograms:
407         if (TMath::Abs(lRapLambda) < 1) {
408           if(!lOnFlyStatus){
409             fHistMassLambdaOff->Fill(lInvMassLambda);
410             fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
411             fHistMassVsRadiusLambdaOff->Fill(lV0Radius,lInvMassLambda);
412             fHistMassVsRadiusAntiLambdaOff->Fill(lV0Radius,lInvMassAntiLambda);
413             fHistPtVsMassLambdaOff->Fill(lInvMassLambda,lPt);
414             fHistPtVsMassAntiLambdaOff->Fill(lInvMassAntiLambda,lPt);
415           }
416           else {
417             fHistMassLambdaOn->Fill(lInvMassLambda);
418             fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
419             fHistMassVsRadiusLambdaOn->Fill(lV0Radius,lInvMassLambda);
420             fHistMassVsRadiusAntiLambdaOn->Fill(lV0Radius,lInvMassAntiLambda);
421             fHistPtVsMassLambdaOn->Fill(lInvMassLambda,lPt);
422             fHistPtVsMassAntiLambdaOn->Fill(lInvMassAntiLambda,lPt);
423           }
424         }
425       } // end V0 loop
426
427   }
428
429   //***********************
430   // AOD loop
431   //***********************
432
433   else if(fAnalysisType == "AOD") {
434
435     const AliAODVertex *primaryVtx = ((AliAODEvent*)lEvent)->GetPrimaryVertex();
436     tPrimaryVtxPosition[0] = primaryVtx->GetX();
437     tPrimaryVtxPosition[1] = primaryVtx->GetY();
438     tPrimaryVtxPosition[2] = primaryVtx->GetZ();
439
440     fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
441     fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
442     fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
443   
444     for (Int_t iV0 = 0; iV0 < nv0s; iV0++) 
445       {// This is the begining of the V0 loop
446         AliAODv0 *myAODv0 = ((AliAODEvent*)lEvent)->GetV0(iV0);
447         if (!myAODv0) continue;
448
449         // common part
450         lV0Radius                = myAODv0->RadiusV0();
451         lDcaPosToPrimVertex      = myAODv0->DcaPosToPrimVertex();
452         lDcaNegToPrimVertex      = myAODv0->DcaNegToPrimVertex();
453         lOnFlyStatus             = myAODv0->GetOnFlyStatus();
454         lChi2V0                  = myAODv0->Chi2V0();
455         lDcaV0Daughters          = myAODv0->DcaV0Daughters();
456         lDcaV0ToPrimVertex       = myAODv0->DcaV0ToPrimVertex();
457         lV0DecayLength           = myAODv0->DecayLengthV0(tPrimaryVtxPosition);
458         lV0CosineOfPointingAngle = myAODv0->CosPointingAngle(tPrimaryVtxPosition);
459
460         lInvMassK0s        = myAODv0->MassK0Short();
461         lInvMassLambda     = myAODv0->MassLambda();
462         lInvMassAntiLambda = myAODv0->MassAntiLambda();
463
464         lPt        = TMath::Sqrt(myAODv0->Pt2V0());
465         lRapK0s    = myAODv0->RapK0Short();
466         lRapLambda = myAODv0->RapLambda();
467         lAlphaV0   = myAODv0->AlphaV0();
468         lPtArmV0   = myAODv0->PtArmV0();
469
470
471         // Selections:
472         if (fUseCut.Contains("yes")) {
473           if ( (lDcaPosToPrimVertex      < 0.036 )||
474                (lDcaNegToPrimVertex      < 0.036 )||
475                (lDcaV0Daughters          > 0.5 )  ||
476                (lV0CosineOfPointingAngle < 0.99)     
477                ) continue;
478         }
479     
480         // Filling histograms
481         fHistDcaPosToPrimVertex->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
482         fHistDcaNegToPrimVertex->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
483         fHistDcaPosToPrimVertexZoom->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
484         fHistDcaNegToPrimVertexZoom->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
485         fHistRadiusV0->Fill(lV0Radius,lOnFlyStatus);
486         fHistDecayLengthV0->Fill(lV0DecayLength,lOnFlyStatus);
487         fHistDcaV0Daughters->Fill(lDcaV0Daughters,lOnFlyStatus);
488         fHistChi2->Fill(lChi2V0,lOnFlyStatus);
489         fHistCosPointAngle->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
490         if (lV0CosineOfPointingAngle >= 0.9) fHistCosPointAngleZoom->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
491         if(!lOnFlyStatus){
492           nv0sOff++;
493           fHistPtVsYK0sOff->Fill(lPt,lRapK0s);
494           fHistPtVsYLambdaOff->Fill(lPt,lRapLambda);
495           fHistPtVsYAntiLambdaOff->Fill(lPt,lRapLambda);
496           fHistArmenterosPodolanskiOff->Fill(lAlphaV0,lPtArmV0);
497         }
498         else {
499           nv0sOn++;
500           fHistPtVsYK0sOn->Fill(lPt,lRapK0s);
501           fHistPtVsYLambdaOn->Fill(lPt,lRapLambda);
502           fHistPtVsYAntiLambdaOn->Fill(lPt,lRapLambda);
503           fHistArmenterosPodolanskiOn->Fill(lAlphaV0,lPtArmV0);
504         }
505         // K0s invariant mass histograms:
506         if (TMath::Abs(lRapK0s) < 1) {  
507           if(!lOnFlyStatus){
508             fHistMassK0sOff->Fill(lInvMassK0s);
509             fHistMassVsRadiusK0sOff->Fill(lV0Radius,lInvMassK0s);
510             fHistPtVsMassK0sOff->Fill(lInvMassK0s,lPt);
511           }
512           else {
513             fHistMassK0sOn->Fill(lInvMassK0s);
514             fHistMassVsRadiusK0sOn->Fill(lV0Radius,lInvMassK0s);
515             fHistPtVsMassK0sOn->Fill(lInvMassK0s,lPt);
516           }
517         }
518         // Lambda and AntiLambda invariant mass histograms:
519         if (TMath::Abs(lRapLambda) < 1) {
520           if(!lOnFlyStatus){
521             fHistMassLambdaOff->Fill(lInvMassLambda);
522             fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
523             fHistMassVsRadiusLambdaOff->Fill(lV0Radius,lInvMassLambda);
524             fHistMassVsRadiusAntiLambdaOff->Fill(lV0Radius,lInvMassAntiLambda);
525             fHistPtVsMassLambdaOff->Fill(lInvMassLambda,lPt);
526             fHistPtVsMassAntiLambdaOff->Fill(lInvMassAntiLambda,lPt);
527           }
528           else {
529             fHistMassLambdaOn->Fill(lInvMassLambda);
530             fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
531             fHistMassVsRadiusLambdaOn->Fill(lV0Radius,lInvMassLambda);
532             fHistMassVsRadiusAntiLambdaOn->Fill(lV0Radius,lInvMassAntiLambda);
533             fHistPtVsMassLambdaOn->Fill(lInvMassLambda,lPt);
534             fHistPtVsMassAntiLambdaOn->Fill(lInvMassAntiLambda,lPt);
535           }
536         }
537       } // end V0 loop
538   }
539
540   fHistV0Multiplicity->Fill(nv0s);
541   fHistV0MultiplicityOff->Fill(nv0sOff);
542   fHistV0MultiplicityOn->Fill(nv0sOn);
543
544   // Post output data.
545   PostData(1, fListHist);
546 }    
547
548 //________________________________________________________________________
549 void AliAnalysisTaskStrange::Terminate(Option_t *) 
550 {
551   // Draw result to the screen
552   // Called once at the end of the query
553 }