CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/GeneratorInterface/CascadeInterface/plugins/Cascade2Hadronizer.cc

Go to the documentation of this file.
00001 #include "GeneratorInterface/CascadeInterface/plugins/Cascade2Hadronizer.h"
00002 
00003 #include "GeneratorInterface/Core/interface/RNDMEngineAccess.h"
00004 
00005 
00006 #include "HepMC/GenEvent.h" 
00007 #include "HepMC/PdfInfo.h"
00008 #include "HepMC/PythiaWrapper6_4.h"  
00009 #include "GeneratorInterface/CascadeInterface/plugins/CascadeWrapper.h"   
00010 
00011 #include "HepMC/HEPEVT_Wrapper.h"
00012 #include "HepMC/IO_HEPEVT.h"
00013 #include "HepMC/IO_GenEvent.h"
00014 
00015 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00016 #include "GeneratorInterface/Core/interface/FortranCallback.h"
00017 
00018 HepMC::IO_HEPEVT hepevtio;
00019 
00020 #include "HepPID/ParticleIDTranslations.hh"
00021 
00022 #include "FWCore/ServiceRegistry/interface/Service.h"
00023 #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
00024 #include "CLHEP/Random/RandFlat.h"
00025 #include "FWCore/Utilities/interface/Exception.h"
00026 
00027 //-- Pythia6 routines and functionalities to pass around Pythia6 params
00028 
00029 #include "GeneratorInterface/Pythia6Interface/interface/Pythia6Service.h"
00030 #include "GeneratorInterface/Pythia6Interface/interface/Pythia6Declarations.h"
00031 
00032 #include "SimDataFormats/GeneratorProducts/interface/GenRunInfoProduct.h"
00033 
00034 using namespace edm;
00035 using namespace std;
00036 
00037 #define debug 0
00038 
00039 CLHEP::RandFlat* fFlat_extern;
00040 
00041 extern "C" {
00042     
00043   double dcasrn_(int *idummy) {
00044 
00045     static int call = 0;
00046 
00047     double rdm_nb = fFlat_extern->fire();
00048 
00049     if(debug && ++call < 100) cout<<"dcasrn from c++, call: "<<call<<" random number: "<<rdm_nb<<endl;
00050 
00051     return rdm_nb;
00052   }
00053 }
00054   
00055 
00056 namespace gen {
00057   
00058   class Pythia6ServiceWithCallback : public Pythia6Service {
00059     
00060   public:
00061     Pythia6ServiceWithCallback(const edm::ParameterSet& pset) : Pythia6Service(pset) {}
00062     
00063   private:
00064     void upInit() {
00065       FortranCallback::getInstance()->fillHeader(); 
00066     }
00067     
00068     void upEvnt() {
00069       FortranCallback::getInstance()->fillEvent(); 
00070     }
00071     
00072     bool upVeto() { 
00073       bool veto = false;
00074       if (!hepeup_.nup) veto = true; //-- LHE Common Blocks
00075       return(veto);
00076     }
00077 
00078   };
00079   
00080   static struct {
00081     int n, npad, k[5][pyjets_maxn];
00082     double p[5][pyjets_maxn], v[5][pyjets_maxn];
00083   } pyjets_local;
00084   
00085   Cascade2Hadronizer::Cascade2Hadronizer(edm::ParameterSet const& pset) 
00086     : BaseHadronizer(pset),
00087       fPy6Service(new Pythia6ServiceWithCallback(pset)), //-- this will store py6 parameters for further settings 
00088       
00089       //-- fRandomEngine(&getEngineReference()),
00090 
00091       //-- defined in GeneratorInterface/Core/src/RNDMEngineAccess.cc
00092       //-- CLHEP::HepRandomEngine& gen::getEngineReference()
00093       //-- { edm::Service<edm::RandomNumberGenerator> rng;
00094       //--  return rng->getEngine(); }
00095 
00096       fComEnergy(pset.getParameter<double>("comEnergy")),
00097       fextCrossSection(pset.getUntrackedParameter<double>("crossSection",-1.)),
00098       fextCrossSectionError(pset.getUntrackedParameter<double>("crossSectionError",-1.)),
00099       fFilterEfficiency(pset.getUntrackedParameter<double>("filterEfficiency",-1.)),
00100       
00101       fMaxEventsToPrint(pset.getUntrackedParameter<int>("maxEventsToPrint",0)),
00102       fHepMCVerbosity(pset.getUntrackedParameter<bool>("pythiaHepMCVerbosity",false)),
00103       fPythiaListVerbosity(pset.getUntrackedParameter<int>("pythiaPylistVerbosity",0)),
00104       
00105       fDisplayPythiaBanner(pset.getUntrackedParameter<bool>("displayPythiaBanner",false)),
00106       fDisplayPythiaCards(pset.getUntrackedParameter<bool>("displayPythiaCards",false)){
00107 
00108     fParameters = pset.getParameter<ParameterSet>("Cascade2Parameters");
00109 
00110     edm::Service<edm::RandomNumberGenerator> rng;
00111     if(debug) cout<<"seed: "<<rng->mySeed()<<endl;
00112     CLHEP::HepRandomEngine& engine = rng->getEngine();
00113     fFlat = new CLHEP::RandFlat(engine);
00114     fFlat_extern = fFlat;
00115 
00116     fConvertToPDG = false;
00117     if(pset.exists("doPDGConvert"))
00118       fConvertToPDG = pset.getParameter<bool>("doPDGConvert");
00119     
00120     //-- silence Pythia6 banner printout, unless display requested
00121    
00122     if(!fDisplayPythiaBanner){
00123       if(!call_pygive("MSTU(12)=12345")){
00124         throw edm::Exception(edm::errors::Configuration,"PythiaError") 
00125           <<" Pythia did not accept MSTU(12)=12345";
00126       }
00127     }
00128     
00129     //-- silence printouts from PYGIVE, unless display requested
00130     
00131     if(!fDisplayPythiaCards){
00132       if(!call_pygive("MSTU(13)=0")){
00133         throw edm::Exception(edm::errors::Configuration,"PythiaError") 
00134           <<" Pythia did not accept MSTU(13)=0";
00135       }
00136     }
00137 
00138     //-- tmp stuff to deal with EvtGen corrupting pyjets
00139     flushTmpStorage();
00140   }
00141 
00142   Cascade2Hadronizer::~Cascade2Hadronizer(){
00143     if(fPy6Service != 0) delete fPy6Service;
00144   }
00145   
00146   void Cascade2Hadronizer::flushTmpStorage(){
00147     
00148     pyjets_local.n = 0 ;
00149     pyjets_local.npad = 0 ;
00150    
00151     for(int ip=0; ip<pyjets_maxn; ip++){
00152       for(int i=0; i<5; i++){
00153         pyjets_local.k[i][ip] = 0;
00154         pyjets_local.p[i][ip] = 0.;
00155         pyjets_local.v[i][ip] = 0.;
00156       }
00157     }
00158     return;
00159   }
00160   
00161   void Cascade2Hadronizer::fillTmpStorage(){
00162 
00163     pyjets_local.n = pyjets.n;
00164     pyjets_local.npad = pyjets.npad;
00165 
00166     for(int ip=0; ip<pyjets_maxn; ip++){
00167       for(int i=0; i<5; i++){
00168         pyjets_local.k[i][ip] = pyjets.k[i][ip];
00169         pyjets_local.p[i][ip] = pyjets.p[i][ip];
00170         pyjets_local.v[i][ip] = pyjets.v[i][ip];
00171       }
00172     }
00173     return ;
00174   }
00175 
00176   void Cascade2Hadronizer::finalizeEvent(){
00177 
00178     HepMC::PdfInfo pdf;
00179 
00180     //-- filling factorization "Q scale" now! pthat moved to binningValues()
00181              
00182     if(event()->signal_process_id() <= 0) event()->set_signal_process_id(pypars.msti[0]);
00183     if(event()->event_scale() <=0 )       event()->set_event_scale(pypars.pari[22]);
00184     if(event()->alphaQED() <= 0)          event()->set_alphaQED(pyint1.vint[56]);
00185     if(event()->alphaQCD() <= 0)          event()->set_alphaQCD(pyint1.vint[57]);
00186    
00187    //-- get pdf info directly from Pythia6 and set it up into HepMC::GenEvent
00188    //-- S. Mrenna: Prefer vint block
00189    
00190    if(pdf.id1() <= 0)      pdf.set_id1(pyint1.mint[14] == 21 ? 0 : pyint1.mint[14]);
00191    if(pdf.id2() <= 0)      pdf.set_id2(pyint1.mint[15] == 21 ? 0 : pyint1.mint[15]);
00192    if(pdf.x1() <= 0)       pdf.set_x1(pyint1.vint[40]);
00193    if(pdf.x2() <= 0)       pdf.set_x2(pyint1.vint[41]);
00194    if(pdf.pdf1() <= 0)     pdf.set_pdf1(pyint1.vint[38] / pyint1.vint[40]);
00195    if(pdf.pdf2() <= 0)     pdf.set_pdf2(pyint1.vint[39] / pyint1.vint[41]);
00196    if(pdf.scalePDF() <= 0) pdf.set_scalePDF(pyint1.vint[50]);   
00197  
00198    event()->set_pdf_info(pdf) ;
00199 
00200    if(debug) {
00201      cout<<"standard Py6 event weight: pyint1.vint[96]: "<<pyint1.vint[96]<<endl;
00202      cout<<"event weight returned by PYEVWT: 1./(pyint1.vint[98]): "<<1./(pyint1.vint[98])<<endl;
00203    }
00204      
00205    //-- this is "standard" Py6 event weight (corresponds to PYINT1/VINT(97))
00206    // event()->weights().push_back(pyint1.vint[96]);
00207  
00208    //-- this is event weight as 1./VINT(99) (PYINT1/VINT(99) is returned by the PYEVWT) 
00209    // event()->weights().push_back(1./(pyint1.vint[98]));
00210 
00211    //-- all cascade events have weight = 1
00212    event()->weights().push_back(1.);
00213 
00214    //-- now create the GenEventInfo product from the GenEvent and fill the missing pieces
00215 
00216    eventInfo().reset(new GenEventInfoProduct(event().get()));
00217 
00218    //-- in Pythia6 pthat is used to subdivide samples into different bins
00219    //-- in LHE mode the binning is done by the external ME generator
00220    //-- which is likely not pthat, so only filling it for Py6 internal mode
00221 
00222    eventInfo()->setBinningValues(vector<double>(1,pypars.pari[16]));
00223   
00224    //-- here we treat long-lived particles
00225   
00226    if (pydat1.mstj[21]==3 || pydat1.mstj[21]==4) imposeProperTime();
00227 
00228    //-- convert particle IDs Py6 -> PDG (if requested)
00229    
00230    if(fConvertToPDG){
00231      for(HepMC::GenEvent::particle_iterator part = event()->particles_begin(); 
00232          part != event()->particles_end(); ++part){
00233        (*part)->set_pdg_id(HepPID::translatePythiatoPDT((*part)->pdg_id()));
00234      }
00235    }
00236    
00237    //-- service printouts, if requested
00238    if(fMaxEventsToPrint > 0){
00239      fMaxEventsToPrint--;
00240      if(fPythiaListVerbosity) call_pylist(fPythiaListVerbosity);
00241      if(fHepMCVerbosity){
00242         cout <<"Event process = "<<pypars.msti[0]<<endl<<"----------------------"<<endl;
00243        event()->print();
00244      }
00245    }
00246    
00247    return;
00248   }
00249 
00250   bool Cascade2Hadronizer::hadronize(){
00251     return false;
00252   }
00253 
00254   bool Cascade2Hadronizer::generatePartonsAndHadronize(){
00255     
00256     //-- grab Py6 instance
00257     Pythia6Service::InstanceWrapper guard(fPy6Service);
00258     
00259     FortranCallback::getInstance()->resetIterationsPerEvent();
00260     
00261     //-- skip event if py6 considers it bad
00262     if(pyint1.mint[50] != 0 ){
00263       event().reset();
00264       return false;
00265     }
00266     
00267     //-- generation of the event with CASCADE
00268     call_event(); 
00269 
00270     //-- pythia pyhepc routine converts common PYJETS in common HEPEVT
00271     call_pyhepc(1);
00272 
00273     //-- delete the created event from memory
00274     event().reset(hepevtio.read_next_event());
00275     
00276     //-- this is to deal with post-gen tools & residualDecay() that may reuse PYJETS
00277     flushTmpStorage();
00278     fillTmpStorage();
00279     
00280     return true;
00281   }
00282 
00283   bool Cascade2Hadronizer::decay(){
00284     return true;
00285   }
00286   
00287   bool Cascade2Hadronizer::residualDecay(){
00288     
00289     //-- grab Py6 instance
00290     Pythia6Service::InstanceWrapper guard(fPy6Service);
00291     
00292     int NPartsBeforeDecays = pyjets_local.n ;
00293     int NPartsAfterDecays = event().get()->particles_size();
00294     int barcode = NPartsAfterDecays;
00295       
00296     // JVY: well, in principle, it's not a 100% fair to go up to NPartsBeforeDecays,
00297     // because Photos will attach gamma's to existing vertices, i.e. in the middle
00298     // of the event rather than at the end; but this will only shift pointers down,
00299     // so we'll be going again over a few "original" particle...
00300     // in the alternative, we may go all the way up to the beginning of the event
00301     // and re-check if anything remains to decay, that's fine even if it'll take
00302     // some extra CPU...
00303     
00304     for(int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- ){
00305       
00306       HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );
00307       int status = part->status();
00308       if ( status != 1 ) continue; // check only "stable" particles, 
00309                                    // as some undecayed may still be marked as such
00310       int pdgid  = part->pdg_id();
00311       int py6ptr = pycomp_( pdgid );
00312       
00313       if(pydat3.mdcy[0][py6ptr-1] != 1) continue; // particle is not expected to decay
00314       
00315       int py6id = HepPID::translatePDTtoPythia(pdgid);
00316       
00317       // first, will need to zero out, then fill up PYJETS
00318       // I better do it directly (by hands) rather than via py1ent 
00319       // - to avoid (additional) mass smearing
00320       
00321       if(part->momentum().t() <= part->generated_mass()){
00322         continue ; // e == m -> 0-momentum, nothing to decay...
00323       }
00324     
00325       else{
00326         pyjets.n = 0;
00327         
00328         for(int i=0; i<5; i++){
00329           pyjets.k[i][0] = 0;
00330           pyjets.p[i][0] = 0.;
00331           pyjets.v[i][0] = 0.;
00332         }
00333         
00334         pyjets.k[0][0] = 1;
00335         pyjets.k[1][0] = py6id;
00336         pyjets.p[4][0] = part->generated_mass();
00337         pyjets.p[3][0] = part->momentum().t();
00338         pyjets.p[0][0] = part->momentum().x();
00339         pyjets.p[1][0] = part->momentum().y();
00340         pyjets.p[2][0] = part->momentum().z();
00341         HepMC::GenVertex* prod_vtx = part->production_vertex();
00342         if ( !prod_vtx ) continue; // in principle, should never happen but...
00343         pyjets.v[0][0] = prod_vtx->position().x();
00344         pyjets.v[1][0] = prod_vtx->position().y();
00345         pyjets.v[2][0] = prod_vtx->position().z();
00346         pyjets.v[3][0] = prod_vtx->position().t();
00347         pyjets.v[4][0] = 0.;
00348         pyjets.n = 1;
00349         pyjets.npad = pyjets_local.npad;
00350       }
00351       
00352       //-- now call Py6 decay routine
00353      
00354       int parent = 1; // since we always pass to Py6 a single particle
00355       pydecy_(parent);
00356       
00357       //-- now attach decay products to mother
00358       
00359       for(int iprt1=1; iprt1<pyjets.n; iprt1++){
00360 
00361         part->set_status( 2 );
00362         
00363         HepMC::GenVertex* DecVtx = new HepMC::GenVertex( HepMC::FourVector(pyjets.v[0][iprt1],
00364                                                                            pyjets.v[1][iprt1],
00365                                                                            pyjets.v[2][iprt1],
00366                                                                            pyjets.v[3][iprt1]) );
00367         DecVtx->add_particle_in( part ); // this will cleanup end_vertex if exists, replace with the new one
00368         // I presume (vtx) barcode will be given automatically
00369         
00370         HepMC::FourVector  pmom(pyjets.p[0][iprt1],pyjets.p[1][iprt1],
00371                                 pyjets.p[2][iprt1],pyjets.p[3][iprt1] );
00372         
00373         int dstatus = 0;
00374         if(pyjets.k[0][iprt1] >= 1 && pyjets.k[0][iprt1] <= 10) dstatus = 1;
00375         else if(pyjets.k[0][iprt1] >= 11 && pyjets.k[0][iprt1] <= 20) dstatus = 2;
00376         else if(pyjets.k[0][iprt1] >= 21 && pyjets.k[0][iprt1] <= 30) dstatus = 3;
00377         else if(pyjets.k[0][iprt1] >= 31 && pyjets.k[0][iprt1] <= 100) dstatus = pyjets.k[0][iprt1];
00378         
00379         HepMC::GenParticle* daughter = new HepMC::GenParticle(pmom,
00380                                                               HepPID::translatePythiatoPDT(pyjets.k[1][iprt1]),
00381                                                               dstatus);
00382         barcode++;
00383         daughter->suggest_barcode( barcode ); 
00384         DecVtx->add_particle_out( daughter );
00385         
00386         int iprt2;
00387         for(iprt2=iprt1+1; iprt2<pyjets.n; iprt2++){ // the pointer is shifted by -1, c++ style
00388          
00389           if(pyjets.k[2][iprt2] != parent){
00390             parent = pyjets.k[2][iprt2];
00391             break; // another parent particle; reset & break the loop
00392           }
00393                
00394           HepMC::FourVector  pmomN(pyjets.p[0][iprt2],pyjets.p[1][iprt2],
00395                                    pyjets.p[2][iprt2],pyjets.p[3][iprt2] );
00396                
00397           dstatus = 0;
00398           if(pyjets.k[0][iprt2] >= 1 && pyjets.k[0][iprt2] <= 10) dstatus = 1;
00399           else if(pyjets.k[0][iprt2] >= 11 && pyjets.k[0][iprt2] <= 20) dstatus = 2;
00400           else if(pyjets.k[0][iprt2] >= 21 && pyjets.k[0][iprt2] <= 30) dstatus = 3;
00401           else if(pyjets.k[0][iprt2] >= 31 && pyjets.k[0][iprt2] <= 100) dstatus = pyjets.k[0][iprt2];
00402           
00403           HepMC::GenParticle* daughterN = new HepMC::GenParticle(pmomN,
00404                                                                  HepPID::translatePythiatoPDT( pyjets.k[1][iprt2] ),
00405                                                                  dstatus);
00406           barcode++;
00407           daughterN->suggest_barcode( barcode ); 
00408           DecVtx->add_particle_out( daughterN );             
00409         } //-- end iprt2 loop
00410          
00411         iprt1 = iprt2-1; // reset counter such that it doesn't go over the same child more than once
00412         // don't forget to offset back into c++ counting, as it's already +1 forward
00413         
00414         event().get()->add_vertex( DecVtx );
00415       } //-- end iprt1 loop
00416     } //-- end loop over decay products
00417    
00418     // now restore the very original Py6 event record 
00419   
00420     if(pyjets_local.n != pyjets.n){
00421    
00422       // restore pyjets to its state as it was before external decays -
00423       // might have been jammed by action above or by py1ent calls in EvtGen
00424 
00425       pyjets.n = pyjets_local.n;
00426       pyjets.npad = pyjets_local.npad;
00427       for(int ip=0; ip<pyjets_local.n; ip++){
00428         for(int i=0; i<5; i++){
00429           pyjets.k[i][ip] = pyjets_local.k[i][ip];
00430           pyjets.p[i][ip] = pyjets_local.p[i][ip];
00431           pyjets.v[i][ip] = pyjets_local.v[i][ip];
00432         }
00433       }
00434     }  
00435          
00436    return true;
00437   }
00438 
00439   bool Cascade2Hadronizer::readSettings(int key) {
00440 
00441     //-- grab Py6 instance
00442     Pythia6Service::InstanceWrapper guard(fPy6Service); 
00443     
00444     fPy6Service->setGeneralParams();
00445     
00446     if (key == 0) fPy6Service->setCSAParams();
00447     fPy6Service->setSLHAParams();
00448     
00449     return true;
00450   }
00451 
00452   bool Cascade2Hadronizer::initializeForExternalPartons(){
00453     return false;
00454   }
00455 
00456   bool Cascade2Hadronizer::initializeForInternalPartons(){
00457 
00458     //-- grab Py6 instance
00459     Pythia6Service::InstanceWrapper guard(fPy6Service);
00460   
00461     //-- change standard parameters of JETSET/PYTHIA - replace call_pytcha()
00462   
00463     fPy6Service->setGeneralParams();   
00464 
00465     fPy6Service->setCSAParams();
00466     fPy6Service->setSLHAParams();
00467     fPy6Service->setPYUPDAParams(false);
00468 
00469     pythia6PrintParameters();
00470     
00471     //-- mstu(8) is set to NMXHEP in this dummy call (version >=6.404)
00472     call_pyhepc(1);
00473 
00474     //-- initialise random number generator: has been changed to be CMSSW compliant    
00475     //-- dcasrn overloaded: call to rluxgo not needed anymore
00476     //-- call_rluxgo(4,314159265,0,0);
00477     
00478     //-- initialise CASCADE parameters (default values)
00479     call_casini();
00480 
00481     //-- Read the parameters and pass them to the common blocks
00482     //-- call_steer();
00483 
00484     //-- initialise CASCADE parameters (user values)
00485 
00486     //-- retrieve all the different sets 
00487     vector<string> AllSets = fParameters.getParameter<vector<string> >("parameterSets");
00488 
00489     //-- loop over the different sets
00490     for(unsigned i=0; i<AllSets.size();++i) {
00491     
00492       string Set = AllSets[i];
00493       vector<string> Para_Set = fParameters.getParameter<vector<string> >(Set);
00494 
00495       //-- loop over all the parameters and stop in case of mistake
00496       for(vector<string>::const_iterator itPara = Para_Set.begin(); itPara != Para_Set.end(); ++itPara) {
00497         
00498         if(!cascadeReadParameters(*itPara)) {
00499           throw edm::Exception(edm::errors::Configuration,"CascadeError")
00500             <<" Cascade did not accept the following parameter: \""<<*itPara<<"\""<<endl;
00501         }
00502       } //-- end loop over all the parameters
00503     } //-- end loop over the different sets
00504 
00505     cainpu.plepin = -fComEnergy/2; 
00506     cainpu.ppin = fComEnergy/2;
00507 
00508     cascadePrintParameters();
00509    
00510     //-- change standard parameters of CASCADE
00511     call_cascha();
00512     
00513     //-- change standard parameters of JETSET/PYTHIA
00514     //-- call_pytcha();
00515 
00516     //-- set up for running CASCADE (integration of the cross-section)
00517     call_cascade(); 
00518 
00519     //-- print cross-section result from integration
00520     call_caend(1);
00521  
00522     fPy6Service->setPYUPDAParams(true);
00523     
00524     fPy6Service->closeSLHA();
00525     
00526     return true;
00527   }
00528 
00529   bool Cascade2Hadronizer::declareStableParticles(vector<int> pdg){
00530    
00531     for(size_t i=0; i<pdg.size(); i++){
00532       int PyID = HepPID::translatePDTtoPythia(pdg[i]);
00533       int pyCode = pycomp_( PyID );
00534 
00535       if(pyCode > 0){
00536         ostringstream pyCard ;
00537         pyCard << "MDCY(" << pyCode << ",1)=0";
00538         //-- cout << "pdg= " << pdg[i] << " " << pyCard.str() << endl; 
00539 
00540         call_pygive(pyCard.str());
00541       }
00542     }
00543       
00544     return true;
00545   }
00546 
00547   void Cascade2Hadronizer::imposeProperTime(){
00548     
00549     // this is practically a copy/paste of the original code by J.Alcaraz, 
00550     // taken directly from PythiaSource
00551     
00552     int dumm=0;
00553     HepMC::GenEvent::vertex_const_iterator vbegin = event()->vertices_begin();
00554     HepMC::GenEvent::vertex_const_iterator vend = event()->vertices_end();
00555     HepMC::GenEvent::vertex_const_iterator vitr = vbegin;
00556     
00557     for(; vitr != vend; ++vitr){
00558       
00559       HepMC::GenVertex::particle_iterator pbegin = (*vitr)->particles_begin(HepMC::children);
00560       HepMC::GenVertex::particle_iterator pend = (*vitr)->particles_end(HepMC::children);
00561       HepMC::GenVertex::particle_iterator pitr = pbegin;
00562       
00563       for(; pitr != pend; ++pitr){
00564         
00565         if((*pitr)->end_vertex()) continue;
00566         if((*pitr)->status()!=1) continue;
00567         
00568         int pdgcode= abs((*pitr)->pdg_id());
00569         // Do nothing if the particle is not expected to decay
00570         if ( pydat3.mdcy[0][pycomp_(pdgcode)-1] !=1 ) continue;
00571       
00572         double ctau = pydat2.pmas[3][pycomp_(pdgcode)-1];
00573         HepMC::FourVector mom = (*pitr)->momentum();
00574         HepMC::FourVector vin = (*vitr)->position();
00575         double x = 0.;
00576         double y = 0.;
00577         double z = 0.;
00578         double t = 0.;
00579         bool decayInRange = false;
00580         while(!decayInRange){
00581 
00582           double unif_rand = fPy6Service->call(pyr_, &dumm);
00583           // Value of 0 is excluded, so following line is OK
00584           double proper_length = - ctau * log(unif_rand);
00585           double factor = proper_length/mom.m();
00586           x = vin.x() + factor * mom.px();
00587           y = vin.y() + factor * mom.py();
00588           z = vin.z() + factor * mom.pz();
00589           t = vin.t() + factor * mom.e();
00590           
00591           // Decay must be happen outside a cylindrical region
00592           if(pydat1.mstj[21]==4){
00593             if (sqrt(x*x+y*y)>pydat1.parj[72] || fabs(z)>pydat1.parj[73]) decayInRange = true;
00594             // Decay must be happen outside a given sphere
00595           } 
00596           else if (pydat1.mstj[21]==3){
00597             if (sqrt(x*x+y*y+z*z)>pydat1.parj[71]) decayInRange = true;
00598           } 
00599           // Decay is always OK otherwise
00600           else {
00601             decayInRange = true;
00602           }
00603         }
00604                   
00605         HepMC::GenVertex* vdec = new HepMC::GenVertex(HepMC::FourVector(x,y,z,t));
00606         event()->add_vertex(vdec);
00607         vdec->add_particle_in((*pitr));
00608       }
00609     }   
00610     
00611     return;
00612   }
00613 
00614   void Cascade2Hadronizer::statistics(){
00615 
00616     runInfo().setExternalXSecLO(GenRunInfoProduct::XSec(fextCrossSection,fextCrossSectionError));
00617     
00618     runInfo().setExternalXSecNLO(GenRunInfoProduct::XSec(-1.,-1.));   
00619 
00620     if(!runInfo().internalXSec()){
00621       double sigma_CMS = 0.001 * caeffic.avgi;
00622       double error_CMS = 0.001 * caeffic.sd;
00623       runInfo().setInternalXSec(GenRunInfoProduct::XSec(sigma_CMS,error_CMS));
00624     }
00625 
00626     //--  print out generated event summary
00627     call_caend(2);
00628 
00629     //-- write out some information from Pythia to the screen
00630     call_pystat(1);
00631   
00632     return;
00633   }
00634 
00635   const char* Cascade2Hadronizer::classname() const {
00636     return "gen::Cascade2Hadronizer";
00637   }
00638 
00639   //-- Read the parameters and pass them to the common blocks
00640  
00641   bool Cascade2Hadronizer::cascadeReadParameters(const string& ParameterString) {
00642 
00643     bool accepted = 1;
00644 
00645     if(!strncmp(ParameterString.c_str(),"KE",2))
00646       caluco.ke = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00647 
00648     else if(!strncmp(ParameterString.c_str(),"IRES(1)",7))
00649       capar6.ires[0] = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00650 
00651     else if(!strncmp(ParameterString.c_str(),"KP",2))
00652       caluco.kp = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00653 
00654     else if(!strncmp(ParameterString.c_str(),"IRES(2)",7))
00655       capar6.ires[1] = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00656 
00657     else if(!strncmp(ParameterString.c_str(),"NFRAG",5))
00658       cainpu.nfrag = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00659 
00660     else if(!strncmp(ParameterString.c_str(),"IPST",4))
00661       cashower.ipst = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00662 
00663     else if(!strncmp(ParameterString.c_str(),"IPSIPOL",7))
00664       jpsi.ipsipol = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00665 
00666     //-- from version 2.2.03 on
00667     // else if(!strncmp(ParameterString.c_str(),"I23S",4))
00668     //  jpsi.i23s = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00669     
00670     else if(!strncmp(ParameterString.c_str(),"IFPS",4))
00671       cainpu.ifps = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00672 
00673     else if(!strncmp(ParameterString.c_str(),"ITIMSHR",7))
00674       casshwr.itimshr = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00675 
00676     else if(!strncmp(ParameterString.c_str(),"IRUNAEM",7))
00677       capar1.irunaem = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00678 
00679     else if(!strncmp(ParameterString.c_str(),"IRUNA",5))
00680       capar1.iruna = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00681 
00682     else if(!strncmp(ParameterString.c_str(),"IQ2",3))
00683       capar1.iq2 = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00684 
00685     else if(!strncmp(ParameterString.c_str(),"IPRO",4))
00686       capar1.ipro = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00687 
00688     else if(!strncmp(ParameterString.c_str(),"NFLAV",5))
00689       caluco.nflav = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00690 
00691     else if(!strncmp(ParameterString.c_str(),"INTER",5))
00692       cainpu.inter = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00693 
00694     else if(!strncmp(ParameterString.c_str(),"IHFLA",5))
00695       cahflav.ihfla = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00696 
00697     else if(!strncmp(ParameterString.c_str(),"IRPA",4))
00698       cascol.irpa = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00699 
00700     else if(!strncmp(ParameterString.c_str(),"IRPB",4))
00701       cascol.irpb = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00702 
00703     else if(!strncmp(ParameterString.c_str(),"IRPC",4))
00704       cascol.irpc = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00705 
00706     else if(!strncmp(ParameterString.c_str(),"ICCFM",5))
00707       casshwr.iccfm = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00708 
00709     else if(!strncmp(ParameterString.c_str(),"IFINAL",6))
00710       cainpu.ifinal = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00711 
00712     else if(!strncmp(ParameterString.c_str(),"IGLU",4))
00713       cagluon.iglu = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00714 
00715     else if(!strncmp(ParameterString.c_str(),"IRspl",5))
00716       casprre.irspl = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00717 
00718     else if(!strncmp(ParameterString.c_str(),"PT2CUT",6))
00719       captcut.pt2cut[capar1.ipro-1] = atof(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00720 
00721     else if(!strncmp(ParameterString.c_str(),"NCB",3))
00722       integr.ncb = atoi(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00723 
00724     else if(!strncmp(ParameterString.c_str(),"ACC1",4))
00725       integr.acc1 = atof(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00726 
00727     else if(!strncmp(ParameterString.c_str(),"ACC2",4))
00728       integr.acc2 = atof(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);      
00729 
00730     else if(!strncmp(ParameterString.c_str(),"SCALFAF",7))
00731       scalf.scalfaf = atof(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00732 
00733     else if(!strncmp(ParameterString.c_str(),"SCALFA",6))
00734       scalf.scalfa = atof(&ParameterString[strcspn(ParameterString.c_str(),"=")+1]);
00735 
00736     else accepted = 0;
00737 
00738     return accepted;
00739   }
00740   
00741   void Cascade2Hadronizer::cascadePrintParameters() {
00742 
00743     cout<<""<<endl;
00744     cout<<"----------------------------------"<<endl;
00745     cout<<"----    Cascade Parameters    ----"<<endl;
00746     cout<<"----------------------------------"<<endl;
00747     cout<<""<<endl;
00748 
00749     cout<<"computation switch: "<<endl;
00750     cout<<"number of calls per iteration for bases: "<<integr.ncb<<endl;
00751     cout<<"relative precision for grid optimisation: "<<integr.acc1<<" %"<<endl;
00752     cout<<"relative precision for integration: "<<integr.acc2<<" %"<<endl;
00753     cout<<""<<endl;
00754 
00755     cout<<"kinematics: "<<endl;
00756     cout<<"flavour code of beam 1: "<<caluco.ke<<endl; 
00757     cout<<"direct or resolved particle 1: "<<capar6.ires[0]<<endl; 
00758     cout<<"pz of incoming beam 1: "<<cainpu.plepin<<" GeV"<<endl;
00759     cout<<"flavour code of beam 2: "<<caluco.kp<<endl;    
00760     cout<<"direct or resolved particle 2: "<<capar6.ires[1]<<endl;  
00761     cout<<"pz of incoming beam 2: "<<cainpu.ppin<<" GeV"<<endl;
00762     cout<<"number of active flavors: "<<caluco.nflav<<endl;
00763     cout<<""<<endl;
00764 
00765     cout<<"hard subprocess selection: "<<endl;
00766     cout<<"hard subprocess number (IPRO): "<<capar1.ipro<<endl;
00767     cout<<"IPRO = 10, switch to select QCD process g* g* -> q qbar: "<<cascol.irpa<<endl;
00768     cout<<"IPRO = 10, switch to select QCD process g* g -> g g: "<<cascol.irpb<<endl;
00769     cout<<"IPRO = 10, switch to select QCD process g* q -> g q: "<<cascol.irpc<<endl;
00770     cout<<"flavor of heavy quark produced (IPRO = 11, 504, 514): "<<cahflav.ihfla<<endl;
00771     // cout<<"select vector meson state: "<<jpsi.i23s<<endl;
00772     cout<<"use matrix element including J/psi (Upsilon) polarisation: "<<jpsi.ipsipol<<endl;
00773     cout<<"pt2 cut in matrix element for massless partons: "<<captcut.pt2cut[capar1.ipro-1]<<" GeV^2"<<endl;
00774     cout<<""<<endl;
00775 
00776     cout<<"parton shower and fragmentation: "<<endl;
00777     cout<<"switch for fragmentation: "<<cainpu.nfrag<<endl;
00778     cout<<"switch for parton shower (0 off - 1 initial state - 2 final state - 3 initial & final state): "<<cainpu.ifps<<endl;
00779     cout<<"switch for time like parton shower in initial state: "<<casshwr.itimshr<<endl;
00780     cout<<"select CCFM (1) or DGLAP (0) evolution: "<<casshwr.iccfm<<endl;
00781     cout<<"scale switch for final state parton shower: "<<cainpu.ifinal<<endl;
00782     cout<<"scale factor for final state parton shower: "<<scalf.scalfaf<<endl;
00783     cout<<"switch for proton remnant treatment: "<<casprre.irspl<<endl;
00784     cout<<"keep track of intermediate state in parton shower: "<<cashower.ipst<<endl;
00785     cout<<"mode of interaction for e p: "<<cainpu.inter<<endl;
00786     cout<<""<<endl;
00787 
00788     cout<<"pdfs, couplings and scales: "<<endl;
00789     cout<<"switch for running alphaem: "<<capar1.irunaem<<endl;
00790     cout<<"switch for running alphas: "<<capar1.iruna<<endl;
00791     cout<<"scale in alphas: "<<capar1.iq2<<endl;
00792     cout<<"scale factor for scale in alphas: "<<scalf.scalfa<<endl;
00793     cout<<"unintegrated pdf: "<<cagluon.iglu<<endl;
00794     cout<<"path where updf are stored: "<<caspdf.pdfpath<<endl;
00795     cout<<""<<endl;
00796 
00797     cout<<"center of mass energy, cross section, filter efficiency: "<<endl;
00798     cout<<"center of mass energy: "<<fComEnergy<<" GeV"<<endl;
00799     cout<<"external LO cross section: "<<fextCrossSection<<" +/- "<<fextCrossSectionError<<" pb"<<endl;
00800     cout<<"filter efficiency: "<<fFilterEfficiency<<endl;
00801     cout<<""<<endl;
00802   }
00803 
00804   void Cascade2Hadronizer::pythia6PrintParameters() {
00805 
00806     cout<<""<<endl;
00807     cout<<"----------------------------------"<<endl;
00808     cout<<"----    Pythia6 Parameters    ----"<<endl;
00809     cout<<"----------------------------------"<<endl;
00810     cout<<""<<endl;
00811 
00812     cout<<"charm mass: "<<pydat2.pmas[0][3]<<" GeV (default value =  1.5 GeV)"<<endl;
00813     cout<<"bottom mass: "<<pydat2.pmas[0][4]<<" GeV (default value = 4.8 GeV)"<<endl;
00814     cout<<"Higgs mass: "<<pydat2.pmas[0][24]<<" GeV (default value = 115 GeV)"<<endl;
00815 
00816     cout<<"lambda QCD: "<<pydat1.paru[111]<<" GeV (default value = 0.25 GeV)"<<endl;
00817 
00818     cout<<"alphas behaviour: "<<pydat1.mstu[110]<<" (default value = 1)"<<endl;
00819     cout<<"nr of flavours wrt lambda QCD: "<<pydat1.mstu[111]<<" (default value = 5)"<<endl;
00820     cout<<"min nr of flavours for alphas: "<<pydat1.mstu[112]<<" (default value = 3)"<<endl;
00821     cout<<"max nr of flavours for alphas: "<<pydat1.mstu[113]<<" (default value = 5)"<<endl;
00822 
00823     cout<<"maximum angle settings: "<<pydat1.mstj[47]<<" (default value = 0)"<<endl;
00824 
00825     cout<<""<<endl;
00826   }
00827 
00828 } //-- namespace gen