StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
genExampleRootFiles.cc
1 // Program to create ROOT files for EvtGen validation plots.
2 // This looks at the 1st generation daughters and stores 4-momenta
3 // info into a ROOT file for further analysis.
4 // Useful for Pythia, Photos and Tauola decay tests.
5 
6 #include "EvtGen/EvtGen.hh"
7 
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"
16 
17 #ifdef EVTGEN_EXTERNAL
18 #include "EvtGenExternal/EvtExternalGenList.hh"
19 #endif
20 
21 //#include "EvtGenBase/EvtHepMCEvent.hh"
22 #include "HepMC/GenEvent.h"
23 
24 #include "TFile.h"
25 #include "TTree.h"
26 
27 #include <iostream>
28 #include <string>
29 #include <vector>
30 
31 using std::cout;
32 using std::endl;
33 using std::string;
34 
35 int main(int argc, char** argv) {
36 
37  string decayFileName("../DECAY_2010.DEC");
38  if (argc > 1) {decayFileName = argv[1];}
39  cout<<"Decay file name is "<<decayFileName<<endl;
40 
41  string rootFileName("evtgenTest.root");
42  if (argc > 2) {rootFileName = argv[2];}
43  cout<<"Root file name is "<<rootFileName<<endl;
44 
45  string parentName("Upsilon(4S)");
46  if (argc > 3) {parentName = argv[3];}
47  cout<<"Parent name is "<<parentName<<endl;
48 
49  int nEvents(10);
50  if (argc > 4) {nEvents = atoi(argv[4]);}
51 
52  bool useXml = false;
53  if(argc > 5) {useXml = (atoi(argv[5])==1);}
54 
55  bool useEvtGenRandom = false;
56  if (argc > 6) {useEvtGenRandom = (atoi(argv[6])==1);}
57 
58  cout<<"Number of events is "<<nEvents<<endl;
59 
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");
64 
65  theTree->SetDirectory(theFile);
66  nDaugTree->SetDirectory(theFile);
67  dalitzTree->SetDirectory(theFile);
68 
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);
75 
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");
90 
91  nDaugTree->Branch("event", &event, "event/I");
92  nDaugTree->Branch("nDaug", &nDaug, "nDaug/I");
93 
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");
100 
101  EvtParticle* baseParticle(0);
102  EvtParticle* theParent(0);
103 
104  // Define the random number generator
105  EvtRandomEngine* myRandomEngine = 0;
106 
107 #ifdef EVTGEN_CPP11
108  // Use the Mersenne-Twister generator (C++11 only)
109  myRandomEngine = new EvtMTRandomEngine();
110 #else
111  myRandomEngine = new EvtSimpleRandomEngine();
112 #endif
113 
114  // Initialize the generator - read in the decay table and particle properties.
115  // For our validation purposes, we just want to read in one decay file and create
116  // plots from that.
117 
118  EvtAbsRadCorr* radCorrEngine = 0;
119  std::list<EvtDecayBase*> extraModels;
120 
121 #ifdef EVTGEN_EXTERNAL
122  EvtExternalGenList genList(true, "", "gamma", useEvtGenRandom);
123  radCorrEngine = genList.getPhotosModel();
124  extraModels = genList.getListOfModels();
125 #endif
126 
127  int mixingType(1);
128 
129  EvtGen myGenerator(decayFileName.c_str(), "../evt.pdl", myRandomEngine,
130  radCorrEngine, &extraModels, mixingType, useXml);
131 
132  // If I wanted a user decay file, I would read it in now, e.g:
133  // myGenerator.readUDecay(otherFileName);
134 
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;
138  return -1;
139  }
140 
141  static EvtId stringId = EvtPDL::getId(std::string("string"));
142  // Loop to create nEvents
143 
144  int i;
145  for (i = 0; i < nEvents; i++) {
146 
147  if (i%1000 == 0) {cout<<"Event number = "<<i+1<<" out of "<<nEvents<<std::endl;}
148 
149  // Set up the parent particle
150  EvtVector4R pInit(EvtPDL::getMass(theId), 0.0, 0.0, 0.0);
151 
152  baseParticle = EvtParticleFactory::particleFactory(theId, pInit);
153  if (baseParticle->getSpinStates() == 3) {baseParticle->setVectorSpinDensity();}
154 
155  // Generate the event
156  myGenerator.generateDecay(baseParticle);
157 
158  // Alternative way to generate decays and print out information:
159  // int PDGId = EvtPDL::getStdHep(theId);
160  // EvtVector4R origin(0.0, 0.0, 0.0, 0.0);
161  // EvtHepMCEvent* theEvent = myGenerator.generateDecay(PDGId, pInit, origin);
162  // HepMC::GenEvent* hepMCEvent = theEvent->getEvent();
163  // hepMCEvent->print();
164  // Extract other info from the HepMC event. Then delete it:
165  // delete theEvent;
166 
167  // Now get the particle decay information, looping through daughter tracks (1st generation only)
168 
169  // Find out if the first daughter is a string.
170  // If so, set this as the new parent particle.
171 
172  EvtId daugEvtId(-1, -1);
173  EvtParticle* baseDaughter = baseParticle->getDaug(0);
174  if (baseDaughter != 0) {daugEvtId = baseDaughter->getId();}
175 
176  if (daugEvtId == stringId) {
177  theParent = baseDaughter;
178  } else {
179  theParent = baseParticle;
180  }
181 
182  nDaug = theParent->getNDaug();
183  int iDaug(0);
184 
185  nDaugTree->Fill();
186 
187  //theParent->printTree();
188 
189  // Loop over the daughter tracks
190  for (iDaug = 0; iDaug < nDaug; iDaug++) {
191 
192  EvtParticle* daug = theParent->getDaug(iDaug);
193 
194  if (daug != 0) {
195 
196  EvtVector4R p4Lab = daug->getP4Lab();
197  EvtVector4R pos4 = daug->get4Pos();
198 
199  // PDG id
200  daugId = EvtPDL::getStdHep(daug->getId());
201 
202  // 4-momenta
203  E = p4Lab.get(0);
204  px = p4Lab.get(1);
205  py = p4Lab.get(2);
206  pz = p4Lab.get(3);
207  p = sqrt(px*px + py*py + pz*pz);
208 
209  // 4-position
210  t = pos4.get(0);
211  x = pos4.get(1);
212  y = pos4.get(2);
213  z = pos4.get(3);
214 
215  mass = daug->mass();
216  lifetime = daug->getLifetime();
217 
218  theTree->Fill();
219 
220  } // Daughter exists
221 
222  } // Number of daughters
223 
224  if (nDaug == 3) {
225 
226  EvtVector4R p4_d1 = theParent->getDaug(0)->getP4Lab();
227  EvtVector4R p4_d2 = theParent->getDaug(1)->getP4Lab();
228  EvtVector4R p4_d3 = theParent->getDaug(2)->getP4Lab();
229 
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);
236 
237  dalitzTree->Fill();
238 
239  }
240 
241  // Cleanup
242  baseParticle->deleteTree();
243 
244  }
245 
246 
247  // Write out the TTree information to the ROOT file
248 
249  theFile->cd();
250  theTree->Write();
251  nDaugTree->Write();
252  dalitzTree->Write();
253 
254  theFile->Close();
255 
256  cout<<"Done."<<endl;
257 
258  delete myRandomEngine;
259  return 0;
260 
261 }
double mass() const
EvtParticle * getDaug(int i)
EvtVector4R getP4Lab() const
size_t getNDaug() const
Definition: EvtId.hh:27
EvtVector4R get4Pos() const
Definition: EvtGen.hh:46
EvtId getId() const
double getLifetime()