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 <AliAnalysisDataSlot.h>
36 #include <AliAnalysisDataContainer.h>
37 #include "AliAnalysisManager.h"
38 #include "AliESDtrack.h"
39 #include "AliVertexerTracks.h"
40 #include "AliAODHandler.h"
41 #include "AliAODEvent.h"
42 #include "AliAODVertex.h"
43 #include "AliAODTrack.h"
44 #include "AliAODMCHeader.h"
45 #include "AliAODMCParticle.h"
46 #include "AliAODRecoDecayHF2Prong.h"
47 #include "AliAODRecoCascadeHF.h"
48 #include "AliAnalysisVertexingHF.h"
49 #include "AliAnalysisTaskSE.h"
50 #include "AliAnalysisTaskSED0Mass.h"
51 #include "AliNormalizationCounter.h"
53 ClassImp(AliAnalysisTaskSED0Mass)
56 //________________________________________________________________________
57 AliAnalysisTaskSED0Mass::AliAnalysisTaskSED0Mass():
74 // Default constructor
77 //________________________________________________________________________
78 AliAnalysisTaskSED0Mass::AliAnalysisTaskSED0Mass(const char *name,AliRDHFCutsD0toKpi* cuts):
79 AliAnalysisTaskSE(name),
94 // Default constructor
96 fNPtBins=cuts->GetNPtBins();
97 // fTotPosPairs=new Int_t[fNPtBins];
98 // fTotNegPairs=new Int_t[fNPtBins];
99 // for(Int_t i=0;i<fNPtBins;i++) {fTotPosPairs[i]=0; fTotNegPairs[i]=0;}
103 // Output slot #1 writes into a TList container (mass with cuts)
104 DefineOutput(1,TList::Class()); //My private output
105 // Output slot #2 writes into a TList container (distributions)
106 DefineOutput(2,TList::Class()); //My private output
107 // Output slot #3 writes into a TH1F container (number of events)
108 DefineOutput(3,TH1F::Class()); //My private output
109 // Output slot #4 writes into a TList container (cuts)
110 DefineOutput(4,AliRDHFCutsD0toKpi::Class()); //My private output
111 // Output slot #5 writes Normalization Counter
112 DefineOutput(5,AliNormalizationCounter::Class());
115 //________________________________________________________________________
116 AliAnalysisTaskSED0Mass::~AliAnalysisTaskSED0Mass()
141 //________________________________________________________________________
142 void AliAnalysisTaskSED0Mass::Init()
146 if(fDebug > 1) printf("AnalysisTaskSED0Mass::Init() \n");
149 AliRDHFCutsD0toKpi* copyfCuts=new AliRDHFCutsD0toKpi(*fCuts);
150 const char* nameoutput=GetOutputSlot(4)->GetContainer()->GetName();
151 copyfCuts->SetName(nameoutput);
153 PostData(4,copyfCuts);
159 //________________________________________________________________________
160 void AliAnalysisTaskSED0Mass::UserCreateOutputObjects()
163 // Create the output container
165 if(fDebug > 1) printf("AnalysisTaskSED0Mass::UserCreateOutputObjects() \n");
167 // Several histograms are more conveniently managed in a TList
168 fOutputMass = new TList();
169 fOutputMass->SetOwner();
170 fOutputMass->SetName("listMass");
172 fDistr = new TList();
174 fDistr->SetName("distributionslist");
176 TString nameMass=" ",nameSgn27=" ",nameSgn=" ", nameBkg=" ", nameRfl=" ",nameMassNocutsS =" ",nameMassNocutsB =" ", namedistr=" ";
178 for(Int_t i=0;i<fCuts->GetNPtBins();i++){
180 nameMass="histMass_";
182 nameSgn27="histSgn27_";
190 nameMassNocutsS="hMassS_";
192 nameMassNocutsB="hMassB_";
195 //histograms of cut variable distributions
200 TH1F *hptpiS = new TH1F(namedistr.Data(), "P_{T} distribution (pions);p_{T} [GeV/c]",200,0.,8.);
204 TH1F *hptKS = new TH1F(namedistr.Data(), "P_{T} distribution (kaons);p_{T} [GeV/c]",200,0.,8.);
208 TH1F *hptB = new TH1F(namedistr.Data(), "P_{T} distribution;p_{T} [GeV/c]",200,0.,8.);
211 // namedistr="hptpiSnoMcut_";
213 // TH1F *hptpiSnoMcut = new TH1F(namedistr.Data(), "P_{T} distribution (pions);p_{T} [GeV/c]",200,0.,8.);
215 // namedistr="hptKSnoMcut_";
217 // TH1F *hptKSnoMcut = new TH1F(namedistr.Data(), "P_{T} distribution (kaons);p_{T} [GeV/c]",200,0.,8.);
219 // namedistr="hptB1prongnoMcut_";
221 // TH1F *hptB1pnoMcut = new TH1F(namedistr.Data(), "P_{T} distribution;p_{T} [GeV/c]",200,0.,8.);
223 // namedistr="hptB2prongsnoMcut_";
225 // TH1F *hptB2pnoMcut = new TH1F(namedistr.Data(), "P_{T} distribution;p_{T} [GeV/c]",200,0.,8.);
227 // fDistr->Add(hptpiSnoMcut);
228 // fDistr->Add(hptKSnoMcut);
229 // fDistr->Add(hptB1pnoMcut);
230 // fDistr->Add(hptB2pnoMcut);
236 TH1F *hdcaS = new TH1F(namedistr.Data(), "DCA distribution;dca [cm]",200,0.,0.1);
239 TH1F *hdcaB = new TH1F(namedistr.Data(), "DCA distribution;dca [cm]",200,0.,0.1);
242 namedistr="hcosthetastarS_";
244 TH1F *hcosthetastarS = new TH1F(namedistr.Data(), "cos#theta* distribution;cos#theta*",200,-1.,1.);
245 namedistr="hcosthetastarB_";
247 TH1F *hcosthetastarB = new TH1F(namedistr.Data(), "cos#theta* distribution;cos#theta*",200,-1.,1.);
252 TH1F *hd0piS = new TH1F(namedistr.Data(), "Impact parameter distribution (pions);d0(#pi) [cm]",200,-0.1,0.1);
256 TH1F *hd0KS = new TH1F(namedistr.Data(), "Impact parameter distribution (kaons);d0(K) [cm]",200,-0.1,0.1);
260 TH1F *hd0B = new TH1F(namedistr.Data(), "Impact parameter distribution (both);d0 [cm]",200,-0.1,0.1);
264 TH1F *hd0p0B = new TH1F(namedistr.Data(), "Impact parameter distribution (prong +);d0 [cm]",200,-0.1,0.1);
268 TH1F *hd0p1B = new TH1F(namedistr.Data(), "Impact parameter distribution (prong -);d0 [cm]",200,-0.1,0.1);
270 namedistr="hd0moresB_";
272 TH1F *hd0moresB = new TH1F(namedistr.Data(), "Impact parameter distribution (both);d0 [cm]",200,-0.1,0.1);
274 namedistr="hd0d0moresB_";
276 TH1F *hd0d0moresB = new TH1F(namedistr.Data(), "Impact parameter distribution (prong +);d0 [cm]",200,-0.001,0.001);
280 TH1F *hd0d0S = new TH1F(namedistr.Data(), "d_{0}#timesd_{0} distribution;d_{0}#timesd_{0} [cm^{2}]",200,-0.001,0.001);
283 TH1F *hd0d0B = new TH1F(namedistr.Data(), "d_{0}#timesd_{0} distribution;d_{0}#timesd_{0} [cm^{2}]",200,-0.001,0.001);
288 TH1F *hdeclengthS = new TH1F(namedistr.Data(), "Decay Length distribution;Decay Length [cm]",200,0,0.15);
292 TH1F *hdeclengthB = new TH1F(namedistr.Data(), "Decay Length distribution;Decay Length [cm]",200,0,0.15);
294 namedistr="hnormdeclS_";
296 TH1F *hnormdeclengthS = new TH1F(namedistr.Data(), "Normalized Decay Length distribution;Decay Length/Err ",200,0,6.);
298 namedistr="hnormdeclB_";
300 TH1F *hnormdeclengthB = new TH1F(namedistr.Data(), "Normalized Decay Length distribution;Decay Length/Err ",200,0,6.);
303 namedistr="hcosthetapointS_";
305 TH1F *hcosthetapointS = new TH1F(namedistr.Data(), "cos#theta_{Point} distribution;cos#theta_{Point}",200,0,1.);
306 namedistr="hcosthetapointB_";
308 TH1F *hcosthetapointB = new TH1F(namedistr.Data(), "cos#theta_{Point} distribution;cos#theta_{Point}",200,0,1.);
310 namedistr="hcosthetapointmoresB_";
312 TH1F *hcosthetapointmoresB = new TH1F(namedistr.Data(), "cos#theta_{Point} distribution;cos#theta_{Point}",200,0,1.);
314 // costhetapoint vs d0 or d0d0
315 namedistr="hcosthpointd0S_";
317 TH2F *hcosthpointd0S= new TH2F(namedistr.Data(),"Correlation cos#theta_{Point}-d_{0};cos#theta_{Point};d_{0} [cm^{2}]",200,0,1.,200,-0.001,0.001);
319 namedistr="hcosthpointd0B_";
321 TH2F *hcosthpointd0B= new TH2F(namedistr.Data(),"Correlation cos#theta_{Point}-d_{0};cos#theta_{Point};d_{0} [cm^{2}]",200,0,1.,200,-0.001,0.001);
323 namedistr="hcosthpointd0d0S_";
325 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);
326 namedistr="hcosthpointd0d0B_";
328 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);
342 fDistr->Add(hd0moresB);
346 fDistr->Add(hd0d0moresB);
348 fDistr->Add(hcosthetastarS);
349 fDistr->Add(hcosthetastarB);
351 fDistr->Add(hcosthetapointS);
352 fDistr->Add(hcosthetapointB);
353 fDistr->Add(hcosthetapointmoresB);
355 fDistr->Add(hdeclengthS);
356 fDistr->Add(hdeclengthB);
358 fDistr->Add(hnormdeclengthS);
359 fDistr->Add(hnormdeclengthB);
361 fDistr->Add(hcosthpointd0S);
362 fDistr->Add(hcosthpointd0B);
364 fDistr->Add(hcosthpointd0d0S);
365 fDistr->Add(hcosthpointd0d0B);
367 //histograms filled only when the secondary vertex is recalculated w/o the daughter tracks (as requested in the cut object)
369 if(fCuts->GetIsPrimaryWithoutDaughters()){
370 namedistr="hd0vmoresB_";
372 TH1F *hd0vmoresB = new TH1F(namedistr.Data(), "Impact parameter distribution (both);d0 [cm]",200,-0.1,0.1);
374 namedistr="hd0d0vmoresB_";
376 TH1F *hd0d0vmoresB = new TH1F(namedistr.Data(), "Impact parameter distribution (prong +);d0 [cm]",200,-0.001,0.001);
379 namedistr="hd0vpiS_";
381 TH1F *hd0vpiS = new TH1F(namedistr.Data(), "Impact parameter distribution (pions)(vtx w/o these tracks);d0(#pi) [cm]",200,-0.1,0.1);
385 TH1F *hd0vKS = new TH1F(namedistr.Data(), "Impact parameter distribution (kaons) (vtx w/o these tracks);d0(K) [cm]",200,-0.1,0.1);
388 TH1F *hd0vB = new TH1F(namedistr.Data(), "Impact parameter distribution (vtx w/o these tracks);d0 [cm]",200,-0.1,0.1);
390 namedistr="hd0vp0B_";
392 TH1F *hd0vp0B = new TH1F(namedistr.Data(), "Impact parameter distribution (prong + ** vtx w/o these tracks);d0 [cm]",200,-0.1,0.1);
394 namedistr="hd0vp1B_";
396 TH1F *hd0vp1B = new TH1F(namedistr.Data(), "Impact parameter distribution (prong - ** vtx w/o these tracks);d0 [cm]",200,-0.1,0.1);
399 namedistr="hd0d0vS_";
401 TH1F *hd0d0vS = new TH1F(namedistr.Data(), "d_{0}#timesd_{0} distribution (vtx w/o these tracks);d_{0}#timesd_{0} [cm^{2}]",200,-0.001,0.001);
402 namedistr="hd0d0vB_";
404 TH1F *hd0d0vB = new TH1F(namedistr.Data(), "d_{0}#timesd_{0} distribution (vtx w/o these tracks);d_{0}#timesd_{0} [cm^{2}]",200,-0.001,0.001);
407 namedistr="hdeclvS_";
409 TH1F *hdeclengthvS = new TH1F(namedistr.Data(), "Decay Length distribution (vtx w/o tracks);Decay Length [cm]",200,0,0.6);
411 namedistr="hdeclvB_";
413 TH1F *hdeclengthvB = new TH1F(namedistr.Data(), "Decay Length distribution (vtx w/o tracks);Decay Length [cm]",200,0,0.6);
415 namedistr="hnormdeclvS_";
417 TH1F *hnormdeclengthvS = new TH1F(namedistr.Data(), "Normalized Decay Length distribution (vtx w/o tracks);Decay Length/Err ",200,0,10.);
419 namedistr="hnormdeclvB_";
421 TH1F *hnormdeclengthvB = new TH1F(namedistr.Data(), "Normalized Decay Length distribution (vtx w/o tracks);Decay Length/Err ",200,0,10.);
425 fDistr->Add(hd0vpiS);
428 fDistr->Add(hd0vp0B);
429 fDistr->Add(hd0vp1B);
430 fDistr->Add(hd0vmoresB);
432 fDistr->Add(hd0d0vS);
433 fDistr->Add(hd0d0vB);
434 fDistr->Add(hd0d0vmoresB);
436 fDistr->Add(hdeclengthvS);
437 fDistr->Add(hdeclengthvB);
439 fDistr->Add(hnormdeclengthvS);
440 fDistr->Add(hnormdeclengthvB);
443 //histograms of invariant mass distributions
445 TH1F* tmpMt = new TH1F(nameMass.Data(),"D^{0} invariant mass; M [GeV]; Entries",200,1.5648,2.1648);
446 TH1F *tmpMl=(TH1F*)tmpMt->Clone();
450 //to compare with AliAnalysisTaskCharmFraction
451 // TH1F* tmpS27t = new TH1F(nameSgn27.Data(),"D^{0} invariant mass in M(D^{0}) +/- 27 MeV - MC; M [GeV]; Entries",200,1.5648,2.1648);
452 // TH1F *tmpS27l=(TH1F*)tmpS27t->Clone();
456 //distribution w/o cuts
457 // TH1F* tmpMS = new TH1F(nameMassNocutsS.Data(),"D^{0} invariant mass; M [GeV]; Entries",300,0.7,3.);
458 TH1F* tmpMS = new TH1F(nameMassNocutsS.Data(),"D^{0} invariant mass; M [GeV]; Entries",300,1.5648,2.1648); //range (MD0-300MeV, mD0 + 300MeV)
459 TH1F *tmpMB=(TH1F*)tmpMS->Clone();
460 tmpMB->SetName(nameMassNocutsB.Data());
464 //MC signal and background
465 TH1F* tmpSt = new TH1F(nameSgn.Data(), "D^{0} invariant mass - MC; M [GeV]; Entries",200,1.5648,2.1648);
466 TH1F *tmpSl=(TH1F*)tmpSt->Clone();
470 TH1F* tmpBt = new TH1F(nameBkg.Data(), "Background invariant mass - MC; M [GeV]; Entries",200,1.5648,2.1648);
471 TH1F *tmpBl=(TH1F*)tmpBt->Clone();
475 //Reflection: histo filled with D0Mass which pass the cut (also) as D0bar and with D0bar which pass (also) the cut as D0
476 TH1F* tmpRt = new TH1F(nameRfl.Data(), "Reflected signal invariant mass - MC; M [GeV]; Entries",200,1.5648,2.1648);
477 TH1F *tmpRl=(TH1F*)tmpRt->Clone();
480 // printf("Created histograms %s\t%s\t%s\t%s\n",tmpM->GetName(),tmpS->GetName(),tmpB->GetName(),tmpR->GetName());
482 fOutputMass->Add(tmpMt);
483 fOutputMass->Add(tmpSt);
484 // fOutputMass->Add(tmpS27t);
485 fOutputMass->Add(tmpBt);
486 fOutputMass->Add(tmpRt);
495 namedistr="hpospair";
496 TH1F* hpospair=new TH1F(namedistr.Data(),"Number of positive pairs",2*fCuts->GetNPtBins(),-0.5,fCuts->GetNPtBins()-0.5);
497 namedistr="hnegpair";
498 TH1F* hnegpair=new TH1F(namedistr.Data(),"Number of negative pairs",fCuts->GetNPtBins(),-0.5,2*fCuts->GetNPtBins()-0.5);
499 fDistr->Add(hpospair);
500 fDistr->Add(hnegpair);
503 const char* nameoutput=GetOutputSlot(3)->GetContainer()->GetName();
505 fNentries=new TH1F(nameoutput, "Integral(1,2) = number of AODs *** Integral(2,3) = number of candidates selected with cuts *** Integral(3,4) = number of D0 selected with cuts *** Integral(4,5) = events with good vertex *** Integral(5,6) = pt out of bounds", 14,-0.5,13.5);
507 fNentries->GetXaxis()->SetBinLabel(1,"nEventsAnal");
508 fNentries->GetXaxis()->SetBinLabel(2,"nCandSel(Cuts)");
509 fNentries->GetXaxis()->SetBinLabel(3,"nD0Selected");
510 fNentries->GetXaxis()->SetBinLabel(4,"nEventsGoodVtxS");
511 fNentries->GetXaxis()->SetBinLabel(5,"ptbin = -1");
512 fNentries->GetXaxis()->SetBinLabel(6,"no daughter");
513 fNentries->GetXaxis()->SetBinLabel(7,"nCandSel(Tr)");
514 fNentries->GetXaxis()->SetBinLabel(8,"PID=0");
515 fNentries->GetXaxis()->SetBinLabel(9,"PID=1");
516 fNentries->GetXaxis()->SetBinLabel(10,"PID=2");
517 fNentries->GetXaxis()->SetBinLabel(11,"PID=3");
518 fNentries->GetXaxis()->SetBinLabel(12,"K");
519 fNentries->GetXaxis()->SetBinLabel(13,"Lambda");
520 fNentries->GetXaxis()->SetBinLabel(14,"Pile-up Rej");
521 fNentries->GetXaxis()->SetNdivisions(1,kFALSE);
523 fCounter = new AliNormalizationCounter(Form("%s",GetOutputSlot(5)->GetContainer()->GetName()));
526 PostData(1,fOutputMass);
528 PostData(3,fNentries);
529 PostData(5,fCounter);
533 //________________________________________________________________________
534 void AliAnalysisTaskSED0Mass::UserExec(Option_t */*option*/)
536 // Execute analysis for current event:
537 // heavy flavor candidates association to MC truth
538 //cout<<"I'm in UserExec"<<endl;
542 // printf(" |M-MD0| [GeV] < %f\n",fD0toKpiCuts[0]);
543 // printf(" dca [cm] < %f\n",fD0toKpiCuts[1]);
544 // printf(" cosThetaStar < %f\n",fD0toKpiCuts[2]);
545 // printf(" pTK [GeV/c] > %f\n",fD0toKpiCuts[3]);
546 // printf(" pTpi [GeV/c] > %f\n",fD0toKpiCuts[4]);
547 // printf(" |d0K| [cm] < %f\n",fD0toKpiCuts[5]);
548 // printf(" |d0pi| [cm] < %f\n",fD0toKpiCuts[6]);
549 // printf(" d0d0 [cm^2] < %f\n",fD0toKpiCuts[7]);
550 // printf(" cosThetaPoint > %f\n",fD0toKpiCuts[8]);
553 AliAODEvent *aod = dynamic_cast<AliAODEvent*> (InputEvent());
556 if(fArray==0){ //D0 candidates
557 // load D0->Kpi candidates
558 //cout<<"D0 candidates"<<endl;
560 } else { //LikeSign candidates
561 // load like sign candidates
562 //cout<<"LS candidates"<<endl;
563 bname="LikeSign2Prong";
566 TClonesArray *inputArray=0;
568 if(!aod && AODEvent() && IsStandardAOD()) {
569 // In case there is an AOD handler writing a standard AOD, use the AOD
570 // event in memory rather than the input (ESD) event.
571 aod = dynamic_cast<AliAODEvent*> (AODEvent());
572 // in this case the braches in the deltaAOD (AliAOD.VertexingHF.root)
573 // have to taken from the AOD event hold by the AliAODExtension
574 AliAODHandler* aodHandler = (AliAODHandler*)
575 ((AliAnalysisManager::GetAnalysisManager())->GetOutputEventHandler());
577 if(aodHandler->GetExtensions()) {
578 AliAODExtension *ext = (AliAODExtension*)aodHandler->GetExtensions()->FindObject("AliAOD.VertexingHF.root");
579 AliAODEvent* aodFromExt = ext->GetAOD();
580 inputArray=(TClonesArray*)aodFromExt->GetList()->FindObject(bname.Data());
583 inputArray=(TClonesArray*)aod->GetList()->FindObject(bname.Data());
588 printf("AliAnalysisTaskSED0Mass::UserExec: input branch not found!\n");
592 // fix for temporary bug in ESDfilter
593 // the AODs with null vertex pointer didn't pass the PhysSel
594 if(!aod->GetPrimaryVertex() || TMath::Abs(aod->GetMagneticField())<0.001) return;
597 TClonesArray *mcArray = 0;
598 AliAODMCHeader *mcHeader = 0;
602 mcArray = (TClonesArray*)aod->GetList()->FindObject(AliAODMCParticle::StdBranchName());
604 printf("AliAnalysisTaskSED0Mass::UserExec: MC particles branch not found!\n");
609 mcHeader = (AliAODMCHeader*)aod->GetList()->FindObject(AliAODMCHeader::StdBranchName());
611 printf("AliAnalysisTaskSED0Mass::UserExec: MC header branch not found!\n");
616 //printf("VERTEX Z %f %f\n",vtx1->GetZ(),mcHeader->GetVtxZ());
618 //histogram filled with 1 for every AOD
620 fCounter->StoreEvent(aod,fReadMC);
621 if(!fCuts->IsEventSelected(aod)) {
622 if(fCuts->GetWhyRejection()==1) // rejected for pileup
627 // AOD primary vertex
628 AliAODVertex *vtx1 = (AliAODVertex*)aod->GetPrimaryVertex();
630 Bool_t isGoodVtx=kFALSE;
633 TString primTitle = vtx1->GetTitle();
634 if(primTitle.Contains("VertexerTracks") && vtx1->GetNContributors()>0) {
639 // loop over candidates
640 Int_t nInD0toKpi = inputArray->GetEntriesFast();
641 if(fDebug>2) printf("Number of D0->Kpi: %d\n",nInD0toKpi);
643 // FILE *f=fopen("4display.txt","a");
644 // fprintf(f,"Number of D0->Kpi: %d\n",nInD0toKpi);
645 Int_t nSelectedloose=0,nSelectedtight=0;
646 for (Int_t iD0toKpi = 0; iD0toKpi < nInD0toKpi; iD0toKpi++) {
647 //Int_t nPosPairs=0, nNegPairs=0;
648 //cout<<"inside the loop"<<endl;
649 AliAODRecoDecayHF2Prong *d = (AliAODRecoDecayHF2Prong*)inputArray->UncheckedAt(iD0toKpi);
652 if(!(d->GetDaughter(0) || d->GetDaughter(1))) {
653 AliDebug(1,"at least one daughter not found!");
658 // Bool_t unsetvtx=kFALSE;
659 // if(!d->GetOwnPrimaryVtx()) {
660 // d->SetOwnPrimaryVtx(vtx1); // needed to compute all variables
665 //check reco daughter in acceptance
667 Double_t eta0=d->EtaProng(0);
668 Double_t eta1=d->EtaProng(1);
670 if ( fCuts->IsInFiducialAcceptance(d->Pt(),d->Y(421)) ) {
671 //if( TMath::Abs(eta0)<0.9 && TMath::Abs(eta1)<0.9 ){
672 //apply cuts on tracks
674 Int_t isSelected = fCuts->IsSelected(d,AliRDHFCuts::kTracks);
675 if(((AliAODTrack*)d->GetDaughter(0))->GetTPCNcls() < 70 || ((AliAODTrack*)d->GetDaughter(1))->GetTPCNcls() < 70) isSelected=kFALSE;
676 if (!isSelected) continue;
679 if(fDebug>2) cout<<"tracks selected"<<endl;
681 Int_t ptbin=fCuts->PtBin(d->Pt());
682 if(ptbin==-1) {fNentries->Fill(4); continue;} //out of bounds
683 FillVarHists(aod,d,mcArray,fCuts,fDistr);
684 FillMassHists(aod,d,mcArray,fCuts,fOutputMass);
687 //if(unsetvtx) d->UnsetOwnPrimaryVtx();
689 fCounter->StoreCandidates(aod,nSelectedloose,kTRUE);
690 fCounter->StoreCandidates(aod,nSelectedtight,kFALSE);
692 PostData(1,fOutputMass);
694 PostData(3,fNentries);
695 PostData(5,fCounter);
699 //____________________________________________________________________________
700 void AliAnalysisTaskSED0Mass::FillVarHists(AliAODEvent* aod,AliAODRecoDecayHF2Prong *part, TClonesArray *arrMC, AliRDHFCutsD0toKpi *cuts, TList *listout){
702 // function used in UserExec to fill variable histograms:
708 Double_t mPDG=TDatabasePDG::Instance()->GetParticle(421)->Mass();
711 Double_t minvD0 = part->InvMassD0();
714 Double_t minvD0bar = part->InvMassD0bar();
715 //cout<<"inside mass cut"<<endl;
717 Int_t pdgDgD0toKpi[2]={321,211};
719 if(fReadMC) lab=part->MatchToMC(421,arrMC,2,pdgDgD0toKpi); //return MC particle label if the array corresponds to a D0, -1 if not (cf. AliAODRecoDecay.cxx)
720 //Double_t pt = d->Pt(); //mother pt
722 Int_t isSelectedPID=99;
724 //isSelectedPID = cuts->IsSelected(part,AliRDHFCuts::kPID); //0 rejected,1 D0,2 Dobar, 3 both
725 isSelectedPID = cuts->IsSelectedPID(part); //0 rejected,1 D0,2 Dobar, 3 both
726 if (isSelectedPID==0)fNentries->Fill(7);
727 if (isSelectedPID==1)fNentries->Fill(8);
728 if (isSelectedPID==2)fNentries->Fill(9);
729 if (isSelectedPID==3)fNentries->Fill(10);
730 //fNentries->Fill(8+isSelectedPID);
733 isSelected = cuts->IsSelected(part,AliRDHFCuts::kCandidate,aod); //cuts with variables recalculated with new vertex (w/o daughters)
735 //cout<<"Not Selected"<<endl;
740 Double_t invmasscut=0.03;
742 TString fillthispi="",fillthisK="",fillthis="";
744 Int_t ptbin=cuts->PtBin(part->Pt());
746 Double_t dz1[2],dz2[2],covar1[3],covar2[3];//,d0xd0proper,errd0xd0proper;
747 dz1[0]=-99; dz2[0]=-99;
749 Double_t decl[2]={-99,-99};
750 Bool_t recalcvtx=kFALSE;
752 if(fCuts->GetIsPrimaryWithoutDaughters()){
755 AliAODVertex *vtxProp=0x0;
756 vtxProp=GetPrimaryVtxSkipped(aod,part);
758 part->SetOwnPrimaryVtx(vtxProp);
759 //Bool_t unsetvtx=kTRUE;
760 //Calculate d0 for daughter tracks with Vtx Skipped
761 AliESDtrack *esdtrack1=new AliESDtrack((AliVTrack*)part->GetDaughter(0));
762 AliESDtrack *esdtrack2=new AliESDtrack((AliVTrack*)part->GetDaughter(1));
763 esdtrack1->PropagateToDCA(vtxProp,aod->GetMagneticField(),1.,dz1,covar1);
764 esdtrack2->PropagateToDCA(vtxProp,aod->GetMagneticField(),1.,dz2,covar2);
765 delete vtxProp; vtxProp=NULL;
773 decl[0]=part->DecayLength();
774 decl[1]=part->NormalizedDecayLength();
775 part->UnsetOwnPrimaryVtx();
780 if(!fCutOnDistr || (fCutOnDistr && isSelected)){ //if no cuts or cuts passed
781 //printf("\nif no cuts or cuts passed\n");
784 if(!fUsePid4Distr) isSelectedPID=0;
785 if((lab>=0 && fReadMC) || (!fReadMC && isSelectedPID)){ //signal (from MC or PID)
787 //check pdg of the prongs
788 AliAODTrack *prong0=(AliAODTrack*)part->GetDaughter(0);
789 AliAODTrack *prong1=(AliAODTrack*)part->GetDaughter(1);
790 if(!prong0 || !prong1) {
795 labprong[0]=prong0->GetLabel();
796 labprong[1]=prong1->GetLabel();
798 AliAODMCParticle *mcprong=0;
799 Int_t pdgProng[2]={0,0};
801 for (Int_t iprong=0;iprong<2;iprong++){
802 if(fReadMC && labprong[iprong]>=0) {
803 mcprong= (AliAODMCParticle*)arrMC->At(labprong[iprong]);
804 pdgProng[iprong]=mcprong->GetPdgCode();
809 //no mass cut ditributions: ptbis
811 // fillthispi="hptpiSnoMcut_";
812 // fillthispi+=ptbin;
814 // fillthisK="hptKSnoMcut_";
817 // if ((TMath::Abs(pdgProng[0]) == 211 && TMath::Abs(pdgProng[1]) == 321)
818 // || (isSelectedPID==1 || isSelectedPID==3)){
819 // ((TH1F*)listout->FindObject(fillthispi))->Fill(part->PtProng(0));
820 // ((TH1F*)listout->FindObject(fillthisK))->Fill(part->PtProng(1));
823 // if ((TMath::Abs(pdgProng[0]) == 321 && TMath::Abs(pdgProng[1]) == 211)
824 // || (isSelectedPID==2 || isSelectedPID==3)){
825 // ((TH1F*)listout->FindObject(fillthisK))->Fill(part->PtProng(0));
826 // ((TH1F*)listout->FindObject(fillthispi))->Fill(part->PtProng(1));
829 //no mass cut ditributions: mass
833 if ((fReadMC && ((AliAODMCParticle*)arrMC->At(lab))->GetPdgCode() == 421)
834 || (!fReadMC && (isSelectedPID==1 || isSelectedPID==3))){//D0
835 ((TH1F*)listout->FindObject(fillthis))->Fill(minvD0);
838 if(fReadMC || (!fReadMC && isSelectedPID > 1))
839 ((TH1F*)listout->FindObject(fillthis))->Fill(minvD0bar);
842 //apply cut on invariant mass on the pair
843 if(TMath::Abs(minvD0-mPDG)<invmasscut || TMath::Abs(minvD0bar-mPDG)<invmasscut){
845 if(fArray==1) cout<<"LS signal: ERROR"<<endl;
846 for (Int_t iprong=0; iprong<2; iprong++){
847 AliAODTrack *prong=(AliAODTrack*)part->GetDaughter(iprong);
848 if (fReadMC) labprong[iprong]=prong->GetLabel();
850 //cout<<"prong name = "<<prong->GetName()<<" label = "<<prong->GetLabel()<<endl;
852 if(fReadMC && labprong[iprong]>=0) {
853 mcprong= (AliAODMCParticle*)arrMC->At(labprong[iprong]);
854 pdgprong=mcprong->GetPdgCode();
857 Bool_t isPionHere[2]={(isSelectedPID==1 || isSelectedPID==3) ? kTRUE : kFALSE,(isSelectedPID==2 || isSelectedPID==3) ? kTRUE : kFALSE};
859 if(TMath::Abs(pdgprong)==211 || isPionHere[iprong]) {
862 fillthispi="hptpiS_";
864 ((TH1F*)listout->FindObject(fillthispi))->Fill(part->PtProng(iprong));
865 fillthispi="hd0piS_";
867 ((TH1F*)listout->FindObject(fillthispi))->Fill(part->Getd0Prong(iprong));
870 fillthispi="hd0vpiS_";
872 ((TH1F*)listout->FindObject(fillthispi))->Fill(d0[iprong]);
877 if(TMath::Abs(pdgprong)==321 || !isPionHere[iprong]) {
878 //cout<<"kappa"<<endl;
882 ((TH1F*)listout->FindObject(fillthisK))->Fill(part->PtProng(iprong));
885 ((TH1F*)listout->FindObject(fillthisK))->Fill(part->Getd0Prong(iprong));
889 ((TH1F*)listout->FindObject(fillthisK))->Fill(d0[iprong]);
893 fillthis="hcosthpointd0S_";
895 ((TH1F*)listout->FindObject(fillthis))->Fill(part->CosPointingAngle(),part->Getd0Prong(iprong));
897 } //end loop on prongs
901 ((TH1F*)listout->FindObject(fillthis))->Fill(part->GetDCA());
903 fillthis="hcosthetapointS_";
905 ((TH1F*)listout->FindObject(fillthis))->Fill(part->CosPointingAngle());
907 fillthis="hcosthpointd0d0S_";
909 ((TH2F*)listout->FindObject(fillthis))->Fill(part->CosPointingAngle(),part->Prodd0d0());
911 fillthis="hcosthetastarS_";
913 if ((fReadMC && ((AliAODMCParticle*)arrMC->At(lab))->GetPdgCode() == 421)) ((TH1F*)listout->FindObject(fillthis))->Fill(part->CosThetaStarD0());
915 if (fReadMC || isSelectedPID>1)((TH1F*)listout->FindObject(fillthis))->Fill(part->CosThetaStarD0bar());
916 if(isSelectedPID==1 || isSelectedPID==3)((TH1F*)listout->FindObject(fillthis))->Fill(part->CosThetaStarD0());
920 ((TH1F*)listout->FindObject(fillthis))->Fill(part->Prodd0d0());
925 ((TH1F*)listout->FindObject(fillthis))->Fill(d0[0]*d0[1]);
930 ((TH1F*)listout->FindObject(fillthis))->Fill(part->DecayLength());
932 fillthis="hnormdeclS_";
934 ((TH1F*)listout->FindObject(fillthis))->Fill(part->NormalizedDecayLength());
939 ((TH1F*)listout->FindObject(fillthis))->Fill(decl[0]);
941 fillthis="hnormdeclvS_";
943 ((TH1F*)listout->FindObject(fillthis))->Fill(decl[1]);
947 } else{ //Background or LS
949 //cout<<"is background"<<endl;
951 //no mass cut distributions: mass, ptbis
955 if (!fCutOnDistr || (fCutOnDistr && (isSelected==1 || isSelected==3))) ((TH1F*)listout->FindObject(fillthis))->Fill(minvD0);
956 if (!fCutOnDistr || (fCutOnDistr && isSelected>1)) ((TH1F*)listout->FindObject(fillthis))->Fill(minvD0bar);
958 // fillthis="hptB1prongnoMcut_";
961 // ((TH1F*)listout->FindObject(fillthis))->Fill(part->PtProng(0));
963 // fillthis="hptB2prongsnoMcut_";
965 // ((TH1F*)listout->FindObject(fillthis))->Fill(part->PtProng(0));
966 // ((TH1F*)listout->FindObject(fillthis))->Fill(part->PtProng(1));
968 //apply cut on invariant mass on the pair
969 if(TMath::Abs(minvD0-mPDG)<invmasscut || TMath::Abs(minvD0bar-mPDG)<invmasscut){
972 AliAODTrack *prongg=(AliAODTrack*)part->GetDaughter(0);
974 if(fDebug>2) cout<<"No daughter found";
979 if(prongg->Charge()==1) {
980 //fTotPosPairs[ptbin]++;
981 ((TH1F*)fDistr->FindObject("hpospair"))->Fill(ptbin);
983 //fTotNegPairs[ptbin]++;
984 ((TH1F*)fDistr->FindObject("hnegpair"))->Fill(ptbin);
989 //normalise pt distr to half afterwards
992 ((TH1F*)listout->FindObject(fillthis))->Fill(part->PtProng(0));
993 ((TH1F*)listout->FindObject(fillthis))->Fill(part->PtProng(1));
997 ((TH1F*)listout->FindObject(fillthis))->Fill(part->Getd0Prong(0));
1000 ((TH1F*)listout->FindObject(fillthis))->Fill(part->Getd0Prong(1));
1003 ((TH1F*)listout->FindObject(fillthis))->Fill(part->Getd0Prong(0));
1004 ((TH1F*)listout->FindObject(fillthis))->Fill(part->Getd0Prong(1));
1007 Int_t pdgMother[2]={0,0};
1008 Double_t factor[2]={1,1};
1010 for(Int_t iprong=0;iprong<2;iprong++){
1011 AliAODTrack *prong=(AliAODTrack*)part->GetDaughter(iprong);
1012 lab=prong->GetLabel();
1014 AliAODMCParticle* mcprong=(AliAODMCParticle*)arrMC->At(lab);
1016 Int_t labmom=mcprong->GetMother();
1018 AliAODMCParticle* mcmother=(AliAODMCParticle*)arrMC->At(labmom);
1019 if(mcmother) pdgMother[iprong]=mcmother->GetPdgCode();
1024 fillthis="hd0moresB_";
1027 if(TMath::Abs(pdgMother[iprong])==310 || TMath::Abs(pdgMother[iprong])==130 || TMath::Abs(pdgMother[iprong])==321){ //K^0_S, K^0_L, K^+-
1028 if(part->PtProng(iprong)<=1)factor[iprong]=1./.7;
1029 else factor[iprong]=1./.6;
1030 fNentries->Fill(11);
1033 if(TMath::Abs(pdgMother[iprong])==3122) { //Lambda
1034 factor[iprong]=1./0.25;
1035 fNentries->Fill(12);
1037 fillthis="hd0moresB_";
1040 ((TH1F*)listout->FindObject(fillthis))->Fill(part->Getd0Prong(iprong),factor[iprong]);
1043 fillthis="hd0vmoresB_";
1045 ((TH1F*)listout->FindObject(fillthis))->Fill(d0[iprong],factor[iprong]);
1050 fillthis="hd0d0moresB_";
1052 ((TH1F*)listout->FindObject(fillthis))->Fill(part->Prodd0d0(),factor[0]*factor[1]);
1055 fillthis="hd0d0vmoresB_";
1057 ((TH1F*)listout->FindObject(fillthis))->Fill(d0[0]*d0[1],factor[0]*factor[1]);
1059 fillthis="hcosthetapointmoresB_";
1061 ((TH1F*)listout->FindObject(fillthis))->Fill(part->CosPointingAngle(),factor[0]*factor[1]);
1065 fillthis="hd0vp0B_";
1067 ((TH1F*)listout->FindObject(fillthis))->Fill(d0[0]);
1068 fillthis="hd0vp1B_";
1070 ((TH1F*)listout->FindObject(fillthis))->Fill(d0[1]);
1074 ((TH1F*)listout->FindObject(fillthis))->Fill(d0[0]);
1075 ((TH1F*)listout->FindObject(fillthis))->Fill(d0[1]);
1078 fillthis="hcosthpointd0B_";
1080 ((TH1F*)listout->FindObject(fillthis))->Fill(part->CosPointingAngle(),part->Getd0Prong(0));
1081 ((TH1F*)listout->FindObject(fillthis))->Fill(part->CosPointingAngle(),part->Getd0Prong(1));
1086 ((TH1F*)listout->FindObject(fillthis))->Fill(part->GetDCA());
1088 fillthis="hcosthetastarB_";
1090 if (!fCutOnDistr || (fCutOnDistr && (isSelected==1 || isSelected==3)))((TH1F*)listout->FindObject(fillthis))->Fill(part->CosThetaStarD0());
1091 if (!fCutOnDistr || (fCutOnDistr && isSelected>1))((TH1F*)listout->FindObject(fillthis))->Fill(part->CosThetaStarD0bar());
1095 ((TH1F*)listout->FindObject(fillthis))->Fill(part->Prodd0d0());
1098 fillthis="hd0d0vB_";
1100 ((TH1F*)listout->FindObject(fillthis))->Fill(d0[0]*d0[1]);
1103 fillthis="hcosthetapointB_";
1105 ((TH1F*)listout->FindObject(fillthis))->Fill(part->CosPointingAngle());
1107 fillthis="hcosthpointd0d0B_";
1109 ((TH2F*)listout->FindObject(fillthis))->Fill(part->CosPointingAngle(),part->Prodd0d0());
1113 ((TH1F*)listout->FindObject(fillthis))->Fill(part->DecayLength());
1115 fillthis="hnormdeclB_";
1117 ((TH1F*)listout->FindObject(fillthis))->Fill(part->NormalizedDecayLength());
1121 fillthis="hdeclvB_";
1123 ((TH1F*)listout->FindObject(fillthis))->Fill(decl[0]);
1125 fillthis="hnormdeclvB_";
1127 ((TH1F*)listout->FindObject(fillthis))->Fill(decl[1]);
1131 }//else (background)
1135 //____________________________________________________________________________
1137 void AliAnalysisTaskSED0Mass::FillMassHists(AliAODEvent* aod,AliAODRecoDecayHF2Prong *part, TClonesArray *arrMC, AliRDHFCutsD0toKpi* cuts, TList *listout){
1139 // function used in UserExec to fill mass histograms:
1143 // Double_t mPDG=TDatabasePDG::Instance()->GetParticle(421)->Mass();
1145 Int_t isSelected=cuts->IsSelected(part,AliRDHFCuts::kCandidate,aod); //selected
1147 //cout<<"is selected = "<<isSelected<<endl;
1149 //cout<<"check cuts = "<<endl;
1152 //cout<<"Not Selected"<<endl;
1153 //cout<<"Rejected because "<<cuts->GetWhy()<<endl;
1158 if(fDebug>2) cout<<"Candidate selected"<<endl;
1160 Double_t invmassD0 = part->InvMassD0(), invmassD0bar = part->InvMassD0bar();
1161 //printf("SELECTED\n");
1162 Int_t ptbin=cuts->PtBin(part->Pt());
1164 AliAODTrack *prong=(AliAODTrack*)part->GetDaughter(0);
1166 AliDebug(2,"No daughter found");
1170 // if(prong->Charge()==1) {
1171 // ((TH1F*)fDistr->FindObject("hpospair"))->Fill(fCuts->GetNPtBins()+ptbin);
1172 // //fTotPosPairs[ptbin]++;
1174 // ((TH1F*)fDistr->FindObject("hnegpair"))->Fill(fCuts->GetNPtBins()+ptbin);
1175 // //fTotNegPairs[ptbin]++;
1179 // for(Int_t it=0;it<2;it++){
1181 // //request on spd points to be addes
1182 // if(/*nSPD==2 && */part->Pt() > 5. && (TMath::Abs(invmassD0-mPDG)<0.01 || TMath::Abs(invmassD0bar-mPDG)<0.01)){
1183 // FILE *f=fopen("4display.txt","a");
1184 // fprintf(f,"pt: %f \n Rapidity: %f \t Period Number: %x \t Run Number: %d \t BunchCrossNumb: %d \t OrbitNumber: %d\n",part->Pt(),part->Y(421),aod->GetPeriodNumber(),aod->GetRunNumber(),aod->GetBunchCrossNumber(),aod->GetOrbitNumber());
1186 // //printf("PrimVtx NContributors: %d \n Prongs Rel Angle: %f \n \n",ncont,relangle);
1190 TString fillthis="";
1191 Int_t pdgDgD0toKpi[2]={321,211};
1193 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)
1195 //count candidates selected by cuts
1197 //count true D0 selected by cuts
1198 if (fReadMC && labD0>=0) fNentries->Fill(2);
1200 if ((isSelected==1 || isSelected==3) && fFillOnlyD0D0bar<2) { //D0
1201 fillthis="histMass_";
1203 //cout<<"Filling "<<fillthis<<endl;
1205 //printf("Fill mass with D0");
1206 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
1209 if(fArray==1) cout<<"LS signal ERROR"<<endl;
1211 AliAODMCParticle *partD0 = (AliAODMCParticle*)arrMC->At(labD0);
1212 Int_t pdgD0 = partD0->GetPdgCode();
1213 //cout<<"pdg = "<<pdgD0<<endl;
1214 if (pdgD0==421){ //D0
1215 //cout<<"Fill S with D0"<<endl;
1216 fillthis="histSgn_";
1218 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
1219 // if(TMath::Abs(invmassD0 - mPDG) < 0.027){
1220 // fillthis="histSgn27_";
1222 // ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
1224 } else{ //it was a D0bar
1225 fillthis="histRfl_";
1227 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
1229 } else {//background
1230 fillthis="histBkg_";
1232 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0);
1236 if (isSelected>1 && (fFillOnlyD0D0bar==0 || fFillOnlyD0D0bar==2)) { //D0bar
1237 fillthis="histMass_";
1239 //printf("Fill mass with D0bar");
1240 ((TH1F*)listout->FindObject(fillthis))->Fill(invmassD0bar);
1243 if(fArray==1) cout<<"LS signal ERROR"<<endl;
1244 AliAODMCParticle *partD0 = (AliAODMCParticle*)arrMC->At(labD0);
1245 Int_t pdgD0 = partD0->GetPdgCode();
1246 //cout<<" pdg = "<<pdgD0<<endl;
1247 if (pdgD0==-421){ //D0bar
1248 fillthis="histSgn_";
1250 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
1251 // if (TMath::Abs(invmassD0bar - mPDG) < 0.027){
1252 // fillthis="histSgn27_";
1254 // ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
1259 fillthis="histRfl_";
1261 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
1263 } else {//background or LS
1264 fillthis="histBkg_";
1266 ((TH1F*)(listout->FindObject(fillthis)))->Fill(invmassD0bar);
1273 //__________________________________________________________________________
1274 AliAODVertex* AliAnalysisTaskSED0Mass::GetPrimaryVtxSkipped(AliAODEvent *aodev,AliAODRecoDecayHF2Prong *d){
1275 //Calculate the primary vertex w/o the daughter tracks of the candidate
1278 Int_t nTrksToSkip=2;
1279 AliAODTrack *dgTrack = (AliAODTrack*)d->GetDaughter(0);
1281 AliDebug(2,"no daughter found!");
1284 skipped[0]=dgTrack->GetID();
1285 dgTrack = (AliAODTrack*)d->GetDaughter(1);
1287 AliDebug(2,"no daughter found!");
1290 skipped[1]=dgTrack->GetID();
1292 AliESDVertex *vertexESD=0x0;
1293 AliAODVertex *vertexAOD=0x0;
1294 AliVertexerTracks *vertexer = new AliVertexerTracks(aodev->GetMagneticField());
1297 vertexer->SetSkipTracks(nTrksToSkip,skipped);
1298 vertexer->SetMinClusters(4);
1299 vertexESD = (AliESDVertex*)vertexer->FindPrimaryVertex(aodev);
1300 if(!vertexESD) return vertexAOD;
1301 if(vertexESD->GetNContributors()<=0) {
1302 AliDebug(2,"vertexing failed");
1303 delete vertexESD; vertexESD=NULL;
1307 delete vertexer; vertexer=NULL;
1310 // convert to AliAODVertex
1311 Double_t pos[3],cov[6],chi2perNDF;
1312 vertexESD->GetXYZ(pos); // position
1313 vertexESD->GetCovMatrix(cov); //covariance matrix
1314 chi2perNDF = vertexESD->GetChi2toNDF();
1315 delete vertexESD; vertexESD=NULL;
1317 vertexAOD = new AliAODVertex(pos,cov,chi2perNDF);
1323 //________________________________________________________________________
1324 void AliAnalysisTaskSED0Mass::Terminate(Option_t */*option*/)
1326 // Terminate analysis
1328 if(fDebug > 1) printf("AnalysisTaskSED0Mass: Terminate() \n");
1331 fOutputMass = dynamic_cast<TList*> (GetOutputData(1));
1333 printf("ERROR: fOutputMass not available\n");
1336 fDistr = dynamic_cast<TList*> (GetOutputData(2));
1338 printf("ERROR: fDistr not available\n");
1342 fNentries = dynamic_cast<TH1F*>(GetOutputData(3));
1345 printf("ERROR: fNEntries not available\n");
1348 fCuts = dynamic_cast<AliRDHFCutsD0toKpi*>(GetOutputData(4));
1350 printf("ERROR: fCuts not available\n");
1353 fCounter = dynamic_cast<AliNormalizationCounter*>(GetOutputData(5));
1355 printf("ERROR: fCounter not available\n");
1358 Int_t nptbins=fCuts->GetNPtBins();
1359 for(Int_t ipt=0;ipt<nptbins;ipt++){
1362 fLsNormalization = 2.*TMath::Sqrt(((TH1F*)fDistr->FindObject("hpospair"))->Integral(nptbins+ipt+1,nptbins+ipt+2)*((TH1F*)fDistr->FindObject("hnegpair"))->Integral(nptbins+ipt+1,nptbins+ipt+2)); //after cuts
1365 if(fLsNormalization>1e-6) {
1367 TString massName="histMass_";
1369 ((TH1F*)fOutputMass->FindObject(massName))->Scale((1/fLsNormalization)*((TH1F*)fOutputMass->FindObject(massName))->GetEntries());
1374 fLsNormalization = 2.*TMath::Sqrt(((TH1F*)fDistr->FindObject("hpospair"))->Integral(ipt+1,ipt+2)*((TH1F*)fDistr->FindObject("hnegpair"))->Integral(ipt+1,ipt+2));
1375 //fLsNormalization = 2.*TMath::Sqrt(fTotPosPairs[4]*fTotNegPairs[4]);
1377 if(fLsNormalization>1e-6) {
1379 TString nameDistr="hptB_";
1381 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
1384 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
1385 nameDistr="hcosthetastarB_";
1387 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
1390 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
1391 nameDistr="hd0d0B_";
1393 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
1394 nameDistr="hcosthetapointB_";
1396 ((TH1F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH1F*)fDistr->FindObject(nameDistr))->GetEntries());
1397 nameDistr="hcosthpointd0d0B_";
1399 ((TH2F*)fDistr->FindObject(nameDistr))->Scale((1/fLsNormalization)*((TH2F*)fDistr->FindObject(nameDistr))->GetEntries());
1404 TString cvname,cstname;
1414 TCanvas *cMass=new TCanvas(cvname,cvname);
1416 ((TH1F*)fOutputMass->FindObject("histMass_3"))->Draw();
1418 TCanvas* cStat=new TCanvas(cstname,Form("Stat%s",fArray ? "LS" : "D0"));
1421 fNentries->Draw("htext0");
1423 // TCanvas *ccheck=new TCanvas(Form("cc%d",fArray),Form("cc%d",fArray));