2 /**************************************************************************
3 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
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 //-----------------------------------------------------------------
28 #include "AliAnalysisTaskSE.h"
30 #include "AliESDVertex.h"
31 #include "AliESDEvent.h"
32 #include "AliAODEvent.h"
35 #include "AliESDtrack.h"
37 #include "AliAODTrack.h"
41 #include "AliAnalysisTaskStrange.h"
43 ClassImp(AliAnalysisTaskStrange)
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)
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)
92 // Define output slots only here
93 // Output slot #1 writes into a TList container
94 DefineOutput(1, TList::Class());
97 //________________________________________________________________________
98 void AliAnalysisTaskStrange::UserCreateOutputObjects()
103 fListHist = new TList();
104 fListHist->SetOwner();
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);
115 if (!fHistTrackMultiplicity) {
116 if (fCollidingSystems)
117 fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 200, 0, 40000);
119 fHistTrackMultiplicity = new TH1F("fHistTrackMultiplicity", "Multiplicity distribution;Number of tracks;Events", 250, 0, 250);
120 fListHist->Add(fHistTrackMultiplicity);
122 if (!fHistV0Multiplicity) {
123 if (fCollidingSystems)
124 fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
126 fHistV0Multiplicity = new TH1F("fHistV0Multiplicity", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
127 fListHist->Add(fHistV0Multiplicity);
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);
152 // bounds of histograms:
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;
157 // V0 offline distributions
158 if (!fHistV0MultiplicityOff) {
159 if (fCollidingSystems)
160 fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 200, 0, 40000);
162 fHistV0MultiplicityOff = new TH1F("fHistV0MultiplicityOff", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
163 fListHist->Add(fHistV0MultiplicityOff);
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);
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);
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);
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);
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);
202 fHistV0MultiplicityOn = new TH1F("fHistV0MultiplicityOn", "Multiplicity distribution;Number of V0s;Events", 50, 0, 50);
203 fListHist->Add(fHistV0MultiplicityOn);
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);
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);
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);
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);
238 PostData(1, fListHist);
241 //________________________________________________________________________
242 void AliAnalysisTaskStrange::UserExec(Option_t *)
245 // Called for each event
246 AliVEvent* lEvent = InputEvent();
248 Printf("ERROR: Event not available");
252 if (!(lEvent->GetNumberOfTracks())) {
253 //Printf("Strange analysis task: There is no track in this event");
256 fHistTrackMultiplicity->Fill(lEvent->GetNumberOfTracks());
258 Double_t tPrimaryVtxPosition[3];
261 nv0s = lEvent->GetNumberOfV0s();
262 //Printf("Strange analysis task: There are %d v0s in this event",nv0s);
264 Int_t lOnFlyStatus = 0, nv0sOn = 0, nv0sOff = 0;
265 Double_t lChi2V0 = 0;
266 Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0;
267 Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0;
268 Double_t lV0CosineOfPointingAngle = 0;
269 Double_t lV0Radius = 0;
270 Double_t lV0DecayLength = 0;
271 Double_t lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0;
272 Double_t lPt = 0, lRapK0s = 0, lRapLambda = 0;
273 Double_t lAlphaV0 = 0, lPtArmV0 = 0;
275 Double_t tV0Position[3];
277 Double_t lMagneticField = 999;
280 //***********************
282 //***********************
284 if(fAnalysisType == "ESD") {
286 const AliESDVertex *primaryVtx = ((AliESDEvent*)lEvent)->GetPrimaryVertex();
287 tPrimaryVtxPosition[0] = primaryVtx->GetXv();
288 tPrimaryVtxPosition[1] = primaryVtx->GetYv();
289 tPrimaryVtxPosition[2] = primaryVtx->GetZv();
291 fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
292 fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
293 fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
295 lMagneticField = ((AliESDEvent*)lEvent)->GetMagneticField();
298 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
299 {// This is the begining of the V0 loop
300 AliESDv0 *v0 = ((AliESDEvent*)lEvent)->GetV0(iV0);
303 // AliAODtrack (V0 Daughters)
304 UInt_t lKeyPos = (UInt_t)TMath::Abs(v0->GetPindex());
305 UInt_t lKeyNeg = (UInt_t)TMath::Abs(v0->GetNindex());
307 AliESDtrack *pTrack = ((AliESDEvent*)lEvent)->GetTrack(lKeyPos);
308 AliESDtrack *nTrack = ((AliESDEvent*)lEvent)->GetTrack(lKeyNeg);
309 if (!pTrack || !nTrack) {
310 Printf("ERROR: Could not retreive one of the daughter track");
315 if ( pTrack->GetSign() == nTrack->GetSign()){
316 //cout<< "like sign, continue"<< endl;
320 // Tracks quality cuts
321 if ( ( (pTrack->GetTPCNcls()) < 80 ) || ( (nTrack->GetTPCNcls()) < 80 ) ) continue;
323 // TPC refit condition (done during reconstruction for Offline but not for On-the-fly)
324 if( !(pTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
325 if( !(nTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue;
327 // DCA between daughter and Primary Vertex:
328 if (pTrack) lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],lMagneticField) );
330 if (nTrack) lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],lMagneticField) );
332 // VO's main characteristics:
333 lOnFlyStatus = v0->GetOnFlyStatus();
334 lChi2V0 = v0->GetChi2V0();
335 lDcaV0Daughters = v0->GetDcaV0Daughters();
336 lDcaV0ToPrimVertex = v0->GetD(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1],tPrimaryVtxPosition[2]);
337 lV0CosineOfPointingAngle = v0->GetV0CosineOfPointingAngle(tPrimaryVtxPosition[0],tPrimaryVtxPosition[1], tPrimaryVtxPosition[2]);
338 v0->GetXYZ(tV0Position[0], tV0Position[1], tV0Position[2]);
339 lV0Radius = TMath::Sqrt(tV0Position[0]*tV0Position[0]+tV0Position[1]*tV0Position[1]);
340 lV0DecayLength = TMath::Sqrt(TMath::Power(tV0Position[0] - tPrimaryVtxPosition[0],2) +
341 TMath::Power(tV0Position[1] - tPrimaryVtxPosition[1],2) +
342 TMath::Power(tV0Position[2] - tPrimaryVtxPosition[2],2 ));
345 v0->ChangeMassHypothesis(310);
346 lInvMassK0s = v0->GetEffMass();
347 v0->ChangeMassHypothesis(3122);
348 lInvMassLambda = v0->GetEffMass();
349 v0->ChangeMassHypothesis(-3122);
350 lInvMassAntiLambda = v0->GetEffMass();
353 lRapK0s = v0->Y(310);
354 lRapLambda = v0->Y(3122);
359 // Armenteros variables: !!
360 lAlphaV0 = v0->AlphaV0();
361 lPtArmV0 = v0->PtArmV0();
364 if (fUseCut.Contains("yes")) {
365 if ( (lDcaPosToPrimVertex < 0.036 )||
366 (lDcaNegToPrimVertex < 0.036 )||
367 (lDcaV0Daughters > 0.5 ) ||
368 (lV0CosineOfPointingAngle < 0.99)
372 // Filling histograms
373 fHistDcaPosToPrimVertex->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
374 fHistDcaNegToPrimVertex->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
375 fHistDcaPosToPrimVertexZoom->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
376 fHistDcaNegToPrimVertexZoom->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
377 fHistRadiusV0->Fill(lV0Radius,lOnFlyStatus);
378 fHistDecayLengthV0->Fill(lV0DecayLength,lOnFlyStatus);
379 fHistDcaV0Daughters->Fill(lDcaV0Daughters,lOnFlyStatus);
380 fHistChi2->Fill(lChi2V0,lOnFlyStatus);
381 fHistCosPointAngle->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
382 if (lV0CosineOfPointingAngle >= 0.9) fHistCosPointAngleZoom->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
385 fHistPtVsYK0sOff->Fill(lPt,lRapK0s);
386 fHistPtVsYLambdaOff->Fill(lPt,lRapLambda);
387 fHistPtVsYAntiLambdaOff->Fill(lPt,lRapLambda);
388 fHistArmenterosPodolanskiOff->Fill(lAlphaV0,lPtArmV0);
392 fHistPtVsYK0sOn->Fill(lPt,lRapK0s);
393 fHistPtVsYLambdaOn->Fill(lPt,lRapLambda);
394 fHistPtVsYAntiLambdaOn->Fill(lPt,lRapLambda);
395 fHistArmenterosPodolanskiOn->Fill(lAlphaV0,lPtArmV0);
397 // K0s invariant mass histograms:
398 if (TMath::Abs(lRapK0s) < 1) {
400 fHistMassK0sOff->Fill(lInvMassK0s);
401 fHistMassVsRadiusK0sOff->Fill(lV0Radius,lInvMassK0s);
402 fHistPtVsMassK0sOff->Fill(lInvMassK0s,lPt);
405 fHistMassK0sOn->Fill(lInvMassK0s);
406 fHistMassVsRadiusK0sOn->Fill(lV0Radius,lInvMassK0s);
407 fHistPtVsMassK0sOn->Fill(lInvMassK0s,lPt);
410 // Lambda and AntiLambda invariant mass histograms:
411 if (TMath::Abs(lRapLambda) < 1) {
413 fHistMassLambdaOff->Fill(lInvMassLambda);
414 fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
415 fHistMassVsRadiusLambdaOff->Fill(lV0Radius,lInvMassLambda);
416 fHistMassVsRadiusAntiLambdaOff->Fill(lV0Radius,lInvMassAntiLambda);
417 fHistPtVsMassLambdaOff->Fill(lInvMassLambda,lPt);
418 fHistPtVsMassAntiLambdaOff->Fill(lInvMassAntiLambda,lPt);
421 fHistMassLambdaOn->Fill(lInvMassLambda);
422 fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
423 fHistMassVsRadiusLambdaOn->Fill(lV0Radius,lInvMassLambda);
424 fHistMassVsRadiusAntiLambdaOn->Fill(lV0Radius,lInvMassAntiLambda);
425 fHistPtVsMassLambdaOn->Fill(lInvMassLambda,lPt);
426 fHistPtVsMassAntiLambdaOn->Fill(lInvMassAntiLambda,lPt);
433 //***********************
435 //***********************
437 else if(fAnalysisType == "AOD") {
439 const AliAODVertex *primaryVtx = ((AliAODEvent*)lEvent)->GetPrimaryVertex();
440 tPrimaryVtxPosition[0] = primaryVtx->GetX();
441 tPrimaryVtxPosition[1] = primaryVtx->GetY();
442 tPrimaryVtxPosition[2] = primaryVtx->GetZ();
444 fHistPrimaryVertexPosX->Fill(tPrimaryVtxPosition[0]);
445 fHistPrimaryVertexPosY->Fill(tPrimaryVtxPosition[1]);
446 fHistPrimaryVertexPosZ->Fill(tPrimaryVtxPosition[2]);
448 for (Int_t iV0 = 0; iV0 < nv0s; iV0++)
449 {// This is the begining of the V0 loop
450 AliAODv0 *myAODv0 = ((AliAODEvent*)lEvent)->GetV0(iV0);
451 if (!myAODv0) continue;
454 lV0Radius = myAODv0->RadiusV0();
455 lDcaPosToPrimVertex = myAODv0->DcaPosToPrimVertex();
456 lDcaNegToPrimVertex = myAODv0->DcaNegToPrimVertex();
457 lOnFlyStatus = myAODv0->GetOnFlyStatus();
458 lChi2V0 = myAODv0->Chi2V0();
459 lDcaV0Daughters = myAODv0->DcaV0Daughters();
460 lDcaV0ToPrimVertex = myAODv0->DcaV0ToPrimVertex();
461 lV0DecayLength = myAODv0->DecayLengthV0(tPrimaryVtxPosition);
462 lV0CosineOfPointingAngle = myAODv0->CosPointingAngle(tPrimaryVtxPosition);
464 lInvMassK0s = myAODv0->MassK0Short();
465 lInvMassLambda = myAODv0->MassLambda();
466 lInvMassAntiLambda = myAODv0->MassAntiLambda();
468 lPt = TMath::Sqrt(myAODv0->Pt2V0());
469 lRapK0s = myAODv0->RapK0Short();
470 lRapLambda = myAODv0->RapLambda();
471 lAlphaV0 = myAODv0->AlphaV0();
472 lPtArmV0 = myAODv0->PtArmV0();
476 if (fUseCut.Contains("yes")) {
477 if ( (lDcaPosToPrimVertex < 0.036 )||
478 (lDcaNegToPrimVertex < 0.036 )||
479 (lDcaV0Daughters > 0.5 ) ||
480 (lV0CosineOfPointingAngle < 0.99)
484 // Filling histograms
485 fHistDcaPosToPrimVertex->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
486 fHistDcaNegToPrimVertex->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
487 fHistDcaPosToPrimVertexZoom->Fill(lDcaPosToPrimVertex,lOnFlyStatus);
488 fHistDcaNegToPrimVertexZoom->Fill(lDcaNegToPrimVertex,lOnFlyStatus);
489 fHistRadiusV0->Fill(lV0Radius,lOnFlyStatus);
490 fHistDecayLengthV0->Fill(lV0DecayLength,lOnFlyStatus);
491 fHistDcaV0Daughters->Fill(lDcaV0Daughters,lOnFlyStatus);
492 fHistChi2->Fill(lChi2V0,lOnFlyStatus);
493 fHistCosPointAngle->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
494 if (lV0CosineOfPointingAngle >= 0.9) fHistCosPointAngleZoom->Fill(lV0CosineOfPointingAngle,lOnFlyStatus);
497 fHistPtVsYK0sOff->Fill(lPt,lRapK0s);
498 fHistPtVsYLambdaOff->Fill(lPt,lRapLambda);
499 fHistPtVsYAntiLambdaOff->Fill(lPt,lRapLambda);
500 fHistArmenterosPodolanskiOff->Fill(lAlphaV0,lPtArmV0);
504 fHistPtVsYK0sOn->Fill(lPt,lRapK0s);
505 fHistPtVsYLambdaOn->Fill(lPt,lRapLambda);
506 fHistPtVsYAntiLambdaOn->Fill(lPt,lRapLambda);
507 fHistArmenterosPodolanskiOn->Fill(lAlphaV0,lPtArmV0);
509 // K0s invariant mass histograms:
510 if (TMath::Abs(lRapK0s) < 1) {
512 fHistMassK0sOff->Fill(lInvMassK0s);
513 fHistMassVsRadiusK0sOff->Fill(lV0Radius,lInvMassK0s);
514 fHistPtVsMassK0sOff->Fill(lInvMassK0s,lPt);
517 fHistMassK0sOn->Fill(lInvMassK0s);
518 fHistMassVsRadiusK0sOn->Fill(lV0Radius,lInvMassK0s);
519 fHistPtVsMassK0sOn->Fill(lInvMassK0s,lPt);
522 // Lambda and AntiLambda invariant mass histograms:
523 if (TMath::Abs(lRapLambda) < 1) {
525 fHistMassLambdaOff->Fill(lInvMassLambda);
526 fHistMassAntiLambdaOff->Fill(lInvMassAntiLambda);
527 fHistMassVsRadiusLambdaOff->Fill(lV0Radius,lInvMassLambda);
528 fHistMassVsRadiusAntiLambdaOff->Fill(lV0Radius,lInvMassAntiLambda);
529 fHistPtVsMassLambdaOff->Fill(lInvMassLambda,lPt);
530 fHistPtVsMassAntiLambdaOff->Fill(lInvMassAntiLambda,lPt);
533 fHistMassLambdaOn->Fill(lInvMassLambda);
534 fHistMassAntiLambdaOn->Fill(lInvMassAntiLambda);
535 fHistMassVsRadiusLambdaOn->Fill(lV0Radius,lInvMassLambda);
536 fHistMassVsRadiusAntiLambdaOn->Fill(lV0Radius,lInvMassAntiLambda);
537 fHistPtVsMassLambdaOn->Fill(lInvMassLambda,lPt);
538 fHistPtVsMassAntiLambdaOn->Fill(lInvMassAntiLambda,lPt);
544 fHistV0Multiplicity->Fill(nv0s);
545 fHistV0MultiplicityOff->Fill(nv0sOff);
546 fHistV0MultiplicityOn->Fill(nv0sOn);
549 PostData(1, fListHist);
552 //________________________________________________________________________
553 void AliAnalysisTaskStrange::Terminate(Option_t *)
555 // Draw result to the screen
556 // Called once at the end of the query