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>
33 #include <TDatabasePDG.h>
35 #include "AliAnalysisManager.h"
36 #include "AliAODHandler.h"
37 #include "AliAODEvent.h"
38 #include "AliAODVertex.h"
39 #include "AliAODTrack.h"
40 #include "AliAODMCHeader.h"
41 #include "AliAODMCParticle.h"
42 #include "AliAODRecoDecayHF2Prong.h"
43 #include "AliAODRecoCascadeHF.h"
44 #include "AliAnalysisVertexingHF.h"
45 #include "AliAnalysisTaskSE.h"
46 #include "AliAnalysisTaskSED0Mass.h"
49 ClassImp(AliAnalysisTaskSED0Mass)
52 //________________________________________________________________________
53 AliAnalysisTaskSED0Mass::AliAnalysisTaskSED0Mass():
66 // Default constructor
67 for(Int_t i=0;i<5;i++) {fTotPosPairs[i]=0; fTotNegPairs[i]=0;}
70 //________________________________________________________________________
71 AliAnalysisTaskSED0Mass::AliAnalysisTaskSED0Mass(const char *name):
72 AliAnalysisTaskSE(name),
83 // Default constructor
84 for(Int_t i=0;i<5;i++) {fTotPosPairs[i]=0; fTotNegPairs[i]=0;}
86 // Output slot #1 writes into a TList container
87 DefineOutput(1,TList::Class()); //My private output
88 // Output slot #2 writes into a TList container
89 DefineOutput(2,TList::Class()); //My private output
90 // Output slot #3 writes into a TH1F container
91 DefineOutput(3,TH1F::Class()); //My private output
92 // Output slot #4 writes into a TList container
93 DefineOutput(4,TList::Class()); //My private output
96 //________________________________________________________________________
97 AliAnalysisTaskSED0Mass::~AliAnalysisTaskSED0Mass()
108 delete fOutputmycuts;
129 //________________________________________________________________________
130 void AliAnalysisTaskSED0Mass::Init()
134 if(fDebug > 1) printf("AnalysisTaskSED0Mass::Init() \n");
136 gROOT->LoadMacro("$ALICE_ROOT/PWG3/vertexingHF/ConfigVertexingHF.C");
138 // 2 sets of dedidcated cuts -- defined in UserExec
139 fVHFPPR = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
140 fVHFmycuts = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
145 //________________________________________________________________________
146 void AliAnalysisTaskSED0Mass::UserCreateOutputObjects()
149 // Create the output container
151 if(fDebug > 1) printf("AnalysisTaskSED0Mass::UserCreateOutputObjects() \n");
153 // Several histograms are more conveniently managed in a TList
154 fOutputPPR = new TList();
155 fOutputPPR->SetOwner();
156 fOutputPPR->SetName("listPPR");
158 fOutputmycuts = new TList();
159 fOutputmycuts->SetOwner();
160 fOutputmycuts->SetName("listloose");
162 fDistr = new TList();
164 fDistr->SetName("distributionslist");
168 TString nameMass=" ",nameSgn27=" ",nameSgn=" ", nameBkg=" ", nameRfl=" ", namedistr=" ";
170 for(Int_t i=0;i<nhist;i++){
171 nameMass="histMass_";
173 nameSgn27="histSgn27_";
182 //histograms of invariant mass distributions
184 TH1F* tmpMt = new TH1F(nameMass.Data(),"D^{0} invariant mass; M [GeV]; Entries",200,1.765,1.965);
185 TH1F *tmpMl=(TH1F*)tmpMt->Clone();
189 //to compare with AliAnalysisTaskCharmFraction
190 TH1F* tmpS27t = new TH1F(nameSgn27.Data(),"D^{0} invariant mass in M(D^{0}) +/- 27 MeV - MC; M [GeV]; Entries",200,1.765,1.965);
191 TH1F *tmpS27l=(TH1F*)tmpS27t->Clone();
195 TH1F* tmpSt = new TH1F(nameSgn.Data(), "D^{0} invariant mass - MC; M [GeV]; Entries",200,1.765,1.965);
196 TH1F *tmpSl=(TH1F*)tmpSt->Clone();
200 TH1F* tmpBt = new TH1F(nameBkg.Data(), "Background invariant mass - MC; M [GeV]; Entries",200,1.765,1.965);
201 TH1F *tmpBl=(TH1F*)tmpBt->Clone();
205 //Reflection: histo filled with D0Mass which pass the cut (also) as D0bar and with D0bar which pass (also) the cut as D0
206 TH1F* tmpRt = new TH1F(nameRfl.Data(), "Reflected signal invariant mass - MC; M [GeV]; Entries",200,1.765,1.965);
207 TH1F *tmpRl=(TH1F*)tmpRt->Clone();
210 // printf("Created histograms %s\t%s\t%s\t%s\n",tmpM->GetName(),tmpS->GetName(),tmpB->GetName(),tmpR->GetName());
212 fOutputPPR->Add(tmpMt);
213 fOutputPPR->Add(tmpSt);
214 fOutputPPR->Add(tmpS27t);
215 fOutputPPR->Add(tmpBt);
216 fOutputPPR->Add(tmpRt);
218 fOutputmycuts->Add(tmpMl);
219 fOutputmycuts->Add(tmpSl);
220 fOutputmycuts->Add(tmpS27l);
221 fOutputmycuts->Add(tmpBl);
222 fOutputmycuts->Add(tmpRl);
226 //histograms of cut variable distributions
229 TH1F *hptpiS = new TH1F(namedistr.Data(), "P_{T} distribution (pions);p_{T} [GeV/c]",200,0.,8.);
232 TH1F *hptKS = new TH1F(namedistr.Data(), "P_{T} distribution (kaons);p_{T} [GeV/c]",200,0.,8.);
235 TH1F *hptB = new TH1F(namedistr.Data(), "P_{T} distribution;p_{T} [GeV/c]",200,0.,8.);
239 TH1F *hdcaS = new TH1F(namedistr.Data(), "DCA distribution;dca [cm]",200,0.,0.1);
241 TH1F *hdcaB = new TH1F(namedistr.Data(), "DCA distribution;dca [cm]",200,0.,0.1);
244 namedistr="hcosthetastarS";
245 TH1F *hcosthetastarS = new TH1F(namedistr.Data(), "cos#theta* distribution;cos#theta*",200,-1.,1.);
246 namedistr="hcosthetastarB";
247 TH1F *hcosthetastarB = new TH1F(namedistr.Data(), "cos#theta* distribution;cos#theta*",200,-1.,1.);
251 TH1F *hd0piS = new TH1F(namedistr.Data(), "Impact parameter distribution (pions);d0(#pi) [cm]",200,-0.1,0.1);
254 TH1F *hd0KS = new TH1F(namedistr.Data(), "Impact parameter distribution (kaons);d0(K) [cm]",200,-0.1,0.1);
256 TH1F *hd0B = new TH1F(namedistr.Data(), "Impact parameter distribution;d0 [cm]",200,-0.1,0.1);
259 TH1F *hd0d0S = new TH1F(namedistr.Data(), "d_{0}#timesd_{0} distribution;d_{0}#timesd_{0} [cm^{2}]",200,-0.001,0.001);
261 TH1F *hd0d0B = new TH1F(namedistr.Data(), "d_{0}#timesd_{0} distribution;d_{0}#timesd_{0} [cm^{2}]",200,-0.001,0.001);
264 namedistr="hcosthetapointS";
265 TH1F *hcosthetapointS = new TH1F(namedistr.Data(), "cos#theta_{Point} distribution;cos#theta_{Point}",200,0,1.);
266 namedistr="hcosthetapointB";
267 TH1F *hcosthetapointB = new TH1F(namedistr.Data(), "cos#theta_{Point} distribution;cos#theta_{Point}",200,0,1.);
269 namedistr="hcosthpointd0d0S";
270 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);
271 namedistr="hcosthpointd0d0B";
272 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);
288 fDistr->Add(hcosthetastarS);
289 fDistr->Add(hcosthetastarB);
291 fDistr->Add(hcosthetapointS);
292 fDistr->Add(hcosthetapointB);
294 fDistr->Add(hcosthpointd0d0S);
295 fDistr->Add(hcosthpointd0d0B);
297 fNentries=new TH1F("nentriesD0", "nentriesD0->Integral(1,2) = number of AODs *** nentriesD0->Integral(3,4) = number of candidates selected with cuts *** nentriesD0->Integral(5,6) = number of D0 selected with cuts", 6,1.,4.);
302 //________________________________________________________________________
303 void AliAnalysisTaskSED0Mass::UserExec(Option_t */*option*/)
305 // Execute analysis for current event:
306 // heavy flavor candidates association to MC truth
307 //cout<<"I'm in UserExec"<<endl;
308 AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
311 if(fArray==0){ //D0 candidates
312 // load D0->Kpi candidates
313 //cout<<"D0 candidates"<<endl;
315 } else { //LikeSign candidates
316 // load like sign candidates
317 //cout<<"LS candidates"<<endl;
318 bname="LikeSign2Prong";
321 TClonesArray *inputArray=0;
323 if(!aod && AODEvent() && IsStandardAOD()) {
324 // In case there is an AOD handler writing a standard AOD, use the AOD
325 // event in memory rather than the input (ESD) event.
326 aod = dynamic_cast<AliAODEvent*> (AODEvent());
327 // in this case the braches in the deltaAOD (AliAOD.VertexingHF.root)
328 // have to taken from the AOD event hold by the AliAODExtension
329 AliAODHandler* aodHandler = (AliAODHandler*)
330 ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
331 if(aodHandler->GetExtensions()) {
332 AliAODExtension *ext = (AliAODExtension*)aodHandler->GetExtensions()->FindObject("AliAOD.VertexingHF.root");
333 AliAODEvent *aodFromExt = ext->GetAOD();
334 inputArray=(TClonesArray*)aodFromExt->GetList()->FindObject(bname.Data());
337 inputArray=(TClonesArray*)aod->GetList()->FindObject(bname.Data());
342 printf("AliAnalysisTaskSED0Mass::UserExec: input branch not found!\n");
346 // AOD primary vertex
347 AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
350 TClonesArray *mcArray = 0;
351 AliAODMCHeader *mcHeader = 0;
355 mcArray = (TClonesArray*)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
357 printf("AliAnalysisTaskSED0Mass::UserExec: MC particles branch not found!\n");
362 mcHeader = (AliAODMCHeader*)aod->GetList()->FindObject(AliAODMCHeader::StdBranchName());
364 printf("AliAnalysisTaskSED0Mass::UserExec: MC header branch not found!\n");
369 //printf("VERTEX Z %f %f\n",vtx1->GetZ(),mcHeader->GetVtxZ());
371 //histogram filled with 1 for every AOD
373 PostData(3,fNentries);
375 // loop over candidates
376 Int_t nInD0toKpi = inputArray->GetEntriesFast();
377 if(fDebug>1) printf("Number of D0->Kpi: %d\n",nInD0toKpi);
379 for (Int_t iD0toKpi = 0; iD0toKpi < nInD0toKpi; iD0toKpi++) {
380 //Int_t nPosPairs=0, nNegPairs=0;
381 //cout<<"inside the loop"<<endl;
382 AliAODRecoDecayHF2Prong *d = (AliAODRecoDecayHF2Prong*)inputArray->UncheckedAt(iD0toKpi);
383 Bool_t unsetvtx=kFALSE;
384 if(!d->GetOwnPrimaryVtx()) {
385 d->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
389 //check reco daughter in acceptance
390 Double_t eta0=d->EtaProng(0);
391 Double_t eta1=d->EtaProng(1);
392 Bool_t prongsinacc=kFALSE;
393 if (TMath::Abs(eta0) < 0.9 && TMath::Abs(eta1) < 0.9) prongsinacc=kTRUE;
398 Double_t mPDG=TDatabasePDG::Instance()->GetParticle(421)->Mass();
401 Double_t minvD0 = d->InvMassD0();
404 Double_t minvD0bar = d->InvMassD0bar();
405 //apply cut on invariant mass on the pair
407 if(TMath::Abs(minvD0-mPDG)<0.03 || TMath::Abs(minvD0bar-mPDG)<0.03){
408 //cout<<"inside mass cut"<<endl;
409 Int_t pdgDgD0toKpi[2]={321,211};
410 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)
413 if(fArray==1) cout<<"LS signal: ERROR"<<endl;
414 for (Int_t iprong=0; iprong<2; iprong++){
415 AliAODTrack *prong=(AliAODTrack*)d->GetDaughter(iprong);
416 Int_t labprong=prong->GetLabel();
418 //cout<<"prong name = "<<prong->GetName()<<" label = "<<prong->GetLabel()<<endl;
419 AliAODMCParticle *mcprong=0;
420 if(labprong>=0) mcprong= (AliAODMCParticle*)mcArray->At(labprong);
421 Int_t pdgprong=mcprong->GetPdgCode();
422 if(TMath::Abs(pdgprong)==211) {
424 ((TH1F*)fDistr->FindObject("hptpiS"))->Fill(d->PtProng(iprong));
425 ((TH1F*)fDistr->FindObject("hd0piS"))->Fill(d->Getd0Prong(iprong));
428 if(TMath::Abs(pdgprong)==321) {
429 //cout<<"kappa"<<endl;
430 ((TH1F*)fDistr->FindObject("hptKS"))->Fill(d->PtProng(iprong));
431 ((TH1F*)fDistr->FindObject("hd0KS"))->Fill(d->Getd0Prong(iprong));
433 ((TH1F*)fDistr->FindObject("hdcaS"))->Fill(d->GetDCA());
437 if (((AliAODMCParticle*)mcArray->At(lab))->GetPdgCode() == 421)
438 ((TH1F*)fDistr->FindObject("hcosthetastarS"))->Fill(d->CosThetaStarD0());
439 else ((TH1F*)fDistr->FindObject("hcosthetastarS"))->Fill(d->CosThetaStarD0bar());
441 ((TH1F*)fDistr->FindObject("hd0d0S"))->Fill(d->Prodd0d0());
443 ((TH1F*)fDistr->FindObject("hcosthetapointS"))->Fill(d->CosPointingAngle());
444 ((TH1F*)fDistr->FindObject("hcosthpointd0d0S"))->Fill(d->CosPointingAngle(),d->Prodd0d0());
446 //cout<<"ok point"<<endl;
449 } else{ //Background or LS
450 //cout<<"is background"<<endl;
451 AliAODTrack *prong=(AliAODTrack*)d->GetDaughter(0);
452 if(!prong) cout<<"No daughter found";
454 if(prong->Charge()==1) {fTotPosPairs[4]++;} else {fTotNegPairs[4]++;}
456 ((TH1F*)fDistr->FindObject("hptB"))->Fill(d->PtProng(0));
457 //cout<<"ptok"<<endl;
458 ((TH1F*)fDistr->FindObject("hd0B"))->Fill(d->Getd0Prong(0));
459 //cout<<"d0ok"<<endl;
460 ((TH1F*)fDistr->FindObject("hdcaB"))->Fill(d->GetDCA());
461 //cout<<"dcaok"<<endl;
462 ((TH1F*)fDistr->FindObject("hcosthetastarB"))->Fill(d->CosThetaStarD0());
463 ((TH1F*)fDistr->FindObject("hcosthetastarB"))->Fill(d->CosThetaStarD0bar());
464 ((TH1F*)fDistr->FindObject("hd0d0B"))->Fill(d->Prodd0d0());
465 //cout<<"d0d0ok"<<endl;
466 ((TH1F*)fDistr->FindObject("hcosthetapointB"))->Fill(d->CosPointingAngle());
467 ((TH1F*)fDistr->FindObject("hcosthpointd0d0B"))->Fill(d->CosPointingAngle(),d->Prodd0d0());
469 //cout<<"pointok"<<endl;
477 // printf(" |M-MD0| [GeV] < %f\n",fD0toKpiCuts[0]);
478 // printf(" dca [cm] < %f\n",fD0toKpiCuts[1]);
479 // printf(" cosThetaStar < %f\n",fD0toKpiCuts[2]);
480 // printf(" pTK [GeV/c] > %f\n",fD0toKpiCuts[3]);
481 // printf(" pTpi [GeV/c] > %f\n",fD0toKpiCuts[4]);
482 // printf(" |d0K| [cm] < %f\n",fD0toKpiCuts[5]);
483 // printf(" |d0pi| [cm] < %f\n",fD0toKpiCuts[6]);
484 // printf(" d0d0 [cm^2] < %f\n",fD0toKpiCuts[7]);
485 // printf(" cosThetaPoint > %f\n",fD0toKpiCuts[8]);
487 Double_t pt = d->Pt();
491 //cout<<"P_t = "<<pt<<endl;
492 if (pt>0. && pt<=1.) {
494 fVHFPPR->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,0.05,0.05,-0.0002,0.5);
495 fVHFmycuts->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,0.05,0.05,-0.00025,0.7);
496 // fVHFPPR->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,0.05,0.05,-0.0002,0.7);
497 // fVHFmycuts->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,1,1,-0.00015,0.5);
498 //printf("I'm in the bin %d\n",ptbin);
502 if(pt>1. && pt<=3.) {
503 if(pt>1. && pt<=2.) ptbin=2;
504 if(pt>2. && pt<=3.) ptbin=3;
505 fVHFPPR->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.0002,0.8);
506 fVHFmycuts->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,1,1,-0.00025,0.8);
507 //printf("I'm in the bin %d\n",ptbin);
512 fVHFPPR->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.0001,0.8);
513 fVHFmycuts->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.00015,0.8);
514 //printf("I'm in the bin %d\n",ptbin);
518 fVHFPPR->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.00005,0.8);
519 fVHFmycuts->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.00015,0.9);
522 //printf("I'm in the bin %d\n",ptbin);
524 //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
526 FillHists(ptbin,d,mcArray,fVHFPPR,fOutputPPR);
527 FillHists(ptbin,d,mcArray,fVHFmycuts,fOutputmycuts);
529 if(unsetvtx) d->UnsetOwnPrimaryVtx();
536 PostData(1,fOutputPPR);
537 PostData(2,fOutputmycuts);
541 //____________________________________________________________________________*
542 void AliAnalysisTaskSED0Mass::FillHists(Int_t ptbin, AliAODRecoDecayHF2Prong *part, TClonesArray *arrMC, AliAnalysisVertexingHF *vhf, TList *listout){
544 // function used in UserExec:
548 Double_t mPDG=TDatabasePDG::Instance()->GetParticle(421)->Mass();
550 Int_t okD0=0,okD0bar=0;
551 //cout<<"inside FillHist"<<endl;
552 if(part->SelectD0(vhf->GetD0toKpiCuts(),okD0,okD0bar)) {//selected
553 Double_t invmassD0 = part->InvMassD0(), invmassD0bar = part->InvMassD0bar();
554 //printf("SELECTED\n");
557 AliAODTrack *prong=(AliAODTrack*)part->GetDaughter(0);
558 if(!prong) cout<<"No daughter found";
560 if(prong->Charge()==1) {fTotPosPairs[ptbin-1]++;} else {fTotNegPairs[ptbin-1]++;}
564 Int_t pdgDgD0toKpi[2]={321,211};
566 if (fReadMC) labD0 = part->MatchToMC(421,arrMC,2,pdgDgD0toKpi); //return MC particle label if the array corresponds to a D0, -1 if not (cf. AliAODRecoDecay.cxx)
568 //count candidates selected by cuts
570 //count true D0 selected by cuts
571 if (labD0>=0) fNentries->Fill(3);
572 PostData(3,fNentries);
575 fillthis="histMass_";
577 //cout<<"Filling "<<fillthis<<endl;
579 //printf("Fill mass with D0");
580 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
583 if(fArray==1) cout<<"LS signal ERROR"<<endl;
585 AliAODMCParticle *partD0 = (AliAODMCParticle*)arrMC->At(labD0);
586 Int_t pdgD0 = partD0->GetPdgCode();
587 //cout<<"pdg = "<<pdgD0<<endl;
588 if (pdgD0==421){ //D0
589 //cout<<"Fill S with D0"<<endl;
592 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
593 if(TMath::Abs(invmassD0 - mPDG) < 0.027){
594 fillthis="histSgn27_";
596 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
598 } else{ //it was a D0bar
601 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
606 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
611 fillthis="histMass_";
613 //printf("Fill mass with D0bar");
614 ((TH1F*)listout->FindObject(fillthis))->Fill(invmassD0bar);
617 if(fArray==1) cout<<"LS signal ERROR"<<endl;
618 AliAODMCParticle *partD0 = (AliAODMCParticle*)arrMC->At(labD0);
619 Int_t pdgD0 = partD0->GetPdgCode();
620 //cout<<" pdg = "<<pdgD0<<endl;
621 if (pdgD0==-421){ //D0bar
624 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
625 if (TMath::Abs(invmassD0bar - mPDG) < 0.027){
626 fillthis="histSgn27_";
628 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
635 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
637 } else {//background or LS
640 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
645 } //else cout<<"NOT SELECTED"<<endl;
648 //________________________________________________________________________
649 void AliAnalysisTaskSED0Mass::Terminate(Option_t */*option*/)
651 // Terminate analysis
653 if(fDebug > 1) printf("AnalysisTaskSED0Mass: Terminate() \n");
655 fOutputPPR = dynamic_cast<TList*> (GetOutputData(1));
657 printf("ERROR: fOutputthight not available\n");
660 fOutputmycuts = dynamic_cast<TList*> (GetOutputData(2));
661 if (!fOutputmycuts) {
662 printf("ERROR: fOutputmycuts not available\n");
665 fDistr = dynamic_cast<TList*> (GetOutputData(4));
667 printf("ERROR: fDistr not available\n");
672 for(Int_t ipt=0;ipt<4;ipt++){
673 fLsNormalization = 2.*TMath::Sqrt(fTotPosPairs[ipt]*fTotNegPairs[ipt]);
676 if(fLsNormalization>0) {
678 TString massName="histMass_";
680 ((TH1F*)fOutputPPR->FindObject(massName))->Scale((1/fLsNormalization)*((TH1F*)fOutputPPR->FindObject(massName))->GetEntries());
681 ((TH1F*)fOutputmycuts->FindObject(massName))->Scale((1/fLsNormalization)*((TH1F*)fOutputmycuts->FindObject(massName))->GetEntries());
685 fLsNormalization = 2.*TMath::Sqrt(fTotPosPairs[4]*fTotNegPairs[4]);
687 if(fLsNormalization>0) {
689 TString nameDistr="hptB";
690 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
692 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
693 nameDistr="hcosthetastarB";
694 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
696 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
698 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
699 nameDistr="hcosthetapointB";
700 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
701 nameDistr="hcosthpointd0d0B";
702 ((TH2F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH2F*)fDistr->FindObject(nameDistr))->GetEntries());
711 } else cvname="LSinvmass";
713 TCanvas *c1=new TCanvas(cvname,cvname);
715 ((TH1F*)fOutputPPR->FindObject("histMass_4"))->Draw();