Main Page | Namespace List | Class Hierarchy | Alphabetical List | Class List | Directories | File List | Namespace Members | Class Members | File Members | Related Pages

EvtGen Class Reference

#include <EvtGen.hh>

List of all members.

Public Member Functions

 EvtGen (const char *const decayName, const char *const pdtTableName, EvtRandomEngine *randomEngine=0, EvtAbsRadCorr *isrEngine=0)
void generateDecay (EvtParticle *p)
void generateDecay (int stdhepid, EvtVector4R P, EvtVector4R D, EvtStdHep *evtStdHep, EvtSpinDensity *spinDensity=0)
void generateEvent (EvtParticle *p, HepLorentzVector D)
void generateEvent (int stdhepid, HepLorentzVector P, HepLorentzVector D)
void readUDecay (const char *const udecay_name)
 ~EvtGen ()

Private Attributes

EvtPDL _pdl


Constructor & Destructor Documentation

EvtGen::EvtGen const char *const   decayName,
const char *const   pdtTableName,
EvtRandomEngine randomEngine = 0,
EvtAbsRadCorr isrEngine = 0
 

00080                                         {
00081 
00082 
00083   report(INFO,"EvtGen") << "Initializing EvtGen"<<endl;
00084 
00085   report(INFO,"EvtGen") << "Storing known decay models"<<endl;
00086   EvtModelReg dummy;
00087 
00088   report(INFO,"EvtGen") << "Main decay file name  :"<<decayName<<endl;
00089   report(INFO,"EvtGen") << "PDT table file name   :"<<pdtTableName<<endl;
00090   
00091   report(INFO,"EvtGen") << "Initializing RadCorr=PHOTOS"<<endl;
00092   if (isrEngine==0){
00093     static EvtPHOTOS defaultRadCorrEngine;
00094     EvtRadCorr::setRadCorrEngine(&defaultRadCorrEngine);
00095     report(INFO,"EvtGen") <<"No RadCorr engine given in "
00096                           <<"EvtGen::EvtGen constructor, "
00097                           <<"will use default EvtPHOTOS."<<endl;
00098   }
00099   else{
00100     EvtRadCorr::setRadCorrEngine(isrEngine);    
00101   }
00102 
00103   _pdl.readPDT(pdtTableName);
00104 
00105 
00106   ifstream indec;
00107 
00108   EvtDecayTable::readDecayFile(decayName);
00109 
00110   if (randomEngine==0){
00111     static EvtRandomEngine defaultRandomEngine;
00112     EvtRandom::setRandomEngine(&defaultRandomEngine);
00113     report(INFO,"EvtGen") <<"No random engine given in "
00114                           <<"EvtGen::EvtGen constructor, "
00115                           <<"will use default EvtRandomEngine."<<endl;
00116   }
00117   else{
00118     EvtRandom::setRandomEngine(randomEngine);    
00119   }
00120 
00121   report(INFO,"EvtGen") << "Done initializing EvtGen"<<endl;
00122 
00123 }

EvtGen::~EvtGen  ) 
 

00064                {
00065 
00066   //This is a bit uggly, should not do anything
00067   //in a destructor. This will fail if EvtGen is made a static
00068   //because then this destructor might be called _after_
00069   //the destructoin of objects that it depends on, e.g., EvtPDL.
00070 
00071   if (getenv("EVTINFO")){
00072     EvtDecayTable::printSummary();
00073   }
00074 
00075 }


Member Function Documentation

void EvtGen::generateDecay EvtParticle p  ) 
 

00183                                         {
00184 
00185   int times=0;
00186   do{
00187     times+=1;
00188     EvtStatus::initRejectFlag();
00189     p->decay();
00190 
00191     //ok then finish.
00192     if ( EvtStatus::getRejectFlag()==0 ) { 
00193       times=0;
00194     }
00195     else{   
00196       
00197       int ii;
00198       for (ii=0;ii<p->getNDaug();ii++){
00199         EvtParticle *temp=p->getDaug(ii);
00200         temp->deleteTree();
00201       }
00202       p->resetFirstOrNot();
00203       p->resetNDaug();
00204       
00205     }
00206 
00207     if ( times==10000) {
00208       report(ERROR,"EvtGen") << "Your event has been rejected 10000 times!"<<endl;
00209       report(ERROR,"EvtGen") << "Will now abort."<<endl;
00210       ::abort();
00211       times=0;
00212     }
00213   } while (times);
00214 
00215 }

void EvtGen::generateDecay int  stdhepid,
EvtVector4R  P,
EvtVector4R  D,
EvtStdHep evtStdHep,
EvtSpinDensity spinDensity = 0
 

00156                                                         {
00157 
00158  
00159   EvtParticle *p;
00160 
00161   if(spinDensity==0){    
00162     p=EvtParticleFactory::particleFactory(EvtPDL::evtIdFromStdHep(stdhepid),P);
00163   }
00164   else{
00165     p=EvtParticleFactory::particleFactory(EvtPDL::evtIdFromStdHep(stdhepid),
00166                                           P,*spinDensity);
00167   }
00168 
00169   generateDecay(p);
00170   //  p->Decay();
00171 
00172   evtStdHep->init();
00173 
00174   p->makeStdHep(*evtStdHep);
00175   
00176   evtStdHep->translate(D);
00177   
00178   p->deleteTree();
00179 
00180 
00181 }

void EvtGen::generateEvent EvtParticle p,
HepLorentzVector  D
 

00241                                                                    {
00242   int i;  
00243   static int nevent=0;
00244   nevent++;  
00245 
00246   static EvtStdHep evtstdhep;
00247   //  static EvtSecondary evtsecondary;
00248 
00249   int j;
00250   int istat;
00251   int partnum;
00252   double px,py,pz,e,m;
00253   double x,y,z,t;
00254 
00255   EvtVector4R p4,x4;
00256   generateDecay(root_part);
00257   //  root_part->Decay();
00258   int          npart=0;
00259   EvtId        list_of_stable[10];
00260   EvtParticle* stable_parent[10];
00261     
00262     
00263   list_of_stable[0]=EvtId(-1,-1);
00264   stable_parent[0]=0;
00265 
00266   evtstdhep.init();
00267   //  evtsecondary.init();
00268   // root_part->makeStdHep(evtstdhep,evtsecondary,list_of_stable);
00269   root_part->makeStdHep(evtstdhep);
00270 
00271   //report(INFO,"EvtGen") << evtstdhep;
00272   //report(INFO,"EvtGen") << evtsecondary;
00273   
00274   npart=evtstdhep.getNPart();
00275   for(i=0;i<evtstdhep.getNPart();i++){
00276 
00277     j=i+1;
00278 
00279     int jmotherfirst=evtstdhep.getFirstMother(i)+1;
00280     int jmotherlast=evtstdhep.getLastMother(i)+1;
00281     int jdaugfirst=evtstdhep.getFirstDaughter(i)+1;
00282     int jdauglast=evtstdhep.getLastDaughter(i)+1;
00283 
00284     partnum=evtstdhep.getStdHepID(i);
00285 
00286     istat=evtstdhep.getIStat(i);
00287 
00288     p4=evtstdhep.getP4(i);
00289     x4=evtstdhep.getX4(i);
00290     px=p4.get(1);
00291     py=p4.get(2);
00292     pz=p4.get(3);
00293     e=p4.get(0);
00294           
00295     x=x4.get(1)+D.x();
00296     y=x4.get(2)+D.y();
00297     z=x4.get(3)+D.z();
00298     t=x4.get(0)+D.t();
00299       
00300     m=p4.mass();
00301 
00302     begevtgenstore_(&j,&nevent,&npart,&istat,
00303                     &partnum,&jmotherfirst,&jmotherlast,
00304                     &jdaugfirst,&jdauglast,
00305                     &px,&py,&pz,&e,
00306                     &m,&x,&y,&z,&t);
00307   }
00308 
00309 }

void EvtGen::generateEvent int  stdhepid,
HepLorentzVector  P,
HepLorentzVector  D
 

00219                                                                             {
00220 
00221   EvtParticle *root_part;
00222   EvtVectorParticle *vector_part;
00223   
00224   vector_part=new EvtVectorParticle;
00225   EvtVector4R p_init;
00226 
00227   p_init.set(P.t(),P.x(),P.y(),P.z());
00228 
00229   vector_part->init(EvtPDL::evtIdFromStdHep(stdhepid),p_init);
00230   
00231   root_part=(EvtParticle *)vector_part;
00232   
00233   root_part->setVectorSpinDensity();      
00234 
00235   generateEvent(root_part,D);
00236 
00237   root_part->deleteTree();  
00238 
00239 }

void EvtGen::readUDecay const char *const   udecay_name  ) 
 

00126                                                    {
00127 
00128   ifstream indec;
00129 
00130   if ( uDecayName[0] == 0) {
00131     report(INFO,"EvtGen") << "Is not reading a user decay file!"<<endl;
00132   }
00133   else{  
00134     indec.open(uDecayName);
00135     if (indec) {
00136       EvtDecayTable::readDecayFile(uDecayName);
00137       
00138       report(INFO,"EvtGen") << "Reading "<<uDecayName
00139                             <<" to override decay table."<<endl;
00140     }    
00141     else{
00142       
00143       report(INFO,"EvtGen") << "Can not find UDECAY file '"
00144                             <<uDecayName<<"'.  Stopping"<<endl;
00145       ::abort();
00146     }
00147   }
00148   
00149 }


Member Data Documentation

EvtPDL EvtGen::_pdl [private]
 


The documentation for this class was generated from the following files:
Generated on Wed Feb 2 16:05:25 2011 for BOSS6.5.5 by  doxygen 1.3.9.1