CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_4/src/GeneratorInterface/ExternalDecays/src/EvtGenInterface.cc

Go to the documentation of this file.
00001 
00002 #include "GeneratorInterface/ExternalDecays/interface/EvtGenInterface.h"
00003 
00004 #include "FWCore/PluginManager/interface/PluginManager.h"
00005 
00006 #include "FWCore/Framework/interface/MakerMacros.h"
00007 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00008 
00009 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00010 
00011 #include "FWCore/Utilities/interface/Exception.h"
00012 #include "FWCore/Framework/interface/Event.h"
00013 #include "FWCore/Framework/interface/Run.h"
00014 #include "Utilities/General/interface/FileInPath.h"
00015 #include "FWCore/ServiceRegistry/interface/Service.h"
00016 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00017 #include "CLHEP/Random/Random.h"
00018 
00019 #include "GeneratorInterface/ExternalDecays/interface/myEvtRandomEngine.h"
00020 #include "GeneratorInterface/Pythia6Interface/interface/Pythia6Service.h"
00021 
00022 #include "HepMC/GenEvent.h"
00023 // #include "HepMC/PythiaWrapper6_2.h"
00024 
00025 #include "DataFormats/GeometryVector/interface/GlobalVector.h"
00026 
00027 //#define PYGIVE pygive_
00028 //extern "C" {
00029 //  void PYGIVE(const char*,int length);
00030 //}
00031 
00032 using namespace gen;
00033 using namespace edm;
00034 
00035 EvtGenInterface::EvtGenInterface( const ParameterSet& pset )
00036 {
00037 
00038   ntotal = 0;
00039   nevent = 0;
00040   std::cout << " EvtGenProducer starting ... " << std::endl;
00041 
00042   // create random engine and initialize seed using Random Number 
00043   // Generator Service 
00044   // as configured in the configuration file
00045 
00046   edm::Service<edm::RandomNumberGenerator> rngen;
00047   
00048   if ( ! rngen.isAvailable()) {             
00049     throw cms::Exception("Configuration")
00050       << "The EvtGenProducer module requires the RandomNumberGeneratorService\n"
00051       "which is not present in the configuration file.  You must add the service\n"
00052       "in the configuration file if you  want to run EvtGenProducer";
00053     }
00054 
00055   CLHEP::HepRandomEngine& m_engine = rngen->getEngine();
00056   m_flat = new CLHEP::RandFlat(m_engine, 0., 1.);
00057   myEvtRandomEngine* the_engine = new myEvtRandomEngine(&m_engine); 
00058 
00059   // Get data from parameter set
00060   edm::FileInPath decay_table = pset.getParameter<edm::FileInPath>("decay_table");
00061   edm::FileInPath pdt = pset.getParameter<edm::FileInPath>("particle_property_file");
00062   bool useDefault = pset.getUntrackedParameter<bool>("use_default_decay",true);
00063   usePythia = pset.getUntrackedParameter<bool>("use_internal_pythia",true);
00064   polarize_ids = pset.getUntrackedParameter<std::vector<int> >("particles_to_polarize",
00065                                                                std::vector<int>());
00066   polarize_pol = pset.getUntrackedParameter<std::vector<double> >("particle_polarizations",
00067                                                                   std::vector<double>());
00068   if (polarize_ids.size() != polarize_pol.size()) {
00069       throw cms::Exception("Configuration")
00070           << "EvtGenProducer requires that the particles_to_polarize and particle_polarization\n"
00071           "vectors be the same size.  Please fix this in your configuration.";
00072   }
00073   for (unsigned int ndx = 0; ndx < polarize_ids.size(); ndx++) {
00074       if (polarize_pol[ndx] < -1. || polarize_pol[ndx] > 1.) {
00075           throw cms::Exception("Configuration")
00076               << "EvtGenProducer error: particle polarizations must be in the range -1 < P < 1";
00077       }
00078       polarizations.insert(std::pair<int, float>(polarize_ids[ndx], polarize_pol[ndx]));
00079   }
00080 
00081   edm::FileInPath user_decay = pset.getParameter<edm::FileInPath>("user_decay_file");
00082   std::string decay_table_s = decay_table.fullPath();
00083   std::string pdt_s = pdt.fullPath();
00084   std::string user_decay_s = user_decay.fullPath();
00085 
00086   //-->pythia_params = pset.getParameter< std::vector<std::string> >("processParameters");
00087   
00088   
00089   // any number of alias names for forced decays can be specified using dynamic std vector of strings 
00090   std::vector<std::string> forced_names = pset.getParameter< std::vector<std::string> >("list_forced_decays");
00091     
00092   m_EvtGen = new EvtGen (decay_table_s.c_str(),pdt_s.c_str(),the_engine);  
00093   // 4th parameter should be rad cor - set to PHOTOS (default)
00094  
00095   if (!useDefault) m_EvtGen->readUDecay( user_decay_s.c_str() );
00096 
00097   std::vector<std::string>::const_iterator i;
00098   nforced=0;
00099 
00100   for (i=forced_names.begin(); i!=forced_names.end(); ++i)   
00101     // i will point to strings containing
00102     // names of particles with forced decay
00103     {
00104       nforced++;
00105       EvtId found = EvtPDL::getId(*i);      
00106       if (found.getId()==-1)
00107         {
00108           throw cms::Exception("Configuration")
00109             << "name in part list for forced decays not found: " << *i; 
00110         }
00111       if (found.getId()==found.getAlias())
00112         {
00113           throw cms::Exception("Configuration")
00114             << "name in part list for forced decays is not an alias: " << *i; 
00115         }
00116       forced_Evt.push_back(found);                      // forced_Evt is the list of EvtId's
00117       forced_Hep.push_back(EvtPDL::getStdHep(found));   // forced_Hep is the list of stdhep codes
00118     }
00119 
00120    // fill up default list of particles to be declared stable in the "master generator"
00121    // these are assumed to be PDG ID's
00122    // in case of combo with Pythia6, translation is done in Pythia6Hadronizer
00123    //
00124    // Note: Pythia6's kc=43, 44, and 84 commented out because they're obsolete (per S.Mrenna)
00125    //
00126    m_PDGs.push_back( 300553 ) ;
00127    m_PDGs.push_back( 511 ) ;
00128    m_PDGs.push_back( 521 ) ;
00129    m_PDGs.push_back( 523 ) ;
00130    m_PDGs.push_back( 513 ) ;
00131    m_PDGs.push_back( 533 ) ;
00132    m_PDGs.push_back( 531 ) ;
00133    
00134    m_PDGs.push_back( 15 ) ;
00135    
00136    m_PDGs.push_back( 413 ) ;
00137    m_PDGs.push_back( 423 ) ;
00138    m_PDGs.push_back( 433 ) ;
00139    m_PDGs.push_back( 411 ) ;
00140    m_PDGs.push_back( 421 ) ;
00141    m_PDGs.push_back( 431 ) ;   
00142    m_PDGs.push_back( 10411 );
00143    m_PDGs.push_back( 10421 );
00144    m_PDGs.push_back( 10413 );
00145    m_PDGs.push_back( 10423 );   
00146    m_PDGs.push_back( 20413 );
00147    m_PDGs.push_back( 20423 );
00148     
00149    m_PDGs.push_back( 415 );
00150    m_PDGs.push_back( 425 );
00151    m_PDGs.push_back( 10431 );
00152    m_PDGs.push_back( 20433 );
00153    m_PDGs.push_back( 10433 );
00154    m_PDGs.push_back( 435 );
00155    
00156    m_PDGs.push_back( 310 );
00157    m_PDGs.push_back( 311 );
00158    m_PDGs.push_back( 313 );
00159    m_PDGs.push_back( 323 );
00160    m_PDGs.push_back( 10321 );
00161    m_PDGs.push_back( 10311 );
00162    m_PDGs.push_back( 10313 );
00163    m_PDGs.push_back( 10323 );
00164    m_PDGs.push_back( 20323 );
00165    m_PDGs.push_back( 20313 );
00166    m_PDGs.push_back( 325 );
00167    m_PDGs.push_back( 315 );
00168    
00169    m_PDGs.push_back( 100313 );
00170    m_PDGs.push_back( 100323 );
00171    m_PDGs.push_back( 30313 );
00172    m_PDGs.push_back( 30323 );
00173    m_PDGs.push_back( 30343 );
00174    m_PDGs.push_back( 30353 );
00175    m_PDGs.push_back( 30363 );
00176 
00177    m_PDGs.push_back( 111 );
00178    m_PDGs.push_back( 221 );
00179    m_PDGs.push_back( 113 );
00180    m_PDGs.push_back( 213 );
00181    m_PDGs.push_back( 223 );
00182    m_PDGs.push_back( 331 );
00183    m_PDGs.push_back( 333 );
00184    m_PDGs.push_back( 20213 );
00185    m_PDGs.push_back( 20113 );
00186    m_PDGs.push_back( 215 );
00187    m_PDGs.push_back( 115 );
00188    m_PDGs.push_back( 10213 );
00189    m_PDGs.push_back( 10113 );
00190    m_PDGs.push_back( 9000111 ); // PDG ID = 9000111, Pythia6 ID = 10111
00191    m_PDGs.push_back( 9000211 ); // PDG ID = 9000211, Pythia6 ID = 10211
00192    m_PDGs.push_back( 9010221 ); // PDG ID = 9010211, Pythia6 ID = ???
00193    m_PDGs.push_back( 10221 );
00194    m_PDGs.push_back( 20223 );
00195    m_PDGs.push_back( 20333 );
00196    m_PDGs.push_back( 225 );
00197    m_PDGs.push_back( 9020221 ); // PDG ID = 9020211, Pythia6 ID = ???
00198    m_PDGs.push_back( 335 );
00199    m_PDGs.push_back( 10223 );
00200    m_PDGs.push_back( 10333 );
00201    m_PDGs.push_back( 100213 );
00202    m_PDGs.push_back( 100113 );
00203    
00204    m_PDGs.push_back( 441 );
00205    m_PDGs.push_back( 100441 );
00206    m_PDGs.push_back( 443 );
00207    m_PDGs.push_back( 100443 );
00208    m_PDGs.push_back( 9000443 );
00209    m_PDGs.push_back( 9010443 );
00210    m_PDGs.push_back( 9020443 );
00211    m_PDGs.push_back( 10441 );
00212    m_PDGs.push_back( 20443 );
00213    m_PDGs.push_back( 445 );
00214 
00215    m_PDGs.push_back( 30443 );
00216    m_PDGs.push_back( 551 );
00217    m_PDGs.push_back( 553 );
00218    m_PDGs.push_back( 100553 );
00219    m_PDGs.push_back( 200553 );
00220    m_PDGs.push_back( 10551 );
00221    m_PDGs.push_back( 20553 );
00222    m_PDGs.push_back( 555 );
00223    m_PDGs.push_back( 10553 );
00224 
00225    m_PDGs.push_back( 110551 );
00226    m_PDGs.push_back( 120553 );
00227    m_PDGs.push_back( 100555 );
00228    m_PDGs.push_back( 210551 );
00229    m_PDGs.push_back( 220553 );
00230    m_PDGs.push_back( 200555 );
00231    m_PDGs.push_back( 30553 );
00232    m_PDGs.push_back( 20555 );
00233 
00234    m_PDGs.push_back( 557 );
00235    m_PDGs.push_back( 130553 ); 
00236    m_PDGs.push_back( 120555 );
00237    m_PDGs.push_back( 100557 );
00238    m_PDGs.push_back( 110553 );
00239    m_PDGs.push_back( 210553 );
00240    m_PDGs.push_back( 10555 );
00241    m_PDGs.push_back( 110555 );
00242 
00243    m_PDGs.push_back( 4122 );
00244    m_PDGs.push_back( 4132 );
00245    // m_PDGs.push_back( 84 ); // obsolete
00246    m_PDGs.push_back( 4112 );
00247    m_PDGs.push_back( 4212 );
00248    m_PDGs.push_back( 4232 );
00249    m_PDGs.push_back( 4222 );
00250    m_PDGs.push_back( 4322 );
00251    m_PDGs.push_back( 4312 );
00252 
00253    m_PDGs.push_back( 13122 );
00254    m_PDGs.push_back( 13124 );
00255    m_PDGs.push_back( 23122 );
00256    m_PDGs.push_back( 33122 );
00257    m_PDGs.push_back( 43122 );
00258    m_PDGs.push_back( 53122 );
00259    m_PDGs.push_back( 13126 );
00260    m_PDGs.push_back( 13212 );
00261    m_PDGs.push_back( 13241 );
00262   
00263    m_PDGs.push_back( 3126 );
00264    m_PDGs.push_back( 3124 );
00265    m_PDGs.push_back( 3122 );
00266    m_PDGs.push_back( 3222 );
00267    m_PDGs.push_back( 2214 );
00268    m_PDGs.push_back( 2224 );
00269    m_PDGs.push_back( 3324 );
00270    m_PDGs.push_back( 2114 );
00271    m_PDGs.push_back( 1114 );
00272    m_PDGs.push_back( 3112 );
00273    m_PDGs.push_back( 3212 );
00274    m_PDGs.push_back( 3114 );
00275    m_PDGs.push_back( 3224 );
00276    m_PDGs.push_back( 3214 );
00277    m_PDGs.push_back( 3216 );
00278    m_PDGs.push_back( 3322 );
00279    m_PDGs.push_back( 3312 );
00280    m_PDGs.push_back( 3314 );
00281    m_PDGs.push_back( 3334 );
00282    
00283    m_PDGs.push_back( 4114 );
00284    m_PDGs.push_back( 4214 );
00285    m_PDGs.push_back( 4224 );
00286    m_PDGs.push_back( 4314 );
00287    m_PDGs.push_back( 4324 );
00288    m_PDGs.push_back( 4332 );
00289    m_PDGs.push_back( 4334 );
00290    //m_PDGs.push_back( 43 ); // obsolete (?)
00291    //m_PDGs.push_back( 44 ); // obsolete (?)
00292    m_PDGs.push_back( 10443 );
00293 
00294    m_PDGs.push_back( 5122 );
00295    m_PDGs.push_back( 5132 );
00296    m_PDGs.push_back( 5232 );
00297    m_PDGs.push_back( 5332 );
00298    m_PDGs.push_back( 5222 );
00299    m_PDGs.push_back( 5112 );
00300    m_PDGs.push_back( 5212 );
00301    m_PDGs.push_back( 541 );
00302    m_PDGs.push_back( 14122 );
00303    m_PDGs.push_back( 14124 );
00304    m_PDGs.push_back( 5312 );
00305    m_PDGs.push_back( 5322 );
00306    m_PDGs.push_back( 10521 );
00307    m_PDGs.push_back( 20523 );
00308    m_PDGs.push_back( 10523 );
00309 
00310    m_PDGs.push_back( 525 );
00311    m_PDGs.push_back( 10511 );
00312    m_PDGs.push_back( 20513 );
00313    m_PDGs.push_back( 10513 );
00314    m_PDGs.push_back( 515 );
00315    m_PDGs.push_back( 10531 );
00316    m_PDGs.push_back( 20533 );
00317    m_PDGs.push_back( 10533 );
00318    m_PDGs.push_back( 535 );
00319    m_PDGs.push_back( 543 );
00320    m_PDGs.push_back( 545 );
00321    m_PDGs.push_back( 5114 );
00322    m_PDGs.push_back( 5224 );
00323    m_PDGs.push_back( 5214 );
00324    m_PDGs.push_back( 5314 );
00325    m_PDGs.push_back( 5324 );
00326    m_PDGs.push_back( 5334 );
00327    m_PDGs.push_back( 10541 );
00328    m_PDGs.push_back( 10543 );
00329    m_PDGs.push_back( 20543 );
00330 
00331    m_PDGs.push_back( 4424 );
00332    m_PDGs.push_back( 4422 );
00333    m_PDGs.push_back( 4414 );
00334    m_PDGs.push_back( 4412 );
00335    m_PDGs.push_back( 4432 );
00336    m_PDGs.push_back( 4434 );
00337 
00338    m_PDGs.push_back( 130 );
00339    
00340    // now check if we need to override default list of particles/IDs
00341    if ( pset.exists("operates_on_particles") )
00342    {
00343       std::vector<int> tmpPIDs = pset.getParameter< std::vector<int> >("operates_on_particles");
00344       if ( tmpPIDs.size() > 0 )
00345       {
00346          if ( tmpPIDs[0] > 0 ) // 0 means default !!!
00347          {
00348             m_PDGs.clear();
00349             m_PDGs = tmpPIDs;
00350          }
00351       }
00352    } 
00353    
00354   m_Py6Service = new Pythia6Service;
00355 } 
00356 
00357 EvtGenInterface::~EvtGenInterface()
00358 {
00359   std::cout << " EvtGenProducer terminating ... " << std::endl; 
00360   delete m_Py6Service;
00361 }
00362 
00363 void EvtGenInterface::init()
00364 {
00365 
00366    Pythia6Service::InstanceWrapper guard(m_Py6Service); // grab Py6 instance
00367    
00368    // Do here initialization of EvtPythia then restore original settings
00369    if (usePythia) EvtPythia::pythiaInit(0);
00370    
00371 // no need - will be done via Pythia6Hadronizer::declareStableParticles
00372 //
00373 //    for( std::vector<std::string>::const_iterator itPar = pythia_params.begin(); itPar != pythia_params.end(); ++itPar ) {
00374 //      call_pygive(*itPar);
00375 //    }
00376 
00377    return ;
00378 
00379 }
00380 
00381 
00382 HepMC::GenEvent* EvtGenInterface::decay( HepMC::GenEvent* evt )
00383 {
00384   Pythia6Service::InstanceWrapper guard(m_Py6Service);  // grab Py6 instance
00385 
00386   nevent++;
00387   npartial = 0;
00388   // std::cout << "nevent = " << nevent << std::endl ;
00389   
00390   int idHep,ipart,status;
00391   EvtId idEvt;
00392 
00393   nPythia = evt->particles_size();
00394   // FIX A MEMORY LEAK (RC)
00395   // HepMC::GenEvent* newEvt = new HepMC::GenEvent( *evt );
00396 
00397   // First pass through undecayed Pythia particles to decay particles known to EvtGen left stable by Pythia
00398   // except candidates to be forced which will be searched later to include EvtGen decay products 
00399   nlist = 0;
00400 
00401   // Notice dynamical use of evt
00402   for (HepMC::GenEvent::particle_const_iterator p= evt->particles_begin(); p != evt->particles_end(); ++p)
00403     {
00404       status = (*p)->status();
00405  
00406       if(status==1) {           // only not decayed (status = 1) particles
00407           
00408 
00409           idHep = (*p)->pdg_id();
00410           int do_force=0;
00411           for(int i=0;i<nforced;i++)           // First check if part with forced decay
00412             {                                  // In that case do not decay immediately 
00413               if(idHep == forced_Hep[i])       // (only 1 per event will be forced)      
00414                 {                              // Fill list
00415                   update_candlist(i,*p);
00416                   do_force=1;
00417                 }
00418             }
00419           if(do_force==0)         // particles with decays not forced are decayed immediately 
00420             {
00421               idEvt = EvtPDL::evtIdFromStdHep(idHep);
00422               ipart = idEvt.getId();
00423               if (ipart==-1) continue;                          // particle not known to EvtGen       
00424               if (EvtDecayTable::getNMode(ipart)==0) continue;  // particles stable for EvtGen
00425               addToHepMC(*p,idEvt,evt,true);                      // generate decay
00426             }
00427         }
00428     }
00429 
00430   if(nlist!=0)   
00431      {
00432       // decide randomly which one to decay as alias
00433       int which = (int)(nlist*m_flat->fire()); 
00434       if (which == nlist) which = nlist-1;
00435   
00436           for(int k=0;k < nlist; k++)
00437             {
00438               if(k == which || !usePythia) {            
00439                 addToHepMC(listp[k],forced_Evt[index[k]],evt,false);  // decay as alias
00440               } 
00441               else
00442                 {
00443                   int id_non_alias = forced_Evt[index[k]].getId();
00444                   EvtId non_alias(id_non_alias,id_non_alias); // create new EvtId with id = alias
00445                   addToHepMC(listp[k],non_alias,evt,false);     // decay as standard (non alias)
00446                 }
00447             }
00448      }
00449 
00450   return evt;
00451   
00452 }
00453 
00454 void EvtGenInterface::addToHepMC(HepMC::GenParticle* partHep, EvtId idEvt, HepMC::GenEvent* theEvent, bool del_daug )
00455 {
00456   // Set spin type
00457   EvtSpinType::spintype stype = EvtPDL::getSpinType(idEvt);
00458   EvtParticle* partEvt;
00459     switch (stype){
00460     case EvtSpinType::SCALAR: 
00461       partEvt = new EvtScalarParticle();
00462       break;
00463     case EvtSpinType::STRING:
00464       partEvt = new EvtStringParticle();    
00465       break;
00466     case EvtSpinType::DIRAC: 
00467       partEvt = new EvtDiracParticle();
00468       break;
00469     case EvtSpinType::VECTOR:
00470       partEvt = new EvtVectorParticle();
00471       break;
00472     case EvtSpinType::RARITASCHWINGER:
00473       partEvt = new EvtRaritaSchwingerParticle();
00474       break;
00475     case EvtSpinType::TENSOR:
00476       partEvt = new EvtTensorParticle();
00477       break;
00478     case EvtSpinType::SPIN5HALF: case EvtSpinType::SPIN3: case EvtSpinType::SPIN7HALF: case EvtSpinType::SPIN4:
00479       partEvt = new EvtHighSpinParticle();
00480       break;
00481     default:
00482       std::cout << "Unknown spintype in EvtSpinType!" << std::endl;   
00483       return;
00484     }
00485 
00486     // Generate decay
00487     EvtVector4R momEvt;
00488     HepMC::FourVector momHep = partHep->momentum();
00489     momEvt.set(momHep.t(),momHep.x(),momHep.y(),momHep.z());
00490     EvtVector4R posEvt;
00491     HepMC::GenVertex* initVert = partHep->production_vertex();
00492     HepMC::FourVector posHep = initVert->position();
00493     posEvt.set(posHep.t(),posHep.x(),posHep.y(),posHep.z());
00494     partEvt->init(idEvt,momEvt);
00495     if (stype == EvtSpinType::DIRAC 
00496         && polarizations.find(partHep->pdg_id()) != polarizations.end()) {
00497          // std::cout << "Polarize particle" << std::endl;
00498         //Particle is spin 1/2, so we can polarize it.
00499         //Check polarizations map for particle, grab its polarization if it exists
00500         // and make the spin density matrix
00501         float pol = polarizations.find(partHep->pdg_id())->second;
00502         GlobalVector pPart(momHep.x(), momHep.y(), momHep.z());
00503         //std::cout << "Polarizing particle with PDG ID "
00504         //  << partHep->pdg_id()
00505         //  << " at " << pol*100 << "%" << std::endl;
00506         GlobalVector zHat(0., 0., 1.);
00507         GlobalVector zCrossP = zHat.cross(pPart);
00508         GlobalVector polVec = pol * zCrossP.unit();
00509 
00510         EvtSpinDensity theSpinDensity;
00511         theSpinDensity.SetDim(2);
00512         theSpinDensity.Set(0, 0, EvtComplex(1./2. + polVec.z()/2., 0.));
00513         theSpinDensity.Set(0, 1, EvtComplex(polVec.x()/2., -polVec.y()/2.));
00514         theSpinDensity.Set(1, 0, EvtComplex(polVec.x()/2., polVec.y()/2.));
00515         theSpinDensity.Set(1, 1, EvtComplex(1./2. - polVec.z()/2., 0.));
00516 
00517         partEvt->setSpinDensityForwardHelicityBasis(theSpinDensity);
00518 
00519     } else {
00520         partEvt->setDiagonalSpinDensity();     
00521     }
00522     partEvt->decay();
00523                        
00524     // extend the search of candidates to be forced to EvtGen decay products and delete their daughters  ** 
00525     // otherwise they wouldn't get their chance to take part in the forced decay lottery                 **
00526     if (del_daug) go_through_daughters(partEvt);    // recursive function go_through_daughters will do   **
00527 
00528     // Change particle in stdHEP format
00529     static EvtStdHep evtstdhep;
00530     
00531     evtstdhep.init();
00532     partEvt->makeStdHep(evtstdhep);
00533 
00534     ntotal++;
00535     partEvt->deleteTree();
00536 
00537     // ********* Now add to the HepMC Event **********
00538 
00539     // Then loop on evtstdhep to add vertexes... 
00540     HepMC::GenVertex* theVerts[200];
00541     for (int ivert = 0; ivert < 200; ivert++) { 
00542       theVerts[ivert] = 0;
00543     }
00544 
00545     for (int ipart = 0; ipart < evtstdhep.getNPart(); ipart++) {
00546       int theMum = evtstdhep.getFirstMother(ipart);
00547       if (theMum != -1 && !theVerts[theMum]) {
00548         EvtVector4R theVpos = evtstdhep.getX4(ipart) + posEvt;
00549         theVerts[theMum] = 
00550           new HepMC::GenVertex(HepMC::FourVector(theVpos.get(1),
00551                                                  theVpos.get(2),
00552                                                  theVpos.get(3),
00553                                                  theVpos.get(0)),0);
00554       }
00555     }
00556 
00557     // ...then particles
00558     partHep->set_status(2);
00559     if (theVerts[0]) theVerts[0]->add_particle_in( partHep );
00560 
00561     for (int ipart2 = 1; ipart2 < evtstdhep.getNPart(); ipart2++) {
00562       int idHep = evtstdhep.getStdHepID(ipart2);
00563       HepMC::GenParticle* thePart = 
00564         new HepMC::GenParticle( HepMC::FourVector(evtstdhep.getP4(ipart2).get(1),
00565                                                   evtstdhep.getP4(ipart2).get(2),
00566                                                   evtstdhep.getP4(ipart2).get(3),
00567                                                   evtstdhep.getP4(ipart2).get(0)),
00568                                 idHep,
00569                                 evtstdhep.getIStat(ipart2));
00570       npartial++;
00571       thePart->suggest_barcode(npartial + nPythia);
00572       int theMum2 = evtstdhep.getFirstMother(ipart2);
00573       if (theMum2 != -1 && theVerts[theMum2]) theVerts[theMum2]->add_particle_out( thePart );
00574       if (theVerts[ipart2]) theVerts[ipart2]->add_particle_in( thePart );
00575        
00576     }
00577     
00578     for (int ipart3 = 0; ipart3 < evtstdhep.getNPart(); ipart3++) {
00579       if (theVerts[ipart3]) theEvent->add_vertex( theVerts[ipart3] );
00580     }
00581     
00582 }        
00583 
00584 /*
00585 void
00586 EvtGenInterface::call_pygive(const std::string& iParm ) {
00587   
00588   //call the fortran routine pygive with a fortran string
00589   PYGIVE( iParm.c_str(), iParm.length() );  
00590  
00591 }
00592 */
00593 
00594 void 
00595 EvtGenInterface::go_through_daughters(EvtParticle* part) {
00596 
00597   int NDaug=part->getNDaug();
00598   if(NDaug)
00599     {
00600       EvtParticle* Daughter;
00601       for (int i=0;i<NDaug;i++)
00602         {
00603           Daughter=part->getDaug(i);
00604           int idHep = EvtPDL::getStdHep(Daughter->getId());
00605           int found=0;
00606           for(int k=0;k<nforced;k++)         
00607             {
00608               if(idHep == forced_Hep[k])
00609                 { 
00610                   found = 1;
00611                   Daughter->deleteDaughters();
00612                 }
00613             }
00614           if (!found) go_through_daughters(Daughter);
00615         }
00616     }
00617 }
00618 
00619 void 
00620 EvtGenInterface::update_candlist( int theIndex, HepMC::GenParticle *thePart )
00621 {
00622   if(nlist<10)                 // not nice ... but is 10 a reasonable maximum?
00623      {
00624        bool isThere = false;
00625        if (nlist) {
00626          for (int check=0; check < nlist; check++) {
00627            if (listp[check]->barcode() == thePart->barcode()) isThere = true;
00628          }
00629        }
00630        if (!isThere) { 
00631          listp[nlist] = thePart;
00632          index[nlist++] = theIndex;
00633        }
00634      }
00635   else
00636     {
00637       throw cms::Exception("runtime")
00638         << "more than 10 candidates to be forced ";
00639       return; 
00640     }  
00641   return;
00642 }
00643