]>
Commit | Line | Data |
---|---|---|
59bd5476 | 1 | /************************************************************************** |
2 | * Copyright(c) 1998-2012, ALICE Experiment at CERN, All rights reserved. * | |
3 | * * | |
4 | * Author: Marek Bombara * | |
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 | /* The task selects candidates for K0s, Lambdas and AntiLambdas (trigger particles) | |
17 | * and calculates correlations with charged unidentified particles (associated particles) in phi and eta. | |
18 | * The task works with AOD events only and containes also mixing for acceptance corrections. | |
19 | * Last update edited by Marek Bombara, October2012, Marek.Bombara@cern.ch | |
20 | */ | |
21 | ||
22 | #include <TCanvas.h> | |
23 | #include <TList.h> | |
24 | #include <TH1.h> | |
25 | #include <TH2.h> | |
26 | #include <TH3.h> | |
27 | #include <THnSparse.h> | |
28 | ||
29 | #include "AliLog.h" | |
30 | #include "AliAnalysisManager.h" | |
31 | #include "AliAODTrack.h" | |
32 | #include "AliAODEvent.h" | |
33 | #include "AliAODv0.h" | |
34 | #include "AliAODcascade.h" | |
35 | #include "AliAODVertex.h" | |
36 | ||
37 | #include "AliAODPid.h" | |
38 | #include "AliPIDResponse.h" | |
39 | #include "AliEventPoolManager.h" | |
40 | #include "AliCentrality.h" | |
41 | ||
42 | #include "AliAnalysisTaskV0ChCorrelations.h" | |
43 | #include "AliPhysicsSelectionTask.h" | |
44 | ||
45 | #include <AliMultiInputEventHandler.h> | |
46 | #include <AliMixInputEventHandler.h> | |
47 | ||
48 | ClassImp(AliAnalysisTaskV0ChCorrelations) | |
49 | ClassImp(AliV0ChBasicParticle) | |
50 | ||
51 | //________________________________________________________________________ | |
52 | AliAnalysisTaskV0ChCorrelations::AliAnalysisTaskV0ChCorrelations() // All data members should be initialised here | |
53 | : AliAnalysisTaskSE(), | |
54 | fFillMixed(kTRUE), | |
55 | fMixingTracks(500), | |
56 | fPoolMgr(0x0), | |
57 | ||
58 | fOutput(0), | |
59 | fPIDResponse(0), | |
60 | fHistCentVtx(0), | |
61 | fHistMultiMain(0), | |
62 | fHistMassK0(0), | |
63 | fHistMassLambda(0), | |
64 | fHistMassAntiLambda(0), | |
65 | ||
66 | fHistdPhidEtaSib(0), | |
67 | fHistdPhidEtaMix(0), | |
68 | fHistTrigSib(0), | |
69 | fHistTrigMix(0), | |
70 | ||
71 | fHistTemp(0)// The last in the above list should not have a comma after it | |
72 | { | |
73 | // Dummy constructor ALWAYS needed for I/O. | |
74 | } | |
75 | ||
76 | //________________________________________________________________________ | |
77 | AliAnalysisTaskV0ChCorrelations::AliAnalysisTaskV0ChCorrelations(const char *name) // All data members should be initialised here | |
78 | : AliAnalysisTaskSE(name), | |
79 | fFillMixed(kTRUE), | |
80 | fMixingTracks(500), | |
81 | fPoolMgr(0x0), | |
82 | ||
83 | fOutput(0), | |
84 | fPIDResponse(0), | |
85 | fHistCentVtx(0), | |
86 | fHistMultiMain(0), | |
87 | fHistMassK0(0), | |
88 | fHistMassLambda(0), | |
89 | fHistMassAntiLambda(0), | |
90 | ||
91 | fHistdPhidEtaSib(0), | |
92 | fHistdPhidEtaMix(0), | |
93 | fHistTrigSib(0), | |
94 | fHistTrigMix(0), | |
95 | ||
96 | fHistTemp(0)// The last in the above list should not have a comma after it | |
97 | { | |
98 | // Constructor | |
99 | // Define input and output slots here (never in the dummy constructor) | |
100 | // Input slot #0 works with a TChain - it is connected to the default input container | |
101 | // Output slot #1 writes into a TH1 container | |
102 | DefineOutput(1, TList::Class()); // for output list | |
103 | } | |
104 | ||
105 | //________________________________________________________________________ | |
106 | AliAnalysisTaskV0ChCorrelations::~AliAnalysisTaskV0ChCorrelations() | |
107 | { | |
108 | // Destructor. Clean-up the output list, but not the histograms that are put inside | |
109 | // (the list is owner and will clean-up these histograms). Protect in PROOF case. | |
110 | if (fOutput && !AliAnalysisManager::GetAnalysisManager()->IsProofMode()) { | |
111 | delete fOutput; | |
112 | } | |
113 | } | |
114 | ||
115 | //________________________________________________________________________ | |
116 | void AliAnalysisTaskV0ChCorrelations::UserCreateOutputObjects() | |
117 | { | |
118 | // Create histograms | |
119 | ||
120 | fOutput = new TList(); | |
121 | fOutput->SetOwner(); // IMPORTANT! | |
122 | // defining bins for centrality | |
123 | Int_t nCentralityBins = 9; | |
124 | Double_t centBins[] = {0.,10.,20.,30.,40.,50.,60.,70.,80.,90.}; | |
125 | const Double_t* centralityBins = centBins; | |
126 | // defining bins for Z vertex | |
127 | Int_t nZvtxBins = 20; | |
128 | Double_t vertexBins[] = {-10.,-9.,-8.,-7.,-6.,-5.,-4.,-3.,-2.,-1.,0.,1.,2.,3.,4.,5.,6.,7.,8.,9.,10.}; | |
129 | const Double_t* zvtxBins = vertexBins; | |
130 | // pt bins of associated particles for the analysis | |
131 | Int_t nPtBins = 4; | |
132 | const Double_t PtBins[5] = {3.0,5.0,7.0,9.0,11.0}; | |
133 | // pt bins of trigger particles for the analysis | |
134 | Int_t nPtBinsV0 = 1; | |
135 | const Double_t PtBinsV0[2] = {8.0,15.0}; | |
136 | // V0 candidate: 1 - K0sig, 2 - Lamsig, 3 - Alamsig, 4 - K0bg, 5 - Lambg, 6 - Alambg | |
137 | Int_t nTrigC = 6; | |
138 | const Double_t TrigC[7] = {0.5,1.5,2.5,3.5,4.5,5.5,6.5}; | |
139 | ||
140 | // Create general histograms | |
141 | fHistCentVtx = new TH2F("fHistCentVtx", "Centrality vs. Z vertex", nCentralityBins, centralityBins, nZvtxBins, zvtxBins); | |
142 | fHistMultiMain = new TH1F("fHistMultiMain", "Multiplicity of main events", 2000, 0, 2000); | |
143 | ||
144 | // defining bins for mass distributions | |
145 | Int_t nBins = 200; | |
146 | Double_t mk0Min = 0.40; | |
147 | Double_t mk0Max = 0.60; | |
148 | Double_t mlaMin = 1.07; | |
149 | Double_t mlaMax = 1.17; | |
150 | Double_t malMin = 1.07; | |
151 | Double_t malMax = 1.17; | |
152 | ||
153 | const Int_t spBins[3] = {nBins, nPtBinsV0, nCentralityBins}; | |
154 | const Double_t spMinK0[3] = {mk0Min, PtBinsV0[0], centralityBins[0]}; | |
155 | const Double_t spMaxK0[3] = {mk0Max, PtBinsV0[1], centralityBins[9]}; | |
156 | const Double_t spMinLa[3] = {mlaMin, PtBinsV0[0], centralityBins[0]}; | |
157 | const Double_t spMaxLa[3] = {mlaMax, PtBinsV0[1], centralityBins[9]}; | |
158 | const Double_t spMinAl[3] = {malMin, PtBinsV0[0], centralityBins[0]}; | |
159 | const Double_t spMaxAl[3] = {malMax, PtBinsV0[1], centralityBins[9]}; | |
160 | // Create mass histograms | |
161 | fHistMassK0 = new THnSparseF("fHistMassK0","V0 mass for K0 hypothesis", 3, spBins, spMinK0, spMaxK0); | |
162 | fHistMassLambda = new THnSparseF("fHistMassLambda","V0 mass for Lambda hypothesis", 3, spBins, spMinLa, spMaxLa); | |
163 | fHistMassAntiLambda = new THnSparseF("fHistMassAntiLambda","V0 mass for AntiLambda hypothesis", 3, spBins, spMinAl, spMaxAl); | |
164 | ||
165 | // defining bins for dPhi distributions | |
166 | const Int_t nbPhiBins = 144; | |
167 | const Double_t kPi = TMath::Pi(); | |
168 | Double_t PhiMin = -kPi/2.; | |
169 | Double_t PhiMax = -kPi/2. + 2*kPi; | |
170 | Double_t PhiBins[nbPhiBins+1] = {0.}; | |
171 | PhiBins[0] = PhiMin; | |
172 | for (Int_t i=0; i<nbPhiBins; i++) { PhiBins[i+1] = PhiBins[i] + (PhiMax - PhiMin)/nbPhiBins; } | |
173 | ||
174 | // defining bins for dEta distributions | |
175 | const Int_t nbEtaBins = 72; | |
176 | Double_t EtaMin = -1.8; | |
177 | Double_t EtaMax = 1.8; | |
178 | Double_t EtaBins[nbEtaBins+1] = {0.}; | |
179 | EtaBins[0] = EtaMin; | |
180 | for (Int_t i=0; i<nbEtaBins; i++) { EtaBins[i+1] = EtaBins[i] + (EtaMax - EtaMin)/nbEtaBins; } | |
181 | ||
182 | const Int_t corBins[7] = {nbPhiBins, nbEtaBins, nPtBinsV0, nPtBins, nCentralityBins, nZvtxBins, nTrigC}; | |
183 | const Double_t corMin[7] = {PhiBins[0], EtaBins[0], PtBinsV0[0], PtBins[0], centralityBins[0], zvtxBins[0], TrigC[0]}; | |
184 | const Double_t corMax[7] = {PhiBins[144], EtaBins[72], PtBinsV0[1], PtBins[4], centralityBins[9], zvtxBins[20], TrigC[6]}; | |
185 | // Create correlation histograms | |
186 | fHistdPhidEtaSib = new THnSparseF("fHistdPhidEtaSib","dPhi vs. dEta siblings", 7, corBins, corMin, corMax); | |
187 | fHistdPhidEtaMix = new THnSparseF("fHistdPhidEtaMix","dPhi vs. dEta mixed", 7, corBins, corMin, corMax); | |
188 | ||
189 | const Int_t corTrigBins[4] = {nPtBinsV0, nCentralityBins, nZvtxBins, nTrigC}; | |
190 | const Double_t corTrigMin[4] = {PtBinsV0[0], centBins[0], vertexBins[0], TrigC[0]}; | |
191 | const Double_t corTrigMax[4] = {PtBinsV0[1], centBins[9], vertexBins[20], TrigC[6]}; | |
192 | // Create histograms for trigger particles | |
193 | fHistTrigSib = new THnSparseF("fHistTrigSib","pt trigger sib", 4, corTrigBins, corTrigMin, corTrigMax); | |
194 | fHistTrigMix = new THnSparseF("fHistTrigMix","pt trigger mix", 4, corTrigBins, corTrigMin, corTrigMax); | |
195 | ||
196 | // Histogram for debugging | |
197 | fHistTemp = new TH1D("fHistTemp", "Temporary", 100, 0., 100.); | |
198 | ||
199 | fOutput->Add(fHistCentVtx); | |
200 | ||
201 | fOutput->Add(fHistMultiMain); | |
202 | fOutput->Add(fHistMassK0); | |
203 | fOutput->Add(fHistMassLambda); | |
204 | fOutput->Add(fHistMassAntiLambda); | |
205 | ||
206 | fOutput->Add(fHistdPhidEtaSib); | |
207 | fOutput->Add(fHistdPhidEtaMix); | |
208 | fOutput->Add(fHistTrigSib); | |
209 | fOutput->Add(fHistTrigMix); | |
210 | ||
211 | fOutput->Add(fHistTemp); | |
212 | ||
213 | PostData(1, fOutput); // Post data for ALL output slots >0 here, to get at least an empty histogram | |
214 | ||
215 | // Settings for event mixing ------------------------------------- | |
216 | Int_t trackDepth = fMixingTracks; | |
217 | Int_t poolSize = 200; // Maximum number of events, ignored in the present implemented of AliEventPoolManager | |
218 | ||
219 | fPoolMgr = new AliEventPoolManager(poolSize, trackDepth, nCentralityBins, centBins, nZvtxBins, vertexBins); | |
220 | //---------------------------------------------- | |
221 | } | |
222 | ||
223 | //________________________________________________________________________ | |
224 | void AliAnalysisTaskV0ChCorrelations::Terminate(Option_t *) | |
225 | { | |
226 | // Draw result to screen, or perform fitting, normalizations | |
227 | // Called once at the end of the query | |
228 | ||
229 | fOutput = dynamic_cast<TList*>(GetOutputData(1)); | |
230 | if (!fOutput) { AliError("Could not retrieve TList fOutput"); return; } | |
231 | ||
232 | // NEW HISTO should be retrieved from the TList container in the above way, | |
233 | // so it is available to draw on a canvas such as below | |
234 | } | |
235 | ||
236 | //_________________________________________________________________________ | |
237 | void AliAnalysisTaskV0ChCorrelations::UserExec(Option_t *) | |
238 | { | |
239 | AliAnalysisManager *mgr = AliAnalysisManager::GetAnalysisManager(); | |
240 | AliInputEventHandler *inEvMain = dynamic_cast<AliInputEventHandler *>(mgr->GetInputEventHandler()); | |
241 | ||
242 | // physics selection | |
243 | UInt_t maskIsSelected = inEvMain->IsEventSelected(); | |
244 | // 2010 data trigger selection | |
245 | //Bool_t isSelected = (maskIsSelected & AliVEvent::kMB); | |
246 | // 2011 data trigger selection | |
247 | Bool_t isSelected = ((maskIsSelected & AliVEvent::kMB) || (maskIsSelected & AliVEvent::kCentral) || (maskIsSelected & AliVEvent::kSemiCentral)); | |
248 | if (!isSelected) return; | |
249 | ||
250 | AliAODEvent* aod = dynamic_cast<AliAODEvent*>(inEvMain->GetEvent()); | |
251 | fPIDResponse = inEvMain->GetPIDResponse(); | |
252 | ||
253 | // pt intervals for trigger particles | |
254 | const Double_t kPi = TMath::Pi(); | |
255 | Double_t PtTrigMin = 8.0; | |
256 | Double_t PtTrigMax = 15.0; | |
257 | // pt interval for associated particles | |
258 | Double_t PtAssocMin = 3.0; | |
259 | fHistMultiMain->Fill(aod->GetNumberOfTracks()); | |
260 | ||
261 | // Vertex cut | |
262 | Double_t cutPrimVertex = 10.0; | |
263 | AliAODVertex *myPrimVertex = aod->GetPrimaryVertex(); | |
264 | if (!myPrimVertex) return; | |
265 | if ( ( TMath::Abs(myPrimVertex->GetZ()) ) >= cutPrimVertex) return ; | |
266 | Double_t lPVx = myPrimVertex->GetX(); | |
267 | Double_t lPVy = myPrimVertex->GetY(); | |
268 | Double_t lPVz = myPrimVertex->GetZ(); | |
269 | ||
270 | if (TMath::Abs(lPVx)<10e-5 && TMath::Abs(lPVy)<10e-5 && TMath::Abs(lPVz)<10e-5) return; | |
271 | ||
272 | // Centrality definition | |
273 | Double_t lCent = 0.0; | |
274 | AliCentrality *centralityObj = 0; | |
275 | centralityObj = aod->GetHeader()->GetCentralityP(); | |
276 | lCent = centralityObj->GetCentralityPercentile("CL1"); | |
277 | ||
278 | if ((lCent < 0.)||(lCent > 90.)) return; | |
279 | fHistCentVtx->Fill(lCent,lPVz); | |
280 | ||
281 | // Track selection loop | |
282 | //-------------------------------- | |
283 | Int_t nTracks = aod->GetNumberOfTracks(); | |
284 | // new tracks array | |
285 | TObjArray * selectedTracks = new TObjArray; | |
286 | selectedTracks->SetOwner(kTRUE); | |
287 | for (Int_t i = 0; i < nTracks; i++) | |
288 | { | |
289 | AliAODTrack* tr = aod->GetTrack(i); | |
290 | if (!(IsMyGoodPrimaryTrack(tr))) continue; | |
291 | if ((tr->Pt())<PtAssocMin) continue; | |
292 | selectedTracks->Add(tr); | |
293 | } | |
294 | //--------------------------------- | |
295 | ||
296 | // V0 loop for reconstructed event | |
297 | TObjArray * selectedV0s = new TObjArray; | |
298 | selectedV0s->SetOwner(kTRUE); | |
299 | ||
300 | Int_t nV0s = aod->GetNumberOfV0s(); | |
301 | ||
302 | for (Int_t i = 0; i < nV0s; i++) | |
303 | { // start of V0 slection loop | |
304 | AliAODv0* aodV0 = dynamic_cast<AliAODv0 *>(aod->GetV0(i)); | |
305 | if (!aodV0) { | |
306 | AliError(Form("ERROR: Could not retrieve aodv0 %d", i)); | |
307 | continue; | |
308 | } | |
309 | if (((aodV0->Pt())<PtTrigMin)||((aodV0->Pt())>PtTrigMax)) continue; | |
310 | // get daughters | |
311 | const AliAODTrack *myTrackPos = NULL; | |
312 | const AliAODTrack *myTrackNeg = NULL; | |
313 | AliAODTrack *myTrackNegTest=dynamic_cast<AliAODTrack *>(aodV0->GetDaughter(1)); | |
314 | AliAODTrack *myTrackPosTest=dynamic_cast<AliAODTrack *>(aodV0->GetDaughter(0)); | |
315 | ||
316 | if (!myTrackPosTest || !myTrackNegTest) { | |
317 | Printf("strange analysis::UserExec:: Error:Could not retreive one of the daughter track\n"); | |
318 | continue; | |
319 | } | |
320 | ||
321 | if( myTrackPosTest->Charge() ==1){ | |
322 | myTrackPos = myTrackPosTest; | |
323 | myTrackNeg = myTrackNegTest; | |
324 | } | |
325 | ||
326 | if( myTrackPosTest->Charge() ==-1){ | |
327 | myTrackPos = myTrackNegTest; | |
328 | myTrackNeg = myTrackPosTest; | |
329 | } | |
330 | ||
331 | if (!IsMyGoodV0(aod,aodV0,myTrackPos,myTrackNeg)) continue; | |
332 | ||
333 | // effective mass calculations for each hypothesis | |
334 | Double_t lInvMassK0 = aodV0->MassK0Short(); | |
335 | Double_t lInvMassAntiLambda = aodV0->MassAntiLambda(); | |
336 | Double_t lInvMassLambda = aodV0->MassLambda(); | |
337 | ||
338 | // calculations for c*tau cut-------------------------------------- | |
339 | Double_t lDVx = aodV0->GetSecVtxX(); | |
340 | Double_t lDVy = aodV0->GetSecVtxY(); | |
341 | Double_t lDVz = aodV0->GetSecVtxZ(); | |
342 | const Double_t kLambdaMass = 1.115683; | |
343 | const Double_t kK0Mass = 0.497648; | |
344 | Double_t cutcTauLam = 3*7.89; | |
345 | Double_t cutcTauK0 = 3*2.68; | |
346 | Double_t lV0DecayLength = TMath::Sqrt(TMath::Power(lDVx - lPVx,2) + TMath::Power(lDVy- lPVy,2) + TMath::Power(lDVz - lPVz,2 )); | |
347 | Double_t lPV0 = TMath::Sqrt((aodV0->Pt())*(aodV0->Pt())+(aodV0->Pz())*(aodV0->Pz())); | |
348 | Double_t lcTauLam = (lV0DecayLength*kLambdaMass)/lPV0; | |
349 | Double_t lcTauK0 = (lV0DecayLength*kK0Mass)/lPV0; | |
350 | // sc - standard cuts | |
351 | Bool_t cutK0sc = (lcTauK0<cutcTauK0); | |
352 | Bool_t cutLambdasc = (lcTauLam<cutcTauLam); | |
353 | Bool_t cutAntiLambdasc = (lcTauLam<cutcTauLam); | |
354 | ||
355 | //------------------------------------------------ | |
356 | ||
357 | // preparations for daughter's PID cut------------------------------ | |
358 | Float_t nSigmaPosPion = 0.; | |
359 | Float_t nSigmaNegPion = 0.; | |
360 | Float_t nSigmaPosProton = 0.; | |
361 | Float_t nSigmaNegProton = 0.; | |
362 | ||
363 | const AliAODPid *pPid = myTrackPos->GetDetPid(); | |
364 | const AliAODPid *nPid = myTrackNeg->GetDetPid(); | |
365 | ||
366 | if (pPid) | |
367 | { | |
368 | Double_t pdMom = pPid->GetTPCmomentum(); | |
369 | if (pdMom<1.) | |
370 | { | |
371 | nSigmaPosPion = fPIDResponse->NumberOfSigmasTPC(myTrackPos, AliPID::kPion); | |
372 | nSigmaPosProton = fPIDResponse->NumberOfSigmasTPC(myTrackPos, AliPID::kProton); | |
373 | } | |
374 | } | |
375 | ||
376 | if (nPid) | |
377 | { | |
378 | Double_t ndMom = nPid->GetTPCmomentum(); | |
379 | if (ndMom<1.) | |
380 | { | |
381 | nSigmaNegPion = fPIDResponse->NumberOfSigmasTPC(myTrackPos, AliPID::kPion); | |
382 | nSigmaNegProton = fPIDResponse->NumberOfSigmasTPC(myTrackPos, AliPID::kProton); | |
383 | } | |
384 | } | |
385 | Bool_t bpPion = TMath::Abs(nSigmaPosPion) <= 3.; | |
386 | Bool_t bpProton = TMath::Abs(nSigmaPosProton) <= 3.; | |
387 | Bool_t bnPion = TMath::Abs(nSigmaNegPion) <= 3.; | |
388 | Bool_t bnProton = TMath::Abs(nSigmaNegProton) <= 3.; | |
389 | ||
390 | Bool_t cutK0Pid = (bpPion && bnPion); | |
391 | Bool_t cutLambdaPid = (bpProton && bnPion); | |
392 | Bool_t cutAntiLambdaPid = (bpPion && bnProton); | |
393 | //-------------------------------------------------- | |
394 | // mass cuts | |
395 | Bool_t cutMassLambda = ((lInvMassLambda>1.10) && (lInvMassLambda<1.13)); | |
396 | Bool_t cutMassAntiLambda = ((lInvMassAntiLambda>1.10) && (lInvMassAntiLambda<1.13)); | |
397 | Bool_t cutMassK0 = (lInvMassK0>0.47) && (lInvMassK0<0.53); | |
398 | ||
399 | cutK0sc = cutK0Pid && (!cutMassLambda) && (!cutMassAntiLambda); | |
400 | cutLambdasc = cutLambdaPid && (!cutMassK0); | |
401 | cutAntiLambdasc = cutAntiLambdaPid && (!cutMassK0); | |
402 | // special cut related to AP diagram for K0s | |
403 | Bool_t k0APcut = (aodV0->Pt()>(TMath::Abs(0.2*aodV0->AlphaV0()))); | |
404 | cutK0sc = cutK0sc && k0APcut; | |
405 | // fill the mass histograms | |
406 | Double_t spK0[3] = {lInvMassK0, aodV0->Pt(), lCent}; | |
407 | Double_t spLa[3] = {lInvMassLambda, aodV0->Pt(), lCent}; | |
408 | Double_t spAl[3] = {lInvMassAntiLambda, aodV0->Pt(), lCent}; | |
409 | if (cutK0sc) fHistMassK0->Fill(spK0); | |
410 | if (cutLambdasc) fHistMassLambda->Fill(spLa); | |
411 | if (cutAntiLambdasc) fHistMassAntiLambda->Fill(spAl); | |
412 | // select final V0s for correlation, selected background to study its contribution to correlation function | |
413 | // the values for signal might change in the future ... | |
414 | Bool_t K0Signal = (lInvMassK0>0.484)&&(lInvMassK0<0.512); | |
415 | Bool_t K0Bckg = ((lInvMassK0>0.450)&&(lInvMassK0<0.462)) || ((lInvMassK0>0.532)&&(lInvMassK0<0.544)); | |
416 | ||
417 | Bool_t LamSignal = (lInvMassLambda>1.110)&&(lInvMassLambda<1.122); | |
418 | Bool_t LamBckg = ((lInvMassLambda>1.097)&&(lInvMassLambda<1.104)) || ((lInvMassLambda>1.127)&&(lInvMassLambda<1.134)); | |
419 | ||
420 | Bool_t ALamSignal = (lInvMassAntiLambda>1.110)&&(lInvMassAntiLambda<1.122); | |
421 | Bool_t ALamBckg = ((lInvMassAntiLambda>1.097)&&(lInvMassAntiLambda<1.104)) || ((lInvMassAntiLambda>1.127)&&(lInvMassAntiLambda<1.134)); | |
422 | ||
423 | // Fill selected V0 particle array | |
424 | if ((cutK0sc)&&(K0Signal)) selectedV0s->Add(new AliV0ChBasicParticle(aodV0->Eta(), aodV0->Phi(), aodV0->Pt(), 1)); | |
425 | if ((cutK0sc)&&(K0Bckg)) selectedV0s->Add(new AliV0ChBasicParticle(aodV0->Eta(), aodV0->Phi(), aodV0->Pt(), 4)); | |
426 | if ((cutLambdasc)&&(LamSignal)) selectedV0s->Add(new AliV0ChBasicParticle(aodV0->Eta(), aodV0->Phi(), aodV0->Pt(), 2)); | |
427 | if ((cutLambdasc)&&(LamBckg)) selectedV0s->Add(new AliV0ChBasicParticle(aodV0->Eta(), aodV0->Phi(), aodV0->Pt(), 5)); | |
428 | if ((cutAntiLambdasc)&&(ALamSignal)) selectedV0s->Add(new AliV0ChBasicParticle(aodV0->Eta(), aodV0->Phi(), aodV0->Pt(), 3)); | |
429 | if ((cutAntiLambdasc)&&(ALamBckg)) selectedV0s->Add(new AliV0ChBasicParticle(aodV0->Eta(), aodV0->Phi(), aodV0->Pt(), 6)); | |
430 | ||
431 | Int_t nSelectedTracks = selectedTracks->GetEntries(); | |
432 | // Correlation part | |
433 | //=================================== | |
434 | for (Int_t j = 0; j < nSelectedTracks; j++) | |
435 | { | |
436 | //const Double_t kPi = TMath::Pi(); | |
437 | AliAODTrack* atr = (AliAODTrack*) selectedTracks->At(j); | |
438 | //if (!(IsMyGoodPrimaryTrack(atr))) continue; | |
439 | if ((atr->Pt())>=(aodV0->Pt())) continue; | |
440 | Double_t dEta = atr->Eta() - aodV0->Eta(); | |
441 | Double_t dPhi = atr->Phi() - aodV0->Phi(); | |
442 | if (dPhi > (1.5*kPi)) dPhi -= 2.0*kPi; | |
443 | if (dPhi < (-0.5*kPi)) dPhi += 2.0*kPi; | |
444 | // removing autocorrelations | |
445 | //---------------------------------- | |
446 | Int_t negID = myTrackNeg->GetID(); | |
447 | Int_t posID = myTrackPos->GetID(); | |
448 | Int_t atrID = atr->GetID(); | |
449 | ||
450 | if ((TMath::Abs(negID)+1)==(TMath::Abs(atrID))) continue; | |
451 | if ((TMath::Abs(posID)+1)==(TMath::Abs(atrID))) continue; | |
452 | //---------------------------------- | |
453 | ||
454 | fHistTrigSib->Sumw2(); | |
455 | fHistdPhidEtaSib->Sumw2(); | |
456 | ||
457 | if ((cutK0sc)&&(K0Signal)) | |
458 | { | |
459 | if (j==0) | |
460 | { | |
461 | Double_t spK0TrigSig[7] = {aodV0->Pt(), lCent, lPVz, 1.}; | |
462 | fHistTrigSib->Fill(spK0TrigSig); | |
463 | } | |
464 | Double_t spK0Sig[7] = {dPhi, dEta, aodV0->Pt(), atr->Pt(), lCent, lPVz, 1.}; | |
465 | fHistdPhidEtaSib->Fill(spK0Sig); | |
466 | } | |
467 | if ((cutK0sc)&&(K0Bckg)) | |
468 | { | |
469 | if (j==0) | |
470 | { | |
471 | Double_t spK0TrigBkg[7] = {aodV0->Pt(), lCent, lPVz, 4.}; | |
472 | fHistTrigSib->Fill(spK0TrigBkg); | |
473 | } | |
474 | Double_t spK0Bkg[7] = {dPhi, dEta, aodV0->Pt(), atr->Pt(), lCent, lPVz, 4.}; | |
475 | fHistdPhidEtaSib->Fill(spK0Bkg); | |
476 | } | |
477 | if ((cutLambdasc)&&(LamSignal)) | |
478 | { | |
479 | if (j==0) | |
480 | { | |
481 | Double_t spLaTrigSig[7] = {aodV0->Pt(), lCent, lPVz, 2.}; | |
482 | fHistTrigSib->Fill(spLaTrigSig); | |
483 | } | |
484 | Double_t spLaSig[7] = {dPhi, dEta, aodV0->Pt(), atr->Pt(), lCent, lPVz, 2.}; | |
485 | fHistdPhidEtaSib->Fill(spLaSig); | |
486 | } | |
487 | if ((cutLambdasc)&&(LamBckg)) | |
488 | { | |
489 | if (j==0) | |
490 | { | |
491 | Double_t spLaTrigBkg[7] = {aodV0->Pt(), lCent, lPVz, 5.}; | |
492 | fHistTrigSib->Fill(spLaTrigBkg); | |
493 | } | |
494 | Double_t spLaBkg[7] = {dPhi, dEta, aodV0->Pt(), atr->Pt(), lCent, lPVz, 5.}; | |
495 | fHistdPhidEtaSib->Fill(spLaBkg); | |
496 | } | |
497 | if ((cutAntiLambdasc)&&(ALamSignal)) | |
498 | { | |
499 | if (j==0) | |
500 | { | |
501 | Double_t spAlTrigSig[7] = {aodV0->Pt(), lCent, lPVz, 3.}; | |
502 | fHistTrigSib->Fill(spAlTrigSig); | |
503 | } | |
504 | Double_t spAlSig[7] = {dPhi, dEta, aodV0->Pt(), atr->Pt(), lCent, lPVz, 3.}; | |
505 | fHistdPhidEtaSib->Fill(spAlSig); | |
506 | } | |
507 | if ((cutAntiLambdasc)&&(ALamBckg)) | |
508 | { | |
509 | if (j==0) | |
510 | { | |
511 | Double_t spAlTrigBkg[7] = {aodV0->Pt(), lCent, lPVz, 6.}; | |
512 | fHistTrigSib->Fill(spAlTrigBkg); | |
513 | } | |
514 | Double_t spAlBkg[7] = {dPhi, dEta, aodV0->Pt(), atr->Pt(), lCent, lPVz, 6.}; | |
515 | fHistdPhidEtaSib->Fill(spAlBkg); | |
516 | } | |
517 | ||
518 | } // end of correlation loop | |
519 | //=================================== | |
520 | ||
521 | ||
522 | } // end of V0 selection loop | |
523 | ||
524 | // Mixing ============================================== | |
525 | ||
526 | fHistTrigMix->Sumw2(); | |
527 | fHistdPhidEtaMix->Sumw2(); | |
528 | AliEventPool* pool = fPoolMgr->GetEventPool(lCent, lPVz); | |
529 | if (!pool) AliFatal(Form("No pool found for centrality = %f, zVtx = %f", lCent, lPVz)); | |
530 | //pool->SetDebug(1); | |
531 | if (pool->IsReady() || pool->NTracksInPool() > fMixingTracks / 10 || pool->GetCurrentNEvents() >= 5) | |
532 | { | |
533 | Int_t nMix = pool->GetCurrentNEvents(); | |
534 | for (Int_t jMix=0; jMix<nMix; jMix++) | |
535 | {// loop through mixing events | |
536 | TObjArray* bgTracks = pool->GetEvent(jMix); | |
537 | for (Int_t i=0; i<selectedV0s->GetEntriesFast(); i++) | |
538 | {// loop through selected V0 particles | |
539 | AliV0ChBasicParticle* trig = (AliV0ChBasicParticle*) selectedV0s->At(i); | |
540 | Double_t trigPhi = trig->Phi(); | |
541 | Double_t trigEta = trig->Eta(); | |
542 | Double_t trigPt = trig->Pt(); | |
543 | Short_t trigC = trig->WhichCandidate(); | |
544 | for (Int_t j=0; j<bgTracks->GetEntriesFast(); j++) | |
545 | { // mixing tracks loop | |
546 | AliVParticle* assoc = (AliVParticle*) bgTracks->At(j); | |
547 | // be careful tracks may have bigger pt than v0s. | |
548 | if ( ( (assoc->Pt())>=trigPt )||( (assoc->Pt())<PtAssocMin ) ) continue; | |
549 | Double_t dEtaMix = assoc->Eta() - trigEta; | |
550 | Double_t dPhiMix = assoc->Phi() - trigPhi; | |
551 | if (dPhiMix > (1.5*kPi)) dPhiMix -= 2.0*kPi; | |
552 | if (dPhiMix < (-0.5*kPi)) dPhiMix += 2.0*kPi; | |
553 | ||
554 | if (j==0) | |
555 | { | |
556 | Double_t spTrigMix[7] = {trigPt, lCent, lPVz, trigC}; | |
557 | fHistTrigMix->Fill(spTrigMix); | |
558 | } | |
559 | Double_t spMix[7] = {dPhiMix, dEtaMix, trigPt, assoc->Pt(), lCent, lPVz, trigC}; | |
560 | fHistdPhidEtaMix->Fill(spMix); | |
561 | ||
562 | } // end of mixing track loop | |
563 | }// end of loop through selected V0 particles | |
564 | }// end of loop of mixing events | |
565 | } | |
566 | ||
567 | TObjArray* tracksClone = (TObjArray*) selectedTracks->Clone(); | |
568 | tracksClone->SetOwner(kTRUE); | |
569 | pool->UpdatePool(tracksClone); | |
570 | //delete selectedtracks; | |
571 | ||
572 | PostData(1, fOutput); | |
573 | ||
574 | } | |
575 | //___________________________________________ | |
576 | Bool_t AliAnalysisTaskV0ChCorrelations::IsMyGoodPrimaryTrack(const AliAODTrack *t) | |
577 | { | |
578 | // Pseudorapidity cut | |
579 | if (TMath::Abs(t->Eta())>0.8) return kFALSE; | |
580 | // Should correspond to set of cuts suitable for correlation analysis | |
581 | if (!t->TestFilterBit(1<<7)) return kFALSE; | |
582 | ||
583 | return kTRUE; | |
584 | } | |
585 | //_____________________________________________ | |
586 | Bool_t AliAnalysisTaskV0ChCorrelations::IsMyGoodDaughterTrack(const AliAODTrack *t) | |
587 | { | |
588 | // TPC refit | |
589 | if (!t->IsOn(AliAODTrack::kTPCrefit)) return kFALSE; | |
590 | // Minimum number of clusters | |
591 | Float_t nCrossedRowsTPC = t->GetTPCClusterInfo(2,1); | |
592 | if (nCrossedRowsTPC < 70) return kFALSE; | |
593 | Int_t findable=t->GetTPCNclsF(); | |
594 | if (findable <= 0) return kFALSE; | |
595 | if (nCrossedRowsTPC/findable < 0.8) return kFALSE; | |
596 | ||
597 | return kTRUE; | |
598 | } | |
599 | //______________________________________________ | |
600 | Bool_t AliAnalysisTaskV0ChCorrelations::IsMyGoodV0(const AliAODEvent* aod, const AliAODv0* aodV0, const AliAODTrack* myTrackPos, const AliAODTrack* myTrackNeg) | |
601 | { | |
602 | if (!aodV0) { | |
603 | AliError(Form("ERROR: Could not retrieve aodV0")); | |
604 | return kFALSE; | |
605 | } | |
606 | ||
607 | // Rapidity cut | |
608 | Double_t lCutRap = 0.75; | |
609 | Double_t lRapK0s = aodV0->Y(310); | |
610 | Double_t lRapLambda = aodV0->Y(3122); | |
611 | Double_t lRapAntiLambda = aodV0->Y(-3122); | |
612 | ||
613 | if (TMath::Abs(lRapK0s)>=lCutRap) return kFALSE; | |
614 | if (TMath::Abs(lRapLambda)>=lCutRap) return kFALSE; | |
615 | if (TMath::Abs(lRapAntiLambda)>=lCutRap) return kFALSE; | |
616 | ||
617 | // Offline reconstructed V0 only | |
618 | if (aodV0->GetOnFlyStatus()) return kFALSE; | |
619 | ||
620 | // DCA of daughter track to Primary Vertex | |
621 | Float_t xyn=aodV0->DcaNegToPrimVertex(); | |
622 | if (TMath::Abs(xyn)<0.1) return kFALSE; | |
623 | Float_t xyp=aodV0->DcaPosToPrimVertex(); | |
624 | if (TMath::Abs(xyp)<0.1) return kFALSE; | |
625 | ||
626 | // DCA of daughter tracks | |
627 | Double_t dca=aodV0->DcaV0Daughters(); | |
628 | if (dca>1.0) return kFALSE; | |
629 | ||
630 | // Cosinus of pointing angle | |
631 | Double_t cpa=aodV0->CosPointingAngle(aod->GetPrimaryVertex()); | |
632 | if (cpa<0.998) return kFALSE; | |
633 | ||
634 | // Fiducial volume cut | |
635 | Double_t xyz[3]; aodV0->GetSecondaryVtx(xyz); | |
636 | Double_t r2=xyz[0]*xyz[0] + xyz[1]*xyz[1]; | |
637 | if (r2<0.9*0.9) return kFALSE; | |
638 | if (r2>100*100) return kFALSE; | |
639 | ||
640 | // c*tau cut - in main V0 loop - depends on particle hypothesis | |
641 | ||
642 | // Get daughters and check them | |
643 | AliAODTrack *myTrackNegTest=dynamic_cast<AliAODTrack *>(aodV0->GetDaughter(1)); | |
644 | AliAODTrack *myTrackPosTest=dynamic_cast<AliAODTrack *>(aodV0->GetDaughter(0)); | |
645 | ||
646 | if (!myTrackPosTest || !myTrackNegTest) { | |
647 | Printf("strange analysis::UserExec:: Error:Could not retreive one of the daughter track\n"); | |
648 | return kFALSE; | |
649 | } | |
650 | ||
651 | if( myTrackPosTest->Charge() ==1){ | |
652 | myTrackPos = myTrackPosTest; | |
653 | myTrackNeg = myTrackNegTest; | |
654 | } | |
655 | ||
656 | if( myTrackPosTest->Charge() ==-1){ | |
657 | myTrackPos = myTrackNegTest; | |
658 | myTrackNeg = myTrackPosTest; | |
659 | } | |
660 | ||
661 | // Track cuts for daugher tracks | |
662 | if ( !(IsMyGoodDaughterTrack(myTrackPos)) || !(IsMyGoodDaughterTrack(myTrackNeg)) ) return kFALSE; | |
663 | ||
664 | // Unlike signs of daughters | |
665 | if (myTrackNegTest->Charge() == myTrackPosTest->Charge()) return kFALSE; | |
666 | ||
667 | // Minimum pt of daughters | |
668 | Double_t lMomPos[3] = {999,999,999}; | |
669 | Double_t lMomNeg[3] = {999,999,999}; | |
670 | ||
671 | lMomPos[0] = aodV0->MomPosX(); | |
672 | lMomPos[1] = aodV0->MomPosY(); | |
673 | lMomPos[2] = aodV0->MomPosZ(); | |
674 | ||
675 | lMomNeg[0] = aodV0->MomNegX(); | |
676 | lMomNeg[1] = aodV0->MomNegY(); | |
677 | lMomNeg[2] = aodV0->MomNegZ(); | |
678 | ||
679 | Double_t lPtPos = TMath::Sqrt(lMomPos[0]*lMomPos[0] + lMomPos[1]*lMomPos[1]); | |
680 | Double_t lPtNeg = TMath::Sqrt(lMomNeg[0]*lMomNeg[0] + lMomNeg[1]*lMomNeg[1]); | |
681 | ||
682 | Double_t cutMinPtDaughter = 0.160; | |
683 | if (lPtPos<cutMinPtDaughter || lPtNeg<cutMinPtDaughter) return kFALSE; | |
684 | ||
685 | // Daughter PID cut - in main V0 loop - depends on particle hypothesis | |
686 | ||
687 | return kTRUE; | |
688 | } | |
689 |