StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
evt_dalitz.cc
1 /*
2  * Generate the kinematics of a Dalitz decay
3  * according to a decay file.
4  */
5 
6 #include <stdio.h>
7 #include <unistd.h>
8 #include <fstream>
9 #include "EvtGenBase/EvtSimpleRandomEngine.hh"
10 #include "EvtGenBase/EvtMTRandomEngine.hh"
11 #include "EvtGenBase/EvtRandom.hh"
12 #include "EvtGenBase/EvtVector4R.hh"
13 #include "EvtGenBase/EvtParticle.hh"
14 #include "EvtGenBase/EvtParticleFactory.hh"
15 #include "EvtGenBase/EvtPDL.hh"
16 #include "EvtGen/EvtGen.hh"
17 #include "EvtGenBase/EvtAbsRadCorr.hh"
18 #include "EvtGenBase/EvtDecayBase.hh"
19 
20 #ifdef EVTGEN_EXTERNAL
21 #include "EvtGenExternal/EvtExternalGenList.hh"
22 #endif
23 
24 #include <cstdlib>
25 #include <list>
26 #include <string>
27 
28 using std::endl;
29 using std::fstream;
30 using std::ifstream;
31 using std::ofstream;
32 
33 void usage(const char *cmd)
34 {
35  printf("Usage: %s -d <decayfile> -m <mother> -n <nevents> -f <datfile> -r <rootfile>\n",cmd);
36 }
37 
38 
39 int main(int argc, char* argv[]) {
40 
41  std::string mother("B+");
42  int N = 10;
43  std::string datfile("output.dat");
44  std::string decfile =("../DECAY_2010.DEC");
45 
46  // Get options
47 
48  extern char* optarg;
49  int c;
50  while ( (c = getopt(argc, argv, "d:m:n:f:")) != EOF ) {
51  switch (c) {
52  case 'm': mother = std::string(optarg); break;
53  case 'n': N = atoi(optarg); break;
54  case 'f': datfile = std::string(optarg); break;
55  case 'd': decfile = std::string(optarg); break;
56  }
57  }
58 
59  // Initialize EvtGen
60 
61  // Define the random number generator
62  EvtRandomEngine* eng = 0;
63 
64 #ifdef EVTGEN_CPP11
65  // Use the Mersenne-Twister generator (C++11 only)
66  eng = new EvtMTRandomEngine();
67 #else
68  eng = new EvtSimpleRandomEngine();
69 #endif
70 
71  EvtRandom::setRandomEngine(eng);
72 
73  EvtAbsRadCorr* radCorrEngine = 0;
74  std::list<EvtDecayBase*> extraModels;
75 
76 #ifdef EVTGEN_EXTERNAL
77  EvtExternalGenList genList;
78  radCorrEngine = genList.getPhotosModel();
79  extraModels = genList.getListOfModels();
80 #endif
81 
82  EvtGen myGenerator(decfile.c_str(),"../evt.pdl",eng,
83  radCorrEngine, &extraModels);
84 
85  // Initialize decay
86 
87  EvtId id = EvtPDL::getId(mother);
88 
89  // generate
90 
91  ofstream df(datfile.c_str());
92 
93  int count = 1;
94  do {
95  EvtVector4R p_init(EvtPDL::getMass(id),0.0,0.0,0.0);
96  EvtParticle *root_part=EvtParticleFactory::particleFactory(id,p_init);
97  root_part->setDiagonalSpinDensity();
98  myGenerator.generateDecay(root_part);
99 
100  int nDaug = root_part->getNDaug();
101  if (nDaug == 3) {
102 
103  EvtVector4R p0 = root_part->getDaug(0)->getP4Lab();
104  EvtVector4R p1 = root_part->getDaug(1)->getP4Lab();
105  EvtVector4R p2 = root_part->getDaug(2)->getP4Lab();
106  double qAB = (p0+p1).mass2();
107  double qBC = (p1+p2).mass2();
108  double qCA = (p2+p0).mass2();
109 
110  df << qAB << " " << qBC << " " << qCA << " " << endl;
111  }
112 
113  root_part->deleteTree();
114 
115  } while(count++ < N);
116 
117  delete eng;
118 
119 }
120 
121 // Macro for displaying resutls
122 /*
123 TTree* readTree(const char* file, int N)
124 {
125  double qAB,qBC,qCA;
126  TTree* t = new TTree();
127  t->Branch("qAB",&qAB,"qAB/D");
128  t->Branch("qBC",&qBC,"qBC/D");
129  t->Branch("qCA",&qCA,"qCA/D");
130 
131  ifstream f(file);
132 
133  int i=0;
134  while(i < N) {
135  i++;
136  f >> qAB >> qBC >> qCA;
137  t->Fill();
138  }
139 
140  return t;
141 }
142 */
143 
144 
145 
146 
EvtParticle * getDaug(int i)
EvtVector4R getP4Lab() const
size_t getNDaug() const
Definition: EvtId.hh:27
Definition: EvtGen.hh:46
void deleteTree()
void setDiagonalSpinDensity()