BOSS 7.0.4
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtGen.cc
Go to the documentation of this file.
1//--------------------------------------------------------------------------
2//
3// Environment:
4// This software is part of the EvtGen package developed jointly
5// for the BaBar and CLEO collaborations. If you use all or part
6// of it, please give an appropriate acknowledgement.
7//
8// Copyright Information: See EvtGen/COPYRIGHT
9// Copyright (C) 1998 Caltech, UCSB
10//
11// Module: EvtGen.cc
12//
13// Description: Main class to provide user interface to EvtGen.
14//
15// Modification history:
16//
17// RYD March 24, 1998 Module created
18//
19//------------------------------------------------------------------------
20//
21
23#include <stdio.h>
24#include <fstream>
25#include <math.h>
27#include <stdlib.h>
28#include "EvtGen.hh"
34#include "EvtGenBase/EvtPDL.hh"
38#include "EvtGenBase/EvtId.hh"
42#include "CLHEP/Vector/LorentzVector.h"
48
49using std::endl;
50using std::fstream;
51using std::ifstream;
52
53extern "C" void begevtgenstore_(int *,int *,int *,int *,
54 int *,int *,int *,int *,int *,
55 double *,double *,double *,
56 double *,double *,double *,
57 double *,double *,double *);
58
59extern "C" {
60extern void evtgen_(float svertex[3],float *e_cms,float *beta_zs,
61 int *mode);
62}
63
65
66 //This is a bit uggly, should not do anything
67 //in a destructor. This will fail if EvtGen is made a static
68 //because then this destructor might be called _after_
69 //the destructoin of objects that it depends on, e.g., EvtPDL.
70
71 if (getenv("EVTINFO")){
73 }
74
75}
76
77EvtGen::EvtGen(const char* const decayName,
78 const char* const pdtTableName,
79 EvtRandomEngine* randomEngine,
80 EvtAbsRadCorr* isrEngine){
81
82
83 report(INFO,"EvtGen") << "Initializing EvtGen"<<endl;
84
85 report(INFO,"EvtGen") << "Storing known decay models"<<endl;
86 EvtModelReg dummy;
87
88 report(INFO,"EvtGen") << "Main decay file name :"<<decayName<<endl;
89 report(INFO,"EvtGen") << "PDT table file name :"<<pdtTableName<<endl;
90
91 report(INFO,"EvtGen") << "Initializing RadCorr=PHOTOS"<<endl;
92 if (isrEngine==0){
93 static EvtPHOTOS defaultRadCorrEngine;
94 EvtRadCorr::setRadCorrEngine(&defaultRadCorrEngine);
95 report(INFO,"EvtGen") <<"No RadCorr engine given in "
96 <<"EvtGen::EvtGen constructor, "
97 <<"will use default EvtPHOTOS."<<endl;
98 }
99 else{
101 }
102
103 _pdl.readPDT(pdtTableName);
104
105
106 ifstream indec;
107
109
110 if (randomEngine==0){
111 static EvtRandomEngine defaultRandomEngine;
112 EvtRandom::setRandomEngine(&defaultRandomEngine);
113 report(INFO,"EvtGen") <<"No random engine given in "
114 <<"EvtGen::EvtGen constructor, "
115 <<"will use default EvtRandomEngine."<<endl;
116 }
117 else{
118 EvtRandom::setRandomEngine(randomEngine);
119 }
120
121 report(INFO,"EvtGen") << "Done initializing EvtGen"<<endl;
122
123}
124
125
126void EvtGen::readUDecay(const char* const uDecayName){
127
128 ifstream indec;
129
130 if ( uDecayName[0] == 0) {
131 report(INFO,"EvtGen") << "Is not reading a user decay file!"<<endl;
132 }
133 else{
134 indec.open(uDecayName);
135 if (indec) {
137
138 report(INFO,"EvtGen") << "Reading "<<uDecayName
139 <<" to override decay table."<<endl;
140 }
141 else{
142
143 report(INFO,"EvtGen") << "Can not find UDECAY file '"
144 <<uDecayName<<"'. Stopping"<<endl;
145 ::abort();
146 }
147 }
148
149}
150
151
152void EvtGen::generateDecay(int stdhepid,
153 EvtVector4R P,
154 EvtVector4R D,
155 EvtStdHep *evtStdHep,
156 EvtSpinDensity *spinDensity ){
157
158
159 EvtParticle *p;
160
161 if(spinDensity==0){
163 }
164 else{
166 P,*spinDensity);
167 }
168
169 generateDecay(p);
170 // p->Decay();
171
172 evtStdHep->init();
173
174 p->makeStdHep(*evtStdHep);
175
176 evtStdHep->translate(D);
177
178 p->deleteTree();
179
180
181}
182
184
185 int times=0;
186 do{
187 times+=1;
189 p->decay();
190 //ok then finish.
191 if ( EvtStatus::getRejectFlag()==0 ) {
192 times=0;
193 }
194 else{
195
196 int ii;
197 for (ii=0;ii<p->getNDaug();ii++){
198 EvtParticle *temp=p->getDaug(ii);
199 temp->deleteTree();
200 }
201 p->resetFirstOrNot();
202 p->resetNDaug();
203
204 }
205
206 if ( times==10000) {
207 report(ERROR,"EvtGen") << "Your event has been rejected 10000 times!"<<endl;
208 report(ERROR,"EvtGen") << "Will now abort."<<endl;
209 ::abort();
210 times=0;
211 }
212 } while (times);
213
214}
215
216
217
218void EvtGen::generateEvent(int stdhepid,HepLorentzVector P,HepLorentzVector D){
219
220 EvtParticle *root_part;
221 EvtVectorParticle *vector_part;
222
223 vector_part=new EvtVectorParticle;
224 EvtVector4R p_init;
225
226 p_init.set(P.t(),P.x(),P.y(),P.z());
227
228 vector_part->init(EvtPDL::evtIdFromStdHep(stdhepid),p_init);
229
230 root_part=(EvtParticle *)vector_part;
231
232 root_part->setVectorSpinDensity();
233
234 generateEvent(root_part,D);
235
236 root_part->deleteTree();
237
238}
239
240void EvtGen::generateEvent(EvtParticle *root_part,HepLorentzVector D){
241 int i;
242 static int nevent=0;
243 nevent++;
244
245 static EvtStdHep evtstdhep;
246 // static EvtSecondary evtsecondary;
247
248 int j;
249 int istat;
250 int partnum;
251 double px,py,pz,e,m;
252 double x,y,z,t;
253
254 EvtVector4R p4,x4;
255 generateDecay(root_part);
256 // root_part->Decay();
257 int npart=0;
258 EvtId list_of_stable[10];
259 EvtParticle* stable_parent[10];
260
261
262 list_of_stable[0]=EvtId(-1,-1);
263 stable_parent[0]=0;
264
265 evtstdhep.init();
266 // evtsecondary.init();
267 // root_part->makeStdHep(evtstdhep,evtsecondary,list_of_stable);
268 root_part->makeStdHep(evtstdhep);
269
270 //report(INFO,"EvtGen") << evtstdhep;
271 //report(INFO,"EvtGen") << evtsecondary;
272
273 npart=evtstdhep.getNPart();
274 for(i=0;i<evtstdhep.getNPart();i++){
275
276 j=i+1;
277
278 int jmotherfirst=evtstdhep.getFirstMother(i)+1;
279 int jmotherlast=evtstdhep.getLastMother(i)+1;
280 int jdaugfirst=evtstdhep.getFirstDaughter(i)+1;
281 int jdauglast=evtstdhep.getLastDaughter(i)+1;
282
283 partnum=evtstdhep.getStdHepID(i);
284
285 istat=evtstdhep.getIStat(i);
286
287 p4=evtstdhep.getP4(i);
288 x4=evtstdhep.getX4(i);
289 px=p4.get(1);
290 py=p4.get(2);
291 pz=p4.get(3);
292 e=p4.get(0);
293
294 x=x4.get(1)+D.x();
295 y=x4.get(2)+D.y();
296 z=x4.get(3)+D.z();
297 t=x4.get(0)+D.t();
298
299 m=p4.mass();
300
301 begevtgenstore_(&j,&nevent,&npart,&istat,
302 &partnum,&jmotherfirst,&jmotherlast,
303 &jdaugfirst,&jdauglast,
304 &px,&py,&pz,&e,
305 &m,&x,&y,&z,&t);
306 }
307
308}
309
310
311
312
313
314
315
316
317
double P(RecMdcKalTrack *trk)
Double_t x[10]
void begevtgenstore_(int *, int *, int *, int *, int *, int *, int *, int *, int *, double *, double *, double *, double *, double *, double *, double *, double *, double *)
void evtgen_(float svertex[3], float *e_cms, float *beta_zs, int *mode)
ostream & report(Severity severity, const char *facility)
Definition: EvtReport.cc:36
@ ERROR
Definition: EvtReport.hh:49
@ INFO
Definition: EvtReport.hh:52
static void readDecayFile(const std::string dec_name)
static void printSummary()
void readUDecay(const char *const udecay_name)
Definition: EvtGen.cc:126
~EvtGen()
Definition: EvtGen.cc:64
void generateDecay(int stdhepid, EvtVector4R P, EvtVector4R D, EvtStdHep *evtStdHep, EvtSpinDensity *spinDensity=0)
Definition: EvtGen.cc:152
void generateEvent(int stdhepid, HepLorentzVector P, HepLorentzVector D)
Definition: EvtGen.cc:218
EvtGen(const char *const decayName, const char *const pdtTableName, EvtRandomEngine *randomEngine=0, EvtAbsRadCorr *isrEngine=0)
Definition: EvtGen.cc:77
Definition: EvtId.hh:27
void readPDT(const std::string fname)
Definition: EvtPDL.cc:63
static EvtId evtIdFromStdHep(int stdhep)
Definition: EvtPDL.cc:244
static EvtParticle * particleFactory(EvtSpinType::spintype spinType)
void resetNDaug()
Definition: EvtParticle.hh:269
void decay()
Definition: EvtParticle.cc:404
void setVectorSpinDensity()
Definition: EvtParticle.cc:138
void resetFirstOrNot()
Definition: EvtParticle.cc:77
int getNDaug() const
Definition: EvtParticle.cc:125
EvtParticle * getDaug(int i)
Definition: EvtParticle.cc:85
void makeStdHep(EvtStdHep &stdhep, EvtSecondary &secondary, EvtId *stable_parent_ihep)
Definition: EvtParticle.cc:759
void deleteTree()
Definition: EvtParticle.cc:557
static void setRadCorrEngine(EvtAbsRadCorr *isrEngine)
Definition: EvtRadCorr.cc:47
static void setRandomEngine(EvtRandomEngine *randomEngine)
Definition: EvtRandom.cc:36
static void initRejectFlag()
Definition: EvtStatus.hh:40
static int getRejectFlag()
Definition: EvtStatus.hh:42
int getIStat(int i)
Definition: EvtStdHep.hh:45
void translate(EvtVector4R d)
Definition: EvtStdHep.cc:69
int getFirstMother(int i)
Definition: EvtStdHep.hh:39
int getLastMother(int i)
Definition: EvtStdHep.hh:40
void init()
Definition: EvtStdHep.cc:33
int getLastDaughter(int i)
Definition: EvtStdHep.hh:42
int getStdHepID(int i)
Definition: EvtStdHep.hh:44
int getFirstDaughter(int i)
Definition: EvtStdHep.hh:41
int getNPart()
Definition: EvtStdHep.cc:37
EvtVector4R getP4(int i)
Definition: EvtStdHep.hh:47
EvtVector4R getX4(int i)
Definition: EvtStdHep.hh:48
double mass() const
Definition: EvtVector4R.cc:39
double get(int i) const
Definition: EvtVector4R.hh:179
void set(int i, double d)
Definition: EvtVector4R.hh:183
void init(EvtId part_n, double e, double px, double py, double pz)
int t()
Definition: t.c:1