1 /**************************************************************************
2 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 /////////////////////////////////////////////////////////////
18 // AliAnalysisTaskSE for D0 candidates invariant mass histogram
19 // and comparison with the MC truth and cut variables distributions.
21 // Authors: A.Dainese, andrea.dainese@lnl.infn.it
22 // Chiara Bianchin, chiara.bianchin@pd.infn.it (invariant mass)
23 // Carmelo Di Giglio, carmelo.digiglio@ba.infn.it (like sign)
24 /////////////////////////////////////////////////////////////
26 #include <Riostream.h>
27 #include <TClonesArray.h>
32 #include <TDatabasePDG.h>
34 #include "AliAODEvent.h"
35 #include "AliAODVertex.h"
36 #include "AliAODTrack.h"
37 #include "AliAODMCHeader.h"
38 #include "AliAODMCParticle.h"
39 #include "AliAODRecoDecayHF2Prong.h"
40 #include "AliAODRecoCascadeHF.h"
41 #include "AliAnalysisVertexingHF.h"
42 #include "AliAnalysisTaskSE.h"
43 #include "AliAnalysisTaskSED0Mass.h"
45 ClassImp(AliAnalysisTaskSED0Mass)
48 //________________________________________________________________________
49 AliAnalysisTaskSED0Mass::AliAnalysisTaskSED0Mass():
61 // Default constructor
64 //________________________________________________________________________
65 AliAnalysisTaskSED0Mass::AliAnalysisTaskSED0Mass(const char *name):
66 AliAnalysisTaskSE(name),
76 // Default constructor
78 // Output slot #1 writes into a TList container
79 DefineOutput(1,TList::Class()); //My private output
80 // Output slot #2 writes into a TList container
81 DefineOutput(2,TList::Class()); //My private output
82 // Output slot #3 writes into a TH1F container
83 DefineOutput(3,TH1F::Class()); //My private output
84 // Output slot #4 writes into a TList container
85 DefineOutput(4,TList::Class()); //My private output
88 //________________________________________________________________________
89 AliAnalysisTaskSED0Mass::~AliAnalysisTaskSED0Mass()
121 //________________________________________________________________________
122 void AliAnalysisTaskSED0Mass::Init()
126 if(fDebug > 1) printf("AnalysisTaskSED0Mass::Init() \n");
128 gROOT->LoadMacro("$ALICE_ROOT/PWG3/vertexingHF/ConfigVertexingHF.C");
130 // 2 sets of dedidcated cuts -- defined in UserExec
131 fVHFtight = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
132 fVHFloose = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
137 //________________________________________________________________________
138 void AliAnalysisTaskSED0Mass::UserCreateOutputObjects()
141 // Create the output container
143 if(fDebug > 1) printf("AnalysisTaskSED0Mass::UserCreateOutputObjects() \n");
145 // Several histograms are more conveniently managed in a TList
146 fOutputtight = new TList();
147 fOutputtight->SetOwner();
148 fOutputtight->SetName("listtight");
150 fOutputloose = new TList();
151 fOutputloose->SetOwner();
152 fOutputloose->SetName("listloose");
154 fDistr = new TList();
156 fDistr->SetName("distributionslist");
160 TString nameMass=" ", nameSgn=" ", nameBkg=" ", nameRfl=" ", namedistr=" ";
162 for(Int_t i=0;i<nhist;i++){
163 nameMass="histMass_";
172 //histograms of invariant mass distributions
174 TH1F* tmpMt = new TH1F(nameMass.Data(),"D^{0} invariant mass; M [GeV]; Entries",200,1.765,1.965);
175 TH1F *tmpMl=(TH1F*)tmpMt->Clone();
179 TH1F* tmpSt = new TH1F(nameSgn.Data(), "D^{0} invariant mass - MC; M [GeV]; Entries",200,1.765,1.965);
180 TH1F *tmpSl=(TH1F*)tmpSt->Clone();
184 TH1F* tmpBt = new TH1F(nameBkg.Data(), "Background invariant mass - MC; M [GeV]; Entries",200,1.765,1.965);
185 TH1F *tmpBl=(TH1F*)tmpBt->Clone();
189 //Reflection: histo filled with D0Mass which pass the cut (also) as D0bar and with D0bar which pass (also) the cut as D0
190 TH1F* tmpRt = new TH1F(nameRfl.Data(), "Reflected signal invariant mass - MC; M [GeV]; Entries",200,1.765,1.965);
191 TH1F *tmpRl=(TH1F*)tmpRt->Clone();
194 // printf("Created histograms %s\t%s\t%s\t%s\n",tmpM->GetName(),tmpS->GetName(),tmpB->GetName(),tmpR->GetName());
196 fOutputtight->Add(tmpMt);
197 fOutputtight->Add(tmpSt);
198 fOutputtight->Add(tmpBt);
199 fOutputtight->Add(tmpRt);
201 fOutputloose->Add(tmpMl);
202 fOutputloose->Add(tmpSl);
203 fOutputloose->Add(tmpBl);
204 fOutputloose->Add(tmpRl);
208 //histograms of cut variable distributions
211 TH1F *hptpiS = new TH1F(namedistr.Data(), "P_{T} distribution (pions);p_{T} [GeV/c]",200,0.,8.);
214 TH1F *hptKS = new TH1F(namedistr.Data(), "P_{T} distribution (kaons);p_{T} [GeV/c]",200,0.,8.);
217 TH1F *hptB = new TH1F(namedistr.Data(), "P_{T} distribution;p_{T} [GeV/c]",200,0.,8.);
221 TH1F *hdcaS = new TH1F(namedistr.Data(), "DCA distribution;dca [cm]",200,0.,0.1);
223 TH1F *hdcaB = new TH1F(namedistr.Data(), "DCA distribution;dca [cm]",200,0.,0.1);
226 namedistr="hcosthetastarS";
227 TH1F *hcosthetastarS = new TH1F(namedistr.Data(), "cos#theta* distribution;cos#theta*",200,-1.,1.);
228 namedistr="hcosthetastarB";
229 TH1F *hcosthetastarB = new TH1F(namedistr.Data(), "cos#theta* distribution;cos#theta*",200,-1.,1.);
233 TH1F *hd0piS = new TH1F(namedistr.Data(), "Impact parameter distribution (pions);d0(#pi) [cm]",200,-0.1,0.1);
236 TH1F *hd0KS = new TH1F(namedistr.Data(), "Impact parameter distribution (kaons);d0(K) [cm]",200,-0.1,0.1);
238 TH1F *hd0B = new TH1F(namedistr.Data(), "Impact parameter distribution;d0 [cm]",200,-0.1,0.1);
241 TH1F *hd0d0S = new TH1F(namedistr.Data(), "d_{0}#timesd_{0} distribution;d_{0}#timesd_{0} [cm^{2}]",200,-0.001,0.001);
243 TH1F *hd0d0B = new TH1F(namedistr.Data(), "d_{0}#timesd_{0} distribution;d_{0}#timesd_{0} [cm^{2}]",200,-0.001,0.001);
246 namedistr="hcosthetapointS";
247 TH1F *hcosthetapointS = new TH1F(namedistr.Data(), "cos#theta_{Point} distribution;cos#theta_{Point}",200,0,1.);
248 namedistr="hcosthetapointB";
249 TH1F *hcosthetapointB = new TH1F(namedistr.Data(), "cos#theta_{Point} distribution;cos#theta_{Point}",200,0,1.);
251 namedistr="hcosthpointd0d0S";
252 TH2F *hcosthpointd0d0S= new TH2F(namedistr.Data(),"Correlation cos#theta_{Point}-d_{0}#timesd_{0};cos#theta_{Point};d_{0}#timesd_{0} [cm^{2}]",200,0,1.,200,-0.001,0.001);
253 namedistr="hcosthpointd0d0B";
254 TH2F *hcosthpointd0d0B= new TH2F(namedistr.Data(),"Correlation cos#theta_{Point}-d_{0}#timesd_{0};cos#theta_{Point};d_{0}#timesd_{0} [cm^{2}]",200,0,1.,200,-0.001,0.001);
270 fDistr->Add(hcosthetastarS);
271 fDistr->Add(hcosthetastarB);
273 fDistr->Add(hcosthetapointS);
274 fDistr->Add(hcosthetapointB);
276 fDistr->Add(hcosthpointd0d0S);
277 fDistr->Add(hcosthpointd0d0B);
279 fNentries=new TH1F("nentriesD0", "Look at the number of entries! it is = to the number of AODs", 2,1.,2.);
284 //________________________________________________________________________
285 void AliAnalysisTaskSED0Mass::UserExec(Option_t */*option*/)
287 // Execute analysis for current event:
288 // heavy flavor candidates association to MC truth
289 //cout<<"I'm in UserExec"<<endl;
290 AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
292 TClonesArray *inputArray=0;
294 if(fArray==0){ //D0 candidates
295 // load D0->Kpi candidates
296 //cout<<"D0 candidates"<<endl;
298 (TClonesArray*)aod->GetList()->FindObject("D0toKpi");
300 printf("AliAnalysisTaskSECompareHFpt::UserExec: D0toKpi branch not found!\n");
304 } else { //LikeSign candidates
306 //cout<<"LS candidates"<<endl;
307 // load like sign candidates
309 (TClonesArray*)aod->GetList()->FindObject("LikeSign2Prong");
311 printf("AliAnalysisTaskSEBkgLikeSignD0::UserExec: LikeSign2Prong branch not found!\n");
317 // AOD primary vertex
318 AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
322 TClonesArray *mcArray =
323 (TClonesArray*)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
325 printf("AliAnalysisTaskSED0Mass::UserExec: MC particles branch not found!\n");
330 AliAODMCHeader *mcHeader =
331 (AliAODMCHeader*)aod->GetList()->FindObject(AliAODMCHeader::StdBranchName());
333 printf("AliAnalysisTaskSED0Mass::UserExec: MC header branch not found!\n");
338 //histogram filled with 1 for every AOD
340 PostData(3,fNentries);
341 // loop over candidates
342 Int_t nInD0toKpi = inputArray->GetEntriesFast();
343 printf("Number of D0->Kpi: %d\n",nInD0toKpi);
345 for (Int_t iD0toKpi = 0; iD0toKpi < nInD0toKpi; iD0toKpi++) {
346 //Int_t nPosPairs=0, nNegPairs=0;
347 //cout<<"inside the loop"<<endl;
348 AliAODRecoDecayHF2Prong *d = (AliAODRecoDecayHF2Prong*)inputArray->UncheckedAt(iD0toKpi);
349 Bool_t unsetvtx=kFALSE;
350 if(!d->GetOwnPrimaryVtx()) {
351 d->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
356 // printf(" |M-MD0| [GeV] < %f\n",fD0toKpiCuts[0]);
357 // printf(" dca [cm] < %f\n",fD0toKpiCuts[1]);
358 // printf(" cosThetaStar < %f\n",fD0toKpiCuts[2]);
359 // printf(" pTK [GeV/c] > %f\n",fD0toKpiCuts[3]);
360 // printf(" pTpi [GeV/c] > %f\n",fD0toKpiCuts[4]);
361 // printf(" |d0K| [cm] < %f\n",fD0toKpiCuts[5]);
362 // printf(" |d0pi| [cm] < %f\n",fD0toKpiCuts[6]);
363 // printf(" d0d0 [cm^2] < %f\n",fD0toKpiCuts[7]);
364 // printf(" cosThetaPoint > %f\n",fD0toKpiCuts[8]);
369 Double_t mPDG=TDatabasePDG::Instance()->GetParticle(421)->Mass();
372 Double_t minvD0 = d->InvMassD0();
375 Double_t minvD0bar = d->InvMassD0bar();
376 //apply cut on invariant mass on the pair
377 if(TMath::Abs(minvD0-mPDG)<0.03 || TMath::Abs(minvD0bar-mPDG)<0.03){
378 //cout<<"inside mass cut"<<endl;
379 Int_t pdgDgD0toKpi[2]={321,211};
380 Int_t lab=d->MatchToMC(421,mcArray,2,pdgDgD0toKpi); //return MC particle label if the array corresponds to a D0, -1 if not (cf. AliAODRecoDecay.cxx)
383 if(fArray==1) cout<<"LS signal: ERROR"<<endl;
384 for (Int_t iprong=0; iprong<2; iprong++){
385 AliAODTrack *prong=(AliAODTrack*)d->GetDaughter(iprong);
386 Int_t labprong=prong->GetLabel();
388 //cout<<"prong name = "<<prong->GetName()<<" label = "<<prong->GetLabel()<<endl;
389 AliAODMCParticle *mcprong=0;
390 if(labprong>=0) mcprong= (AliAODMCParticle*)mcArray->At(labprong);
391 Int_t pdgprong=mcprong->GetPdgCode();
392 if(TMath::Abs(pdgprong)==211) {
394 ((TH1F*)fDistr->FindObject("hptpiS"))->Fill(d->PtProng(iprong));
395 ((TH1F*)fDistr->FindObject("hd0piS"))->Fill(d->Getd0Prong(iprong));
398 if(TMath::Abs(pdgprong)==321) {
399 //cout<<"kappa"<<endl;
400 ((TH1F*)fDistr->FindObject("hptKS"))->Fill(d->PtProng(iprong));
401 ((TH1F*)fDistr->FindObject("hd0KS"))->Fill(d->Getd0Prong(iprong));
403 ((TH1F*)fDistr->FindObject("hdcaS"))->Fill(d->GetDCA());
407 if (((AliAODMCParticle*)mcArray->At(lab))->GetPdgCode() == 421)
408 ((TH1F*)fDistr->FindObject("hcosthetastarS"))->Fill(d->CosThetaStarD0());
409 else ((TH1F*)fDistr->FindObject("hcosthetastarS"))->Fill(d->CosThetaStarD0bar());
411 ((TH1F*)fDistr->FindObject("hd0d0S"))->Fill(d->Prodd0d0());
413 ((TH1F*)fDistr->FindObject("hcosthetapointS"))->Fill(d->CosPointingAngle());
414 ((TH1F*)fDistr->FindObject("hcosthpointd0d0S"))->Fill(d->CosPointingAngle(),d->Prodd0d0());
416 //cout<<"ok point"<<endl;
419 } else{ //Background or LS
420 //cout<<"is background"<<endl;
421 AliAODTrack *prong=(AliAODTrack*)d->GetDaughter(0);
422 if(!prong) cout<<"No daughter found";
424 if(prong->Charge()==1) fTotPosPairs[4]++; else fTotNegPairs[4]++;
426 ((TH1F*)fDistr->FindObject("hptB"))->Fill(d->PtProng(0));
427 //cout<<"ptok"<<endl;
428 ((TH1F*)fDistr->FindObject("hd0B"))->Fill(d->Getd0Prong(0));
429 //cout<<"d0ok"<<endl;
430 ((TH1F*)fDistr->FindObject("hdcaB"))->Fill(d->GetDCA());
431 //cout<<"dcaok"<<endl;
432 ((TH1F*)fDistr->FindObject("hcosthetastarB"))->Fill(d->CosThetaStarD0());
433 ((TH1F*)fDistr->FindObject("hcosthetastarB"))->Fill(d->CosThetaStarD0bar());
434 ((TH1F*)fDistr->FindObject("hd0d0B"))->Fill(d->Prodd0d0());
435 //cout<<"d0d0ok"<<endl;
436 ((TH1F*)fDistr->FindObject("hcosthetapointB"))->Fill(d->CosPointingAngle());
437 ((TH1F*)fDistr->FindObject("hcosthpointd0d0B"))->Fill(d->CosPointingAngle(),d->Prodd0d0());
439 //cout<<"pointok"<<endl;
446 Double_t pt = d->Pt();
450 //cout<<"P_t = "<<pt<<endl;
451 if (pt>0. && pt<=1.) {
453 fVHFtight->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,0.05,0.05,-0.0003,0.7);
454 fVHFloose->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,0.05,0.05,-0.00025,0.7);
455 // fVHFtight->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,0.05,0.05,-0.0002,0.7);
456 // fVHFloose->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,1,1,-0.00015,0.5);
457 //printf("I'm in the bin %d\n",ptbin);
461 if(pt>1. && pt<=3.) {
463 fVHFtight->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.0003,0.9);
464 fVHFloose->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,1,1,-0.00025,0.8);
465 //printf("I'm in the bin %d\n",ptbin);
470 fVHFtight->SetD0toKpiCuts(0.7,0.015,0.8,0.7,0.7,0.05,0.05,-0.0002,0.9);
471 fVHFloose->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.00015,0.8);
472 //printf("I'm in the bin %d\n",ptbin);
476 fVHFtight->SetD0toKpiCuts(0.7,0.015,0.8,0.7,0.7,0.05,0.05,-0.0002,0.95);
477 fVHFloose->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.00015,0.9);
480 //printf("I'm in the bin %d\n",ptbin);
482 //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
484 FillHists(ptbin,d,mcArray,fVHFtight,fOutputtight);
485 FillHists(ptbin,d,mcArray,fVHFloose,fOutputloose);
486 if(unsetvtx) d->UnsetOwnPrimaryVtx();
493 PostData(1,fOutputtight);
494 PostData(2,fOutputloose);
498 //____________________________________________________________________________*
499 void AliAnalysisTaskSED0Mass::FillHists(Int_t ptbin, AliAODRecoDecayHF2Prong *part, TClonesArray *arrMC, AliAnalysisVertexingHF *vhf, TList *listout){
501 // function used in UserExec:
503 Int_t okD0=0,okD0bar=0;
504 //cout<<"inside FillHist"<<endl;
505 if(part->SelectD0(vhf->GetD0toKpiCuts(),okD0,okD0bar)) {//selected
506 Double_t invmassD0 = part->InvMassD0(), invmassD0bar = part->InvMassD0bar();
507 //printf("SELECTED\n");
509 AliAODTrack *prong=(AliAODTrack*)part->GetDaughter(0);
510 if(!prong) cout<<"No daughter found";
512 if(prong->Charge()==1) fTotPosPairs[ptbin-1]++; else fTotNegPairs[ptbin-1]++;
516 Int_t pdgDgD0toKpi[2]={321,211};
517 Int_t labD0 = part->MatchToMC(421,arrMC,2,pdgDgD0toKpi); //return MC particle label if the array corresponds to a D0, -1 if not (cf. AliAODRecoDecay.cxx)
520 fillthis="histMass_";
522 //cout<<"Filling "<<fillthis<<endl;
524 //printf("Fill mass with D0");
525 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
528 if(fArray==1) cout<<"LS signal ERROR"<<endl;
530 AliAODMCParticle *partD0 = (AliAODMCParticle*)arrMC->At(labD0);
531 Int_t pdgD0 = partD0->GetPdgCode();
532 //cout<<"pdg = "<<pdgD0<<endl;
533 if (pdgD0==421){ //D0
534 //cout<<"Fill S with D0"<<endl;
537 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
539 } else{ //it was a D0bar
542 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
547 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
552 fillthis="histMass_";
554 //printf("Fill mass with D0bar");
555 ((TH1F*)listout->FindObject(fillthis))->Fill(invmassD0bar);
558 if(fArray==1) cout<<"LS signal ERROR"<<endl;
559 AliAODMCParticle *partD0 = (AliAODMCParticle*)arrMC->At(labD0);
560 Int_t pdgD0 = partD0->GetPdgCode();
561 //cout<<" pdg = "<<pdgD0<<endl;
562 if (pdgD0==-421){ //D0bar
565 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
570 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
572 } else {//background or LS
575 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
580 } //else cout<<"NOT SELECTED"<<endl;
583 //________________________________________________________________________
584 void AliAnalysisTaskSED0Mass::Terminate(Option_t */*option*/)
586 // Terminate analysis
588 if(fDebug > 1) printf("AnalysisTaskSED0Mass: Terminate() \n");
590 fOutputtight = dynamic_cast<TList*> (GetOutputData(1));
592 printf("ERROR: fOutputthight not available\n");
595 fOutputloose = dynamic_cast<TList*> (GetOutputData(2));
597 printf("ERROR: fOutputloose not available\n");
600 fDistr = dynamic_cast<TList*> (GetOutputData(3));
602 printf("ERROR: fDistr not available\n");
607 for(Int_t ipt=0;ipt<4;ipt++){
608 fLsNormalization = 2.*TMath::Sqrt(fTotPosPairs[ipt]*fTotNegPairs[ipt]);
611 if(fLsNormalization>0) {
613 TString massName="histMass_";
615 ((TH1F*)fOutputtight->FindObject(massName))->Scale((1/fLsNormalization)*((TH1F*)fOutputtight->FindObject(massName))->GetEntries());
616 ((TH1F*)fOutputloose->FindObject(massName))->Scale((1/fLsNormalization)*((TH1F*)fOutputloose->FindObject(massName))->GetEntries());
620 fLsNormalization = 2.*TMath::Sqrt(fTotPosPairs[4]*fTotNegPairs[4]);
622 if(fLsNormalization>0) {
624 TString nameDistr="hptB";
625 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
627 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
628 nameDistr="hcosthetastarB";
629 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
631 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
633 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
634 nameDistr="hcosthetapointB";
635 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
636 nameDistr="hcosthpointd0d0B";
637 ((TH2F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH2F*)fDistr->FindObject(nameDistr))->GetEntries());