1 #if !defined(__CINT__) || defined(__MAKECINT__)
2 #include <TDatabasePDG.h>
5 #include <TLorentzVector.h>
8 #include <TClonesArray.h>
17 #include <TPythia6Decayer.h>
18 #include <TPaveStats.h>
45 gSystem->Load(
"libEGPythia6.so");
49 TPythia6Decayer* pdec=TPythia6Decayer::Instance();
50 if(fDecayTableFileName.CompareTo(
"")!=0){
51 if(fDecayTableFileName.Contains(
"ALICE_PHYSICS")){
52 gSystem->Exec(Form(
"cp %s .",fDecayTableFileName.Data()));
53 fDecayTableFileName.ReplaceAll(
"$ALICE_PHYSICS/../src/PWGHF/vertexingHF/macros/",
"./");
55 pdec->SetDecayTableFile(fDecayTableFileName.Data());
56 pdec->ReadDecayTable();
64 TString outFileName=
"Acceptance_Toy_";
70 outFileName.Append(
"D0Kpi_");
76 outFileName.Append(
"DplusKpipi_");
82 outFileName.Append(
"DStarD0pi_");
88 outFileName.Append(
"DsKKpi_");
94 outFileName.Append(
"LcpKpi_");
100 outFileName.Append(
"LcK0Sp_");
102 printf(
"ERROR: Wrong decay selected\n");
106 else outFileName.Append(
"yfidPtDep_");
109 TDatabasePDG* db=TDatabasePDG::Instance();
111 TClonesArray *array =
new TClonesArray(
"TParticle",100);
113 TH2D* hPtVsYGen=
new TH2D(
"hPtVsYGen",
"",400,0.,40.,20.,-1.,1.);
114 hPtVsYGen->GetXaxis()->SetTitle(
"p_{T} (GeV/c)");
115 hPtVsYGen->GetYaxis()->SetTitle(
"y");
116 TH2D* hPtVsYGenLimAcc=
new TH2D(
"hPtVsYGenLimAcc",
"",400,0.,40.,20.,-1.,1.);
117 hPtVsYGenLimAcc->GetXaxis()->SetTitle(
"p_{T} (GeV/c)");
118 hPtVsYGenLimAcc->GetYaxis()->SetTitle(
"y");
119 TH2D* hPtVsYGenAcc=
new TH2D(
"hPtVsYGenAcc",
"",400,0.,40.,20.,-1.,1.);
120 hPtVsYGenAcc->GetXaxis()->SetTitle(
"p_{T} (GeV/c)");
121 hPtVsYGenAcc->GetYaxis()->SetTitle(
"y");
125 funcPt=
new TF1(
"fFONLL",
"[0]*x/TMath::Power((1+TMath::Power(x/[1],[3])),[2])",0.,40.);
126 funcPt->SetParameters(0.518046,3.01138,3.38914,1.75899);
127 outFileName.Append(
"FONLL8ptshape.root");
129 funcPt=
new TF1(
"fFONLL",
"[0]*x/TMath::Power((1+TMath::Power(x/[1],[3])),[2])",0.,40.);
130 funcPt->SetParameters(0.398252, 3.9603, 3.915, 1.51853);
131 outFileName.Append(
"FONLL8ptshapeFeedDown.root");
133 funcPt=
new TF1(
"fFONLL",
"[0]*x/TMath::Power((1+TMath::Power(x/[1],[3])),[2])",0.,40.);
134 funcPt->SetParameters(0.322643,2.96275,2.30301,2.5);
135 outFileName.Append(
"FONLL7ptshape.root");
137 funcPt=
new TF1(
"fFONLL",
"[0]*x/TMath::Power((1+TMath::Power(x/[1],[3])),[2])",0.,40.);
138 funcPt->SetParameters(0.302879,2.9750,3.68139,1.68855);
139 outFileName.Append(
"FONLL5ptshape.root");
141 funcPt=
new TF1(
"fFONLL",
"[0]*x/TMath::Power((1+TMath::Power(x/[1],[3])),[2])",0.,40.);
142 funcPt->SetParameters(0.322643,1.94635,1.40463,2.5);
143 outFileName.Append(
"PYTHIA7ptshape.root");
145 funcPt=
new TF1(
"fFlat",
"pol0",0.,40.);
146 funcPt->SetParameter(0,1.);
147 outFileName.Append(
"flatpt.root");
149 TRandom3* gener=
new TRandom3(0);
152 if(itry%10000==0) printf(
"Event %d\n",itry);
153 Float_t ptD=funcPt->GetRandom();
154 Float_t phiD=gener->Rndm()*2*TMath::Pi();
155 Float_t yD=gener->Rndm()*2.-1.;
156 Float_t px=ptD*TMath::Cos(phiD);
157 Float_t py=ptD*TMath::Sin(phiD);
158 Float_t mt=TMath::Sqrt(massD*massD+ptD*ptD);
160 Float_t E=TMath::Sqrt(massD*massD+px*px+py*py+pz*pz);
162 TLorentzVector* vec=
new TLorentzVector(px,py,pz,E);
163 pdec->Decay(pdgCode,vec);
165 Int_t nentries = pdec->ImportParticles(array);
166 TParticle* dmes=(TParticle*)array->At(0);
167 Int_t nDaughters=dmes->GetNDaughters();
171 Int_t nProtonsInAcc=0;
176 Bool_t isOk=
CountPKpi(array,nentries,nPions,nKaons,nProtons,nPionsInAcc,nKaonsInAcc,nProtonsInAcc);
179 if(nPions==nPionDau && nKaons==nKaonDau && nProtons==nProtonDau){
180 hPtVsYGen->Fill(ptD,yD);
181 if(TMath::Abs(yD)<0.5){
182 hPtVsYGenLimAcc->Fill(ptD,yD);
185 if(nPionsInAcc==nPionDau && nKaonsInAcc==nKaonDau && nProtonsInAcc==nProtonDau){
186 hPtVsYGenAcc->Fill(ptD,yD);
194 TH1D* hPtGenAcc=(
TH1D*)hPtVsYGenAcc->ProjectionX(
"hPtGenAcc");
195 hPtGenAcc->GetYaxis()->SetTitle(
"Entries");
196 TH1D* hPtGenLimAcc=(
TH1D*)hPtVsYGenLimAcc->ProjectionX(
"hPtGenLimAcc");
197 hPtGenLimAcc->GetYaxis()->SetTitle(
"Entries");
199 hPtGenLimAcc->Sumw2();
200 TH1D* hAccVsPt=(
TH1D*)hPtGenAcc->Clone(
"hAccVsPt");
201 hAccVsPt->Divide(hPtGenAcc,hPtGenLimAcc,1,1,
"B");
202 hAccVsPt->GetYaxis()->SetTitle(
"Acceptance");
203 hAccVsPt->SetStats(0);
205 TCanvas* c2d=
new TCanvas(
"c2d",
"Pt vs y",1200,600);
208 hPtVsYGen->Draw(
"colz");
210 hPtVsYGenLimAcc->Draw(
"colz");
212 hPtVsYGenAcc->Draw(
"colz");
214 TCanvas* c1d=
new TCanvas(
"c1d",
"Acceptance",1200,600);
217 hPtGenLimAcc->Draw(
"");
218 Double_t ymax=1.2*TMath::Max(hPtGenLimAcc->GetMaximum(),hPtGenAcc->GetMaximum());
219 hPtGenLimAcc->SetMaximum(ymax);
221 TPaveStats *st1=(TPaveStats*)hPtGenLimAcc->GetListOfFunctions()->FindObject(
"stats");
224 hPtGenAcc->SetLineColor(kRed+1);
225 hPtGenAcc->Draw(
"sames");
227 TPaveStats *st2=(TPaveStats*)hPtGenAcc->GetListOfFunctions()->FindObject(
"stats");
230 st2->SetTextColor(hPtGenAcc->GetLineColor());
236 TFile* outfil=
new TFile(outFileName.Data(),
"recreate");
238 hPtVsYGenLimAcc->Write();
239 hPtVsYGenAcc->Write();
240 hPtGenLimAcc->Write();
257 if (TMath::Abs(y) > 0.8)
return kFALSE;
259 Double_t maxFiducialY = -0.2/15*pt*pt+1.9/15*pt+0.5;
260 Double_t minFiducialY = 0.2/15*pt*pt-1.9/15*pt-0.5;
261 if (y < minFiducialY || y > maxFiducialY)
return kFALSE;
270 TParticle* dmes=(TParticle*)array->At(0);
274 for(
int j=0; j<nentries; j++){
275 TParticle * o = (TParticle*)array->At(j);
276 Int_t pdgdau=TMath::Abs(o->GetPdgCode());
279 if(
fDebugLevel>0) printf(
"K0 dacaying into K0L\n");
282 Float_t ptdau=TMath::Sqrt(o->Px()*o->Px()+o->Py()*o->Py());
297 if(pdgdau==211) nPionsInAcc++;
298 if(pdgdau==321) nKaonsInAcc++;
302 if(TMath::Abs(sumPx-dmes->Px())>0.001 ||
303 TMath::Abs(sumPy-dmes->Py())>0.001 ||
304 TMath::Abs(sumPz-dmes->Pz())>0.001){
305 printf(
"Momentum conservation violation\n");
316 TParticle* dmes=(TParticle*)array->At(0);
321 for(
int j=0; j<nentries; j++){
322 TParticle * o = (TParticle*)array->At(j);
323 Int_t pdgdau=TMath::Abs(o->GetPdgCode());
326 if(
fDebugLevel>0) printf(
"K0 dacaying into K0L\n");
329 Float_t ptdau=TMath::Sqrt(o->Px()*o->Px()+o->Py()*o->Py());
350 if(pdgdau==211) nPionsInAcc++;
351 if(pdgdau==321) nKaonsInAcc++;
352 if(pdgdau==2212) nProtonsInAcc++;
356 if(TMath::Abs(sumPx-dmes->Px())>0.001 ||
357 TMath::Abs(sumPy-dmes->Py())>0.001 ||
358 TMath::Abs(sumPz-dmes->Pz())>0.001){
359 printf(
"Momentum conservation violation\n");
Bool_t CountKpi(TClonesArray *array, Int_t nentries, Int_t &nPions, Int_t &nKaons, Int_t &nPionsInAcc, Int_t &nKaonsInAcc)
TString fDecayTableFileName
Bool_t CountPKpi(TClonesArray *array, Int_t nentries, Int_t &nPions, Int_t &nKaons, Int_t &nProtons, Int_t &nPionsInAcc, Int_t &nKaonsInAcc, Int_t &nProtonsInAcc)
Bool_t IsInFiducialAcceptance(Double_t pt, Double_t y)