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