]>
Commit | Line | Data |
---|---|---|
397cdc9f | 1 | /************************************************************************** |
2 | * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: The ALICE Off-line Project. * | |
5 | * Contributors are mentioned in the code where appropriate. * | |
6 | * * | |
7 | * Permission to use, copy, modify and distribute this software and its * | |
8 | * documentation strictly for non-commercial purposes is hereby granted * | |
9 | * without fee, provided that the above copyright notice appears in all * | |
10 | * copies and that both the copyright notice and this permission notice * | |
11 | * appear in the supporting documentation. The authors make no claims * | |
12 | * about the suitability of this software for any purpose. It is * | |
13 | * provided "as is" without express or implied warranty. * | |
14 | **************************************************************************/ | |
15 | ||
16 | // +-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+ | |
17 | // | |
18 | // Modified version of AliAnalysisTaskCheckCascade.cxx. | |
19 | // This is a 'hybrid' output version, in that it uses a classic TTree | |
20 | // ROOT object to store the candidates, plus a couple of histograms filled on | |
21 | // a per-event basis for storing variables too numerous to put in a tree. | |
22 | // | |
23 | // --- Added bits of code for checking V0s | |
24 | // (from AliAnalysisTaskCheckStrange.cxx) | |
25 | // | |
26 | // --- Algorithm Description | |
27 | // 1. Perform Physics Selection | |
28 | // 2. Perform Primary Vertex |z|<10cm selection | |
29 | // 3. Perform Primary Vertex NoTPCOnly vertexing selection (>0 contrib.) | |
30 | // 4. Perform Pileup Rejection | |
31 | // 5. Analysis Loops: | |
32 | // 5a. Fill TTree object with V0 information, candidates | |
33 | // | |
34 | // Please Report Any Bugs! | |
35 | // | |
36 | // --- David Dobrigkeit Chinellato | |
37 | // (david.chinellato@gmail.com) | |
38 | // | |
39 | // +-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+-+ | |
40 | ||
41 | class TTree; | |
42 | class TParticle; | |
43 | class TVector3; | |
44 | ||
45 | //class AliMCEventHandler; | |
46 | //class AliMCEvent; | |
47 | //class AliStack; | |
48 | ||
49 | class AliESDVertex; | |
50 | class AliAODVertex; | |
51 | class AliESDv0; | |
52 | class AliAODv0; | |
53 | ||
54 | #include <Riostream.h> | |
55 | #include "TList.h" | |
56 | #include "TH1.h" | |
57 | #include "TH2.h" | |
58 | #include "TH3.h" | |
59 | #include "THnSparse.h" | |
60 | #include "TVector3.h" | |
61 | #include "TCanvas.h" | |
62 | #include "TMath.h" | |
63 | #include "TLegend.h" | |
64 | #include "AliLog.h" | |
65 | #include "AliCentrality.h" | |
66 | #include "AliESDEvent.h" | |
67 | #include "AliAODEvent.h" | |
68 | #include "AliV0vertexer.h" | |
69 | #include "AliCascadeVertexer.h" | |
70 | #include "AliESDpid.h" | |
71 | #include "AliESDtrack.h" | |
72 | ||
73 | #include "AliESDtrackCuts.h" | |
74 | #include "AliInputEventHandler.h" | |
75 | #include "AliAnalysisManager.h" | |
76 | #include "AliMCEventHandler.h" | |
77 | ||
78 | #include "AliCFContainer.h" | |
79 | #include "AliMultiplicity.h" | |
80 | ||
81 | #include "AliESDcascade.h" | |
82 | #include "AliAODcascade.h" | |
83 | #include "AliESDUtils.h" | |
84 | #include "AliESDHeader.h" | |
85 | #include "AliAODTrack.h" | |
86 | #include "AliAnalysisTaskQAV0AOD.h" | |
87 | ||
88 | //debugging purposes | |
89 | #include "TObjectTable.h" | |
90 | ||
91 | ClassImp(AliAnalysisTaskQAV0AOD) | |
92 | ||
93 | AliAnalysisTaskQAV0AOD::AliAnalysisTaskQAV0AOD() | |
94 | : AliAnalysisTaskSE(), | |
95 | //Output lists | |
96 | fOutput(0), | |
97 | ||
98 | //Histos | |
99 | fHistEvent(0), | |
100 | fHistTopDCANegToPV(0), | |
101 | fHistTopDCAPosToPV(0), | |
102 | fHistTopDCAV0Daughters(0), | |
103 | fHistTopCosinePA(0), | |
104 | fHistTopV0Radius(0), | |
105 | fHistSelectedTopDCANegToPV(0), | |
106 | fHistSelectedTopDCAPosToPV(0), | |
107 | fHistSelectedTopDCAV0Daughters(0), | |
108 | fHistSelectedTopCosinePA(0), | |
109 | fHistSelectedTopV0Radius(0), | |
110 | ||
111 | f2dHistInvMassK0Short(0), | |
112 | f2dHistInvMassLambda(0), | |
113 | f2dHistInvMassAntiLambda(0), | |
114 | ||
115 | f2dHistInvMassWithdEdxK0Short(0), | |
116 | f2dHistInvMassWithdEdxLambda(0), | |
117 | f2dHistInvMassWithdEdxAntiLambda(0), | |
118 | ||
119 | f2dHistResponseNegativeAsPion(0), | |
120 | f2dHistResponseNegativeAsProton(0), | |
121 | f2dHistResponsePositiveAsPion(0), | |
122 | f2dHistResponsePositiveAsProton(0), | |
123 | ||
124 | f2dHistdEdxSignalPionFromLambda(0), | |
125 | f2dHistdEdxSignalProtonFromLambda(0), | |
126 | f2dHistResponsePionFromLambda(0), | |
127 | f2dHistResponseProtonFromLambda(0), | |
128 | ||
129 | //Task Control / Utils | |
130 | fPIDResponse(0) | |
131 | { | |
132 | // Dummy Constructor | |
133 | for(Int_t iV0selIdx = 0; iV0selIdx < 7; iV0selIdx++ ) { fV0Sels [iV0selIdx ] = -1.; } | |
134 | } | |
135 | ||
136 | AliAnalysisTaskQAV0AOD::AliAnalysisTaskQAV0AOD(const char *name) | |
137 | : AliAnalysisTaskSE(name), | |
138 | //Output lists | |
139 | fOutput(0), | |
140 | ||
141 | //Histos | |
142 | fHistEvent(0), | |
143 | fHistTopDCANegToPV(0), | |
144 | fHistTopDCAPosToPV(0), | |
145 | fHistTopDCAV0Daughters(0), | |
146 | fHistTopCosinePA(0), | |
147 | fHistTopV0Radius(0), | |
148 | fHistSelectedTopDCANegToPV(0), | |
149 | fHistSelectedTopDCAPosToPV(0), | |
150 | fHistSelectedTopDCAV0Daughters(0), | |
151 | fHistSelectedTopCosinePA(0), | |
152 | fHistSelectedTopV0Radius(0), | |
153 | ||
154 | f2dHistInvMassK0Short(0), | |
155 | f2dHistInvMassLambda(0), | |
156 | f2dHistInvMassAntiLambda(0), | |
157 | ||
158 | f2dHistInvMassWithdEdxK0Short(0), | |
159 | f2dHistInvMassWithdEdxLambda(0), | |
160 | f2dHistInvMassWithdEdxAntiLambda(0), | |
161 | ||
162 | f2dHistResponseNegativeAsPion(0), | |
163 | f2dHistResponseNegativeAsProton(0), | |
164 | f2dHistResponsePositiveAsPion(0), | |
165 | f2dHistResponsePositiveAsProton(0), | |
166 | ||
167 | f2dHistdEdxSignalPionFromLambda(0), | |
168 | f2dHistdEdxSignalProtonFromLambda(0), | |
169 | f2dHistResponsePionFromLambda(0), | |
170 | f2dHistResponseProtonFromLambda(0), | |
171 | ||
172 | //Task Control / Utils | |
173 | fPIDResponse(0) | |
174 | { | |
175 | // Constructor | |
176 | //SELECTION CUTS | |
177 | // REALLY LOOSE? Be careful when attempting to run over PbPb if fkRunV0Vertexer is set! | |
178 | fV0Sels[0] = 33. ; // max allowed chi2 | |
179 | fV0Sels[1] = 0.02; // min allowed impact parameter for the 1st daughter (LHC09a4 : 0.05) | |
180 | fV0Sels[2] = 0.02; // min allowed impact parameter for the 2nd daughter (LHC09a4 : 0.05) | |
181 | fV0Sels[3] = 2.0 ; // max allowed DCA between the daughter tracks (LHC09a4 : 0.5) | |
182 | fV0Sels[4] = 0.95; // min allowed cosine of V0's pointing angle (LHC09a4 : 0.99) | |
183 | fV0Sels[5] = 0.5 ; // min radius of the fiducial volume (LHC09a4 : 0.2) | |
184 | fV0Sels[6] = 200. ; // max radius of the fiducial volume (LHC09a4 : 100.0) | |
185 | ||
186 | // Output slot #0 writes into a TList container (Lambda Histos and fTree) | |
187 | DefineOutput(1, TList::Class()); | |
188 | } | |
189 | ||
190 | ||
191 | AliAnalysisTaskQAV0AOD::~AliAnalysisTaskQAV0AOD() | |
192 | { | |
193 | //------------------------------------------------ | |
194 | // DESTRUCTOR | |
195 | //------------------------------------------------ | |
196 | ||
197 | if (fOutput){ | |
198 | delete fOutput; | |
199 | fOutput = 0x0; | |
200 | } | |
201 | } | |
202 | ||
203 | ||
204 | ||
205 | //________________________________________________________________________ | |
206 | void AliAnalysisTaskQAV0AOD::UserCreateOutputObjects() | |
207 | { | |
208 | ||
209 | //Define Output Lists | |
210 | fOutput = new TList(); | |
211 | fOutput->SetOwner(); | |
212 | ||
213 | ||
214 | //Histogram Output: Event-by-Event | |
215 | fHistEvent = new TH1D( "fHistEvent", ";Evt. Sel. Step;Count",4,0,4); | |
216 | fHistEvent->GetXaxis()->SetBinLabel(1, "Processed"); | |
217 | fHistEvent->GetXaxis()->SetBinLabel(2, "Phys-Sel"); | |
218 | fHistEvent->GetXaxis()->SetBinLabel(3, "Has Vtx"); | |
219 | fHistEvent->GetXaxis()->SetBinLabel(4, "Vtx |z|<10cm"); | |
220 | fOutput->Add(fHistEvent); | |
221 | ||
222 | //Topological Selection Histograms, 1D | |
223 | fHistTopDCANegToPV = new TH1D( "fHistTopDCANegToPV",";DCA Neg. Daughter to PV (cm);Counts",200,0,1); | |
224 | fHistTopDCAPosToPV = new TH1D( "fHistTopDCAPosToPV",";DCA Pos. Daughter to PV (cm);Counts",200,0,1); | |
225 | fHistTopDCAV0Daughters = new TH1D( "fHistTopDCAV0Daughters",";DCA V0 Daughters (#sigma);Counts",200,0,2); | |
226 | fHistTopCosinePA = new TH1D( "fHistTopCosinePA",";Cosine of PA;Counts",200,-1,1); | |
227 | fHistTopV0Radius = new TH1D( "fHistTopV0Radius",";Decay Radius (cm);Counts",200,0.,10); | |
228 | ||
229 | fOutput->Add( fHistTopDCANegToPV ); | |
230 | fOutput->Add( fHistTopDCAPosToPV ); | |
231 | fOutput->Add( fHistTopDCAV0Daughters ); | |
232 | fOutput->Add( fHistTopCosinePA ); | |
233 | fOutput->Add( fHistTopV0Radius ); | |
234 | ||
235 | //Zoomed In | |
236 | fHistSelectedTopDCANegToPV = new TH1D( "fHistSelectedTopDCANegToPV",";DCA Neg. Daughter to PV (cm);Counts",200,fV0Sels[1],1); | |
237 | fHistSelectedTopDCAPosToPV = new TH1D( "fHistSelectedTopDCAPosToPV",";DCA Pos. Daughter to PV (cm);Counts",200,fV0Sels[2],1); | |
238 | fHistSelectedTopDCAV0Daughters = new TH1D( "fHistSelectedTopDCAV0Daughters",";DCA V0 Daughters (#sigma);Counts",200,0,fV0Sels[3]); | |
239 | fHistSelectedTopCosinePA = new TH1D( "fHistSelectedTopCosinePA",";Cosine of PA;Counts",200,fV0Sels[4],1); | |
240 | fHistSelectedTopV0Radius = new TH1D( "fHistSelectedTopV0Radius",";Decay Radius (cm);Counts",200,fV0Sels[5],10); | |
241 | ||
242 | fOutput->Add( fHistSelectedTopDCANegToPV ); | |
243 | fOutput->Add( fHistSelectedTopDCAPosToPV ); | |
244 | fOutput->Add( fHistSelectedTopDCAV0Daughters ); | |
245 | fOutput->Add( fHistSelectedTopCosinePA ); | |
246 | fOutput->Add( fHistSelectedTopV0Radius ); | |
247 | ||
248 | //Invariant Mass Plots | |
249 | f2dHistInvMassK0Short = new TH2D( "f2dHistInvMassK0Short" , ";p_{T};M(#pi^{+},#pi^{-})" ,250,0,25,500,0.25,0.75); | |
250 | f2dHistInvMassLambda = new TH2D( "f2dHistInvMassLambda" , ";p_{T};M(p,#pi^{-})" ,250,0,25,300,1.07,1.115+0.255); | |
251 | f2dHistInvMassAntiLambda = new TH2D( "f2dHistInvMassAntiLambda", ";p_{T};M(#pi^{+},#bar{p})" ,250,0,25,300,1.07,1.115+0.255); | |
252 | ||
253 | fOutput->Add( f2dHistInvMassK0Short ); | |
254 | fOutput->Add( f2dHistInvMassLambda ); | |
255 | fOutput->Add( f2dHistInvMassAntiLambda ); | |
256 | ||
257 | //Invariant Mass Plots, with dEdx | |
258 | f2dHistInvMassWithdEdxK0Short = new TH2D( "f2dHistInvMassWithdEdxK0Short" , ";p_{T};M(#pi^{+},#pi^{-})" ,250,0,25,500,0.25,0.75); | |
259 | f2dHistInvMassWithdEdxLambda = new TH2D( "f2dHistInvMassWithdEdxLambda" , ";p_{T};M(p,#pi^{-})" ,250,0,25,300,1.07,1.115+0.255); | |
260 | f2dHistInvMassWithdEdxAntiLambda = new TH2D( "f2dHistInvMassWithdEdxAntiLambda", ";p_{T};M(#pi^{+},#bar{p})" ,250,0,25,300,1.07,1.115+0.255); | |
261 | ||
262 | fOutput->Add( f2dHistInvMassWithdEdxK0Short ); | |
263 | fOutput->Add( f2dHistInvMassWithdEdxLambda ); | |
264 | fOutput->Add( f2dHistInvMassWithdEdxAntiLambda ); | |
265 | ||
266 | //dE/dx QA for main analysis and for calibration check | |
267 | f2dHistResponseNegativeAsPion = new TH2D( "f2dHistResponseNegativeAsPion", ";p_{T}^{V0};N#sigma",500,0,5,400,-20,20); | |
268 | f2dHistResponseNegativeAsProton = new TH2D( "f2dHistResponseNegativeAsProton", ";p_{T}^{V0};N#sigma",500,0,5,400,-20,20); | |
269 | f2dHistResponsePositiveAsPion = new TH2D( "f2dHistResponsePositiveAsPion", ";p_{T}^{V0};N#sigma",500,0,5,400,-20,20); | |
270 | f2dHistResponsePositiveAsProton = new TH2D( "f2dHistResponsePositiveAsProton", ";p_{T}^{V0};N#sigma",500,0,5,400,-20,20); | |
271 | ||
272 | //Clean Signal Check from Lambdas: stricter cuts, raw signal check | |
273 | f2dHistdEdxSignalPionFromLambda = new TH2D( "f2dHistdEdxSignalPionFromLambda", ";p_{T}^{V0};TPC Signal",500,0,5,8000,0,800); | |
274 | f2dHistdEdxSignalProtonFromLambda = new TH2D( "f2dHistdEdxSignalProtonFromLambda", ";p_{T}^{V0};TPC Signal",500,0,5,8000,0,800); | |
275 | f2dHistResponsePionFromLambda = new TH2D( "f2dHistResponsePionFromLambda", ";p_{T}^{V0};N#sigma",500,0,5,400,-20,20); | |
276 | f2dHistResponseProtonFromLambda = new TH2D( "f2dHistResponseProtonFromLambda", ";p_{T}^{V0};N#sigma",500,0,5,400,-20,20); | |
277 | ||
278 | ||
279 | fOutput->Add( f2dHistResponseNegativeAsPion ); | |
280 | fOutput->Add( f2dHistResponseNegativeAsProton ); | |
281 | fOutput->Add( f2dHistResponsePositiveAsPion ); | |
282 | fOutput->Add( f2dHistResponsePositiveAsProton ); | |
283 | ||
284 | fOutput->Add( f2dHistdEdxSignalPionFromLambda ); | |
285 | fOutput->Add( f2dHistdEdxSignalProtonFromLambda ); | |
286 | fOutput->Add( f2dHistResponsePionFromLambda ); | |
287 | fOutput->Add( f2dHistResponseProtonFromLambda ); | |
288 | ||
289 | //------------------------------------------------ | |
290 | // Particle Identification Setup | |
291 | //------------------------------------------------ | |
292 | ||
293 | AliAnalysisManager *man=AliAnalysisManager::GetAnalysisManager(); | |
294 | AliInputEventHandler* inputHandler = (AliInputEventHandler*) (man->GetInputEventHandler()); | |
295 | fPIDResponse = inputHandler->GetPIDResponse(); | |
296 | ||
297 | //Regular output: Histograms | |
298 | PostData(1, fOutput); | |
299 | }// end UserCreateOutputObjects | |
300 | ||
301 | ||
302 | //________________________________________________________________________ | |
303 | void AliAnalysisTaskQAV0AOD::UserExec(Option_t *) | |
304 | { | |
305 | ||
306 | // Main loop | |
307 | // Called for each event | |
308 | //gObjectTable->Print(); | |
309 | AliAODEvent *lAODevent = 0x0; | |
310 | ||
311 | //AliAODEvent *lAODevent = 0x0; | |
312 | Int_t nV0s = -1; | |
313 | ||
314 | Double_t lBestPrimaryVtxPos[3] = {-100.0, -100.0, -100.0}; | |
315 | Double_t lMagneticField = -10.; | |
316 | ||
317 | // Connect to the InputEvent | |
318 | // After these lines, we should have an ESD/AOD event + the number of cascades in it. | |
319 | ||
320 | // Appropriate for ESD analysis! | |
321 | ||
322 | lAODevent = dynamic_cast<AliAODEvent*>( InputEvent() ); | |
323 | if (!lAODevent) { | |
324 | AliWarning("ERROR: lAODevent not available \n"); | |
325 | return; | |
326 | } | |
327 | fHistEvent->Fill(0.5); | |
328 | ||
329 | //------------------------------------------------ | |
330 | // Physics Selection | |
331 | //------------------------------------------------ | |
332 | ||
333 | // new method | |
334 | UInt_t maskIsSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected(); | |
335 | Bool_t isSelected = 0; | |
336 | //kMB: default selection, also if fTriggerMask is something not understood... | |
337 | isSelected = (maskIsSelected & AliVEvent::kMB) == AliVEvent::kMB; | |
338 | ||
339 | //Standard Min-Bias Selection | |
340 | if ( ! isSelected ) { | |
341 | PostData(1, fOutput); | |
342 | } | |
343 | fHistEvent->Fill(1.5); | |
344 | ||
345 | //------------------------------------------------ | |
346 | // After Trigger Selection | |
347 | //------------------------------------------------ | |
348 | ||
349 | nV0s = lAODevent->GetNumberOfV0s(); | |
350 | ||
351 | //------------------------------------------------ | |
352 | // Getting: Primary Vertex + MagField Info | |
353 | //------------------------------------------------ | |
354 | ||
355 | const AliAODVertex *lPrimaryBestAODVtx = lAODevent->GetPrimaryVertex(); | |
356 | // get the best primary vertex available for the event | |
357 | // As done in AliCascadeVertexer, we keep the one which is the best one available. | |
358 | // between : Tracking vertex > SPD vertex > TPC vertex > default SPD vertex | |
359 | // This one will be used for next calculations (DCA essentially) | |
360 | lPrimaryBestAODVtx->GetXYZ( lBestPrimaryVtxPos ); | |
361 | ||
362 | Double_t tPrimaryVtxPosition[3]; | |
363 | const AliVVertex *primaryVtx = lAODevent->GetPrimaryVertex(); | |
364 | tPrimaryVtxPosition[0] = primaryVtx->GetX(); | |
365 | tPrimaryVtxPosition[1] = primaryVtx->GetY(); | |
366 | tPrimaryVtxPosition[2] = primaryVtx->GetZ(); | |
367 | ||
368 | //------------------------------------------------ | |
369 | // Only look at events with well-established PV | |
370 | //------------------------------------------------ | |
371 | ||
372 | const AliAODVertex *lPrimaryTrackingAODVtxCheck = lAODevent->GetPrimaryVertex(); | |
373 | const AliAODVertex *lPrimarySPDVtx = lAODevent->GetPrimaryVertexSPD(); | |
374 | if (!lPrimarySPDVtx && !lPrimaryTrackingAODVtxCheck ){ | |
375 | AliWarning("Pb / No SPD prim. vertex nor prim. Tracking vertex ... return !"); | |
376 | PostData(1, fOutput); | |
377 | return; | |
378 | } | |
379 | fHistEvent->Fill(2.5); | |
380 | ||
381 | //------------------------------------------------ | |
382 | // Primary Vertex Z position: SKIP | |
383 | //------------------------------------------------ | |
384 | ||
385 | if(TMath::Abs(lBestPrimaryVtxPos[2]) > 10.0 ) { | |
386 | AliWarning("Pb / | Z position of Best Prim Vtx | > 10.0 cm ... return !"); | |
387 | PostData(1, fOutput); | |
388 | return; | |
389 | } | |
390 | ||
391 | lMagneticField = lAODevent->GetMagneticField( ); | |
392 | fHistEvent->Fill(3.5); | |
393 | ||
394 | //------------------------------------------------ | |
395 | // MAIN LAMBDA LOOP STARTS HERE | |
396 | //------------------------------------------------ | |
397 | ||
398 | //Variable definition | |
399 | Int_t lOnFlyStatus = 0;// nv0sOn = 0, nv0sOff = 0; | |
400 | Double_t lChi2V0 = 0; | |
401 | Double_t lDcaV0Daughters = 0, lDcaV0ToPrimVertex = 0; | |
402 | Double_t lDcaPosToPrimVertex = 0, lDcaNegToPrimVertex = 0; | |
403 | Double_t lV0CosineOfPointingAngle = 0; | |
404 | Double_t lV0Radius = 0, lPt = 0; | |
405 | Double_t lRapK0Short = 0, lRapLambda = 0; | |
406 | Double_t lInvMassK0s = 0, lInvMassLambda = 0, lInvMassAntiLambda = 0; | |
407 | Double_t lAlphaV0 = 0, lPtArmV0 = 0; | |
408 | ||
409 | Double_t fMinV0Pt = 0; | |
410 | Double_t fMaxV0Pt = 100; | |
411 | ||
412 | Int_t nv0s = 0; | |
413 | nv0s = lAODevent->GetNumberOfV0s(); | |
414 | ||
415 | //for (Int_t iV0 = 0; iV0 < nv0s; iV0++) | |
416 | for (Int_t iV0 = 0; iV0 < nv0s; iV0++) //extra-crazy test | |
417 | {// This is the begining of the V0 loop | |
418 | AliAODv0 *v0 = lAODevent->GetV0(iV0); | |
419 | if (!v0) continue; | |
420 | ||
421 | //Only use Offline Candidates for QA | |
422 | lOnFlyStatus = v0->GetOnFlyStatus(); | |
423 | if( lOnFlyStatus == kTRUE ) continue; | |
424 | ||
425 | Double_t tDecayVertexV0[3]; v0->GetXYZ(tDecayVertexV0); | |
426 | Double_t tV0mom[3]; | |
427 | v0->GetPxPyPz( tV0mom ); | |
428 | Double_t lV0TotalMomentum = TMath::Sqrt( | |
429 | tV0mom[0]*tV0mom[0]+tV0mom[1]*tV0mom[1]+tV0mom[2]*tV0mom[2] ); | |
430 | ||
431 | lV0Radius = TMath::Sqrt(tDecayVertexV0[0]*tDecayVertexV0[0]+tDecayVertexV0[1]*tDecayVertexV0[1]); | |
432 | lPt = v0->Pt(); | |
433 | lRapK0Short = v0->RapK0Short(); | |
434 | lRapLambda = v0->RapLambda(); | |
435 | if ((lPt<fMinV0Pt)||(fMaxV0Pt<lPt)) continue; | |
436 | ||
437 | Double_t lMomPos[3]; //v0->GetPPxPyPz(lMomPos[0],lMomPos[1],lMomPos[2]); | |
438 | Double_t lMomNeg[3]; //v0->GetNPxPyPz(lMomNeg[0],lMomNeg[1],lMomNeg[2]); | |
439 | lMomPos[0] = v0->MomPosX(); | |
440 | lMomPos[1] = v0->MomPosY(); | |
441 | lMomPos[2] = v0->MomPosZ(); | |
442 | lMomNeg[0] = v0->MomNegX(); | |
443 | lMomNeg[1] = v0->MomNegY(); | |
444 | lMomNeg[2] = v0->MomNegZ(); | |
445 | ||
446 | AliAODTrack *pTrack=(AliAODTrack *)v0->GetDaughter(0); //0->Positive Daughter | |
447 | AliAODTrack *nTrack=(AliAODTrack *)v0->GetDaughter(1); //1->Negative Daughter | |
448 | if (!pTrack || !nTrack) { | |
449 | Printf("ERROR: Could not retreive one of the daughter track"); | |
450 | continue; | |
451 | } | |
452 | ||
453 | //Daughter Eta for Eta selection, afterwards | |
454 | Float_t lNegEta = nTrack->Eta(); | |
455 | Float_t lPosEta = pTrack->Eta(); | |
456 | ||
457 | // Filter like-sign V0 (next: add counter and distribution) | |
458 | if ( pTrack->Charge() == nTrack->Charge()){ | |
459 | continue; | |
460 | } | |
461 | ||
462 | //Quick test this far! | |
463 | ||
464 | //________________________________________________________________________ | |
465 | // Track quality cuts | |
466 | Float_t lPosTrackCrossedRows = pTrack->GetTPCClusterInfo(2,1); | |
467 | Float_t lNegTrackCrossedRows = nTrack->GetTPCClusterInfo(2,1); | |
468 | Int_t lLeastNbrCrossedRows = (Int_t) lPosTrackCrossedRows; | |
469 | if( lNegTrackCrossedRows < lLeastNbrCrossedRows ) | |
470 | lLeastNbrCrossedRows = (Int_t) lNegTrackCrossedRows; | |
471 | ||
472 | // TPC refit condition (done during reconstruction for Offline but not for On-the-fly) | |
473 | if( !(pTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue; | |
474 | if( !(nTrack->GetStatus() & AliESDtrack::kTPCrefit)) continue; | |
475 | ||
476 | if ( ( ( pTrack->GetTPCClusterInfo(2,1) ) < 70 ) || ( ( nTrack->GetTPCClusterInfo(2,1) ) < 70 ) ) continue; | |
477 | ||
478 | //GetKinkIndex condition | |
479 | //if( pTrack->GetKinkIndex(0)>0 || nTrack->GetKinkIndex(0)>0 ) continue; | |
480 | ||
481 | //Findable clusters > 0 condition | |
482 | if( pTrack->GetTPCNclsF()<=0 || nTrack->GetTPCNclsF()<=0 ) continue; | |
483 | ||
484 | //Compute ratio Crossed Rows / Findable clusters | |
485 | //Note: above test avoids division by zero! | |
486 | Float_t lPosTrackCrossedRowsOverFindable = lPosTrackCrossedRows / ((double)(pTrack->GetTPCNclsF())); | |
487 | Float_t lNegTrackCrossedRowsOverFindable = lNegTrackCrossedRows / ((double)(nTrack->GetTPCNclsF())); | |
488 | ||
489 | Float_t lLeastRatioCrossedRowsOverFindable = lPosTrackCrossedRowsOverFindable; | |
490 | if( lNegTrackCrossedRowsOverFindable < lLeastRatioCrossedRowsOverFindable ) | |
491 | lLeastRatioCrossedRowsOverFindable = lNegTrackCrossedRowsOverFindable; | |
492 | ||
493 | //Lowest Cut Level for Ratio Crossed Rows / Findable = 0.8, set here | |
494 | if ( lLeastRatioCrossedRowsOverFindable < 0.8 ) continue; | |
495 | ||
496 | //End track Quality Cuts | |
497 | //________________________________________________________________________ | |
498 | ||
499 | ||
500 | //ESD code: obsolete | |
501 | //lDcaPosToPrimVertex = TMath::Abs(pTrack->GetD(tPrimaryVtxPosition[0], | |
502 | // tPrimaryVtxPosition[1], | |
503 | // lMagneticField) ); | |
504 | ||
505 | //lDcaNegToPrimVertex = TMath::Abs(nTrack->GetD(tPrimaryVtxPosition[0], | |
506 | // tPrimaryVtxPosition[1], | |
507 | // lMagneticField) ); | |
508 | ||
509 | lDcaPosToPrimVertex = v0->DcaPosToPrimVertex(); | |
510 | lDcaNegToPrimVertex = v0->DcaNegToPrimVertex(); | |
511 | ||
512 | ||
513 | lOnFlyStatus = v0->GetOnFlyStatus(); | |
514 | lChi2V0 = v0->Chi2V0(); | |
515 | lDcaV0Daughters = v0->DcaV0Daughters(); | |
516 | lDcaV0ToPrimVertex = v0->DcaV0ToPrimVertex(); | |
517 | lV0CosineOfPointingAngle = v0->CosPointingAngle(tPrimaryVtxPosition); | |
518 | ||
519 | // Getting invariant mass infos directly from ESD | |
520 | lInvMassK0s = v0->MassK0Short(); | |
521 | lInvMassLambda = v0->MassLambda(); | |
522 | lInvMassAntiLambda = v0->MassAntiLambda(); | |
523 | lAlphaV0 = v0->AlphaV0(); | |
524 | lPtArmV0 = v0->PtArmV0(); | |
525 | ||
526 | //Official means of acquiring N-sigmas | |
527 | Float_t lNSigmasPosProton = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kProton ); | |
528 | Float_t lNSigmasPosPion = fPIDResponse->NumberOfSigmasTPC( pTrack, AliPID::kPion ); | |
529 | Float_t lNSigmasNegProton = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kProton ); | |
530 | Float_t lNSigmasNegPion = fPIDResponse->NumberOfSigmasTPC( nTrack, AliPID::kPion ); | |
531 | ||
532 | //This requires an Invariant Mass Hypothesis afterwards | |
533 | Float_t lDistOverTotMom = TMath::Sqrt( | |
534 | TMath::Power( tDecayVertexV0[0] - lBestPrimaryVtxPos[0] , 2) + | |
535 | TMath::Power( tDecayVertexV0[1] - lBestPrimaryVtxPos[1] , 2) + | |
536 | TMath::Power( tDecayVertexV0[2] - lBestPrimaryVtxPos[2] , 2) | |
537 | ); | |
538 | lDistOverTotMom /= (lV0TotalMomentum+1e-10); //avoid division by zero, to be sure | |
539 | ||
540 | //------------------------------------------------ | |
541 | // Fill Main Output Histograms | |
542 | //------------------------------------------------ | |
543 | ||
544 | //Topological Variable Checks, One-Dimensional | |
545 | fHistTopDCANegToPV -> Fill( lDcaNegToPrimVertex ) ; | |
546 | fHistTopDCAPosToPV -> Fill( lDcaPosToPrimVertex ) ; | |
547 | fHistTopDCAV0Daughters -> Fill( lDcaV0Daughters ) ; | |
548 | fHistTopCosinePA -> Fill( lV0CosineOfPointingAngle ) ; | |
549 | fHistTopV0Radius -> Fill( lV0Radius ) ; | |
550 | ||
551 | ||
552 | if( lDcaNegToPrimVertex > fV0Sels[1] && lDcaPosToPrimVertex > fV0Sels[2] && | |
553 | lDcaV0Daughters < fV0Sels[3] && lV0CosineOfPointingAngle > fV0Sels[4] && | |
554 | lV0Radius > fV0Sels[5] && lV0Radius < fV0Sels [6] ){ | |
555 | ||
556 | //Topological Variables zoomed in at selection level (whatever that may be) | |
557 | //May be slightly redundant if no specific extra configuration was done | |
558 | fHistSelectedTopDCANegToPV -> Fill( lDcaNegToPrimVertex ) ; | |
559 | fHistSelectedTopDCAPosToPV -> Fill( lDcaPosToPrimVertex ) ; | |
560 | fHistSelectedTopDCAV0Daughters -> Fill( lDcaV0Daughters ) ; | |
561 | fHistSelectedTopCosinePA -> Fill( lV0CosineOfPointingAngle ) ; | |
562 | fHistSelectedTopV0Radius -> Fill( lV0Radius ) ; | |
563 | ||
564 | //Specific fV0Sel selection level, but no dEdx applied | |
565 | f2dHistInvMassK0Short -> Fill ( lPt , lInvMassK0s ) ; | |
566 | f2dHistInvMassLambda -> Fill ( lPt , lInvMassLambda ) ; | |
567 | f2dHistInvMassAntiLambda -> Fill ( lPt , lInvMassAntiLambda ) ; | |
568 | ||
569 | //General dE/dx QA | |
570 | f2dHistResponseNegativeAsPion -> Fill( lPt, lNSigmasNegPion ); | |
571 | f2dHistResponseNegativeAsProton -> Fill( lPt, lNSigmasNegProton ); | |
572 | f2dHistResponsePositiveAsPion -> Fill( lPt, lNSigmasPosPion ); | |
573 | f2dHistResponsePositiveAsProton -> Fill( lPt, lNSigmasPosProton ); | |
574 | ||
575 | //Clean Sample From Lambdas | |
576 | //Very strict cuts to ensure dealing with good Lambdas | |
577 | if ( lDcaV0Daughters < 1.0 && lV0CosineOfPointingAngle > 0.999 && TMath::Abs( lInvMassK0s - 0.497614 ) > 0.012 | |
578 | && TMath::Abs( lInvMassAntiLambda - 1.115683) > 0.08 && TMath::Abs( lInvMassLambda - 1.115683) < 0.002 ) { | |
579 | ||
580 | f2dHistdEdxSignalPionFromLambda -> Fill( lPt, nTrack-> GetTPCsignal() ); | |
581 | f2dHistdEdxSignalProtonFromLambda -> Fill( lPt, pTrack-> GetTPCsignal() ); | |
582 | f2dHistResponsePionFromLambda -> Fill( lPt, lNSigmasNegPion ); | |
583 | f2dHistResponseProtonFromLambda -> Fill( lPt, lNSigmasPosProton ); | |
584 | } | |
585 | ||
586 | //Specific fV0Sel selection level, dE/dx applied | |
587 | if ( TMath::Abs(lNSigmasPosPion) < 5 && TMath::Abs(lNSigmasNegPion) < 5 ) f2dHistInvMassWithdEdxK0Short -> Fill ( lPt , lInvMassK0s ) ; | |
588 | if ( TMath::Abs(lNSigmasPosProton) < 5 && TMath::Abs(lNSigmasNegPion) < 5 ) f2dHistInvMassWithdEdxLambda -> Fill ( lPt , lInvMassLambda ) ; | |
589 | if ( TMath::Abs(lNSigmasPosPion) < 5 && TMath::Abs(lNSigmasNegProton) < 5 ) f2dHistInvMassWithdEdxAntiLambda -> Fill ( lPt , lInvMassAntiLambda ) ; | |
590 | } | |
591 | ||
592 | ||
593 | }// This is the end of the V0 loop | |
594 | ||
595 | // Post output data. | |
596 | PostData(1, fOutput); | |
597 | ||
598 | ||
599 | } | |
600 | ||
601 | //________________________________________________________________________ | |
602 | void AliAnalysisTaskQAV0AOD::Terminate(Option_t *) | |
603 | { | |
604 | // Draw result to the screen | |
605 | // Called once at the end of the query | |
606 | // This will draw the V0 candidate multiplicity, whose | |
607 | // number of entries corresponds to the number of triggered events. | |
608 | TList *cRetrievedList = 0x0; | |
609 | cRetrievedList = (TList*)GetOutputData(1); | |
610 | if(!cRetrievedList){ | |
611 | Printf("ERROR - AliAnalysisTaskQAV0 : ouput data container list not available\n"); | |
612 | return; | |
613 | } | |
614 | fHistEvent = dynamic_cast<TH1D*> ( cRetrievedList->FindObject("fHistEvent") ); | |
615 | if (!fHistEvent) { | |
616 | Printf("ERROR - AliAnalysisTaskQAV0 : fHistEvent not available"); | |
617 | return; | |
618 | } | |
619 | TCanvas *canCheck = new TCanvas("AliAnalysisTaskQAV0","V0 Multiplicity",10,10,510,510); | |
620 | canCheck->cd(1)->SetLogy(); | |
621 | fHistEvent->SetMarkerStyle(22); | |
622 | fHistEvent->DrawCopy("E"); | |
623 | } |