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