CMS 3D CMS Logo

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