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/EvtDecayBase.hh"
17 #ifdef EVTGEN_EXTERNAL
18 #include "EvtGenExternal/EvtExternalGenList.hh"
22 #include "HepMC/GenEvent.h"
35 int main(
int argc,
char** argv) {
37 string decayFileName(
"../DECAY_2010.DEC");
38 if (argc > 1) {decayFileName = argv[1];}
39 cout<<
"Decay file name is "<<decayFileName<<endl;
41 string rootFileName(
"evtgenTest.root");
42 if (argc > 2) {rootFileName = argv[2];}
43 cout<<
"Root file name is "<<rootFileName<<endl;
45 string parentName(
"Upsilon(4S)");
46 if (argc > 3) {parentName = argv[3];}
47 cout<<
"Parent name is "<<parentName<<endl;
50 if (argc > 4) {nEvents = atoi(argv[4]);}
53 if(argc > 5) {useXml = (atoi(argv[5])==1);}
55 bool useEvtGenRandom =
false;
56 if (argc > 6) {useEvtGenRandom = (atoi(argv[6])==1);}
58 cout<<
"Number of events is "<<nEvents<<endl;
60 TFile* theFile =
new TFile(rootFileName.c_str(),
"recreate");
61 TTree* theTree =
new TTree(
"Data",
"Data");
62 TTree* nDaugTree =
new TTree(
"nDaugTree",
"nDaugTree");
63 TTree* dalitzTree =
new TTree(
"dalitzTree",
"dalitzTree");
65 theTree->SetDirectory(theFile);
66 nDaugTree->SetDirectory(theFile);
67 dalitzTree->SetDirectory(theFile);
69 int event(0), nDaug(0), daugId(0);
70 double E(0.0), p(0.0), px(0.0), py(0.0), pz(0.0);
71 double t(0.0), x(0.0), y(0.0), z(0.0);
72 double mass(0.0), lifetime(0.0);
73 double inv12(0.0), inv13(0.0), inv23(0.0);
74 double inv12Sq(0.0), inv13Sq(0.0), inv23Sq(0.0);
76 theTree->Branch(
"event", &event,
"event/I");
77 theTree->Branch(
"nDaug", &nDaug,
"nDaug/I");
78 theTree->Branch(
"id", &daugId,
"id/I");
79 theTree->Branch(
"E", &E,
"E/D");
80 theTree->Branch(
"p", &p,
"p/D");
81 theTree->Branch(
"px", &px,
"px/D");
82 theTree->Branch(
"py", &py,
"py/D");
83 theTree->Branch(
"pz", &pz,
"pz/D");
84 theTree->Branch(
"t", &t,
"t/D");
85 theTree->Branch(
"x", &x,
"x/D");
86 theTree->Branch(
"y", &y,
"x/D");
87 theTree->Branch(
"z", &z,
"x/D");
88 theTree->Branch(
"mass", &mass,
"mass/D");
89 theTree->Branch(
"lifetime", &lifetime,
"lifetime/D");
91 nDaugTree->Branch(
"event", &event,
"event/I");
92 nDaugTree->Branch(
"nDaug", &nDaug,
"nDaug/I");
94 dalitzTree->Branch(
"invMass12", &inv12,
"invMass12/D");
95 dalitzTree->Branch(
"invMass13", &inv13,
"invMass13/D");
96 dalitzTree->Branch(
"invMass23", &inv23,
"invMass23/D");
97 dalitzTree->Branch(
"invMass12Sq", &inv12Sq,
"invMass12Sq/D");
98 dalitzTree->Branch(
"invMass13Sq", &inv13Sq,
"invMass13Sq/D");
99 dalitzTree->Branch(
"invMass23Sq", &inv23Sq,
"invMass23Sq/D");
109 myRandomEngine =
new EvtMTRandomEngine();
119 std::list<EvtDecayBase*> extraModels;
121 #ifdef EVTGEN_EXTERNAL
123 radCorrEngine = genList.getPhotosModel();
124 extraModels = genList.getListOfModels();
129 EvtGen myGenerator(decayFileName.c_str(),
"../evt.pdl", myRandomEngine,
130 radCorrEngine, &extraModels, mixingType, useXml);
135 EvtId theId = EvtPDL::getId(parentName);
136 if (theId.getId() == -1 && theId.getAlias() == -1) {
137 cout<<
"Error. Could not find valid EvtId for "<<parentName<<endl;
141 static EvtId stringId = EvtPDL::getId(std::string(
"string"));
145 for (i = 0; i < nEvents; i++) {
147 if (i%1000 == 0) {cout<<
"Event number = "<<i+1<<
" out of "<<nEvents<<std::endl;}
150 EvtVector4R pInit(EvtPDL::getMass(theId), 0.0, 0.0, 0.0);
152 baseParticle = EvtParticleFactory::particleFactory(theId, pInit);
153 if (baseParticle->getSpinStates() == 3) {baseParticle->setVectorSpinDensity();}
156 myGenerator.generateDecay(baseParticle);
172 EvtId daugEvtId(-1, -1);
174 if (baseDaughter != 0) {daugEvtId = baseDaughter->
getId();}
176 if (daugEvtId == stringId) {
177 theParent = baseDaughter;
179 theParent = baseParticle;
190 for (iDaug = 0; iDaug < nDaug; iDaug++) {
200 daugId = EvtPDL::getStdHep(daug->
getId());
207 p = sqrt(px*px + py*py + pz*pz);
226 EvtVector4R p4_d1 = theParent->getDaug(0)->getP4Lab();
227 EvtVector4R p4_d2 = theParent->getDaug(1)->getP4Lab();
228 EvtVector4R p4_d3 = theParent->getDaug(2)->getP4Lab();
230 inv12Sq = (p4_d1+p4_d2)*(p4_d1+p4_d2);
231 inv13Sq = (p4_d1+p4_d3)*(p4_d1+p4_d3);
232 inv23Sq = (p4_d2+p4_d3)*(p4_d2+p4_d3);
233 inv12 = sqrt(inv12Sq);
234 inv13 = sqrt(inv13Sq);
235 inv23 = sqrt(inv23Sq);
242 baseParticle->deleteTree();
258 delete myRandomEngine;
EvtParticle * getDaug(int i)
EvtVector4R getP4Lab() const
EvtVector4R get4Pos() const