]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG3/vertexingHF/AliAnalysisTaskSEBkgLikeSignJPSI.cxx
New class to analyze unlike-sign and like-sign background for Jpsi from B (Carmelo)
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliAnalysisTaskSEBkgLikeSignJPSI.cxx
CommitLineData
adef5d4e 1/**************************************************************************
2 * Copyright(c) 1998-2009, 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// AliAnalysisTaskSE for reading both reconstructed JPSI -> ee candidates
19// and like sign pairs and for drawing corresponding distributions
20//
21// Author: C.Di Giglio, carmelo.digiglio@ba.infn.it
22///////////////////////////////////////////////////////////////////////////
23
24#include <TSystem.h>
25#include <TROOT.h>
26#include <TClonesArray.h>
27
28#include "AliAODEvent.h"
29#include "AliAODVertex.h"
30#include "AliAODTrack.h"
31#include "AliAODRecoDecayHF2Prong.h"
32#include "AliAnalysisTaskSE.h"
33#include "AliAnalysisTaskSEBkgLikeSignJPSI.h"
34
35ClassImp(AliAnalysisTaskSEBkgLikeSignJPSI)
36
37//________________________________________________________________________
38AliAnalysisTaskSEBkgLikeSignJPSI::AliAnalysisTaskSEBkgLikeSignJPSI():
39AliAnalysisTaskSE(),
40fOutput(0),
41fHistMassJPSI(0),
42fHistMassLS(0),
43fHistCtsJPSI(0),
44fHistCtsLS(0),
45fHistCtsLSpos(0),
46fHistCtsLSneg(0),
47fHistCPtaJPSI(0),
48fHistCPtaLS(0),
49fHistd0d0JPSI(0),
50fHistd0d0LS(0),
51fHistDCAJPSI(0),
52fHistDCALS(0),
53fVHF(0),
54fTotPosPairs(0),
55fTotNegPairs(0),
56fLsNormalization(1.)
57{
58 //
59 // Default constructor
60 //
61 // Output slot #1 writes into a TList container
62 DefineOutput(1,TList::Class()); //My private output
63}
64
65//________________________________________________________________________
66AliAnalysisTaskSEBkgLikeSignJPSI::AliAnalysisTaskSEBkgLikeSignJPSI(const char *name):
67AliAnalysisTaskSE(name),
68fOutput(0),
69fHistMassJPSI(0),
70fHistMassLS(0),
71fHistCtsJPSI(0),
72fHistCtsLS(0),
73fHistCtsLSpos(0),
74fHistCtsLSneg(0),
75fHistCPtaJPSI(0),
76fHistCPtaLS(0),
77fHistd0d0JPSI(0),
78fHistd0d0LS(0),
79fHistDCAJPSI(0),
80fHistDCALS(0),
81fVHF(0),
82fTotPosPairs(0),
83fTotNegPairs(0),
84fLsNormalization(1.)
85{
86 //
87 // Default constructor
88 //
89 // Output slot #1 writes into a TList container
90 DefineOutput(1,TList::Class()); //My private output
91}
92
93//________________________________________________________________________
94AliAnalysisTaskSEBkgLikeSignJPSI::~AliAnalysisTaskSEBkgLikeSignJPSI()
95{
96 // Destructor
97 if (fOutput) {
98 delete fOutput;
99 fOutput = 0;
100 }
101
102 if (fVHF) {
103 delete fVHF;
104 fVHF = 0;
105 }
106
107}
108//________________________________________________________________________
109void AliAnalysisTaskSEBkgLikeSignJPSI::Init()
110{
111 // Initialization
112
113 if(fDebug > 1) printf("AnalysisTaskSEBkgLikeSignJPSI::Init() \n");
114
115 gROOT->LoadMacro("ConfigVertexingHF.C");
116
117 fVHF = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
118 fVHF->PrintStatus();
119
120 return;
121
122 return;
123}
124
125//________________________________________________________________________
126void AliAnalysisTaskSEBkgLikeSignJPSI::UserCreateOutputObjects()
127{
128 // Create the output container
129 //
130 if(fDebug > 1) printf("AnalysisTaskSEBkgLikeSignJPSI::UserCreateOutputObjects() \n");
131
132 // Several histograms are more conveniently managed in a TList
133 fOutput = new TList();
134 fOutput->SetOwner();
135
136 fHistMassJPSI = new TH1F("fHistMassJPSI", "J/#Psi invariant mass; M [GeV]; Entries",200,2.8,3.25);
137 fHistMassJPSI->Sumw2();
138 fHistMassJPSI->SetMinimum(0);
139 fOutput->Add(fHistMassJPSI);
140
141 fHistMassLS = new TH1F("fHistMassLS", "Like sign pairs invariant mass; M [GeV]; Entries",200,2.8,3.25);
142 fHistMassLS->Sumw2();
143 fHistMassLS->SetMinimum(0);
144 fOutput->Add(fHistMassLS);
145
146 fHistCtsJPSI = new TH1F("fHistCtsJPSI", "J/#Psi cosine of decay angle; Cos#Theta*; Entries",200,-1.,1.);
147 fHistCtsJPSI->Sumw2();
148 fHistCtsJPSI->SetMinimum(0);
149 fOutput->Add(fHistCtsJPSI);
150
151 fHistCtsLS = new TH1F("fHistCtsLS", "Like sign pairs cosine of decay angle; Cos#Theta*; Entries",200,-1.,1.);
152 fHistCtsLS->Sumw2();
153 fHistCtsLS->SetMinimum(0);
154 fOutput->Add(fHistCtsLS);
155
156 fHistCtsLSpos = new TH1F("fHistCtsLSpos", "Like sign ++ pairs cosine of decay angle; Cos#Theta*; Entries",200,-1.,1.);
157 fHistCtsLSpos->Sumw2();
158 fHistCtsLSpos->SetMinimum(0);
159 fOutput->Add(fHistCtsLSpos);
160
161 fHistCtsLSneg = new TH1F("fHistCtsLSneg", "Like sign -- pairs cosine of decay angle; Cos#Theta*; Entries",200,-1.,1.);
162 fHistCtsLSneg->Sumw2();
163 fHistCtsLSneg->SetMinimum(0);
164 fOutput->Add(fHistCtsLSneg);
165
166 fHistCPtaJPSI = new TH1F("fHistCPtaJPSI", "J/#Psi cosine of pointing angle; Cos#Theta_{point}; Entries",200,-1.,1.);
167 fHistCPtaJPSI->Sumw2();
168 fHistCPtaJPSI->SetMinimum(0);
169 fOutput->Add(fHistCPtaJPSI);
170
171 fHistCPtaLS = new TH1F("fHistCPtaLS", "Like sign pairs cosine of pointing angle; Cos#Theta_{point}; Entries",200,-1.,1.);
172 fHistCPtaLS->Sumw2();
173 fHistCPtaLS->SetMinimum(0);
174 fOutput->Add(fHistCPtaLS);
175
176 fHistd0d0JPSI = new TH1F("fHistd0d0JPSI", "J/#Psi product of impact parameters; d0xd0 [#mu m^{2}]; Entries",200,-100000.,100000.);
177 fHistd0d0JPSI->Sumw2();
178 fHistd0d0JPSI->SetMinimum(0);
179 fOutput->Add(fHistd0d0JPSI);
180
181 fHistd0d0LS = new TH1F("fHistd0d0LS", "Like sign pairs product of impact parameters; d0xd0 [#mu m^{2}]; Entries",200,-100000.,100000.);
182 fHistd0d0LS->Sumw2();
183 fHistd0d0LS->SetMinimum(0);
184 fOutput->Add(fHistd0d0LS);
185
186 fHistDCAJPSI = new TH1F("fHistDCAJPSI", "J/#Psi distance of closest approach; dca [10^{2}#mu m]; Entries",100,0.,5.);
187 fHistDCAJPSI->Sumw2();
188 fHistDCAJPSI->SetMinimum(0);
189 fOutput->Add(fHistDCAJPSI);
190
191 fHistDCALS = new TH1F("fHistDCALS", "Like sign pairs distance of closest approach; dca [10^{2}#mu m]; Entries",100,0.,5.);
192 fHistDCALS->Sumw2();
193 fHistDCALS->SetMinimum(0);
194 fOutput->Add(fHistDCALS);
195
196 return;
197}
198
199//________________________________________________________________________
200void AliAnalysisTaskSEBkgLikeSignJPSI::UserExec(Option_t */*option*/)
201{
202 // Execute analysis for current event:
203 // heavy flavor candidates association to MC truth
204
205 AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
206
207 // load heavy flavour vertices
208 TClonesArray *arrayVerticesHF =
209 (TClonesArray*)aod->GetList()->FindObject("VerticesHF");
210 if(!arrayVerticesHF) {
211 printf("AliAnalysisTaskSEBkgLikeSignJPSI::UserExec: VerticesHF branch not found!\n");
212 return;
213 }
214
215 // load JPSI->ee candidates
216 TClonesArray *arrayJPSItoEle =
217 (TClonesArray*)aod->GetList()->FindObject("JPSItoEle");
218 if(!arrayJPSItoEle) {
219 printf("AliAnalysisTaskSEBkgLikeSignJPSI::UserExec: JPSItoEle branch not found!\n");
220 return;
221 }
222
223 // load like sign candidates
224 TClonesArray *arrayLikeSign =
225 (TClonesArray*)aod->GetList()->FindObject("LikeSign2Prong");
226 if(!arrayLikeSign) {
227 printf("AliAnalysisTaskSEBkgLikeSignJPSI::UserExec: LikeSign2Prong branch not found!\n");
228 return;
229 }
230
231 // AOD primary vertex
232 AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
233
234 // make trkIDtoEntry register (temporary)
235 Int_t trkIDtoEntry[100000];
236 for(Int_t it=0;it<aod->GetNumberOfTracks();it++) {
237 AliAODTrack *track = aod->GetTrack(it);
238 trkIDtoEntry[track->GetID()]=it;
239 }
240
241 // loop over Like sign candidates
242 Int_t nPosPairs=0,nNegPairs=0;
243 Int_t nLikeSign = arrayLikeSign->GetEntriesFast();
244 printf("+++\n+++Number of like sign pairs ---> %d \n+++\n", nLikeSign);
245
246 for(Int_t iLikeSign = 0; iLikeSign < nLikeSign; iLikeSign++) {
247 AliAODRecoDecayHF2Prong *d = (AliAODRecoDecayHF2Prong*)arrayLikeSign->UncheckedAt(iLikeSign);
248 Bool_t unsetvtx=kFALSE;
249 if(!d->GetOwnPrimaryVtx()) {
250 d->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
251 unsetvtx=kTRUE;
252 }
253 Int_t okBtoJPSIls=0;
254 if(d->SelectBtoJPSI(fVHF->GetBtoJPSICuts(),okBtoJPSIls)) {
255 fHistMassLS->Fill(d->InvMassJPSIee());
256 fHistCPtaLS->Fill(d->CosPointingAngle());
257 fHistd0d0LS->Fill(1e8*d->Prodd0d0());
258 fHistCtsLS->Fill(d->CosThetaStarJPSI());
259 fHistDCALS->Fill(100*d->GetDCA());
260 PostData(1,fOutput);
261 AliAODTrack *trk0 = (AliAODTrack*)d->GetDaughter(0);
262 if(!trk0) {
263 trk0=aod->GetTrack(trkIDtoEntry[d->GetProngID(0)]);
264 printf("references to standard AOD not available \n");
265 }
266 if((trk0->Charge())==1) {
267 nPosPairs++;
268 fHistCtsLSpos->Fill(d->CosThetaStarJPSI());
269 PostData(1,fOutput);
270 } else {
271 nNegPairs++;
272 fHistCtsLSneg->Fill(d->CosThetaStarJPSI());
273 PostData(1,fOutput);
274 }
275 PostData(1,fOutput);
276 }
277 if(unsetvtx) d->UnsetOwnPrimaryVtx();
278 }
279
280 printf("------------ N. of positive pairs in Event ----- %d \n", nPosPairs);
281 printf("------------ N. of negative pairs in Event ----- %d \n", nNegPairs);
282
283 fTotPosPairs += nPosPairs;
284 fTotNegPairs += nNegPairs;
285
286 // loop over JPSI candidates
287 Int_t nBtoJpsiToEle = arrayJPSItoEle->GetEntriesFast();
288 printf("Number of like JPSI -> ee candidates ---> %d \n", nBtoJpsiToEle);
289
290 for (Int_t iBtoJpsiToEle = 0; iBtoJpsiToEle < nBtoJpsiToEle; iBtoJpsiToEle++) {
291 AliAODRecoDecayHF2Prong *d = (AliAODRecoDecayHF2Prong*)arrayJPSItoEle->UncheckedAt(iBtoJpsiToEle);
292 Bool_t unsetvtx=kFALSE;
293 if(!d->GetOwnPrimaryVtx()) {
294 d->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
295 unsetvtx=kTRUE;
296 }
297 Int_t okBtoJPSI=0;
298 if(d->SelectBtoJPSI(fVHF->GetBtoJPSICuts(),okBtoJPSI)) {
299 fHistMassJPSI->Fill(d->InvMassJPSIee());
300 fHistCtsJPSI->Fill(d->CosThetaStarJPSI());
301 fHistd0d0JPSI->Fill(1e8*d->Prodd0d0());
302 fHistCPtaJPSI->Fill(d->CosPointingAngle());
303 fHistDCAJPSI->Fill(100*d->GetDCA());
304 PostData(1,fOutput);
305 }
306 if(unsetvtx) d->UnsetOwnPrimaryVtx();
307 }
308
309 return;
310}
311
312//________________________________________________________________________
313void AliAnalysisTaskSEBkgLikeSignJPSI::Terminate(Option_t */*option*/)
314{
315 // Terminate analysis
316 //
317 if(fDebug > 1) printf("AnalysisTaskSEBkgLikeSignJPSI: Terminate() \n");
318
319 fOutput = dynamic_cast<TList*> (GetOutputData(1));
320 if (!fOutput) {
321 printf("ERROR: fOutput not available\n");
322 return;
323 }
324
325 fLsNormalization = 2.*TMath::Sqrt(fTotPosPairs*fTotNegPairs);
326
327 fHistMassJPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistMassJPSI"));
328 fHistMassLS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistMassLS"));
329 fHistCtsJPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCtsJPSI"));
330 fHistCtsLS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCtsLS"));
331 fHistCtsLSpos = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCtsLSpos"));
332 fHistCtsLSneg = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCtsLSneg"));
333 fHistCPtaJPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCPtaJPSI"));
334 fHistCPtaLS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistCPtaLS"));
335 fHistd0d0JPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistd0d0JPSI"));
336 fHistd0d0LS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistd0d0LS"));
337 fHistDCAJPSI = dynamic_cast<TH1F*>(fOutput->FindObject("fHistDCAJPSI"));
338 fHistDCALS = dynamic_cast<TH1F*>(fOutput->FindObject("fHistDCALS"));
339
340 fHistMassLS->Scale((1/fLsNormalization)*fHistMassLS->GetEntries());
341 fHistCtsLS->Scale((1/fLsNormalization)*fHistCtsLS->GetEntries());
342 fHistCtsLSpos->Scale((1/fLsNormalization)*fHistCtsLSpos->GetEntries());
343 fHistCtsLSneg->Scale((1/fLsNormalization)*fHistCtsLSneg->GetEntries());
344 fHistCPtaLS->Scale((1/fLsNormalization)*fHistCPtaLS->GetEntries());
345 fHistd0d0LS->Scale((1/fLsNormalization)*fHistd0d0LS->GetEntries());
346 fHistDCALS->Scale((1/fLsNormalization)*fHistDCALS->GetEntries());
347
348 return;
349}