]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/SPECTRA/AliAnalysisTaskStrange.cxx
Updates of the macro that dispays the results (efficiency plots) - Efficiency calcula...
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliAnalysisTaskStrange.cxx
CommitLineData
9ea746fc 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#include "AliAODVertex.h"
34
35#include "AliESDv0.h"
36#include "AliESDtrack.h"
37#include "AliAODv0.h"
38#include "AliAODTrack.h"
39
40#include "AliLog.h"
41
42#include "AliAnalysisTaskStrange.h"
43
44ClassImp(AliAnalysisTaskStrange)
45
46//________________________________________________________________________
47AliAnalysisTaskStrange::AliAnalysisTaskStrange()
48 : AliAnalysisTaskSE(), fAnalysisType("ESD"), fCollidingSystems(0), fOption("no"), fListHist(),
49 fHistPrimaryVertexPosX(0), fHistPrimaryVertexPosY(0), fHistPrimaryVertexPosZ(0),
50 fHistTrackMultiplicity(0), fHistV0Multiplicity(0),
51 fHistDcaPosToPrimVertex(0), fHistDcaNegToPrimVertex(0),
52 fHistDcaPosToPrimVertexZoom(0), fHistDcaNegToPrimVertexZoom(0),
53 fHistRadiusV0(0), fHistDecayLengthV0(0), fHistDcaV0Daughters(0), fHistChi2(0),
54 fHistCosPointAngle(0), fHistCosPointAngleZoom(0),
55 fHistV0MultiplicityOff(0),
56 fHistPtVsYK0sOff(0), fHistPtVsYLambdaOff(0), fHistPtVsYAntiLambdaOff(0),
57 fHistMassK0sOff(0), fHistMassLambdaOff(0), fHistMassAntiLambdaOff(0),
58 fHistMassVsRadiusK0sOff(0), fHistMassVsRadiusLambdaOff(0), fHistMassVsRadiusAntiLambdaOff(0),
59 fHistPtVsMassK0sOff(0), fHistPtVsMassLambdaOff(0), fHistPtVsMassAntiLambdaOff(0),
60 fHistArmenterosPodolanskiOff(0),
61 fHistV0MultiplicityOn(0),
62 fHistPtVsYK0sOn(0), fHistPtVsYLambdaOn(0), fHistPtVsYAntiLambdaOn(0),
63 fHistMassK0sOn(0), fHistMassLambdaOn(0), fHistMassAntiLambdaOn(0),
64 fHistMassVsRadiusK0sOn(0), fHistMassVsRadiusLambdaOn(0), fHistMassVsRadiusAntiLambdaOn(0),
65 fHistPtVsMassK0sOn(0), fHistPtVsMassLambdaOn(0), fHistPtVsMassAntiLambdaOn(0),
66 fHistArmenterosPodolanskiOn(0)
67{
68 // Dummy constructor
69}
70//________________________________________________________________________
71AliAnalysisTaskStrange::AliAnalysisTaskStrange(const char *name, const char *optCuts)
72 : AliAnalysisTaskSE(name), fAnalysisType("ESD"), fCollidingSystems(0), fOption(optCuts), fListHist(),
73 fHistPrimaryVertexPosX(0), fHistPrimaryVertexPosY(0), fHistPrimaryVertexPosZ(0),
74 fHistTrackMultiplicity(0), fHistV0Multiplicity(0),
75 fHistDcaPosToPrimVertex(0), fHistDcaNegToPrimVertex(0),
76 fHistDcaPosToPrimVertexZoom(0), fHistDcaNegToPrimVertexZoom(0),
77 fHistRadiusV0(0), fHistDecayLengthV0(0), fHistDcaV0Daughters(0), fHistChi2(0),
78 fHistCosPointAngle(0), fHistCosPointAngleZoom(0),
79 fHistV0MultiplicityOff(0),
80 fHistPtVsYK0sOff(0), fHistPtVsYLambdaOff(0), fHistPtVsYAntiLambdaOff(0),
81 fHistMassK0sOff(0), fHistMassLambdaOff(0), fHistMassAntiLambdaOff(0),
82 fHistMassVsRadiusK0sOff(0), fHistMassVsRadiusLambdaOff(0), fHistMassVsRadiusAntiLambdaOff(0),
83 fHistPtVsMassK0sOff(0), fHistPtVsMassLambdaOff(0), fHistPtVsMassAntiLambdaOff(0),
84 fHistArmenterosPodolanskiOff(0),
85 fHistV0MultiplicityOn(0),
86 fHistPtVsYK0sOn(0), fHistPtVsYLambdaOn(0), fHistPtVsYAntiLambdaOn(0),
87 fHistMassK0sOn(0), fHistMassLambdaOn(0), fHistMassAntiLambdaOn(0),
88 fHistMassVsRadiusK0sOn(0), fHistMassVsRadiusLambdaOn(0), fHistMassVsRadiusAntiLambdaOn(0),
89 fHistPtVsMassK0sOn(0), fHistPtVsMassLambdaOn(0), fHistPtVsMassAntiLambdaOn(0),
90 fHistArmenterosPodolanskiOn(0)
91{
92 // Constructor
93 // Define output slots only here
94 // Output slot #1 writes into a TList container
95 DefineOutput(1, TList::Class());
96}
97
98//________________________________________________________________________
99void AliAnalysisTaskStrange::UserCreateOutputObjects()
100{
101 // Create histograms
102 // Called once
103
104 fListHist = new TList();
105
106 // Primary Vertex:
107 fHistPrimaryVertexPosX = new TH1F("h1PrimaryVertexPosX", "Primary Vertex Position X;Primary Vertex Position X (cm);Events",100,-0.5,0.5);
108 fListHist->Add(fHistPrimaryVertexPosX);
109 fHistPrimaryVertexPosY = new TH1F("h1PrimaryVertexPosY", "Primary Vertex Position Y;Primary Vertex Position Y (cm);Events",100,-0.5,0.5);
110 fListHist->Add(fHistPrimaryVertexPosY);
111 fHistPrimaryVertexPosZ = new TH1F("h1PrimaryVertexPosZ", "Primary Vertex Position Z;Primary Vertex Position Z (cm);Events",200,-2.0,2.0);
112 fListHist->Add(fHistPrimaryVertexPosZ);
113
114 // Multiplicity:
115 if (!fHistTrackMultiplicity) {
116 if (fCollidingSystems)
117 fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 200, 0, 40000);
118 else
119 fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 250, 0, 250);
120 fListHist->Add(fHistTrackMultiplicity);
121 }
122 if (!fHistV0Multiplicity) {
123 if (fCollidingSystems)
124 fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
125 else
126 fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
127 fListHist->Add(fHistV0Multiplicity);
128 }
129
130 // Selection checks:
131 fHistDcaPosToPrimVertex = new TH2F("h2DcaPosToPrimVertex", "Positive V0 daughter;dca(cm);Status",500,0,5,2,-0.5,1.5);
132 fListHist->Add(fHistDcaPosToPrimVertex);
133 fHistDcaNegToPrimVertex = new TH2F("h2DcaNegToPrimVertex", "Negative V0 daughter;dca(cm);Status",500,0,5,2,-0.5,1.5);
134 fListHist->Add(fHistDcaNegToPrimVertex);
135 fHistDcaPosToPrimVertexZoom = new TH2F("h2DcaPosToPrimVertexZoom", "Positive V0 daughter;dca(cm);Status",100,0,0.1,2,-0.5,1.5);
136 fListHist->Add(fHistDcaPosToPrimVertexZoom);
137 fHistDcaNegToPrimVertexZoom = new TH2F("h2DcaNegToPrimVertexZoom", "Negative V0 daughter;dca(cm);Status",100,0,0.1,2,-0.5,1.5);
138 fListHist->Add(fHistDcaNegToPrimVertexZoom);
139 fHistRadiusV0 = new TH2F("h2RadiusV0", "Radius;Radius(cm);Status",1000,0,100,2,-0.5,1.5);
140 fListHist->Add(fHistRadiusV0);
141 fHistDecayLengthV0 = new TH2F("h2DecayLengthV0", "V0s decay Length;decay length(cm);Status", 200, 0, 100,2,-0.5,1.5);
142 fListHist->Add(fHistDecayLengthV0);
143 fHistDcaV0Daughters = new TH2F("h2DcaV0Daughters", "DCA between daughters;dca(cm);Status", 160, 0, 4,2,-0.5,1.5);
144 fListHist->Add(fHistDcaV0Daughters);
145 fHistChi2 = new TH2F("h2Chi2", "V0s chi2;chi2;Status", 33, 0, 33,2,-0.5,1.5);
146 fListHist->Add(fHistChi2);
147 fHistCosPointAngle = new TH2F("h2CosPointAngle", "Cosine of V0's pointing angle", 100,0,1,2,-0.5,1.5);
148 fListHist->Add(fHistCosPointAngle);
149 fHistCosPointAngleZoom = new TH2F("h2CosPointAngleZoom", "Cosine of V0's pointing angle", 100,0.9,1,2,-0.5,1.5);
150 fListHist->Add(fHistCosPointAngleZoom);
151
152 // bounds of histograms:
153 // Radius
154 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};
155 Int_t NbinRadius = 9;
156
157 // V0 offline distributions
158 if (!fHistV0MultiplicityOff) {
159 if (fCollidingSystems)
160 fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
161 else
162 fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
163 fListHist->Add(fHistV0MultiplicityOff);
164 }
165 // Pt vs rapidity:
166 fHistPtVsYK0sOff = new TH2F("h2PtVsYK0sOff", "K^{0} Offline candidates;p_{t} (GeV/c);rapidity",30,0,15,30,-1.5,1.5);
167 fListHist->Add(fHistPtVsYK0sOff);
168 fHistPtVsYLambdaOff = new TH2F("h2PtVsYLambdaOff", "#Lambda^{0} Offline candidates;p_{t} (GeV/c);rapidity",30,0,15,30,-1.5,1.5);
169 fListHist->Add(fHistPtVsYLambdaOff);
170 fHistPtVsYAntiLambdaOff = new TH2F("h2PtVsYAntiLambdaOff", "#bar{#Lambda}^{0} Offline candidates;p_{t} (GeV/c);rapidity",30,0,15,30,-1.5,1.5);
171 fListHist->Add(fHistPtVsYAntiLambdaOff);
172 // Mass:
173 fHistMassK0sOff = new TH1F("h1MassK0sOff", "K^{0} Offline candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});Counts", 100, 0.4, 0.6);
174 fListHist->Add(fHistMassK0sOff);
175 fHistMassLambdaOff = new TH1F("h1MassLambdaOff", "#Lambda^{0} Offline candidates;M(p#pi^{-}) (GeV/c^{2});Counts", 75, 1.05, 1.2);
176 fListHist->Add(fHistMassLambdaOff);
177 fHistMassAntiLambdaOff = new TH1F("h1MassAntiLambdaOff", "#bar{#Lambda}^{0} Offline candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts", 75, 1.05, 1.2);
178 fListHist->Add(fHistMassAntiLambdaOff);
179 // Mass vs radius:
180 fHistMassVsRadiusK0sOff = new TH2F("h2MassVsRadiusK0sOff", "K^{0} Offline candidates;radius (cm);M(#pi^{+}#pi^{-}) (GeV/c^{2})",NbinRadius,radius, 200, 0.4, 0.6);
181 fListHist->Add(fHistMassVsRadiusK0sOff);
182 fHistMassVsRadiusLambdaOff = new TH2F("h2MassVsRadiusLambdaOff", "#Lambda Offline candidates;radius (cm);M(p#pi^{-}) (GeV/c^{2})",NbinRadius,radius, 140, 1.06, 1.2);
183 fListHist->Add(fHistMassVsRadiusLambdaOff);
184 fHistMassVsRadiusAntiLambdaOff = new TH2F("h2MassVsRadiusAntiLambdaOff", "#bar{#Lambda} Offline candidates;radius (cm);M(#bar{p}#pi^{+}) (GeV/c^{2})",NbinRadius,radius, 140, 1.06, 1.2);
185 fListHist->Add(fHistMassVsRadiusAntiLambdaOff);
186 // Pt Vs Mass:
187 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);
188 fListHist->Add(fHistPtVsMassK0sOff);
189 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);
190 fListHist->Add(fHistPtVsMassLambdaOff);
191 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);
192 fListHist->Add(fHistPtVsMassAntiLambdaOff);
193 //ArmenterosPodolanski:
194 fHistArmenterosPodolanskiOff = new TH2F("h2ArmenterosPodolanskiOff","Armenteros-Podolanski Offline phase space;#alpha;p_{t} arm",100,-1.0,1.0,50,0,0.5);
195 fListHist->Add(fHistArmenterosPodolanskiOff);
196
197 // V0 on-the-fly distributions
198 if (!fHistV0MultiplicityOn) {
199 if (fCollidingSystems)
200 fHistV0MultiplicityOn = new TH1F("fHistV0MultiplicityOn", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
201 else
202 fHistV0MultiplicityOn = new TH1F("fHistV0MultiplicityOn", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
203 fListHist->Add(fHistV0MultiplicityOn);
204 }
205 // Pt vs rapidity:
206 fHistPtVsYK0sOn = new TH2F("h2PtVsYK0sOn", "K^{0} Onthefly candidates;p_{t} (GeV/c);rapidity",30,0,15,30,-1.5,1.5);
207 fListHist->Add(fHistPtVsYK0sOn);
208 fHistPtVsYLambdaOn = new TH2F("h2PtVsYLambdaOn", "#Lambda^{0} Onthefly candidates;p_{t} (GeV/c);rapidity",30,0,15,30,-1.5,1.5);
209 fListHist->Add(fHistPtVsYLambdaOn);
210 fHistPtVsYAntiLambdaOn = new TH2F("h2PtVsYAntiLambdaOn", "#bar{#Lambda}^{0} Onthefly candidates;p_{t} (GeV/c);rapidity",30,0,15,30,-1.5,1.5);
211 fListHist->Add(fHistPtVsYAntiLambdaOn);
212 // Mass:
213 fHistMassK0sOn = new TH1F("h1MassK0sOn", "K^{0} Onthefly candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});Counts", 100, 0.4, 0.6);
214 fListHist->Add(fHistMassK0sOn);
215 fHistMassLambdaOn = new TH1F("h1MassLambdaOn", "#Lambda^{0} Onthefly candidates;M(p#pi^{-}) (GeV/c^{2});Counts", 75, 1.05, 1.2);
216 fListHist->Add(fHistMassLambdaOn);
217 fHistMassAntiLambdaOn = new TH1F("h1MassAntiLambdaOn", "#bar{#Lambda}^{0} Onthefly candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts", 75, 1.05, 1.2);
218 fListHist->Add(fHistMassAntiLambdaOn);
219 // Mass vs radius:
220 fHistMassVsRadiusK0sOn = new TH2F("h2MassVsRadiusK0sOn", "K^{0} Onthefly candidates;radius (cm);M(#pi^{+}#pi^{-}) (GeV/c^{2})",NbinRadius,radius, 200, 0.4, 0.6);
221 fListHist->Add(fHistMassVsRadiusK0sOn);
222 fHistMassVsRadiusLambdaOn = new TH2F("h2MassVsRadiusLambdaOn", "#Lambda Onthefly candidates;radius (cm);M(p#pi^{-}) (GeV/c^{2})",NbinRadius,radius, 140, 1.06, 1.2);
223 fListHist->Add(fHistMassVsRadiusLambdaOn);
224 fHistMassVsRadiusAntiLambdaOn = new TH2F("h2MassVsRadiusAntiLambdaOn", "#bar{#Lambda} Onthefly candidates;radius (cm);M(#bar{p}#pi^{+}) (GeV/c^{2})",NbinRadius,radius, 140, 1.06, 1.2);
225 fListHist->Add(fHistMassVsRadiusAntiLambdaOn);
226 // Pt Vs Mass:
227 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);
228 fListHist->Add(fHistPtVsMassK0sOn);
229 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);
230 fListHist->Add(fHistPtVsMassLambdaOn);
231 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);
232 fListHist->Add(fHistPtVsMassAntiLambdaOn);
233 //ArmenterosPodolanski:
234 fHistArmenterosPodolanskiOn = new TH2F("h2ArmenterosPodolanskiOn","Armenteros-Podolanski Onthefly phase space;#alpha;p_{t} arm",100,-1.0,1.0,50,0,0.5);
235 fListHist->Add(fHistArmenterosPodolanskiOn);
236}
237
238//________________________________________________________________________
239void AliAnalysisTaskStrange::UserExec(Option_t *)
240{
241 // Main loop
242 // Called for each event
243 AliVEvent* lEvent = InputEvent();
244 if (!lEvent) {
245 Printf("ERROR: Event not available");
246 return;
247 }
248
249 if (!(lEvent->GetNumberOfTracks())) {
250 Printf("Strange analysis task: There is no track in this event");
251 return;
252 }
253 fHistTrackMultiplicity->Fill(lEvent->GetNumberOfTracks());
254
255 Double_t tPrimaryVtxPosition[3];
256 Double_t tPrimaryVtxCov[3];
257 Int_t nv0s = 0;
258 nv0s = lEvent->GetNumberOfV0s();
259 Printf("Strange analysis task: There are %d v0s in this event",nv0s);
260
261 Int_t lOnFlyStatus = 0, nv0sOn = 0, nv0sOff = 0;
262 Double_t lChi2V0 = 0;
263 Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0;
264 Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0;
265 Double_t lV0CosineOfPointingAngle = 0;
266 Double_t lV0Radius = 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 if(fAnalysisType == "ESD") {
272
273 const AliESDVertex *primaryVtx = ((AliESDEvent*)lEvent)->GetPrimaryVertex();
274 tPrimaryVtxPosition[0] = primaryVtx->GetXv();
275 tPrimaryVtxPosition[1] = primaryVtx->GetYv();
276 tPrimaryVtxPosition[2] = primaryVtx->GetZv();
277
278 fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
279 fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
280 fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
281
282 primaryVtx->GetCovMatrix(tPrimaryVtxCov);
283 AliAODVertex *primary = new AliAODVertex(tPrimaryVtxPosition, tPrimaryVtxCov, primaryVtx->GetChi2toNDF(), NULL, -1, AliAODVertex::kPrimary);
284
285 // V0 variables:
286 // to get info from ESD files and fill AliAODVertex:
287 Float_t tdcaPosToPrimVertexXYZ[2], tdcaNegToPrimVertexXYZ[2]; // ..[0] = Impact parameter in XY plane and ..[1] = Impact parameter in Z
288 Double_t tdcaDaughterToPrimVertex[2]; // ..[0] = Pos and ..[1] = Neg
289 Double_t tMomPos[3];
290 Double_t tMomNeg[3];
291 Double_t V0Position[3];
292 Double_t V0Cov[6];
293
294 // to fill AliAODtrack:
295 Double_t TrackP[3];
296 Double_t TrackPosition[3];
297 Double_t TrackcovTr[21];
298 Double_t Trackpid[10];
299
300 AliAODTrack *myPosAodTrack = new AliAODTrack();
301 AliAODTrack *myNegAodTrack = new AliAODTrack();
302 AliAODVertex *myAODVertex = new AliAODVertex();
303 AliAODv0 *myAODv0 = new AliAODv0();
304
305 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
306 {// This is the begining of the V0 loop
307 AliESDv0 *v0 = ((AliESDEvent*)lEvent)->GetV0(iV0);
308 if (!v0) continue;
309
310 // AliAODVertex
311 v0->GetXYZ(V0Position[0], V0Position[1], V0Position[2]);
312 v0->GetPosCov(V0Cov);
313 myAODVertex->SetPosition(V0Position[0],V0Position[1],V0Position[2]);
314 myAODVertex->SetCovMatrix(V0Cov);
315 myAODVertex->SetChi2perNDF(v0->GetChi2V0());
316 myAODVertex->SetID((Short_t)iV0);
317 myAODVertex->SetParent(primary);
318 myAODVertex->SetType(AliAODVertex::kV0);
319
320 // AliAODtrack (V0 Daughters)
321 UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
322 UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
323
324 AliESDtrack *TrackPos = ((AliESDEvent*)lEvent)->GetTrack(lKeyPos);
325 AliESDtrack *TrackNeg = ((AliESDEvent*)lEvent)->GetTrack(lKeyNeg);
326 if (!TrackPos || !TrackNeg) {
327 Printf("ERROR: Could not retreive one of the daughter track");
328 continue;
329 }
330 UInt_t lLabelPos = (UInt_t)TMath::Abs(TrackPos->GetLabel());
331 UInt_t lLabelNeg = (UInt_t)TMath::Abs(TrackNeg->GetLabel());
332
333 myPosAodTrack->SetID((Short_t)(TrackPos->GetID()));
334 myPosAodTrack->SetLabel(lLabelPos);
335 TrackPos->GetPxPyPz(TrackP);
336 myPosAodTrack->SetP(TrackP);
337 TrackPos->GetXYZ(TrackPosition);
338 myPosAodTrack->SetPosition(TrackPosition,kFALSE);
339 TrackPos->GetCovarianceXYZPxPyPz(TrackcovTr);
340 myPosAodTrack->SetCovMatrix(TrackcovTr);
341 TrackPos->GetESDpid(Trackpid);
342 myPosAodTrack->SetPID(Trackpid);
343 myPosAodTrack->SetCharge((Short_t)(TrackPos->Charge()));
344 myPosAodTrack->SetITSClusterMap(TrackPos->GetITSClusterMap());
345 myPosAodTrack->SetProdVertex(myAODVertex);
346 myPosAodTrack->SetUsedForVtxFit(kTRUE);
347 myPosAodTrack->SetUsedForPrimVtxFit(kFALSE);
348 myPosAodTrack->SetType(AliAODTrack::kSecondary);
349 myPosAodTrack->ConvertAliPIDtoAODPID();
350
351 myNegAodTrack->SetID((Short_t)(TrackNeg->GetID()));
352 myNegAodTrack->SetLabel(lLabelNeg);
353 TrackNeg->GetPxPyPz(TrackP);
354 myNegAodTrack->SetP(TrackP);
355 TrackNeg->GetXYZ(TrackPosition);
356 myNegAodTrack->SetPosition(TrackPosition,kFALSE);
357 TrackNeg->GetCovarianceXYZPxPyPz(TrackcovTr);
358 myNegAodTrack->SetCovMatrix(TrackcovTr);
359 TrackNeg->GetESDpid(Trackpid);
360 myNegAodTrack->SetPID(Trackpid);
361 myNegAodTrack->SetCharge((Short_t)(TrackNeg->Charge()));
362 myNegAodTrack->SetITSClusterMap(TrackPos->GetITSClusterMap());
363 myNegAodTrack->SetProdVertex(myAODVertex);
364 myNegAodTrack->SetUsedForVtxFit(kTRUE);
365 myNegAodTrack->SetUsedForPrimVtxFit(kFALSE);
366 myNegAodTrack->SetType(AliAODTrack::kSecondary);
367 myNegAodTrack->ConvertAliPIDtoAODPID();
368
369 myAODVertex->AddDaughter(myPosAodTrack);
370 myAODVertex->AddDaughter(myNegAodTrack);
371
372 // filling myAODv0
373 lDcaV0Daughters = v0->GetDcaV0Daughters();
374 lDcaV0ToPrimVertex = v0->GetD();
375
376 if (TrackPos) TrackPos->GetImpactParameters(tdcaPosToPrimVertexXYZ[0],tdcaPosToPrimVertexXYZ[1]);
377 if (TrackNeg) TrackNeg->GetImpactParameters(tdcaNegToPrimVertexXYZ[0],tdcaNegToPrimVertexXYZ[1]);
378 tdcaDaughterToPrimVertex[0] = TMath::Sqrt(tdcaPosToPrimVertexXYZ[0]*tdcaPosToPrimVertexXYZ[0]+tdcaPosToPrimVertexXYZ[1]*tdcaPosToPrimVertexXYZ[1]);
379 tdcaDaughterToPrimVertex[1] = TMath::Sqrt(tdcaNegToPrimVertexXYZ[0]*tdcaNegToPrimVertexXYZ[0]+tdcaNegToPrimVertexXYZ[1]*tdcaNegToPrimVertexXYZ[1]);
380
381 v0->GetPPxPyPz(tMomPos[0],tMomPos[1],tMomPos[2]);
382 v0->GetNPxPyPz(tMomNeg[0],tMomNeg[1],tMomNeg[2]);
383
384 myAODv0->Fill(myAODVertex, lDcaV0Daughters, lDcaV0ToPrimVertex, tMomPos, tMomNeg, tdcaDaughterToPrimVertex);
385 myAODv0->SetOnFlyStatus(v0->GetOnFlyStatus());
386
387 // common part
388 lV0Radius = myAODv0->RadiusV0();
389 lDcaPosToPrimVertex = myAODv0->DcaPosToPrimVertex();
390 lDcaNegToPrimVertex = myAODv0->DcaNegToPrimVertex();
391
392 lOnFlyStatus = myAODv0->GetOnFlyStatus();
393 lChi2V0 = myAODv0->Chi2V0();
394 lDcaV0Daughters = myAODv0->DcaV0Daughters();
395 lDcaV0ToPrimVertex = myAODv0->DcaV0ToPrimVertex();
396 lV0CosineOfPointingAngle = myAODv0->CosPointingAngle(tPrimaryVtxPosition);
397
398 lInvMassK0s = myAODv0->MassK0Short();
399 lInvMassLambda = myAODv0->MassLambda();
400 lInvMassAntiLambda = myAODv0->MassAntiLambda();
401
402 lPt = TMath::Sqrt(myAODv0->Pt2V0());
403 lRapK0s = myAODv0->RapK0Short();
404 lRapLambda = myAODv0->RapLambda();
405 lAlphaV0 = myAODv0->AlphaV0();
406 lPtArmV0 = myAODv0->PtArmV0();
407
408
409 // Selections:
410 if (fOption.Contains("yes")) {
411 if ( (lDcaPosToPrimVertex < 0.036 )||
412 (lDcaPosToPrimVertex < 0.036 )||
413 (lDcaV0Daughters > 0.5 ) ||
414 (lV0CosineOfPointingAngle < 0.99)
415 ) continue;
416 }
417
418 // Filling histograms
419 fHistDcaPosToPrimVertex->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
420 fHistDcaNegToPrimVertex->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
421 fHistDcaPosToPrimVertexZoom->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
422 fHistDcaNegToPrimVertexZoom->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
423 fHistRadiusV0->Fill(lV0Radius,lOnFlyStatus);
424 fHistDecayLengthV0->Fill(myAODv0->DecayLengthV0(tPrimaryVtxPosition),lOnFlyStatus);
425 fHistDcaV0Daughters->Fill(lDcaV0Daughters,lOnFlyStatus);
426 fHistChi2->Fill(lChi2V0,lOnFlyStatus);
427 fHistCosPointAngle->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
428 if (lV0CosineOfPointingAngle >= 0.9) fHistCosPointAngleZoom->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
429 if(!lOnFlyStatus){
430 nv0sOff++;
431 fHistPtVsYK0sOff->Fill(lPt,lRapK0s);
432 fHistPtVsYLambdaOff->Fill(lPt,lRapLambda);
433 fHistPtVsYAntiLambdaOff->Fill(lPt,lRapLambda);
434 fHistArmenterosPodolanskiOff->Fill(lAlphaV0,lPtArmV0);
435 }
436 else {
437 nv0sOn++;
438 fHistPtVsYK0sOn->Fill(lPt,lRapK0s);
439 fHistPtVsYLambdaOn->Fill(lPt,lRapLambda);
440 fHistPtVsYAntiLambdaOn->Fill(lPt,lRapLambda);
441 fHistArmenterosPodolanskiOn->Fill(lAlphaV0,lPtArmV0);
442 }
443 // K0s invariant mass histograms:
444 if (TMath::Abs(lRapK0s) < 1) {
445 if(!lOnFlyStatus){
446 fHistMassK0sOff->Fill(lInvMassK0s);
447 fHistMassVsRadiusK0sOff->Fill(lV0Radius,lInvMassK0s);
448 fHistPtVsMassK0sOff->Fill(lInvMassK0s,lPt);
449 }
450 else {
451 fHistMassK0sOn->Fill(lInvMassK0s);
452 fHistMassVsRadiusK0sOn->Fill(lV0Radius,lInvMassK0s);
453 fHistPtVsMassK0sOn->Fill(lInvMassK0s,lPt);
454 }
455 }
456 // Lambda and AntiLambda invariant mass histograms:
457 if (TMath::Abs(lRapLambda) < 1) {
458 if(!lOnFlyStatus){
459 fHistMassLambdaOff->Fill(lInvMassLambda);
460 fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
461 fHistMassVsRadiusLambdaOff->Fill(lV0Radius,lInvMassLambda);
462 fHistMassVsRadiusAntiLambdaOff->Fill(lV0Radius,lInvMassAntiLambda);
463 fHistPtVsMassLambdaOff->Fill(lInvMassLambda,lPt);
464 fHistPtVsMassAntiLambdaOff->Fill(lInvMassAntiLambda,lPt);
465 }
466 else {
467 fHistMassLambdaOn->Fill(lInvMassLambda);
468 fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
469 fHistMassVsRadiusLambdaOn->Fill(lV0Radius,lInvMassLambda);
470 fHistMassVsRadiusAntiLambdaOn->Fill(lV0Radius,lInvMassAntiLambda);
471 fHistPtVsMassLambdaOn->Fill(lInvMassLambda,lPt);
472 fHistPtVsMassAntiLambdaOn->Fill(lInvMassAntiLambda,lPt);
473 }
474 }
475 } // end V0 loop
476
477 if (primary) delete primary;
478 if (myPosAodTrack) delete myPosAodTrack;
479 if (myNegAodTrack) delete myNegAodTrack;
480 if (myAODVertex) delete myAODVertex;
481 if (myAODv0) delete myAODv0;
482 }
483 else if(fAnalysisType == "AOD") {
484
485 const AliAODVertex *primaryVtx = ((AliAODEvent*)lEvent)->GetPrimaryVertex();
486 tPrimaryVtxPosition[0] = primaryVtx->GetX();
487 tPrimaryVtxPosition[1] = primaryVtx->GetY();
488 tPrimaryVtxPosition[2] = primaryVtx->GetZ();
489
490 fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
491 fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
492 fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
493
494 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
495 {// This is the begining of the V0 loop
496 AliAODv0 *v0 = ((AliAODEvent*)lEvent)->GetV0(iV0);
497 if (!v0) continue;
498
499 // common part
500 lV0Radius = v0->RadiusV0();
501 lDcaPosToPrimVertex = v0->DcaPosToPrimVertex();
502 lDcaNegToPrimVertex = v0->DcaNegToPrimVertex();
503
504 lOnFlyStatus = v0->GetOnFlyStatus();
505 lChi2V0 = v0->Chi2V0();
506 lDcaV0Daughters = v0->DcaV0Daughters();
507 lDcaV0ToPrimVertex = v0->DcaV0ToPrimVertex();
508 lV0CosineOfPointingAngle = v0->CosPointingAngle(tPrimaryVtxPosition);
509
510 lInvMassK0s = v0->MassK0Short();
511 lInvMassLambda = v0->MassLambda();
512 lInvMassAntiLambda = v0->MassAntiLambda();
513
514 lPt = TMath::Sqrt(v0->Pt2V0());
515 lRapK0s = v0->RapK0Short();
516 lRapLambda = v0->RapLambda();
517 lAlphaV0 = v0->AlphaV0();
518 lPtArmV0 = v0->PtArmV0();
519
520
521 // Selections:
522 if (fOption.Contains("yes")) {
523 if ( (lDcaPosToPrimVertex < 0.036 )||
524 (lDcaPosToPrimVertex < 0.036 )||
525 (lDcaV0Daughters > 0.5 ) ||
526 (lV0CosineOfPointingAngle < 0.99)
527 ) continue;
528 }
529
530 // Filling histograms
531 fHistDcaPosToPrimVertex->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
532 fHistDcaNegToPrimVertex->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
533 fHistDcaPosToPrimVertexZoom->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
534 fHistDcaNegToPrimVertexZoom->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
535 fHistRadiusV0->Fill(lV0Radius,lOnFlyStatus);
536 fHistDecayLengthV0->Fill(v0->DecayLengthV0(tPrimaryVtxPosition),lOnFlyStatus);
537 fHistDcaV0Daughters->Fill(lDcaV0Daughters,lOnFlyStatus);
538 fHistChi2->Fill(lChi2V0,lOnFlyStatus);
539 fHistCosPointAngle->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
540 if (lV0CosineOfPointingAngle >= 0.9) fHistCosPointAngleZoom->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
541 if(!lOnFlyStatus){
542 nv0sOff++;
543 fHistPtVsYK0sOff->Fill(lPt,lRapK0s);
544 fHistPtVsYLambdaOff->Fill(lPt,lRapLambda);
545 fHistPtVsYAntiLambdaOff->Fill(lPt,lRapLambda);
546 fHistArmenterosPodolanskiOff->Fill(lAlphaV0,lPtArmV0);
547 }
548 else {
549 nv0sOn++;
550 fHistPtVsYK0sOn->Fill(lPt,lRapK0s);
551 fHistPtVsYLambdaOn->Fill(lPt,lRapLambda);
552 fHistPtVsYAntiLambdaOn->Fill(lPt,lRapLambda);
553 fHistArmenterosPodolanskiOn->Fill(lAlphaV0,lPtArmV0);
554 }
555 // K0s invariant mass histograms:
556 if (TMath::Abs(lRapK0s) < 1) {
557 if(!lOnFlyStatus){
558 fHistMassK0sOff->Fill(lInvMassK0s);
559 fHistMassVsRadiusK0sOff->Fill(lV0Radius,lInvMassK0s);
560 fHistPtVsMassK0sOff->Fill(lInvMassK0s,lPt);
561 }
562 else {
563 fHistMassK0sOn->Fill(lInvMassK0s);
564 fHistMassVsRadiusK0sOn->Fill(lV0Radius,lInvMassK0s);
565 fHistPtVsMassK0sOn->Fill(lInvMassK0s,lPt);
566 }
567 }
568 // Lambda and AntiLambda invariant mass histograms:
569 if (TMath::Abs(lRapLambda) < 1) {
570 if(!lOnFlyStatus){
571 fHistMassLambdaOff->Fill(lInvMassLambda);
572 fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
573 fHistMassVsRadiusLambdaOff->Fill(lV0Radius,lInvMassLambda);
574 fHistMassVsRadiusAntiLambdaOff->Fill(lV0Radius,lInvMassAntiLambda);
575 fHistPtVsMassLambdaOff->Fill(lInvMassLambda,lPt);
576 fHistPtVsMassAntiLambdaOff->Fill(lInvMassAntiLambda,lPt);
577 }
578 else {
579 fHistMassLambdaOn->Fill(lInvMassLambda);
580 fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
581 fHistMassVsRadiusLambdaOn->Fill(lV0Radius,lInvMassLambda);
582 fHistMassVsRadiusAntiLambdaOn->Fill(lV0Radius,lInvMassAntiLambda);
583 fHistPtVsMassLambdaOn->Fill(lInvMassLambda,lPt);
584 fHistPtVsMassAntiLambdaOn->Fill(lInvMassAntiLambda,lPt);
585 }
586 }
587 } // end V0 loop
588 }
589
590 fHistV0Multiplicity->Fill(nv0s);
591 fHistV0MultiplicityOff->Fill(nv0sOff);
592 fHistV0MultiplicityOn->Fill(nv0sOn);
593
594 // Post output data.
595 PostData(1, fListHist);
596}
597
598//________________________________________________________________________
599void AliAnalysisTaskStrange::Terminate(Option_t *)
600{
601 // Draw result to the screen
602 // Called once at the end of the query
603}