2 /**************************************************************************
3 * Copyright(c) 1998-2009, ALICE Experiment at CERN, All rights reserved. *
5 * Author: The ALICE Off-line Project. *
6 * Contributors are mentioned in the code where appropriate. *
8 * Permission to use, copy, modify and distribute this software and its *
9 * documentation strictly for non-commercial purposes is hereby granted *
10 * without fee, provided that the above copyright notice appears in all *
11 * copies and that both the copyright notice and this permission notice *
12 * appear in the supporting documentation. The authors make no claims *
13 * about the suitability of this software for any purpose. It is *
14 * provided "as is" without express or implied warranty. *
15 **************************************************************************/
17 /////////////////////////////////////////////////////////////
19 // AliAnalysisTaskSE for D0 candidates invariant mass histogram
20 // and comparison with the MC truth and cut variables distributions.
22 // Authors: A.Dainese, andrea.dainese@lnl.infn.it
23 // Chiara Bianchin, chiara.bianchin@pd.infn.it (invariant mass)
24 // Carmelo Di Giglio, carmelo.digiglio@ba.infn.it (like sign)
25 /////////////////////////////////////////////////////////////
27 #include <Riostream.h>
28 #include <TClonesArray.h>
34 #include <TDatabasePDG.h>
36 #include "AliAnalysisManager.h"
37 #include "AliAODHandler.h"
38 #include "AliAODEvent.h"
39 #include "AliAODVertex.h"
40 #include "AliAODTrack.h"
41 #include "AliAODMCHeader.h"
42 #include "AliAODMCParticle.h"
43 #include "AliAODRecoDecayHF2Prong.h"
44 #include "AliAODRecoCascadeHF.h"
45 #include "AliAnalysisVertexingHF.h"
46 #include "AliAnalysisTaskSE.h"
47 #include "AliAnalysisTaskSED0Mass.h"
50 ClassImp(AliAnalysisTaskSED0Mass)
53 //________________________________________________________________________
54 AliAnalysisTaskSED0Mass::AliAnalysisTaskSED0Mass():
67 // Default constructor
68 for(Int_t i=0;i<5;i++) {fTotPosPairs[i]=0; fTotNegPairs[i]=0;}
71 //________________________________________________________________________
72 AliAnalysisTaskSED0Mass::AliAnalysisTaskSED0Mass(const char *name):
73 AliAnalysisTaskSE(name),
84 // Default constructor
85 for(Int_t i=0;i<5;i++) {fTotPosPairs[i]=0; fTotNegPairs[i]=0;}
87 // Output slot #1 writes into a TList container
88 DefineOutput(1,TList::Class()); //My private output
89 // Output slot #2 writes into a TList container
90 DefineOutput(2,TList::Class()); //My private output
91 // Output slot #3 writes into a TH1F container
92 DefineOutput(3,TH1F::Class()); //My private output
93 // Output slot #4 writes into a TList container
94 DefineOutput(4,TList::Class()); //My private output
97 //________________________________________________________________________
98 AliAnalysisTaskSED0Mass::~AliAnalysisTaskSED0Mass()
109 delete fOutputmycuts;
130 //________________________________________________________________________
131 void AliAnalysisTaskSED0Mass::Init()
135 if(fDebug > 1) printf("AnalysisTaskSED0Mass::Init() \n");
137 gROOT->LoadMacro("$ALICE_ROOT/PWG3/vertexingHF/ConfigVertexingHF.C");
139 // 2 sets of dedidcated cuts -- defined in UserExec
140 fVHFPPR = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
141 fVHFmycuts = (AliAnalysisVertexingHF*)gROOT->ProcessLine("ConfigVertexingHF()");
146 //________________________________________________________________________
147 void AliAnalysisTaskSED0Mass::UserCreateOutputObjects()
150 // Create the output container
152 if(fDebug > 1) printf("AnalysisTaskSED0Mass::UserCreateOutputObjects() \n");
154 // Several histograms are more conveniently managed in a TList
155 fOutputPPR = new TList();
156 fOutputPPR->SetOwner();
157 fOutputPPR->SetName("listPPR");
159 fOutputmycuts = new TList();
160 fOutputmycuts->SetOwner();
161 fOutputmycuts->SetName("listloose");
163 fDistr = new TList();
165 fDistr->SetName("distributionslist");
169 TString nameMass=" ",nameSgn27=" ",nameSgn=" ", nameBkg=" ", nameRfl=" ",nameMassNocutsS =" ",nameMassNocutsB =" ", namedistr=" ";
171 for(Int_t i=0;i<nhist;i++){
172 nameMass="histMass_";
174 nameSgn27="histSgn27_";
182 nameMassNocutsS="hMassS_";
183 nameMassNocutsS+=i+1;
184 nameMassNocutsB="hMassB_";
185 nameMassNocutsB+=i+1;
187 //histograms of invariant mass distributions
189 TH1F* tmpMt = new TH1F(nameMass.Data(),"D^{0} invariant mass; M [GeV]; Entries",200,1.765,1.965);
190 TH1F *tmpMl=(TH1F*)tmpMt->Clone();
194 //to compare with AliAnalysisTaskCharmFraction
195 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);
196 TH1F *tmpS27l=(TH1F*)tmpS27t->Clone();
200 //distribution w/o cuts
201 TH1F* tmpMS = new TH1F(nameMassNocutsS.Data(),"D^{0} invariant mass; M [GeV]; Entries",200,1.765,1.965);
202 TH1F *tmpMB=(TH1F*)tmpMt->Clone();
203 tmpMB->SetName(nameMassNocutsB.Data());
207 //MC signal and background
208 TH1F* tmpSt = new TH1F(nameSgn.Data(), "D^{0} invariant mass - MC; M [GeV]; Entries",200,1.765,1.965);
209 TH1F *tmpSl=(TH1F*)tmpSt->Clone();
213 TH1F* tmpBt = new TH1F(nameBkg.Data(), "Background invariant mass - MC; M [GeV]; Entries",200,1.765,1.965);
214 TH1F *tmpBl=(TH1F*)tmpBt->Clone();
218 //Reflection: histo filled with D0Mass which pass the cut (also) as D0bar and with D0bar which pass (also) the cut as D0
219 TH1F* tmpRt = new TH1F(nameRfl.Data(), "Reflected signal invariant mass - MC; M [GeV]; Entries",200,1.765,1.965);
220 TH1F *tmpRl=(TH1F*)tmpRt->Clone();
223 // printf("Created histograms %s\t%s\t%s\t%s\n",tmpM->GetName(),tmpS->GetName(),tmpB->GetName(),tmpR->GetName());
225 fOutputPPR->Add(tmpMt);
226 fOutputPPR->Add(tmpSt);
227 fOutputPPR->Add(tmpS27t);
228 fOutputPPR->Add(tmpBt);
229 fOutputPPR->Add(tmpRt);
231 fOutputmycuts->Add(tmpMl);
232 fOutputmycuts->Add(tmpSl);
233 fOutputmycuts->Add(tmpS27l);
234 fOutputmycuts->Add(tmpBl);
235 fOutputmycuts->Add(tmpRl);
242 //histograms of cut variable distributions
245 TH1F *hptpiS = new TH1F(namedistr.Data(), "P_{T} distribution (pions);p_{T} [GeV/c]",200,0.,8.);
248 TH1F *hptKS = new TH1F(namedistr.Data(), "P_{T} distribution (kaons);p_{T} [GeV/c]",200,0.,8.);
251 TH1F *hptB = new TH1F(namedistr.Data(), "P_{T} distribution;p_{T} [GeV/c]",200,0.,8.);
255 TH1F *hdcaS = new TH1F(namedistr.Data(), "DCA distribution;dca [cm]",200,0.,0.1);
257 TH1F *hdcaB = new TH1F(namedistr.Data(), "DCA distribution;dca [cm]",200,0.,0.1);
260 namedistr="hcosthetastarS";
261 TH1F *hcosthetastarS = new TH1F(namedistr.Data(), "cos#theta* distribution;cos#theta*",200,-1.,1.);
262 namedistr="hcosthetastarB";
263 TH1F *hcosthetastarB = new TH1F(namedistr.Data(), "cos#theta* distribution;cos#theta*",200,-1.,1.);
267 TH1F *hd0piS = new TH1F(namedistr.Data(), "Impact parameter distribution (pions);d0(#pi) [cm]",200,-0.1,0.1);
270 TH1F *hd0KS = new TH1F(namedistr.Data(), "Impact parameter distribution (kaons);d0(K) [cm]",200,-0.1,0.1);
272 TH1F *hd0B = new TH1F(namedistr.Data(), "Impact parameter distribution;d0 [cm]",200,-0.1,0.1);
275 TH1F *hd0d0S = new TH1F(namedistr.Data(), "d_{0}#timesd_{0} distribution;d_{0}#timesd_{0} [cm^{2}]",200,-0.001,0.001);
277 TH1F *hd0d0B = new TH1F(namedistr.Data(), "d_{0}#timesd_{0} distribution;d_{0}#timesd_{0} [cm^{2}]",200,-0.001,0.001);
280 namedistr="hcosthetapointS";
281 TH1F *hcosthetapointS = new TH1F(namedistr.Data(), "cos#theta_{Point} distribution;cos#theta_{Point}",200,0,1.);
282 namedistr="hcosthetapointB";
283 TH1F *hcosthetapointB = new TH1F(namedistr.Data(), "cos#theta_{Point} distribution;cos#theta_{Point}",200,0,1.);
285 namedistr="hcosthpointd0d0S";
286 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);
287 namedistr="hcosthpointd0d0B";
288 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);
304 fDistr->Add(hcosthetastarS);
305 fDistr->Add(hcosthetastarB);
307 fDistr->Add(hcosthetapointS);
308 fDistr->Add(hcosthetapointB);
310 fDistr->Add(hcosthpointd0d0S);
311 fDistr->Add(hcosthpointd0d0B);
313 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.);
318 //________________________________________________________________________
319 void AliAnalysisTaskSED0Mass::UserExec(Option_t */*option*/)
321 // Execute analysis for current event:
322 // heavy flavor candidates association to MC truth
323 //cout<<"I'm in UserExec"<<endl;
324 AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
327 if(fArray==0){ //D0 candidates
328 // load D0->Kpi candidates
329 //cout<<"D0 candidates"<<endl;
331 } else { //LikeSign candidates
332 // load like sign candidates
333 //cout<<"LS candidates"<<endl;
334 bname="LikeSign2Prong";
337 TClonesArray *inputArray=0;
339 if(!aod && AODEvent() && IsStandardAOD()) {
340 // In case there is an AOD handler writing a standard AOD, use the AOD
341 // event in memory rather than the input (ESD) event.
342 aod = dynamic_cast<AliAODEvent*> (AODEvent());
343 // in this case the braches in the deltaAOD (AliAOD.VertexingHF.root)
344 // have to taken from the AOD event hold by the AliAODExtension
345 AliAODHandler* aodHandler = (AliAODHandler*)
346 ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
347 if(aodHandler->GetExtensions()) {
348 AliAODExtension *ext = (AliAODExtension*)aodHandler->GetExtensions()->FindObject("AliAOD.VertexingHF.root");
349 AliAODEvent *aodFromExt = ext->GetAOD();
350 inputArray=(TClonesArray*)aodFromExt->GetList()->FindObject(bname.Data());
353 inputArray=(TClonesArray*)aod->GetList()->FindObject(bname.Data());
358 printf("AliAnalysisTaskSED0Mass::UserExec: input branch not found!\n");
362 // AOD primary vertex
363 AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
366 TClonesArray *mcArray = 0;
367 AliAODMCHeader *mcHeader = 0;
371 mcArray = (TClonesArray*)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
373 printf("AliAnalysisTaskSED0Mass::UserExec: MC particles branch not found!\n");
378 mcHeader = (AliAODMCHeader*)aod->GetList()->FindObject(AliAODMCHeader::StdBranchName());
380 printf("AliAnalysisTaskSED0Mass::UserExec: MC header branch not found!\n");
385 //printf("VERTEX Z %f %f\n",vtx1->GetZ(),mcHeader->GetVtxZ());
387 //histogram filled with 1 for every AOD
389 PostData(3,fNentries);
391 // loop over candidates
392 Int_t nInD0toKpi = inputArray->GetEntriesFast();
393 if(fDebug>1) printf("Number of D0->Kpi: %d\n",nInD0toKpi);
395 for (Int_t iD0toKpi = 0; iD0toKpi < nInD0toKpi; iD0toKpi++) {
396 //Int_t nPosPairs=0, nNegPairs=0;
397 //cout<<"inside the loop"<<endl;
398 AliAODRecoDecayHF2Prong *d = (AliAODRecoDecayHF2Prong*)inputArray->UncheckedAt(iD0toKpi);
399 Bool_t unsetvtx=kFALSE;
400 if(!d->GetOwnPrimaryVtx()) {
401 d->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
405 //check reco daughter in acceptance
406 Double_t eta0=d->EtaProng(0);
407 Double_t eta1=d->EtaProng(1);
408 Bool_t prongsinacc=kFALSE;
409 if (TMath::Abs(eta0) < 0.9 && TMath::Abs(eta1) < 0.9) prongsinacc=kTRUE;
414 Double_t mPDG=TDatabasePDG::Instance()->GetParticle(421)->Mass();
417 Double_t minvD0 = d->InvMassD0();
420 Double_t minvD0bar = d->InvMassD0bar();
421 //cout<<"inside mass cut"<<endl;
423 Int_t pdgDgD0toKpi[2]={321,211};
425 if(fReadMC) lab=d->MatchToMC(421,mcArray,2,pdgDgD0toKpi); //return MC particle label if the array corresponds to a D0, -1 if not (cf. AliAODRecoDecay.cxx)
426 Double_t pt = d->Pt();
428 if(lab>=0 && fReadMC){ //signal
429 if (((AliAODMCParticle*)mcArray->At(lab))->GetPdgCode() == 421){//D0
430 if(TMath::Abs(minvD0-mPDG)<0.03 || TMath::Abs(minvD0bar-mPDG)<0.03){
431 ((TH1F*)fDistr->FindObject("hcosthetastarS"))->Fill(d->CosThetaStarD0());
434 //no mass cut on mass distribution
435 if(pt>0. && pt<=1.) ((TH1F*)fDistr->FindObject("hMassS_1"))->Fill(minvD0);
436 if(pt>1. && pt<=2.) ((TH1F*)fDistr->FindObject("hMassS_2"))->Fill(minvD0);
437 if(pt>2. && pt<=3.) ((TH1F*)fDistr->FindObject("hMassS_3"))->Fill(minvD0);
438 if(pt>3. && pt<=5.) ((TH1F*)fDistr->FindObject("hMassS_4"))->Fill(minvD0);
439 if(pt>5.) ((TH1F*)fDistr->FindObject("hMassS_5"))->Fill(minvD0);
443 if (((AliAODMCParticle*)mcArray->At(lab))->GetPdgCode() == 421){
444 ((TH1F*)fDistr->FindObject("hcosthetastarS"))->Fill(d->CosThetaStarD0bar());
447 //no mass cut on mass distribution
448 if(pt>0. && pt<=1.) ((TH1F*)fDistr->FindObject("hMassS_1"))->Fill(minvD0bar);
449 if(pt>1. && pt<=2.) ((TH1F*)fDistr->FindObject("hMassS_2"))->Fill(minvD0bar);
450 if(pt>2. && pt<=3.) ((TH1F*)fDistr->FindObject("hMassS_3"))->Fill(minvD0bar);
451 if(pt>3. && pt<=5.) ((TH1F*)fDistr->FindObject("hMassS_4"))->Fill(minvD0bar);
452 if(pt>5.) ((TH1F*)fDistr->FindObject("hMassS_5"))->Fill(minvD0bar);
456 //apply cut on invariant mass on the pair
457 if(TMath::Abs(minvD0-mPDG)<0.03 || TMath::Abs(minvD0bar-mPDG)<0.03){
459 if(fArray==1) cout<<"LS signal: ERROR"<<endl;
460 for (Int_t iprong=0; iprong<2; iprong++){
461 AliAODTrack *prong=(AliAODTrack*)d->GetDaughter(iprong);
462 Int_t labprong=prong->GetLabel();
464 //cout<<"prong name = "<<prong->GetName()<<" label = "<<prong->GetLabel()<<endl;
465 AliAODMCParticle *mcprong=0;
466 if(labprong>=0) mcprong= (AliAODMCParticle*)mcArray->At(labprong);
467 Int_t pdgprong=mcprong->GetPdgCode();
468 if(TMath::Abs(pdgprong)==211) {
470 ((TH1F*)fDistr->FindObject("hptpiS"))->Fill(d->PtProng(iprong));
471 ((TH1F*)fDistr->FindObject("hd0piS"))->Fill(d->Getd0Prong(iprong));
474 if(TMath::Abs(pdgprong)==321) {
475 //cout<<"kappa"<<endl;
476 ((TH1F*)fDistr->FindObject("hptKS"))->Fill(d->PtProng(iprong));
477 ((TH1F*)fDistr->FindObject("hd0KS"))->Fill(d->Getd0Prong(iprong));
479 ((TH1F*)fDistr->FindObject("hdcaS"))->Fill(d->GetDCA());
484 ((TH1F*)fDistr->FindObject("hd0d0S"))->Fill(d->Prodd0d0());
486 ((TH1F*)fDistr->FindObject("hcosthetapointS"))->Fill(d->CosPointingAngle());
487 ((TH1F*)fDistr->FindObject("hcosthpointd0d0S"))->Fill(d->CosPointingAngle(),d->Prodd0d0());
489 //cout<<"ok point"<<endl;
493 } else{ //Background or LS
494 //cout<<"is background"<<endl;
495 Double_t pt = d->Pt();
496 if(pt>0. && pt<=1.) {
497 ((TH1F*)fDistr->FindObject("hMassB_1"))->Fill(minvD0);
498 ((TH1F*)fDistr->FindObject("hMassB_1"))->Fill(minvD0bar);
500 if(pt>1. && pt<=2.) {
501 ((TH1F*)fDistr->FindObject("hMassB_2"))->Fill(minvD0);
502 ((TH1F*)fDistr->FindObject("hMassB_2"))->Fill(minvD0bar);
504 if(pt>2. && pt<=3.) {
505 ((TH1F*)fDistr->FindObject("hMassB_3"))->Fill(minvD0);
506 ((TH1F*)fDistr->FindObject("hMassB_3"))->Fill(minvD0bar);
508 if(pt>3. && pt<=5.) {
509 ((TH1F*)fDistr->FindObject("hMassB_4"))->Fill(minvD0);
510 ((TH1F*)fDistr->FindObject("hMassB_4"))->Fill(minvD0bar);
513 ((TH1F*)fDistr->FindObject("hMassB_5"))->Fill(minvD0);
514 ((TH1F*)fDistr->FindObject("hMassB_5"))->Fill(minvD0bar);
517 //apply cut on invariant mass on the pair
518 if(TMath::Abs(minvD0-mPDG)<0.03 || TMath::Abs(minvD0bar-mPDG)<0.03){
521 AliAODTrack *prong=(AliAODTrack*)d->GetDaughter(0);
522 if(!prong) cout<<"No daughter found";
524 if(prong->Charge()==1) {fTotPosPairs[4]++;} else {fTotNegPairs[4]++;}
526 ((TH1F*)fDistr->FindObject("hptB"))->Fill(d->PtProng(0));
527 //cout<<"ptok"<<endl;
528 ((TH1F*)fDistr->FindObject("hd0B"))->Fill(d->Getd0Prong(0));
529 //cout<<"d0ok"<<endl;
530 ((TH1F*)fDistr->FindObject("hdcaB"))->Fill(d->GetDCA());
531 //cout<<"dcaok"<<endl;
532 ((TH1F*)fDistr->FindObject("hcosthetastarB"))->Fill(d->CosThetaStarD0());
533 ((TH1F*)fDistr->FindObject("hcosthetastarB"))->Fill(d->CosThetaStarD0bar());
534 ((TH1F*)fDistr->FindObject("hd0d0B"))->Fill(d->Prodd0d0());
535 //cout<<"d0d0ok"<<endl;
536 ((TH1F*)fDistr->FindObject("hcosthetapointB"))->Fill(d->CosPointingAngle());
537 ((TH1F*)fDistr->FindObject("hcosthpointd0d0B"))->Fill(d->CosPointingAngle(),d->Prodd0d0());
540 }// end if inv mass cut
546 // printf(" |M-MD0| [GeV] < %f\n",fD0toKpiCuts[0]);
547 // printf(" dca [cm] < %f\n",fD0toKpiCuts[1]);
548 // printf(" cosThetaStar < %f\n",fD0toKpiCuts[2]);
549 // printf(" pTK [GeV/c] > %f\n",fD0toKpiCuts[3]);
550 // printf(" pTpi [GeV/c] > %f\n",fD0toKpiCuts[4]);
551 // printf(" |d0K| [cm] < %f\n",fD0toKpiCuts[5]);
552 // printf(" |d0pi| [cm] < %f\n",fD0toKpiCuts[6]);
553 // printf(" d0d0 [cm^2] < %f\n",fD0toKpiCuts[7]);
554 // printf(" cosThetaPoint > %f\n",fD0toKpiCuts[8]);
558 //cout<<"P_t = "<<pt<<endl;
559 if (pt>0. && pt<=1.) {
561 fVHFPPR->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,0.05,0.05,-0.0002,0.5);
562 fVHFmycuts->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,0.05,0.05,-0.00025,0.7);
563 // fVHFPPR->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,0.05,0.05,-0.0002,0.7);
564 // fVHFmycuts->SetD0toKpiCuts(0.7,0.04,0.8,0.5,0.5,1,1,-0.00015,0.5);
565 //printf("I'm in the bin %d\n",ptbin);
569 //determination of pt bin
570 if(pt>1. && pt<=3.) {
571 if(pt>1. && pt<=2.) ptbin=2;
572 if(pt>2. && pt<=3.) ptbin=3;
573 fVHFPPR->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.0002,0.8);
574 fVHFmycuts->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,1,1,-0.00025,0.8);
575 //printf("I'm in the bin %d\n",ptbin);
580 fVHFPPR->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.0001,0.8);
581 fVHFmycuts->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.00015,0.8);
582 //printf("I'm in the bin %d\n",ptbin);
586 fVHFPPR->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.00005,0.8);
587 fVHFmycuts->SetD0toKpiCuts(0.7,0.02,0.8,0.7,0.7,0.05,0.05,-0.00015,0.9);
590 //printf("I'm in the bin %d\n",ptbin);
592 //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
594 FillHists(ptbin,d,mcArray,fVHFPPR,fOutputPPR);
595 FillHists(ptbin,d,mcArray,fVHFmycuts,fOutputmycuts);
597 if(unsetvtx) d->UnsetOwnPrimaryVtx();
604 PostData(1,fOutputPPR);
605 PostData(2,fOutputmycuts);
609 //____________________________________________________________________________*
610 void AliAnalysisTaskSED0Mass::FillHists(Int_t ptbin, AliAODRecoDecayHF2Prong *part, TClonesArray *arrMC, AliAnalysisVertexingHF *vhf, TList *listout){
612 // function used in UserExec:
616 Double_t mPDG=TDatabasePDG::Instance()->GetParticle(421)->Mass();
618 Int_t okD0=0,okD0bar=0;
619 //cout<<"inside FillHist"<<endl;
620 if(part->SelectD0(vhf->GetD0toKpiCuts(),okD0,okD0bar)) {//selected
621 Double_t invmassD0 = part->InvMassD0(), invmassD0bar = part->InvMassD0bar();
622 //printf("SELECTED\n");
625 AliAODTrack *prong=(AliAODTrack*)part->GetDaughter(0);
626 if(!prong) cout<<"No daughter found";
628 if(prong->Charge()==1) {fTotPosPairs[ptbin-1]++;} else {fTotNegPairs[ptbin-1]++;}
632 Int_t pdgDgD0toKpi[2]={321,211};
634 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)
636 //count candidates selected by cuts
638 //count true D0 selected by cuts
639 if (labD0>=0) fNentries->Fill(3);
640 PostData(3,fNentries);
643 fillthis="histMass_";
645 //cout<<"Filling "<<fillthis<<endl;
647 //printf("Fill mass with D0");
648 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
651 if(fArray==1) cout<<"LS signal ERROR"<<endl;
653 AliAODMCParticle *partD0 = (AliAODMCParticle*)arrMC->At(labD0);
654 Int_t pdgD0 = partD0->GetPdgCode();
655 //cout<<"pdg = "<<pdgD0<<endl;
656 if (pdgD0==421){ //D0
657 //cout<<"Fill S with D0"<<endl;
660 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
661 if(TMath::Abs(invmassD0 - mPDG) < 0.027){
662 fillthis="histSgn27_";
664 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
666 } else{ //it was a D0bar
669 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
674 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
679 fillthis="histMass_";
681 //printf("Fill mass with D0bar");
682 ((TH1F*)listout->FindObject(fillthis))->Fill(invmassD0bar);
685 if(fArray==1) cout<<"LS signal ERROR"<<endl;
686 AliAODMCParticle *partD0 = (AliAODMCParticle*)arrMC->At(labD0);
687 Int_t pdgD0 = partD0->GetPdgCode();
688 //cout<<" pdg = "<<pdgD0<<endl;
689 if (pdgD0==-421){ //D0bar
692 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
693 if (TMath::Abs(invmassD0bar - mPDG) < 0.027){
694 fillthis="histSgn27_";
696 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
703 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
705 } else {//background or LS
708 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
713 } //else cout<<"NOT SELECTED"<<endl;
716 //________________________________________________________________________
717 void AliAnalysisTaskSED0Mass::Terminate(Option_t */*option*/)
719 // Terminate analysis
721 if(fDebug > 1) printf("AnalysisTaskSED0Mass: Terminate() \n");
723 fOutputPPR = dynamic_cast<TList*> (GetOutputData(1));
725 printf("ERROR: fOutputthight not available\n");
728 fOutputmycuts = dynamic_cast<TList*> (GetOutputData(2));
729 if (!fOutputmycuts) {
730 printf("ERROR: fOutputmycuts not available\n");
733 fDistr = dynamic_cast<TList*> (GetOutputData(4));
735 printf("ERROR: fDistr not available\n");
740 for(Int_t ipt=0;ipt<4;ipt++){
741 fLsNormalization = 2.*TMath::Sqrt(fTotPosPairs[ipt]*fTotNegPairs[ipt]);
744 if(fLsNormalization>0) {
746 TString massName="histMass_";
748 ((TH1F*)fOutputPPR->FindObject(massName))->Scale((1/fLsNormalization)*((TH1F*)fOutputPPR->FindObject(massName))->GetEntries());
749 ((TH1F*)fOutputmycuts->FindObject(massName))->Scale((1/fLsNormalization)*((TH1F*)fOutputmycuts->FindObject(massName))->GetEntries());
753 fLsNormalization = 2.*TMath::Sqrt(fTotPosPairs[4]*fTotNegPairs[4]);
755 if(fLsNormalization>0) {
757 TString nameDistr="hptB";
758 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
760 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
761 nameDistr="hcosthetastarB";
762 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
764 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
766 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
767 nameDistr="hcosthetapointB";
768 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
769 nameDistr="hcosthpointd0d0B";
770 ((TH2F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH2F*)fDistr->FindObject(nameDistr))->GetEntries());
779 } else cvname="LSinvmass";
781 TCanvas *c1=new TCanvas(cvname,cvname);
783 ((TH1F*)fOutputPPR->FindObject("histMass_4"))->Draw();