]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG3/vertexingHF/AliAnalysisTaskSED0Mass.cxx
Reduced default debug level
[u/mrichter/AliRoot.git] / PWG3 / vertexingHF / AliAnalysisTaskSED0Mass.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 D0 candidates invariant mass histogram
19 // and comparison with the MC truth.
20 //
21 // Authors: A.Dainese, andrea.dainese@lnl.infn.it
22 // and Chiara Bianchin, chiara.bianchin@pd.infn.it
23 /////////////////////////////////////////////////////////////
24
25 #include <Riostream.h>
26 #include <TClonesArray.h>
27 #include <TNtuple.h>
28 #include <TList.h>
29 #include <TH1F.h>
30
31
32 #include "AliAODEvent.h"
33 #include "AliAODVertex.h"
34 #include "AliAODTrack.h"
35 #include "AliAODMCHeader.h"
36 #include "AliAODMCParticle.h"
37 #include "AliAODRecoDecayHF2Prong.h"
38 #include "AliAODRecoCascadeHF.h"
39 #include "AliAnalysisVertexingHF.h"
40 #include "AliAnalysisTaskSE.h"
41 #include "AliAnalysisTaskSED0Mass.h"
42
43 ClassImp(AliAnalysisTaskSED0Mass)
44
45
46 //________________________________________________________________________
47 AliAnalysisTaskSED0Mass::AliAnalysisTaskSED0Mass():
48 AliAnalysisTaskSE(),
49 fOutputtight(0),
50 fOutputloose(0), 
51 fVHFtight(0),
52 fVHFloose(0),
53 fNentries(0)
54 {
55   // Default constructor
56 }
57
58 //________________________________________________________________________
59 AliAnalysisTaskSED0Mass::AliAnalysisTaskSED0Mass(const char *name):
60 AliAnalysisTaskSE(name),
61 fOutputtight(0), 
62 fOutputloose(0), 
63 fVHFtight(0),
64 fVHFloose(0),
65 fNentries(0)
66 {
67   // Default constructor
68
69   // Output slot #1 writes into a TList container
70   DefineOutput(1,TList::Class());  //My private output
71   // Output slot #2 writes into a TList container
72   DefineOutput(2,TList::Class());  //My private output
73   // Output slot #3 writes into a TH1F container
74   DefineOutput(3,TH1F::Class());  //My private output
75 }
76
77 //________________________________________________________________________
78 AliAnalysisTaskSED0Mass::~AliAnalysisTaskSED0Mass()
79 {
80    if (fOutputtight) {
81     delete fOutputtight;
82     fOutputtight = 0;
83   }
84   if (fVHFtight) {
85     delete fVHFtight;
86     fVHFtight = 0;
87   }
88   if (fOutputloose) {
89     delete fOutputloose;
90     fOutputloose = 0;
91   }
92   if (fVHFloose) {
93     delete fVHFloose;
94     fVHFloose = 0;
95   }
96   if (fNentries){
97     delete fNentries;
98     fNentries = 0;
99   }
100 }  
101
102 //________________________________________________________________________
103 void AliAnalysisTaskSED0Mass::Init()
104 {
105   // Initialization
106
107   if(fDebug > 1) printf("AnalysisTaskSED0Mass::Init() \n");
108
109   gROOT->LoadMacro("$ALICE_ROOT/PWG3/vertexingHF/ConfigVertexingHF.C");
110
111   // 2 sets of dedidcated cuts -- defined in UserExec
112   fVHFtight = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
113   fVHFloose = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");  
114   
115   return;
116 }
117
118 //________________________________________________________________________
119 void AliAnalysisTaskSED0Mass::UserCreateOutputObjects()
120 {
121
122   // Create the output container
123   //
124   if(fDebug > 1) printf("AnalysisTaskSED0Mass::UserCreateOutputObjects() \n");
125
126   // Several histograms are more conveniently managed in a TList
127   fOutputtight = new TList();
128   fOutputtight->SetOwner();
129   fOutputtight->SetName("listtight");
130
131   fOutputloose = new TList();
132   fOutputloose->SetOwner();
133   fOutputloose->SetName("listloose");
134
135   const Int_t nhist=4;
136
137   TString nameMass=" ", nameSgn=" ", nameBkg=" ", nameRfl=" ";
138
139   for(Int_t i=0;i<nhist;i++){
140     nameMass="histMass_";
141     nameMass+=i+1;
142     nameSgn="histSgn_";
143     nameSgn+=i+1;
144     nameBkg="histBkg_";
145     nameBkg+=i+1;
146     nameRfl="histRfl_";
147     nameRfl+=i+1;
148
149     TH1F* tmpMt = new TH1F(nameMass.Data(),"D^{0} invariant mass; M [GeV]; Entries",200,1.765,1.965);
150     TH1F *tmpMl=(TH1F*)tmpMt->Clone();
151     tmpMt->Sumw2();
152     tmpMl->Sumw2();
153
154     TH1F* tmpSt = new TH1F(nameSgn.Data(), "D^{0} invariant mass - MC; M [GeV]; Entries",200,1.765,1.965);
155     TH1F *tmpSl=(TH1F*)tmpSt->Clone();
156     tmpSt->Sumw2();
157     tmpSl->Sumw2();
158
159     TH1F* tmpBt = new TH1F(nameBkg.Data(), "Background invariant mass - MC; M [GeV]; Entries",200,1.765,1.965);
160     TH1F *tmpBl=(TH1F*)tmpBt->Clone();
161     tmpBt->Sumw2();
162     tmpBl->Sumw2();
163
164     //Reflection: histo filled with D0 which pass the cut (also) as D0bar and with D0bar which pass (also) the cut as D0
165     TH1F* tmpRt = new TH1F(nameRfl.Data(), "Reflected signal invariant mass - MC; M [GeV]; Entries",200,1.765,1.965);
166     TH1F *tmpRl=(TH1F*)tmpRt->Clone();
167     tmpRt->Sumw2();
168     tmpRl->Sumw2();
169   //  printf("Created histograms %s\t%s\t%s\t%s\n",tmpM->GetName(),tmpS->GetName(),tmpB->GetName(),tmpR->GetName());
170
171     fOutputtight->Add(tmpMt);
172     fOutputtight->Add(tmpSt);
173     fOutputtight->Add(tmpBt);
174     fOutputtight->Add(tmpRt);
175
176     fOutputloose->Add(tmpMl);
177     fOutputloose->Add(tmpSl);
178     fOutputloose->Add(tmpBl);
179     fOutputloose->Add(tmpRl);
180     
181   }
182
183
184   //fOutputloose->ls();
185   //fOutputtight->ls();
186
187   fNentries=new TH1F("nentries", "Look at the number of entries! = number of AODs", 2,1.,2.);
188
189   return;
190 }
191
192 //________________________________________________________________________
193 void AliAnalysisTaskSED0Mass::UserExec(Option_t */*option*/)
194 {
195   // Execute analysis for current event:
196   // heavy flavor candidates association to MC truth
197   //cout<<"I'm in UserExec"<<endl;
198   AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
199  
200   // load D0->Kpi candidates                                                   
201   TClonesArray *inputArrayD0toKpi =
202     (TClonesArray*)aod->GetList()->FindObject("D0toKpi");
203   if(!inputArrayD0toKpi) {
204     printf("AliAnalysisTaskSECompareHFpt::UserExec: D0toKpi branch not found!\n");
205     return;
206   }
207   
208   // AOD primary vertex
209   AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
210   //vtx1->Print();
211   
212   // load MC particles
213   TClonesArray *mcArray = 
214     (TClonesArray*)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
215   if(!mcArray) {
216     printf("AliAnalysisTaskSED0Mass::UserExec: MC particles branch not found!\n");
217     return;
218   }
219   
220   // load MC header
221   AliAODMCHeader *mcHeader = 
222     (AliAODMCHeader*)aod->GetList()->FindObject(AliAODMCHeader::StdBranchName());
223   if(!mcHeader) {
224     printf("AliAnalysisTaskSED0Mass::UserExec: MC header branch not found!\n");
225     return;
226   }
227   
228
229   //histogram filled with 1 for every AOD
230   fNentries->Fill(1);
231   PostData(3,fNentries);
232   // loop over D0->Kpi candidates
233   Int_t nInD0toKpi = inputArrayD0toKpi->GetEntriesFast();
234   if(fDebug>1) printf("Number of D0->Kpi: %d\n",nInD0toKpi);
235   
236   for (Int_t iD0toKpi = 0; iD0toKpi < nInD0toKpi; iD0toKpi++) {
237     
238     AliAODRecoDecayHF2Prong *d = (AliAODRecoDecayHF2Prong*)inputArrayD0toKpi->UncheckedAt(iD0toKpi);
239     Bool_t unsetvtx=kFALSE;
240     if(!d->GetOwnPrimaryVtx()) {
241       d->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
242       unsetvtx=kTRUE;
243     }
244     
245      //cuts order
246 //       printf("    |M-MD0| [GeV]    < %f\n",fD0toKpiCuts[0]);
247 //     printf("    dca    [cm]  < %f\n",fD0toKpiCuts[1]);
248 //     printf("    cosThetaStar     < %f\n",fD0toKpiCuts[2]);
249 //     printf("    pTK     [GeV/c]    > %f\n",fD0toKpiCuts[3]);
250 //     printf("    pTpi    [GeV/c]    > %f\n",fD0toKpiCuts[4]);
251 //     printf("    |d0K|  [cm]  < %f\n",fD0toKpiCuts[5]);
252 //     printf("    |d0pi| [cm]  < %f\n",fD0toKpiCuts[6]);
253 //     printf("    d0d0  [cm^2] < %f\n",fD0toKpiCuts[7]);
254 //     printf("    cosThetaPoint    > %f\n",fD0toKpiCuts[8]);
255   
256     
257     Double_t pt = d->Pt();
258     Int_t ptbin=0;
259     
260     //cout<<"P_t = "<<pt<<endl;
261     if (pt>0. && pt<=1.) {
262       ptbin=1;
263       fVHFtight->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,0.05,0.05,-0.0003,0.7);
264       fVHFloose->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,0.05,0.05,-0.00025,0.7);
265 //       fVHFtight->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,0.05,0.05,-0.0002,0.7);
266 //       fVHFloose->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,1,1,-0.00015,0.5);
267       //printf("I'm in the bin %d\n",ptbin);
268       
269     }
270     
271     if(pt>1. && pt<=3.) {
272       ptbin=2;  
273       fVHFtight->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.0003,0.9);
274       fVHFloose->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,1,1,-0.00025,0.8);
275       //printf("I'm in the bin %d\n",ptbin);
276     }
277  
278     if(pt>3. && pt<=5.){
279         ptbin=3;  
280         fVHFtight->SetD0toKpiCuts(0.7,0.015,0.8,0.7,0.7,0.05,0.05,-0.0002,0.9);
281         fVHFloose->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.00015,0.8);
282         //printf("I'm in the bin %d\n",ptbin);
283     }
284     if(pt>5.){
285       ptbin=4;
286       fVHFtight->SetD0toKpiCuts(0.7,0.015,0.8,0.7,0.7,0.05,0.05,-0.0002,0.95);
287       fVHFloose->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.00015,0.9);
288     }//if(pt>5)
289   
290     //printf("I'm in the bin %d\n",ptbin);
291     //old
292     //fVHF->SetD0toKpiCuts(0.7,0.03,0.8,0.06,0.06,0.05,0.05,-0.0002,0.6); //2.p-p vertex reconstructed    
293
294     FillHists(ptbin,d,mcArray,fVHFtight,fOutputtight);
295     FillHists(ptbin,d,mcArray,fVHFloose,fOutputloose);
296     if(unsetvtx) d->UnsetOwnPrimaryVtx();
297   }
298   
299      
300   
301    
302   // Post the data
303   PostData(1,fOutputtight);
304   PostData(2,fOutputloose);
305
306   return;
307 }
308 //____________________________________________________________________________*
309 void AliAnalysisTaskSED0Mass::FillHists(Int_t ptbin, AliAODRecoDecayHF2Prong *part, TClonesArray *arrMC, AliAnalysisVertexingHF *vhf, TList *listout){
310   //
311   // function used in UserExec:
312   //
313   Int_t okD0=0,okD0bar=0;
314
315   if(part->SelectD0(vhf->GetD0toKpiCuts(),okD0,okD0bar)) {//selected
316     Double_t invmassD0 = part->InvMassD0(), invmassD0bar = part->InvMassD0bar();
317     //printf("SELECTED\n");
318     Int_t labD0 = part->MatchToMC(421,arrMC); //return MC particle label if the array corresponds to a D0, -1 if not (cf. AliAODRecoDecay.cxx)
319     //printf("labD0 %d",labD0);
320     
321     TString fillthis="";
322     
323     if (okD0==1) {
324       fillthis="histMass_";
325       fillthis+=ptbin;
326       //cout<<"Filling "<<fillthis<<endl;
327
328       //printf("Fill mass with D0");
329       ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
330       if(labD0>=0) {
331         AliAODMCParticle *partD0 = (AliAODMCParticle*)arrMC->At(labD0);
332         Int_t pdgD0 = partD0->GetPdgCode();
333         //cout<<"pdg = "<<pdgD0<<endl;
334         if (pdgD0==421){ //D0
335           //cout<<"Fill S with D0"<<endl;
336           fillthis="histSgn_";
337           fillthis+=ptbin;
338           ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
339
340         } else{ //it was a D0bar
341           fillthis="histRfl_";
342           fillthis+=ptbin;
343           ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
344         }
345       } else {//background
346         fillthis="histBkg_";
347         fillthis+=ptbin;
348         ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
349       }
350     }
351     if (okD0bar==1) {
352       fillthis="histMass_";
353       fillthis+=ptbin;
354       //printf("Fill mass with D0bar");
355       ((TH1F*)listout->FindObject(fillthis))->Fill(invmassD0bar);
356       
357       if(labD0>=0) {
358         AliAODMCParticle *partD0 = (AliAODMCParticle*)arrMC->At(labD0);
359         Int_t pdgD0 = partD0->GetPdgCode();
360         //cout<<" pdg = "<<pdgD0<<endl;
361         if (pdgD0==-421){ //D0bar
362           fillthis="histSgn_";
363           fillthis+=ptbin;
364           ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
365           
366         } else{
367           fillthis="histRfl_";
368           fillthis+=ptbin;
369           ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
370         }
371       } else {//background
372         fillthis="histBkg_";
373         fillthis+=ptbin;
374         ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
375       }
376     }
377
378
379   } //else cout<<"NOT SELECTED"<<endl;
380
381 }
382 //________________________________________________________________________
383 void AliAnalysisTaskSED0Mass::Terminate(Option_t */*option*/)
384 {
385   // Terminate analysis
386   //
387   if(fDebug > 1) printf("AnalysisTaskSED0Mass: Terminate() \n");
388
389   fOutputtight = dynamic_cast<TList*> (GetOutputData(1));
390   if (!fOutputtight) {     
391     printf("ERROR: fOutput not available\n");
392     return;
393   }
394   fOutputloose = dynamic_cast<TList*> (GetOutputData(2));
395   if (!fOutputloose) {     
396     printf("ERROR: fOutput not available\n");
397     return;
398   }
399
400
401   return;
402 }
403