]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/SPECTRA/AliAnalysisTaskCheckV0.cxx
Important update for the task CheckV0 (B.Hippolyte)
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliAnalysisTaskCheckV0.cxx
CommitLineData
96cad837 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//-----------------------------------------------------------------
c8eada4f 15// AliAnalysisTaskCheckV0 class
16// This task is for QAing the V0s from ESD/AOD
96cad837 17// Origin: B.H. Nov2007, hippolyt@in2p3.fr
18//-----------------------------------------------------------------
96cad837 19#include "TList.h"
20#include "TH1F.h"
50b713dc 21#include "TH2F.h"
22
23#include "TStyle.h"
96cad837 24#include "TCanvas.h"
c8eada4f 25#include "TLegend.h"
96cad837 26
c8eada4f 27#include "AliAnalysisTaskSE.h"
7285f641 28#include "AliAnalysisManager.h"
29#include "AliInputEventHandler.h"
96cad837 30
31#include "AliESDEvent.h"
b638d7f5 32#include "AliESDVertex.h"
c8eada4f 33#include "AliAODEvent.h"
96cad837 34
35#include "AliESDv0.h"
36
c8eada4f 37#include "AliAnalysisTaskCheckV0.h"
96cad837 38
c8eada4f 39ClassImp(AliAnalysisTaskCheckV0)
96cad837 40
b638d7f5 41//________________________________________________________________________
42AliAnalysisTaskCheckV0::AliAnalysisTaskCheckV0()
7285f641 43 : AliAnalysisTaskSE(), fAnalysisType("ESD"), fCollidingSystems(0), fUsePhysicsSelection(0),
44 fMaxPrimaryVtxPosZ(100.), fMinV0Pt(0.), fMaxV0Pt(100.), fMaxV0Rapidity(1.), fMinDaughterTpcClusters(80),
45 fListHist(),
b638d7f5 46 fHistPrimaryVertexPosX(0), fHistPrimaryVertexPosY(0), fHistPrimaryVertexPosZ(0),
7285f641 47 fHistKeptPrimaryVertexPosX(0), fHistKeptPrimaryVertexPosY(0), fHistKeptPrimaryVertexPosZ(0),
b638d7f5 48 fHistTrackMultiplicity(0), fHistV0Multiplicity(0), fHistV0OnFlyStatus(0),
50b713dc 49 fHistV0MultiplicityOff(0),
50 fHistV0Chi2Off(0),
b638d7f5 51 fHistDcaV0DaughtersOff(0), fHistV0CosineOfPointingAngleOff(0),
52 fHistV0RadiusOff(0),fHistDcaV0ToPrimVertexOff(0),
53 fHistDcaPosToPrimVertexOff(0),fHistDcaNegToPrimVertexOff(0),
9ea746fc 54 fHistMassK0sOff(0),fHistMassLambdaOff(0),fHistMassAntiLambdaOff(0),
7285f641 55 fHistMassK0sOffVsPt(0),fHistMassLambdaOffVsPt(0),fHistMassAntiLambdaOffVsPt(0),
50b713dc 56 fHistArmenterosPodolanskiOff(0),
57 fHistV0MultiplicityOn(0),
58 fHistV0Chi2On(0),
b638d7f5 59 fHistDcaV0DaughtersOn(0), fHistV0CosineOfPointingAngleOn(0),
60 fHistV0RadiusOn(0),fHistDcaV0ToPrimVertexOn(0),
61 fHistDcaPosToPrimVertexOn(0),fHistDcaNegToPrimVertexOn(0),
50b713dc 62 fHistMassK0sOn(0),fHistMassLambdaOn(0),fHistMassAntiLambdaOn(0),
7285f641 63 fHistMassK0sOnVsPt(0),fHistMassLambdaOnVsPt(0),fHistMassAntiLambdaOnVsPt(0),
50b713dc 64 fHistArmenterosPodolanskiOn(0)
b638d7f5 65{
66 // Dummy constructor
67}
96cad837 68//________________________________________________________________________
c8eada4f 69AliAnalysisTaskCheckV0::AliAnalysisTaskCheckV0(const char *name)
7285f641 70 : AliAnalysisTaskSE(name), fAnalysisType("ESD"), fCollidingSystems(0), fUsePhysicsSelection(0),
71 fMaxPrimaryVtxPosZ(100.), fMinV0Pt(0.), fMaxV0Pt(100.), fMaxV0Rapidity(1.), fMinDaughterTpcClusters(80),
72 fListHist(),
b638d7f5 73 fHistPrimaryVertexPosX(0), fHistPrimaryVertexPosY(0), fHistPrimaryVertexPosZ(0),
7285f641 74 fHistKeptPrimaryVertexPosX(0), fHistKeptPrimaryVertexPosY(0), fHistKeptPrimaryVertexPosZ(0),
96cad837 75 fHistTrackMultiplicity(0), fHistV0Multiplicity(0), fHistV0OnFlyStatus(0),
50b713dc 76 fHistV0MultiplicityOff(0),
77 fHistV0Chi2Off(0),
96cad837 78 fHistDcaV0DaughtersOff(0), fHistV0CosineOfPointingAngleOff(0),
79 fHistV0RadiusOff(0),fHistDcaV0ToPrimVertexOff(0),
80 fHistDcaPosToPrimVertexOff(0),fHistDcaNegToPrimVertexOff(0),
9ea746fc 81 fHistMassK0sOff(0),fHistMassLambdaOff(0),fHistMassAntiLambdaOff(0),
7285f641 82 fHistMassK0sOffVsPt(0),fHistMassLambdaOffVsPt(0),fHistMassAntiLambdaOffVsPt(0),
50b713dc 83 fHistArmenterosPodolanskiOff(0),
84 fHistV0MultiplicityOn(0),
85 fHistV0Chi2On(0),
96cad837 86 fHistDcaV0DaughtersOn(0), fHistV0CosineOfPointingAngleOn(0),
87 fHistV0RadiusOn(0),fHistDcaV0ToPrimVertexOn(0),
88 fHistDcaPosToPrimVertexOn(0),fHistDcaNegToPrimVertexOn(0),
50b713dc 89 fHistMassK0sOn(0),fHistMassLambdaOn(0),fHistMassAntiLambdaOn(0),
7285f641 90 fHistMassK0sOnVsPt(0),fHistMassLambdaOnVsPt(0),fHistMassAntiLambdaOnVsPt(0),
50b713dc 91 fHistArmenterosPodolanskiOn(0)
96cad837 92{
93 // Constructor
976dc218 94 // Define output slots only here
95 // Output slot #1 writes into a TList container
96 DefineOutput(1, TList::Class());
96cad837 97}
98
99//________________________________________________________________________
976dc218 100void AliAnalysisTaskCheckV0::UserCreateOutputObjects()
96cad837 101{
102 // Create histograms
103 // Called once
104
105 // Distinguish Track and V0 Multiplicity !
106
107 fListHist = new TList();
b638d7f5 108
109 if (!fHistPrimaryVertexPosX) {
110 fHistPrimaryVertexPosX = new TH1F("fHistPrimaryVertexPosX", "Primary vertex position in x;Position in x (cm);Events;",100,-1,1);
111 fListHist->Add(fHistPrimaryVertexPosX);
112 }
113 if (!fHistPrimaryVertexPosY) {
114 fHistPrimaryVertexPosY = new TH1F("fHistPrimaryVertexPosY", "Primary vertex position in y;Position in y (cm);Events;",100,-1,1);
115 fListHist->Add(fHistPrimaryVertexPosY);
116 }
117 if (!fHistPrimaryVertexPosZ) {
118 fHistPrimaryVertexPosZ = new TH1F("fHistPrimaryVertexPosZ", "Primary vertex position in z;Position in z (cm);Events;",100,-1,1);
119 fListHist->Add(fHistPrimaryVertexPosZ);
120 }
7285f641 121 if (!fHistKeptPrimaryVertexPosX) {
122 fHistKeptPrimaryVertexPosX = new TH1F("fHistKeptPrimaryVertexPosX", "Kept Primary vertex position in x;Position in x (cm);Events;",100,-1,1);
123 fListHist->Add(fHistKeptPrimaryVertexPosX);
124 }
125 if (!fHistKeptPrimaryVertexPosY) {
126 fHistKeptPrimaryVertexPosY = new TH1F("fHistKeptPrimaryVertexPosY", "Kept Primary vertex position in y;Position in y (cm);Events;",100,-1,1);
127 fListHist->Add(fHistKeptPrimaryVertexPosY);
128 }
129 if (!fHistKeptPrimaryVertexPosZ) {
130 fHistKeptPrimaryVertexPosZ = new TH1F("fHistKeptPrimaryVertexPosZ", "Kept Primary vertex position in z;Position in z (cm);Events;",100,-1,1);
131 fListHist->Add(fHistKeptPrimaryVertexPosZ);
132 }
96cad837 133
134 if (!fHistTrackMultiplicity) {
c8eada4f 135 if (fCollidingSystems)
136 fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 200, 0, 40000);
137 else
138 fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 250, 0, 250);
96cad837 139 fListHist->Add(fHistTrackMultiplicity);
140 }
141 if (!fHistV0Multiplicity) {
c8eada4f 142 if (fCollidingSystems)
143 fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
144 else
145 fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
96cad837 146 fListHist->Add(fHistV0Multiplicity);
147 }
148 if (!fHistV0OnFlyStatus) {
149 fHistV0OnFlyStatus = new TH1F("fHistV0OnFlyStatus", "V0 On fly status;status;Number of V0s", 3, 0, 3);
150 fListHist->Add(fHistV0OnFlyStatus);
151 }
152
153 // V0 offline distributions
154 if (!fHistV0MultiplicityOff) {
c8eada4f 155 if (fCollidingSystems)
156 fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
157 else
158 fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
96cad837 159 fListHist->Add(fHistV0MultiplicityOff);
160 }
161 if (!fHistV0Chi2Off) {
162 fHistV0Chi2Off = new TH1F("fHistV0Chi2Off", "V0 chi2;chi2;Number of V0s", 33, 0, 33);
163 fListHist->Add(fHistV0Chi2Off);
164 }
165 if (!fHistDcaV0DaughtersOff) {
166 fHistDcaV0DaughtersOff = new TH1F("fHistDcaV0DaughtersOff", "DCA between V0 daughters;DCA (cm);Number of V0s", 300, 0, 3);
167 fListHist->Add(fHistDcaV0DaughtersOff);
168 }
169 if (!fHistV0CosineOfPointingAngleOff) {
170 fHistV0CosineOfPointingAngleOff = new TH1F("fHistV0CosineOfPointingAngleOff", "V0 Cosine of Pointing Angle;Number of V0s", 200, 0, 1);
171 fListHist->Add(fHistV0CosineOfPointingAngleOff);
172 }
173 if (!fHistV0RadiusOff) {
174 fHistV0RadiusOff = new TH1F("fHistV0RadiusOff", "V0 decay radius;Radius (cm);Number of V0s", 33, 0, 33);
175 fListHist->Add(fHistV0RadiusOff);
176 }
177 if (!fHistDcaV0ToPrimVertexOff) {
178 fHistDcaV0ToPrimVertexOff = new TH1F("fHistDcaV0ToPrimVertexOff", "DCA of V0 to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
179 fListHist->Add(fHistDcaV0ToPrimVertexOff);
180 }
181 if (!fHistDcaPosToPrimVertexOff) {
182 fHistDcaPosToPrimVertexOff = new TH1F("fHistDcaPosToPrimVertexOff", "DCA of V0 neg daughter to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
183 fListHist->Add(fHistDcaPosToPrimVertexOff);
184 }
185 if (!fHistDcaNegToPrimVertexOff) {
186 fHistDcaNegToPrimVertexOff = new TH1F("fHistDcaNegToPrimVertexOff", "DCA of V0 pos daughter to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
187 fListHist->Add(fHistDcaNegToPrimVertexOff);
188 }
189
9ea746fc 190 if (!fHistMassK0sOff) {
191 fHistMassK0sOff = new TH1F("fHistMassK0sOff","K^{0} candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});Counts",100,0.4,0.6);
192 fListHist->Add(fHistMassK0sOff);
96cad837 193 }
194 if (!fHistMassLambdaOff) {
c8eada4f 195 fHistMassLambdaOff = new TH1F("fHistMassLambdaOff","#Lambda^{0} candidates;M(p#pi^{-}) (GeV/c^{2});Counts",75,1.05,1.2);
96cad837 196 fListHist->Add(fHistMassLambdaOff);
197 }
198 if (!fHistMassAntiLambdaOff) {
c8eada4f 199 fHistMassAntiLambdaOff = new TH1F("fHistMassAntiLambdaOff","#bar{#Lambda}^{0} candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts",75,1.05,1.2);
96cad837 200 fListHist->Add(fHistMassAntiLambdaOff);
201 }
7285f641 202 if (!fHistMassK0sOffVsPt) {
203 fHistMassK0sOffVsPt = new TH2F("fHistMassK0sOffVsPt","K^{0} candidates;p_{T} (GeV/c);M(#pi^{+}#pi^{-}) (GeV/c^{2})",200,0,10,100,0.4,0.6);
204 fListHist->Add(fHistMassK0sOffVsPt);
205 }
206 if (!fHistMassLambdaOffVsPt) {
207 fHistMassLambdaOffVsPt = new TH2F("fHistMassLambdaOffVsPt","#Lambda^{0} candidates;p_{T} (GeV/c);M(p#pi^{-}) (GeV/c^{2})",200,0,10,75,1.05,1.2);
208 fListHist->Add(fHistMassLambdaOffVsPt);
209 }
210 if (!fHistMassAntiLambdaOffVsPt) {
211 fHistMassAntiLambdaOffVsPt = new TH2F("fHistMassAntiLambdaOffVsPt","#bar{#Lambda}^{0} candidates;p_{T} (GeV/c);M(#bar{p}#pi^{+}) (GeV/c^{2})",200,0,10,75,1.05,1.2);
212 fListHist->Add(fHistMassAntiLambdaOffVsPt);
213 }
50b713dc 214 if (!fHistArmenterosPodolanskiOff) {
abcd19bd 215 fHistArmenterosPodolanskiOff = new TH2F("fHistArmenterosPodolanskiOff","Armenteros-Podolanski Offline phase space;#alpha;p_{t} arm",200,-1.0,1.0,150,0,0.3);
50b713dc 216 fListHist->Add(fHistArmenterosPodolanskiOff);
217 }
96cad837 218
219 // V0 on-the-fly distributions
220 if (!fHistV0MultiplicityOn) {
c8eada4f 221 if (fCollidingSystems)
222 fHistV0MultiplicityOn = new TH1F("fHistV0MultiplicityOn", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
223 else
224 fHistV0MultiplicityOn = new TH1F("fHistV0MultiplicityOn", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
96cad837 225 fListHist->Add(fHistV0MultiplicityOn);
226 }
227 if (!fHistV0Chi2On) {
228 fHistV0Chi2On = new TH1F("fHistV0Chi2On", "V0 chi2;chi2;Number of V0s", 33, 0, 33);
229 fListHist->Add(fHistV0Chi2On);
230 }
231 if (!fHistDcaV0DaughtersOn) {
232 fHistDcaV0DaughtersOn = new TH1F("fHistDcaV0DaughtersOn", "DCA between V0 daughters;DCA (cm);Number of V0s", 300, 0, 3);
233 fListHist->Add(fHistDcaV0DaughtersOn);
234 }
235 if (!fHistV0CosineOfPointingAngleOn) {
236 fHistV0CosineOfPointingAngleOn = new TH1F("fHistV0CosineOfPointingAngleOn", "V0 Cosine of Pointing Angle;Number of V0s", 200, 0, 1);
237 fListHist->Add(fHistV0CosineOfPointingAngleOn);
238 }
239 if (!fHistV0RadiusOn) {
240 fHistV0RadiusOn = new TH1F("fHistV0RadiusOn", "V0 decay radius;Radius (cm);Number of V0s", 33, 0, 33);
241 fListHist->Add(fHistV0RadiusOn);
242 }
243 if (!fHistDcaV0ToPrimVertexOn) {
244 fHistDcaV0ToPrimVertexOn = new TH1F("fHistDcaV0ToPrimVertexOn", "DCA of V0 to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
245 fListHist->Add(fHistDcaV0ToPrimVertexOn);
246 }
247 if (!fHistDcaPosToPrimVertexOn) {
248 fHistDcaPosToPrimVertexOn = new TH1F("fHistDcaPosToPrimVertexOn", "DCA of V0 neg daughter to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
249 fListHist->Add(fHistDcaPosToPrimVertexOn);
250 }
251 if (!fHistDcaNegToPrimVertexOn) {
252 fHistDcaNegToPrimVertexOn = new TH1F("fHistDcaNegToPrimVertexOn", "DCA of V0 pos daughter to Prim. Vertex;DCA (cm);Number of V0s", 300, 0, 3);
253 fListHist->Add(fHistDcaNegToPrimVertexOn);
254 }
255
9ea746fc 256 if (!fHistMassK0sOn) {
257 fHistMassK0sOn = new TH1F("fHistMassK0sOn","K^{0} candidates;M(#pi^{+}#pi^{-}) (GeV/c^{2});Counts",100,0.4,0.6);
258 fListHist->Add(fHistMassK0sOn);
96cad837 259 }
260 if (!fHistMassLambdaOn) {
c8eada4f 261 fHistMassLambdaOn = new TH1F("fHistMassLambdaOn","#Lambda^{0} candidates;M(p#pi^{-}) (GeV/c^{2});Counts",75,1.05,1.2);
96cad837 262 fListHist->Add(fHistMassLambdaOn);
263 }
264 if (!fHistMassAntiLambdaOn) {
c8eada4f 265 fHistMassAntiLambdaOn = new TH1F("fHistMassAntiLambdaOn","#bar{#Lambda}^{0} candidates;M(#bar{p}#pi^{+}) (GeV/c^{2});Counts",75,1.05,1.2);
96cad837 266 fListHist->Add(fHistMassAntiLambdaOn);
267 }
7285f641 268 if (!fHistMassK0sOnVsPt) {
269 fHistMassK0sOnVsPt = new TH2F("fHistMassK0sOnVsPt","K^{0} candidates;p_{T} (GeV/c);M(#pi^{+}#pi^{-}) (GeV/c^{2})",200,0,10,100,0.4,0.6);
270 fListHist->Add(fHistMassK0sOnVsPt);
271 }
272 if (!fHistMassLambdaOnVsPt) {
273 fHistMassLambdaOnVsPt = new TH2F("fHistMassLambdaOnVsPt","#Lambda^{0} candidates;p_{T} (GeV/c);M(p#pi^{-}) (GeV/c^{2})",200,0,10,75,1.05,1.2);
274 fListHist->Add(fHistMassLambdaOnVsPt);
275 }
276 if (!fHistMassAntiLambdaOnVsPt) {
277 fHistMassAntiLambdaOnVsPt = new TH2F("fHistMassAntiLambdaOnVsPt","#bar{#Lambda}^{0} candidates;p_{T} (GeV/c);M(#bar{p}#pi^{+}) (GeV/c^{2})",200,0,10,75,1.05,1.2);
278 fListHist->Add(fHistMassAntiLambdaOnVsPt);
279 }
50b713dc 280 if (!fHistArmenterosPodolanskiOn) {
abcd19bd 281 fHistArmenterosPodolanskiOn = new TH2F("fHistArmenterosPodolanskiOn","Armenteros-Podolanski Onthefly phase space;#alpha;p_{t} arm",200,-1.0,1.0,150,0,0.3);
50b713dc 282 fListHist->Add(fHistArmenterosPodolanskiOn);
283 }
96cad837 284
285}
286
287//________________________________________________________________________
976dc218 288void AliAnalysisTaskCheckV0::UserExec(Option_t *)
96cad837 289{
290 // Main loop
291 // Called for each event
976dc218 292 AliVEvent* lEvent = InputEvent();
293 if (!lEvent) {
294 Printf("ERROR: Event not available");
295 return;
296 }
7285f641 297
298 Bool_t isSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
299 if ((fUsePhysicsSelection)&&(!isSelected)) return;
300
976dc218 301 fHistTrackMultiplicity->Fill(lEvent->GetNumberOfTracks());
302
303 Double_t tPrimaryVtxPosition[3];
304 Int_t nv0s = 0;
305 nv0s = lEvent->GetNumberOfV0s();
83db4a5d 306 // Printf("CheckV0 analysis task: There are %d v0s in this event",nv0s);
976dc218 307
308 Int_t lOnFlyStatus = 0, nv0sOn = 0, nv0sOff = 0;
309 Double_t lChi2V0 = 0;
310 Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0;
311 Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0;
312 Double_t lV0CosineOfPointingAngle = 0;
7285f641 313 Double_t lV0Radius = 0, lPt = 0;
314 Double_t lRapK0Short = 0, lRapLambda = 0;
9ea746fc 315 Double_t lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
50b713dc 316 Double_t lAlphaV0 = 0, lPtArmV0 = 0;
96cad837 317
420b11e4 318 const AliVVertex *primaryVtx = lEvent->GetPrimaryVertex();
f80d1a51 319 tPrimaryVtxPosition[0] = primaryVtx->GetX();
320 tPrimaryVtxPosition[1] = primaryVtx->GetY();
321 tPrimaryVtxPosition[2] = primaryVtx->GetZ();
322 fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
323 fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
324 fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
b638d7f5 325
7285f641 326 if (TMath::Abs(tPrimaryVtxPosition[2])<fMaxPrimaryVtxPosZ){// event selections
327
328 fHistKeptPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
329 fHistKeptPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
330 fHistKeptPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
331
332 if(fAnalysisType == "ESD") {
333
334 Double_t lMagneticField = ((AliESDEvent*)lEvent)->GetMagneticField();
335
336 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
337 {// This is the begining of the V0 loop
338 AliESDv0 *v0 = ((AliESDEvent*)lEvent)->GetV0(iV0);
339 if (!v0) continue;
340
341 Double_t tDecayVertexV0[3]; v0->GetXYZ(tDecayVertexV0[0],tDecayVertexV0[1],tDecayVertexV0[2]);
342
343 lV0Radius = TMath::Sqrt(tDecayVertexV0[0]*tDecayVertexV0[0]+tDecayVertexV0[1]*tDecayVertexV0[1]);
344 lPt = v0->Pt();
345 lRapK0Short = v0->RapK0Short();
346 lRapLambda = v0->RapLambda();
347 if ((lPt<fMinV0Pt)||(fMaxV0Pt<lPt)) continue;
348
349 UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
350 UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
351
352 Double_t lMomPos[3]; v0->GetPPxPyPz(lMomPos[0],lMomPos[1],lMomPos[2]);
353 Double_t lMomNeg[3]; v0->GetNPxPyPz(lMomNeg[0],lMomNeg[1],lMomNeg[2]);
354
355 AliESDtrack *pTrack=((AliESDEvent*)lEvent)->GetTrack(lKeyPos);
356 AliESDtrack *nTrack=((AliESDEvent*)lEvent)->GetTrack(lKeyNeg);
357 if (!pTrack || !nTrack) {
358 Printf("ERROR: Could not retreive one of the daughter track");
359 continue;
360 }
361
362 // Filter like-sign V0 (next: add counter and distribution)
363 if ( pTrack->GetSign() == nTrack->GetSign()){
364 continue;
365 }
366
367 // WARNING: the following selections cannot be done for AOD yet...
368
369 // Tracks quality cuts
370 if ( ( (pTrack->GetTPCNcls()) < fMinDaughterTpcClusters ) || ( (nTrack->GetTPCNcls()) < fMinDaughterTpcClusters ) ) continue;
371
372 // TPC refit condition (done during reconstruction for Offline but not for On-the-fly)
373 if( !(pTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
374 if( !(nTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
375
376
377 lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(tPrimaryVtxPosition[0],
378 tPrimaryVtxPosition[1],
379 lMagneticField) );
380
381 lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(tPrimaryVtxPosition[0],
382 tPrimaryVtxPosition[1],
383 lMagneticField) );
384
385 lOnFlyStatus = v0->GetOnFlyStatus();
386 lChi2V0 = v0->GetChi2V0();
387 lDcaV0Daughters = v0->GetDcaV0Daughters();
388 lDcaV0ToPrimVertex = v0->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
389 lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
390
391 // Getting invariant mass infos directly from ESD
392 v0->ChangeMassHypothesis(310);
393 lInvMassK0s = v0->GetEffMass();
394 v0->ChangeMassHypothesis(3122);
395 lInvMassLambda = v0->GetEffMass();
396 v0->ChangeMassHypothesis(-3122);
397 lInvMassAntiLambda = v0->GetEffMass();
398 lAlphaV0 = v0->AlphaV0();
399 lPtArmV0 = v0->PtArmV0();
400
401 fHistV0OnFlyStatus->Fill(lOnFlyStatus);
402 if(!lOnFlyStatus){
403 nv0sOff++;
404 fHistV0Chi2Off->Fill(lChi2V0);
405 fHistDcaV0ToPrimVertexOff->Fill(lDcaV0ToPrimVertex);
406 fHistDcaV0DaughtersOff->Fill(lDcaV0Daughters);
407 fHistV0CosineOfPointingAngleOff->Fill(lV0CosineOfPointingAngle);
408
409 fHistV0RadiusOff->Fill(lV0Radius);
410 fHistDcaPosToPrimVertexOff->Fill(lDcaPosToPrimVertex);
411 fHistDcaNegToPrimVertexOff->Fill(lDcaNegToPrimVertex);
412
413 // Filling invariant mass histos for all candidates
414 fHistMassK0sOff->Fill(lInvMassK0s);
415 fHistMassLambdaOff->Fill(lInvMassLambda);
416 fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
417 if (TMath::Abs(lRapK0Short)<fMaxV0Rapidity)
418 fHistMassK0sOffVsPt->Fill(lPt,lInvMassK0s);
419 if (TMath::Abs(lRapLambda)<fMaxV0Rapidity){
420 fHistMassLambdaOffVsPt->Fill(lPt,lInvMassLambda);
421 fHistMassAntiLambdaOffVsPt->Fill(lPt,lInvMassAntiLambda);
422 }
423 fHistArmenterosPodolanskiOff->Fill(lAlphaV0,lPtArmV0);
424 }
425 else {
426 nv0sOn++;
427 fHistV0Chi2On->Fill(lChi2V0);
428 fHistDcaV0ToPrimVertexOn->Fill(lDcaV0ToPrimVertex);
429 fHistDcaV0DaughtersOn->Fill(lDcaV0Daughters);
430 fHistV0CosineOfPointingAngleOn->Fill(lV0CosineOfPointingAngle);
431
432 fHistV0RadiusOn->Fill(lV0Radius);
433 fHistDcaPosToPrimVertexOn->Fill(lDcaPosToPrimVertex);
434 fHistDcaNegToPrimVertexOn->Fill(lDcaNegToPrimVertex);
435
436 // Filling invariant mass histos for all candidates
437 fHistMassK0sOn->Fill(lInvMassK0s);
438 fHistMassLambdaOn->Fill(lInvMassLambda);
439 fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
440 if (TMath::Abs(lRapK0Short)<fMaxV0Rapidity)
441 fHistMassK0sOnVsPt->Fill(lPt,lInvMassK0s);
442 if (TMath::Abs(lRapLambda)<fMaxV0Rapidity){
443 fHistMassLambdaOnVsPt->Fill(lPt,lInvMassLambda);
444 fHistMassAntiLambdaOnVsPt->Fill(lPt,lInvMassAntiLambda);
445 }
446 fHistArmenterosPodolanskiOn->Fill(lAlphaV0,lPtArmV0);
447 }
448 }// This is the end of the V0 loop
449 } // end of "ESD" analysis
450
451 else if(fAnalysisType == "AOD") {
452
453 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
454 {// This is the begining of the V0 loop
455 AliAODv0 *v0 = ((AliAODEvent*)lEvent)->GetV0(iV0);
456 if (!v0) continue;
457
458 lV0Radius = v0->RadiusV0();
459 lPt = v0->Pt();
460 lRapK0Short = v0->RapK0Short();
461 lRapLambda = v0->RapLambda();
462 if ((lPt<fMinV0Pt)||(fMaxV0Pt<lPt)) continue;
463 lDcaPosToPrimVertex = v0->DcaPosToPrimVertex();
464 lDcaNegToPrimVertex = v0->DcaNegToPrimVertex();
465
466 // propose getOnFlyStatus() = -1 for like-sign
467
468 lOnFlyStatus = v0->GetOnFlyStatus();
469 lChi2V0 = v0->Chi2V0();
470 lDcaV0Daughters = v0->DcaV0Daughters();
471 lDcaV0ToPrimVertex = v0->DcaV0ToPrimVertex();
472 lV0CosineOfPointingAngle = v0->CosPointingAngle(tPrimaryVtxPosition);
473
474 lInvMassK0s = v0->MassK0Short();
475 lInvMassLambda = v0->MassLambda();
476 lInvMassAntiLambda = v0->MassAntiLambda();
477 lAlphaV0 = v0->AlphaV0();
478 lPtArmV0 = v0->PtArmV0();
479
480 fHistV0OnFlyStatus->Fill(lOnFlyStatus);
481 if(!lOnFlyStatus){
482 nv0sOff++;
483 fHistV0Chi2Off->Fill(lChi2V0);
484 fHistDcaV0ToPrimVertexOff->Fill(lDcaV0ToPrimVertex);
485 fHistDcaV0DaughtersOff->Fill(lDcaV0Daughters);
486 fHistV0CosineOfPointingAngleOff->Fill(lV0CosineOfPointingAngle);
487
488 fHistV0RadiusOff->Fill(lV0Radius);
489 fHistDcaPosToPrimVertexOff->Fill(lDcaPosToPrimVertex);
490 fHistDcaNegToPrimVertexOff->Fill(lDcaNegToPrimVertex);
491
492 // Filling invariant mass histos for all candidates
493 fHistMassK0sOff->Fill(lInvMassK0s);
494 fHistMassLambdaOff->Fill(lInvMassLambda);
495 fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
496 if (TMath::Abs(lRapK0Short)<fMaxV0Rapidity)
497 fHistMassK0sOffVsPt->Fill(lPt,lInvMassK0s);
498 if (TMath::Abs(lRapLambda)<fMaxV0Rapidity){
499 fHistMassLambdaOffVsPt->Fill(lPt,lInvMassLambda);
500 fHistMassAntiLambdaOffVsPt->Fill(lPt,lInvMassAntiLambda);
501 }
502 fHistArmenterosPodolanskiOff->Fill(lAlphaV0,lPtArmV0);
503 }
504 else {
505 nv0sOn++;
506 fHistV0Chi2On->Fill(lChi2V0);
507 fHistDcaV0ToPrimVertexOn->Fill(lDcaV0ToPrimVertex);
508 fHistDcaV0DaughtersOn->Fill(lDcaV0Daughters);
509 fHistV0CosineOfPointingAngleOn->Fill(lV0CosineOfPointingAngle);
510
511 fHistV0RadiusOn->Fill(lV0Radius);
512 fHistDcaPosToPrimVertexOn->Fill(lDcaPosToPrimVertex);
513 fHistDcaNegToPrimVertexOn->Fill(lDcaNegToPrimVertex);
514
515 // Filling invariant mass histos for all candidates
516 fHistMassK0sOn->Fill(lInvMassK0s);
517 fHistMassLambdaOn->Fill(lInvMassLambda);
518 fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
519 if (TMath::Abs(lRapK0Short)<fMaxV0Rapidity)
520 fHistMassK0sOnVsPt->Fill(lPt,lInvMassK0s);
521 if (TMath::Abs(lRapLambda)<fMaxV0Rapidity){
522 fHistMassLambdaOnVsPt->Fill(lPt,lInvMassLambda);
523 fHistMassAntiLambdaOnVsPt->Fill(lPt,lInvMassAntiLambda);
524 }
525 fHistArmenterosPodolanskiOn->Fill(lAlphaV0,lPtArmV0);
526 }
527 }// This is the end of the V0 loop
528 } // end of "AOD" analysis
529
530 fHistV0Multiplicity->Fill(nv0s);
531 fHistV0MultiplicityOff->Fill(nv0sOff);
532 fHistV0MultiplicityOn->Fill(nv0sOn);
533
534 }// end of event selection
976dc218 535
96cad837 536 // Post output data.
976dc218 537 PostData(1, fListHist);
7285f641 538}
96cad837 539
540//________________________________________________________________________
c8eada4f 541void AliAnalysisTaskCheckV0::Terminate(Option_t *)
96cad837 542{
543 // Draw result to the screen
544 // Called once at the end of the query
545
f00218b3 546 TList *cRetrievedList = 0x0;
547 cRetrievedList = (TList*)GetOutputData(1);
548
549 if(!cRetrievedList){
550 AliWarning("ERROR - AliAnalysisTaskPerformanceStrange: output data container list not available\n"); return;
551 }
552
50b713dc 553 // Implement a decent style
554 TStyle *myStyle = new TStyle("myStyle","my style");
555 Int_t font = 42;
556 myStyle->SetCanvasColor(10);
557 myStyle->SetStatColor(10);
558 myStyle->SetPadColor(10);
559 myStyle->SetDrawBorder(0);
560 myStyle->SetCanvasBorderMode(0);
561 myStyle->SetPadBorderMode(0);
562 myStyle->SetTextFont(font);
563 myStyle->SetStatFont(font);
564 myStyle->SetLabelFont(font,"xyz");
565 myStyle->SetTitleFont(font);
566 myStyle->SetPalette(1,0);
567 myStyle->cd();
568
f00218b3 569 fHistTrackMultiplicity = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistTrackMultiplicity"));
96cad837 570 if (!fHistTrackMultiplicity) {
571 Printf("ERROR: fHistTrackMultiplicity not available");
572 return;
573 }
f00218b3 574 fHistV0Multiplicity = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistV0Multiplicity"));
96cad837 575 if (!fHistV0Multiplicity) {
576 Printf("ERROR: fHistV0Multiplicity not available");
577 return;
578 }
f00218b3 579 fHistV0MultiplicityOff = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistV0MultiplicityOff"));
96cad837 580 if (!fHistV0MultiplicityOff) {
581 Printf("ERROR: fHistV0MultiplicityOff not available");
582 return;
583 }
f00218b3 584 fHistV0MultiplicityOn = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistV0MultiplicityOn"));
96cad837 585 if (!fHistV0MultiplicityOn) {
586 Printf("ERROR: fHistV0MultiplicityOn not available");
587 return;
588 }
50b713dc 589
590 TCanvas *canCheckV0 = new TCanvas("AliAnalysisTaskCheckV0","Check V0",10,10,510,700);
591 canCheckV0->Divide(2,3);
976dc218 592 if (fHistTrackMultiplicity->GetMaximum() > 0.) canCheckV0->cd(1)->SetLogy();
b638d7f5 593 fHistTrackMultiplicity->SetMarkerStyle(26);
96cad837 594 fHistTrackMultiplicity->DrawCopy("E");
b638d7f5 595 fHistV0Multiplicity->SetMarkerStyle(25);
96cad837 596 fHistV0Multiplicity->DrawCopy("ESAME");
597 fHistV0MultiplicityOff->SetMarkerStyle(24);
598 fHistV0MultiplicityOff->DrawCopy("ESAME");
599 fHistV0MultiplicityOn->SetMarkerStyle(20);
600 fHistV0MultiplicityOn->DrawCopy("ESAME");
c8eada4f 601
a95976c8 602 TLegend *legendMultiplicity = new TLegend(0.5,0.5,0.75,0.75);
c8eada4f 603 legendMultiplicity->AddEntry(fHistTrackMultiplicity,"tracks");
604 legendMultiplicity->AddEntry(fHistV0Multiplicity,"all V^{0}");
605 legendMultiplicity->AddEntry(fHistV0MultiplicityOff,"offline V^{0}");
606 legendMultiplicity->AddEntry(fHistV0MultiplicityOn,"onthefly V^{0}");
607 legendMultiplicity->Draw();
608
f00218b3 609 fHistMassK0sOff = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistMassK0sOff"));
9ea746fc 610 if (!fHistMassK0sOff) {
611 Printf("ERROR: fHistMassK0sOff not available");
c8eada4f 612 return;
613 }
f00218b3 614 fHistMassK0sOn = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistMassK0sOn"));
9ea746fc 615 if (!fHistMassK0sOn) {
616 Printf("ERROR: fHistMassK0sOn not available");
c8eada4f 617 return;
618 }
f00218b3 619 fHistMassLambdaOff = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistMassLambdaOff"));
c8eada4f 620 if (!fHistMassLambdaOff) {
621 Printf("ERROR: fHistMassLambdaOff not available");
622 return;
623 }
f00218b3 624 fHistMassLambdaOn = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistMassLambdaOn"));
c8eada4f 625 if (!fHistMassLambdaOn) {
626 Printf("ERROR: fHistMassLambdaOn not available");
627 return;
628 }
f00218b3 629 fHistMassAntiLambdaOff = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistMassAntiLambdaOff"));
c8eada4f 630 if (!fHistMassAntiLambdaOff) {
631 Printf("ERROR: fHistMassAntiLambdaOff not available");
632 return;
633 }
f00218b3 634 fHistMassAntiLambdaOn = dynamic_cast<TH1F*> (cRetrievedList->FindObject("fHistMassAntiLambdaOn"));
c8eada4f 635 if (!fHistMassAntiLambdaOn) {
636 Printf("ERROR: fHistMassAntiLambdaOn not available");
637 return;
638 }
f00218b3 639 fHistArmenterosPodolanskiOff = dynamic_cast<TH2F*> (cRetrievedList->FindObject("fHistArmenterosPodolanskiOff"));
50b713dc 640 if (!fHistArmenterosPodolanskiOff) {
641 Printf("ERROR: fHistArmenterosPodolanskiOff not available");
642 return;
643 }
f00218b3 644 fHistArmenterosPodolanskiOn = dynamic_cast<TH2F*> (cRetrievedList->FindObject("fHistArmenterosPodolanskiOn"));
50b713dc 645 if (!fHistArmenterosPodolanskiOn) {
646 Printf("ERROR: fHistArmenterosPodolanskiOn not available");
647 return;
648 }
c8eada4f 649
650 canCheckV0->cd(2);
9ea746fc 651 fHistMassK0sOn->SetMarkerStyle(20);
652 fHistMassK0sOn->DrawCopy("E");
653 fHistMassK0sOff->SetMarkerStyle(24);
654 fHistMassK0sOff->DrawCopy("ESAME");
c8eada4f 655
656 canCheckV0->cd(3);
c8eada4f 657 fHistMassLambdaOn->SetMarkerStyle(20);
a95976c8 658 fHistMassLambdaOn->DrawCopy("E");
659 fHistMassLambdaOff->SetMarkerStyle(24);
660 fHistMassLambdaOff->DrawCopy("ESAME");
c8eada4f 661
662 canCheckV0->cd(4);
c8eada4f 663 fHistMassAntiLambdaOn->SetMarkerStyle(20);
a95976c8 664 fHistMassAntiLambdaOn->DrawCopy("E");
665 fHistMassAntiLambdaOff->SetMarkerStyle(24);
666 fHistMassAntiLambdaOff->DrawCopy("ESAME");
976dc218 667
50b713dc 668 canCheckV0->cd(5);
669 fHistArmenterosPodolanskiOff->DrawCopy("COL2Z");
670 canCheckV0->cd(6);
671 fHistArmenterosPodolanskiOn->DrawCopy("COL2Z");
672
96cad837 673}