]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/SPECTRA/AliAnalysisTaskStrange.cxx
Adding more QA histograms: eta vs phi vs Nclusters
[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"
9ea746fc 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
43ClassImp(AliAnalysisTaskStrange)
44
45//________________________________________________________________________
46AliAnalysisTaskStrange::AliAnalysisTaskStrange()
47 : AliAnalysisTaskSE(), fAnalysisType("ESD"), fCollidingSystems(0), fOption("no"), 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//________________________________________________________________________
70AliAnalysisTaskStrange::AliAnalysisTaskStrange(const char *name, const char *optCuts)
71 : AliAnalysisTaskSE(name), fAnalysisType("ESD"), fCollidingSystems(0), fOption(optCuts), 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//________________________________________________________________________
98void 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};
c2c631b3 154 Int_t nBinRadius = 9;
9ea746fc 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:
c2c631b3 179 fHistMassVsRadiusK0sOff = new TH2F("h2MassVsRadiusK0sOff", "K^{0} Offline candidates;radius (cm);M(#pi^{+}#pi^{-}) (GeV/c^{2})",nBinRadius,radius, 200, 0.4, 0.6);
9ea746fc 180 fListHist->Add(fHistMassVsRadiusK0sOff);
c2c631b3 181 fHistMassVsRadiusLambdaOff = new TH2F("h2MassVsRadiusLambdaOff", "#Lambda Offline candidates;radius (cm);M(p#pi^{-}) (GeV/c^{2})",nBinRadius,radius, 140, 1.06, 1.2);
9ea746fc 182 fListHist->Add(fHistMassVsRadiusLambdaOff);
c2c631b3 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);
9ea746fc 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:
c2c631b3 219 fHistMassVsRadiusK0sOn = new TH2F("h2MassVsRadiusK0sOn", "K^{0} Onthefly candidates;radius (cm);M(#pi^{+}#pi^{-}) (GeV/c^{2})",nBinRadius,radius, 200, 0.4, 0.6);
9ea746fc 220 fListHist->Add(fHistMassVsRadiusK0sOn);
c2c631b3 221 fHistMassVsRadiusLambdaOn = new TH2F("h2MassVsRadiusLambdaOn", "#Lambda Onthefly candidates;radius (cm);M(p#pi^{-}) (GeV/c^{2})",nBinRadius,radius, 140, 1.06, 1.2);
9ea746fc 222 fListHist->Add(fHistMassVsRadiusLambdaOn);
c2c631b3 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);
9ea746fc 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//________________________________________________________________________
238void 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 Double_t tPrimaryVtxCov[3];
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 lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
267 Double_t lPt = 0, lRapK0s = 0, lRapLambda = 0;
268 Double_t lAlphaV0 = 0, lPtArmV0 = 0;
269
270 if(fAnalysisType == "ESD") {
271
272 const AliESDVertex *primaryVtx = ((AliESDEvent*)lEvent)->GetPrimaryVertex();
273 tPrimaryVtxPosition[0] = primaryVtx->GetXv();
274 tPrimaryVtxPosition[1] = primaryVtx->GetYv();
275 tPrimaryVtxPosition[2] = primaryVtx->GetZv();
276
277 fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
278 fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
279 fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
280
281 primaryVtx->GetCovMatrix(tPrimaryVtxCov);
282 AliAODVertex *primary = new AliAODVertex(tPrimaryVtxPosition, tPrimaryVtxCov, primaryVtx->GetChi2toNDF(), NULL, -1, AliAODVertex::kPrimary);
283
284 // V0 variables:
285 // to get info from ESD files and fill AliAODVertex:
286 Float_t tdcaPosToPrimVertexXYZ[2], tdcaNegToPrimVertexXYZ[2]; // ..[0] = Impact parameter in XY plane and ..[1] = Impact parameter in Z
287 Double_t tdcaDaughterToPrimVertex[2]; // ..[0] = Pos and ..[1] = Neg
288 Double_t tMomPos[3];
289 Double_t tMomNeg[3];
c2c631b3 290 Double_t tV0Position[3];
291 Double_t tV0Cov[6];
9ea746fc 292
293 // to fill AliAODtrack:
c2c631b3 294 Double_t tTrackP[3];
295 Double_t tTrackPosition[3];
296 Double_t tTrackCovTr[21];
297 Double_t tTrackPid[10];
9ea746fc 298
299 AliAODTrack *myPosAodTrack = new AliAODTrack();
300 AliAODTrack *myNegAodTrack = new AliAODTrack();
301 AliAODVertex *myAODVertex = new AliAODVertex();
302 AliAODv0 *myAODv0 = new AliAODv0();
303
304 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
305 {// This is the begining of the V0 loop
306 AliESDv0 *v0 = ((AliESDEvent*)lEvent)->GetV0(iV0);
307 if (!v0) continue;
308
309 // AliAODVertex
c2c631b3 310 v0->GetXYZ(tV0Position[0], tV0Position[1], tV0Position[2]);
311 v0->GetPosCov(tV0Cov);
312 myAODVertex->SetPosition(tV0Position[0],tV0Position[1],tV0Position[2]);
313 myAODVertex->SetCovMatrix(tV0Cov);
9ea746fc 314 myAODVertex->SetChi2perNDF(v0->GetChi2V0());
315 myAODVertex->SetID((Short_t)iV0);
316 myAODVertex->SetParent(primary);
317 myAODVertex->SetType(AliAODVertex::kV0);
318
319 // AliAODtrack (V0 Daughters)
320 UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
321 UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
322
c2c631b3 323 AliESDtrack *pTrack = ((AliESDEvent*)lEvent)->GetTrack(lKeyPos);
324 AliESDtrack *nTrack = ((AliESDEvent*)lEvent)->GetTrack(lKeyNeg);
325 if (!pTrack || !nTrack) {
9ea746fc 326 Printf("ERROR: Could not retreive one of the daughter track");
327 continue;
328 }
c2c631b3 329 UInt_t lLabelPos = (UInt_t)TMath::Abs(pTrack->GetLabel());
330 UInt_t lLabelNeg = (UInt_t)TMath::Abs(nTrack->GetLabel());
9ea746fc 331
c2c631b3 332 myPosAodTrack->SetID((Short_t)(pTrack->GetID()));
9ea746fc 333 myPosAodTrack->SetLabel(lLabelPos);
c2c631b3 334 pTrack->GetPxPyPz(tTrackP);
335 myPosAodTrack->SetP(tTrackP);
336 pTrack->GetXYZ(tTrackPosition);
337 myPosAodTrack->SetPosition(tTrackPosition,kFALSE);
338 pTrack->GetCovarianceXYZPxPyPz(tTrackCovTr);
339 myPosAodTrack->SetCovMatrix(tTrackCovTr);
340 pTrack->GetESDpid(tTrackPid);
341 myPosAodTrack->SetPID(tTrackPid);
342 myPosAodTrack->SetCharge((Short_t)(pTrack->Charge()));
343 myPosAodTrack->SetITSClusterMap(pTrack->GetITSClusterMap());
9ea746fc 344 myPosAodTrack->SetProdVertex(myAODVertex);
345 myPosAodTrack->SetUsedForVtxFit(kTRUE);
346 myPosAodTrack->SetUsedForPrimVtxFit(kFALSE);
347 myPosAodTrack->SetType(AliAODTrack::kSecondary);
348 myPosAodTrack->ConvertAliPIDtoAODPID();
349
c2c631b3 350 myNegAodTrack->SetID((Short_t)(nTrack->GetID()));
9ea746fc 351 myNegAodTrack->SetLabel(lLabelNeg);
c2c631b3 352 nTrack->GetPxPyPz(tTrackP);
353 myNegAodTrack->SetP(tTrackP);
354 nTrack->GetXYZ(tTrackPosition);
355 myNegAodTrack->SetPosition(tTrackPosition,kFALSE);
356 nTrack->GetCovarianceXYZPxPyPz(tTrackCovTr);
357 myNegAodTrack->SetCovMatrix(tTrackCovTr);
358 nTrack->GetESDpid(tTrackPid);
359 myNegAodTrack->SetPID(tTrackPid);
360 myNegAodTrack->SetCharge((Short_t)(nTrack->Charge()));
361 myNegAodTrack->SetITSClusterMap(nTrack->GetITSClusterMap());
9ea746fc 362 myNegAodTrack->SetProdVertex(myAODVertex);
363 myNegAodTrack->SetUsedForVtxFit(kTRUE);
364 myNegAodTrack->SetUsedForPrimVtxFit(kFALSE);
365 myNegAodTrack->SetType(AliAODTrack::kSecondary);
366 myNegAodTrack->ConvertAliPIDtoAODPID();
367
368 myAODVertex->AddDaughter(myPosAodTrack);
369 myAODVertex->AddDaughter(myNegAodTrack);
370
371 // filling myAODv0
372 lDcaV0Daughters = v0->GetDcaV0Daughters();
77f2c47e 373 lDcaV0ToPrimVertex = v0->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
9ea746fc 374
c2c631b3 375 if (pTrack) pTrack->GetImpactParameters(tdcaPosToPrimVertexXYZ[0],tdcaPosToPrimVertexXYZ[1]);
376 if (nTrack) nTrack->GetImpactParameters(tdcaNegToPrimVertexXYZ[0],tdcaNegToPrimVertexXYZ[1]);
9ea746fc 377 tdcaDaughterToPrimVertex[0] = TMath::Sqrt(tdcaPosToPrimVertexXYZ[0]*tdcaPosToPrimVertexXYZ[0]+tdcaPosToPrimVertexXYZ[1]*tdcaPosToPrimVertexXYZ[1]);
378 tdcaDaughterToPrimVertex[1] = TMath::Sqrt(tdcaNegToPrimVertexXYZ[0]*tdcaNegToPrimVertexXYZ[0]+tdcaNegToPrimVertexXYZ[1]*tdcaNegToPrimVertexXYZ[1]);
379
380 v0->GetPPxPyPz(tMomPos[0],tMomPos[1],tMomPos[2]);
381 v0->GetNPxPyPz(tMomNeg[0],tMomNeg[1],tMomNeg[2]);
382
383 myAODv0->Fill(myAODVertex, lDcaV0Daughters, lDcaV0ToPrimVertex, tMomPos, tMomNeg, tdcaDaughterToPrimVertex);
384 myAODv0->SetOnFlyStatus(v0->GetOnFlyStatus());
385
386 // common part
387 lV0Radius = myAODv0->RadiusV0();
388 lDcaPosToPrimVertex = myAODv0->DcaPosToPrimVertex();
389 lDcaNegToPrimVertex = myAODv0->DcaNegToPrimVertex();
390
391 lOnFlyStatus = myAODv0->GetOnFlyStatus();
392 lChi2V0 = myAODv0->Chi2V0();
393 lDcaV0Daughters = myAODv0->DcaV0Daughters();
394 lDcaV0ToPrimVertex = myAODv0->DcaV0ToPrimVertex();
395 lV0CosineOfPointingAngle = myAODv0->CosPointingAngle(tPrimaryVtxPosition);
396
397 lInvMassK0s = myAODv0->MassK0Short();
398 lInvMassLambda = myAODv0->MassLambda();
399 lInvMassAntiLambda = myAODv0->MassAntiLambda();
400
401 lPt = TMath::Sqrt(myAODv0->Pt2V0());
402 lRapK0s = myAODv0->RapK0Short();
403 lRapLambda = myAODv0->RapLambda();
404 lAlphaV0 = myAODv0->AlphaV0();
405 lPtArmV0 = myAODv0->PtArmV0();
406
407
408 // Selections:
409 if (fOption.Contains("yes")) {
410 if ( (lDcaPosToPrimVertex < 0.036 )||
411 (lDcaPosToPrimVertex < 0.036 )||
412 (lDcaV0Daughters > 0.5 ) ||
413 (lV0CosineOfPointingAngle < 0.99)
414 ) continue;
415 }
416
417 // Filling histograms
418 fHistDcaPosToPrimVertex->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
419 fHistDcaNegToPrimVertex->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
420 fHistDcaPosToPrimVertexZoom->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
421 fHistDcaNegToPrimVertexZoom->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
422 fHistRadiusV0->Fill(lV0Radius,lOnFlyStatus);
423 fHistDecayLengthV0->Fill(myAODv0->DecayLengthV0(tPrimaryVtxPosition),lOnFlyStatus);
424 fHistDcaV0Daughters->Fill(lDcaV0Daughters,lOnFlyStatus);
425 fHistChi2->Fill(lChi2V0,lOnFlyStatus);
426 fHistCosPointAngle->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
427 if (lV0CosineOfPointingAngle >= 0.9) fHistCosPointAngleZoom->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
428 if(!lOnFlyStatus){
429 nv0sOff++;
430 fHistPtVsYK0sOff->Fill(lPt,lRapK0s);
431 fHistPtVsYLambdaOff->Fill(lPt,lRapLambda);
432 fHistPtVsYAntiLambdaOff->Fill(lPt,lRapLambda);
433 fHistArmenterosPodolanskiOff->Fill(lAlphaV0,lPtArmV0);
434 }
435 else {
436 nv0sOn++;
437 fHistPtVsYK0sOn->Fill(lPt,lRapK0s);
438 fHistPtVsYLambdaOn->Fill(lPt,lRapLambda);
439 fHistPtVsYAntiLambdaOn->Fill(lPt,lRapLambda);
440 fHistArmenterosPodolanskiOn->Fill(lAlphaV0,lPtArmV0);
441 }
442 // K0s invariant mass histograms:
443 if (TMath::Abs(lRapK0s) < 1) {
444 if(!lOnFlyStatus){
445 fHistMassK0sOff->Fill(lInvMassK0s);
446 fHistMassVsRadiusK0sOff->Fill(lV0Radius,lInvMassK0s);
447 fHistPtVsMassK0sOff->Fill(lInvMassK0s,lPt);
448 }
449 else {
450 fHistMassK0sOn->Fill(lInvMassK0s);
451 fHistMassVsRadiusK0sOn->Fill(lV0Radius,lInvMassK0s);
452 fHistPtVsMassK0sOn->Fill(lInvMassK0s,lPt);
453 }
454 }
455 // Lambda and AntiLambda invariant mass histograms:
456 if (TMath::Abs(lRapLambda) < 1) {
457 if(!lOnFlyStatus){
458 fHistMassLambdaOff->Fill(lInvMassLambda);
459 fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
460 fHistMassVsRadiusLambdaOff->Fill(lV0Radius,lInvMassLambda);
461 fHistMassVsRadiusAntiLambdaOff->Fill(lV0Radius,lInvMassAntiLambda);
462 fHistPtVsMassLambdaOff->Fill(lInvMassLambda,lPt);
463 fHistPtVsMassAntiLambdaOff->Fill(lInvMassAntiLambda,lPt);
464 }
465 else {
466 fHistMassLambdaOn->Fill(lInvMassLambda);
467 fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
468 fHistMassVsRadiusLambdaOn->Fill(lV0Radius,lInvMassLambda);
469 fHistMassVsRadiusAntiLambdaOn->Fill(lV0Radius,lInvMassAntiLambda);
470 fHistPtVsMassLambdaOn->Fill(lInvMassLambda,lPt);
471 fHistPtVsMassAntiLambdaOn->Fill(lInvMassAntiLambda,lPt);
472 }
473 }
474 } // end V0 loop
475
476 if (primary) delete primary;
477 if (myPosAodTrack) delete myPosAodTrack;
478 if (myNegAodTrack) delete myNegAodTrack;
479 if (myAODVertex) delete myAODVertex;
480 if (myAODv0) delete myAODv0;
481 }
482 else if(fAnalysisType == "AOD") {
483
484 const AliAODVertex *primaryVtx = ((AliAODEvent*)lEvent)->GetPrimaryVertex();
485 tPrimaryVtxPosition[0] = primaryVtx->GetX();
486 tPrimaryVtxPosition[1] = primaryVtx->GetY();
487 tPrimaryVtxPosition[2] = primaryVtx->GetZ();
488
489 fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
490 fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
491 fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
492
493 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
494 {// This is the begining of the V0 loop
495 AliAODv0 *v0 = ((AliAODEvent*)lEvent)->GetV0(iV0);
496 if (!v0) continue;
497
498 // common part
499 lV0Radius = v0->RadiusV0();
500 lDcaPosToPrimVertex = v0->DcaPosToPrimVertex();
501 lDcaNegToPrimVertex = v0->DcaNegToPrimVertex();
502
503 lOnFlyStatus = v0->GetOnFlyStatus();
504 lChi2V0 = v0->Chi2V0();
505 lDcaV0Daughters = v0->DcaV0Daughters();
506 lDcaV0ToPrimVertex = v0->DcaV0ToPrimVertex();
507 lV0CosineOfPointingAngle = v0->CosPointingAngle(tPrimaryVtxPosition);
508
509 lInvMassK0s = v0->MassK0Short();
510 lInvMassLambda = v0->MassLambda();
511 lInvMassAntiLambda = v0->MassAntiLambda();
512
513 lPt = TMath::Sqrt(v0->Pt2V0());
514 lRapK0s = v0->RapK0Short();
515 lRapLambda = v0->RapLambda();
516 lAlphaV0 = v0->AlphaV0();
517 lPtArmV0 = v0->PtArmV0();
518
519
520 // Selections:
521 if (fOption.Contains("yes")) {
522 if ( (lDcaPosToPrimVertex < 0.036 )||
523 (lDcaPosToPrimVertex < 0.036 )||
524 (lDcaV0Daughters > 0.5 ) ||
525 (lV0CosineOfPointingAngle < 0.99)
526 ) continue;
527 }
528
529 // Filling histograms
530 fHistDcaPosToPrimVertex->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
531 fHistDcaNegToPrimVertex->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
532 fHistDcaPosToPrimVertexZoom->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
533 fHistDcaNegToPrimVertexZoom->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
534 fHistRadiusV0->Fill(lV0Radius,lOnFlyStatus);
535 fHistDecayLengthV0->Fill(v0->DecayLengthV0(tPrimaryVtxPosition),lOnFlyStatus);
536 fHistDcaV0Daughters->Fill(lDcaV0Daughters,lOnFlyStatus);
537 fHistChi2->Fill(lChi2V0,lOnFlyStatus);
538 fHistCosPointAngle->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
539 if (lV0CosineOfPointingAngle >= 0.9) fHistCosPointAngleZoom->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
540 if(!lOnFlyStatus){
541 nv0sOff++;
542 fHistPtVsYK0sOff->Fill(lPt,lRapK0s);
543 fHistPtVsYLambdaOff->Fill(lPt,lRapLambda);
544 fHistPtVsYAntiLambdaOff->Fill(lPt,lRapLambda);
545 fHistArmenterosPodolanskiOff->Fill(lAlphaV0,lPtArmV0);
546 }
547 else {
548 nv0sOn++;
549 fHistPtVsYK0sOn->Fill(lPt,lRapK0s);
550 fHistPtVsYLambdaOn->Fill(lPt,lRapLambda);
551 fHistPtVsYAntiLambdaOn->Fill(lPt,lRapLambda);
552 fHistArmenterosPodolanskiOn->Fill(lAlphaV0,lPtArmV0);
553 }
554 // K0s invariant mass histograms:
555 if (TMath::Abs(lRapK0s) < 1) {
556 if(!lOnFlyStatus){
557 fHistMassK0sOff->Fill(lInvMassK0s);
558 fHistMassVsRadiusK0sOff->Fill(lV0Radius,lInvMassK0s);
559 fHistPtVsMassK0sOff->Fill(lInvMassK0s,lPt);
560 }
561 else {
562 fHistMassK0sOn->Fill(lInvMassK0s);
563 fHistMassVsRadiusK0sOn->Fill(lV0Radius,lInvMassK0s);
564 fHistPtVsMassK0sOn->Fill(lInvMassK0s,lPt);
565 }
566 }
567 // Lambda and AntiLambda invariant mass histograms:
568 if (TMath::Abs(lRapLambda) < 1) {
569 if(!lOnFlyStatus){
570 fHistMassLambdaOff->Fill(lInvMassLambda);
571 fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
572 fHistMassVsRadiusLambdaOff->Fill(lV0Radius,lInvMassLambda);
573 fHistMassVsRadiusAntiLambdaOff->Fill(lV0Radius,lInvMassAntiLambda);
574 fHistPtVsMassLambdaOff->Fill(lInvMassLambda,lPt);
575 fHistPtVsMassAntiLambdaOff->Fill(lInvMassAntiLambda,lPt);
576 }
577 else {
578 fHistMassLambdaOn->Fill(lInvMassLambda);
579 fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
580 fHistMassVsRadiusLambdaOn->Fill(lV0Radius,lInvMassLambda);
581 fHistMassVsRadiusAntiLambdaOn->Fill(lV0Radius,lInvMassAntiLambda);
582 fHistPtVsMassLambdaOn->Fill(lInvMassLambda,lPt);
583 fHistPtVsMassAntiLambdaOn->Fill(lInvMassAntiLambda,lPt);
584 }
585 }
586 } // end V0 loop
587 }
588
589 fHistV0Multiplicity->Fill(nv0s);
590 fHistV0MultiplicityOff->Fill(nv0sOff);
591 fHistV0MultiplicityOn->Fill(nv0sOn);
592
593 // Post output data.
594 PostData(1, fListHist);
595}
596
597//________________________________________________________________________
598void AliAnalysisTaskStrange::Terminate(Option_t *)
599{
600 // Draw result to the screen
601 // Called once at the end of the query
602}