6 #include "EvtGen/EvtGen.hh"
8 #include "EvtGenBase/EvtParticle.hh"
9 #include "EvtGenBase/EvtParticleFactory.hh"
10 #include "EvtGenBase/EvtPatches.hh"
11 #include "EvtGenBase/EvtPDL.hh"
12 #include "EvtGenBase/EvtRandom.hh"
13 #include "EvtGenBase/EvtSimpleRandomEngine.hh"
14 #include "EvtGenBase/EvtMTRandomEngine.hh"
15 #include "EvtGenBase/EvtHepMCEvent.hh"
16 #include "EvtGenBase/EvtSpinDensity.hh"
17 #include "EvtGenBase/EvtSpinType.hh"
18 #include "EvtGenBase/EvtComplex.hh"
19 #include "EvtGenBase/EvtCPUtil.hh"
21 #include "EvtGenBase/EvtAbsRadCorr.hh"
22 #include "EvtGenBase/EvtDecayBase.hh"
24 #include "HepMC/GenEvent.h"
25 #include "HepMC/GenVertex.h"
26 #include "HepMC/GenParticle.h"
27 #include "HepMC/SimpleVector.h"
29 #ifdef EVTGEN_EXTERNAL
30 #include "EvtGenExternal/EvtExternalGenList.hh"
52 TH1F* H_total =
new TH1F(
"Total",
"", 300, 0.0, 12.0);
53 TH1F* H_B0 =
new TH1F(
"B0",
"", 300, 0.0, 12.0);
54 TH1F* H_B0bar =
new TH1F(
"B0bar",
"", 300, 0.0, 12.0);
56 int B0Id(511), B0barId(-511);
59 bool checkSignal(std::vector<int>& daugIdVect);
61 double sineFitFun(
double* x,
double* p);
62 double timeFitFun(
double* x,
double* p);
64 int main(
int argc,
char** argv) {
66 string decayFileName(
"CPVDecayFiles/Bd_JpsiKSeeCPV.dec");
67 if (argc > 1) {decayFileName = argv[1];}
68 cout<<
"Decay file name is "<<decayFileName<<endl;
70 string rootFileName(
"rootFiles/CPVDecayTest.root");
71 if (argc > 2) {rootFileName = argv[2];}
72 cout<<
"Root file name is "<<rootFileName<<endl;
74 string parentName(
"Upsilon(4S)");
75 if (argc > 3) {parentName = argv[3];}
76 cout<<
"Parent name is "<<parentName<<endl;
79 if (argc > 4) {nEvents = atoi(argv[4]);}
81 double sin2Beta = sin(0.775);
82 if (argc > 5) {sin2Beta = atof(argv[5]);}
84 cout<<
"Number of events is "<<nEvents<<endl;
85 cout<<
"sin2Beta = "<<sin2Beta<<
" (used to draw oscillation maxima lines)"<<endl;
93 myRandomEngine =
new EvtMTRandomEngine();
99 std::list<EvtDecayBase*> extraModels;
101 #ifdef EVTGEN_EXTERNAL
103 radCorrEngine = genList.getPhotosModel();
104 extraModels = genList.getListOfModels();
107 int mixingType = EvtCPUtil::Incoherent;
110 EvtGen myGenerator(
"../DECAY_2010.DEC",
"../evt.pdl", myRandomEngine,
111 radCorrEngine, &extraModels, mixingType);
113 myGenerator.readUDecay(decayFileName.c_str());
115 EvtId theId = EvtPDL::getId(parentName);
116 if (theId.getId() == -1 && theId.getAlias() == -1) {
117 cout<<
"Error. Could not find valid EvtId for "<<parentName<<endl;
125 EvtSpinType::spintype baseSpin = EvtPDL::getSpinType(theId);
127 if (baseSpin == EvtSpinType::VECTOR) {
128 cout<<
"Setting spin density for vector particle "<<parentName<<endl;
130 spinDensity->setDiag(EvtSpinType::getSpinStates(EvtSpinType::VECTOR));
137 for (i = 0; i < nEvents; i++) {
139 if (i%1000 == 0) {cout<<
"Event number = "<<i+1<<
" out of "<<nEvents<<std::endl;}
142 int PDGId = EvtPDL::getStdHep(theId);
143 EvtVector4R pInit(EvtPDL::getMass(theId), 0.0, 0.0, 0.0);
147 EvtHepMCEvent* theEvent = myGenerator.generateDecay(PDGId, pInit, origin, spinDensity);
154 storeBFlightTimes(hepMCEvent);
165 TH1F* H_Diff =
dynamic_cast<TH1F*
>(H_B0->Clone(
"H_Diff"));
166 H_Diff->Add(H_B0bar, -1.0);
168 TH1F* H_DiffSum =
dynamic_cast<TH1F*
>(H_Diff->Clone(
"H_DiffSum"));
169 H_DiffSum->Divide(H_total);
171 TF1* sineFit =
new TF1(
"sineFit", sineFitFun, 0.0, 12.0, 3);
173 sineFit->SetParName(0,
"N");
174 sineFit->SetParName(1,
"a");
175 sineFit->SetParName(2,
"#phi");
176 sineFit->SetParameter(0, 0.5);
177 sineFit->SetParameter(1, 0.7);
178 sineFit->SetParameter(2, -0.7);
179 H_DiffSum->Fit(sineFit );
182 TF1* timeFit =
new TF1(
"timeFit", timeFitFun, 0.0, 12.0, 2);
183 timeFit->SetParName(0,
"N");
184 timeFit->SetParName(1,
"#Gamma");
185 timeFit->SetParameter(0, 500);
186 timeFit->SetParameter(1, 0.6);
187 timeFit->SetParLimits(1, 0.0, 1.0);
188 H_total->Fit(timeFit);
190 gROOT->SetStyle(
"Plain");
191 gStyle->SetOptFit(1111);
192 TCanvas* theCanvas =
new TCanvas(
"theCanvas",
"", 900, 700);
193 theCanvas->UseCurrentStyle();
195 H_DiffSum->SetXTitle(
"t (ps)");
199 TLine line1(0.0, sin2Beta, 12.0, sin2Beta); line1.Draw();
200 TLine line2(0.0, -sin2Beta, 12.0, -sin2Beta); line2.Draw();
201 theCanvas->Print(
"BCPVSinFit.gif");
203 H_total->SetXTitle(
"t (ps)");
205 theCanvas->Print(
"BTimeFit.gif");
207 TFile* theFile =
new TFile(rootFileName.c_str(),
"recreate");
219 delete myRandomEngine;
229 std::list<HepMC::GenVertex*> allVertices;
239 if (theVertex == 0) {
continue;}
246 bool gotB0(
false), gotB0bar(
false);
255 if (inParticle == 0) {
continue;}
257 int inPDGId = inParticle->
pdg_id();
258 if (inPDGId == B0Id) {
260 }
else if (inPDGId == B0barId) {
264 if (gotB0 ==
true || gotB0bar ==
true) {
270 if (gotB0 ==
true || gotB0bar ==
true) {
273 std::vector<int> daugIdVect;
280 if (outParticle != 0) {
281 int outPDGId = outParticle->
pdg_id();
282 daugIdVect.push_back(outPDGId);
288 bool gotSignal = checkSignal(daugIdVect);
291 if (gotSignal ==
true) {
294 double flightTime = calcFlightTime(BDecayVtx, B4mtm);
298 H_B0->Fill(flightTime);
299 H_total->Fill(flightTime);
303 H_B0bar->Fill(flightTime);
304 H_total->Fill(flightTime);
318 double flightTime(0.0);
320 double distance = BDecayVtx.
rho()*1e-3;
321 double momentum = B4mtm.
rho();
322 double BMass = 5.2795;
323 double c0 = 299792458.0;
325 if (momentum > 0.0) {
327 flightTime = 1.0e12*distance*BMass/(momentum*c0);
335 bool checkSignal(std::vector<int>& daugIdVect) {
337 bool gotSignal(
false);
339 int nDaug = daugIdVect.size();
344 int daug1Id = daugIdVect[0];
345 int daug2Id = daugIdVect[1];
346 if ((daug1Id == 443 && daug2Id == 310) || (daug1Id == 310 && daug2Id == 443)) {
356 double sineFitFun(
double* x,
double* p) {
363 double funcVal = N*sin(a*t + phi);
369 double timeFitFun(
double* x,
double* p) {
374 double funcVal = N0*(exp(-gamma*t));
int pdg_id() const
particle ID
std::vector< HepMC::GenParticle * >::const_iterator particles_in_const_iterator
const iterator for incoming particles
particles_out_const_iterator particles_out_const_end() const
end iteration of outgoing particles
std::vector< HepMC::GenParticle * >::const_iterator particles_out_const_iterator
const iterator for outgoing particles
GenVertex contains information about decay vertices.
double rho() const
spatial vector component magnitude
The GenEvent class is the core of HepMC.
particles_in_const_iterator particles_in_const_end() const
end iteration of incoming particles
vertex_const_iterator vertices_end() const
end vertex iteration
particles_out_const_iterator particles_out_const_begin() const
begin iteration of outgoing particles
const FourVector & position() const
vertex position and time
const FourVector & momentum() const
standard 4 momentum
FourVector is a simple representation of a physics 4 vector.
vertex_const_iterator vertices_begin() const
begin vertex iteration
non-const vertex iterator
particles_in_const_iterator particles_in_const_begin() const
begin iteration of incoming particles
The GenParticle class contains information about generated particles.