CMS 3D CMS Logo

Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes

gen::Py8GunBase Class Reference

#include <Py8GunBase.h>

Inheritance diagram for gen::Py8GunBase:
gen::Py8InterfaceBase gen::Py8EGun gen::Py8JetGun gen::Py8PtGun

List of all members.

Public Member Functions

void finalizeEvent ()
edm::EventgetEDMEvent () const
HepMC::GenEvent * getGenEvent ()
GenEventInfoProductgetGenEventInfo ()
GenRunInfoProductgetGenRunInfo ()
bool initializeForInternalPartons ()
 Py8GunBase (edm::ParameterSet const &ps)
void resetEvent (HepMC::GenEvent *event)
void resetEventInfo (GenEventInfoProduct *eventInfo)
virtual bool residualDecay ()
virtual bool select (HepMC::GenEvent *) const
void setEDMEvent (edm::Event &event)
void statistics ()
 ~Py8GunBase ()

Protected Member Functions

std::auto_ptr< HepMC::GenEvent > & event ()
std::auto_ptr
< GenEventInfoProduct > & 
eventInfo ()
GenRunInfoProductrunInfo ()

Protected Attributes

double fMaxPhi
double fMinPhi
std::vector< int > fPartIDs

Private Attributes

edm::EventedmEvent_
std::auto_ptr< HepMC::GenEvent > genEvent_
std::auto_ptr
< GenEventInfoProduct
genEventInfo_
GenRunInfoProduct genRunInfo_

Detailed Description

Definition at line 25 of file Py8GunBase.h.


Constructor & Destructor Documentation

gen::Py8GunBase::Py8GunBase ( edm::ParameterSet const &  ps)

Definition at line 9 of file Py8GunBase.cc.

References fMaxPhi, fMinPhi, fPartIDs, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), runInfo(), GenRunInfoProduct::setExternalXSecLO(), GenRunInfoProduct::setExternalXSecNLO(), and GenRunInfoProduct::setFilterEfficiency().

   : Py8InterfaceBase(ps)
{

   runInfo().setFilterEfficiency(
      ps.getUntrackedParameter<double>("filterEfficiency", -1.) );
   runInfo().setExternalXSecLO(
      GenRunInfoProduct::XSec(ps.getUntrackedParameter<double>("crossSection", -1.)) );
   runInfo().setExternalXSecNLO(
       GenRunInfoProduct::XSec(ps.getUntrackedParameter<double>("crossSectionNLO", -1.)) );

  
  // PGun specs
  //
   edm::ParameterSet pgun_params = 
      ps.getParameter<edm::ParameterSet>("PGunParameters");
      
   // although there's the method ParameterSet::empty(),  
   // it looks like it's NOT even necessary to check if it is,
   // before trying to extract parameters - if it is empty,
   // the default values seem to be taken
   //
   fPartIDs    = pgun_params.getParameter< std::vector<int> >("ParticleID");
   fMinPhi     = pgun_params.getParameter<double>("MinPhi"); // ,-3.14159265358979323846);
   fMaxPhi     = pgun_params.getParameter<double>("MaxPhi"); // , 3.14159265358979323846);
   
}
gen::Py8GunBase::~Py8GunBase ( ) [inline]

Definition at line 28 of file Py8GunBase.h.

{}

Member Function Documentation

std::auto_ptr<HepMC::GenEvent>& gen::Py8GunBase::event ( ) [inline, protected]
std::auto_ptr<GenEventInfoProduct>& gen::Py8GunBase::eventInfo ( ) [inline, protected]

Definition at line 51 of file Py8GunBase.h.

References genEventInfo_.

{ return genEventInfo_; }
void gen::Py8GunBase::finalizeEvent ( ) [virtual]

Implements gen::Py8InterfaceBase.

Definition at line 138 of file Py8GunBase.cc.

References gather_cfg::cout, event(), gen::Py8InterfaceBase::fMasterGen, gen::Py8InterfaceBase::maxEventsToPrint, gen::Py8InterfaceBase::pythiaHepMCVerbosity, and gen::Py8InterfaceBase::pythiaPylistVerbosity.

{
   

// FIXME !!!

  //******** Verbosity ********

   if (maxEventsToPrint > 0 &&
      (pythiaPylistVerbosity || pythiaHepMCVerbosity)) 
   {
      maxEventsToPrint--;
      if (pythiaPylistVerbosity) 
      {
         fMasterGen->info.list(std::cout); 
         fMasterGen->event.list(std::cout);
      } 

      if (pythiaHepMCVerbosity) 
      {
         std::cout << "Event process = "
                   << fMasterGen->info.code() << "\n"
                   << "----------------------" << std::endl;
         event()->print();
      }
   }
      
   return;
}
edm::Event& gen::Py8GunBase::getEDMEvent ( ) const [inline]

Definition at line 40 of file Py8GunBase.h.

References edmEvent_.

{ return *edmEvent_; }
HepMC::GenEvent* gen::Py8GunBase::getGenEvent ( ) [inline]

Definition at line 32 of file Py8GunBase.h.

References genEvent_.

{ return genEvent_.release(); }
GenEventInfoProduct* gen::Py8GunBase::getGenEventInfo ( ) [inline]

Definition at line 33 of file Py8GunBase.h.

References genEventInfo_.

{ return genEventInfo_.release(); }
GenRunInfoProduct& gen::Py8GunBase::getGenRunInfo ( ) [inline]

Definition at line 31 of file Py8GunBase.h.

References genRunInfo_.

{ return genRunInfo_; }
bool gen::Py8GunBase::initializeForInternalPartons ( ) [virtual]

Implements gen::Py8InterfaceBase.

Definition at line 39 of file Py8GunBase.cc.

References gen::Py8InterfaceBase::fDecayer, and gen::Py8InterfaceBase::fMasterGen.

{

   // NO MATTER what was this setting below, override it before init 
   // - this is essencial for the PGun mode 
   
   // Key requirement: switch off ProcessLevel, and thereby also PartonLevel.
   fMasterGen->readString("ProcessLevel:all = off");
   fMasterGen->readString("Standalone:allowResDec=on");
   // pythia->readString("ProcessLevel::resonanceDecays=on");
   fMasterGen->init();
   
   // init decayer
   fDecayer->readString("ProcessLevel:all = off"); // Same trick!
   fDecayer->readString("Standalone:allowResDec=on");
   // pythia->readString("ProcessLevel::resonanceDecays=on");
   fDecayer->init();
  
   return true;

}
void gen::Py8GunBase::resetEvent ( HepMC::GenEvent *  event) [inline]

Definition at line 35 of file Py8GunBase.h.

References genEvent_.

{ genEvent_.reset(event); }
void gen::Py8GunBase::resetEventInfo ( GenEventInfoProduct eventInfo) [inline]

Definition at line 36 of file Py8GunBase.h.

References genEventInfo_.

{ genEventInfo_.reset(eventInfo); }
bool gen::Py8GunBase::residualDecay ( ) [virtual]

Definition at line 61 of file Py8GunBase.cc.

References event(), gen::Py8InterfaceBase::fDecayer, gen::Py8InterfaceBase::fMasterGen, and configurableAnalysis::GenParticle.

{

   Event* pythiaEvent = &(fMasterGen->event);
  
   int NPartsBeforeDecays = pythiaEvent->size()-1; // do NOT count the very 1st "system" particle 
                                                   // in Pythia8::Event record; it does NOT even
                                                   // get translated by the HepMCInterface to the
                                                   // HepMC::GenEvent record !!!
   int NPartsAfterDecays = event().get()->particles_size();
   int NewBarcode = NPartsAfterDecays;
   
   for ( int ipart=NPartsAfterDecays; ipart>NPartsBeforeDecays; ipart-- )
   {

      HepMC::GenParticle* part = event().get()->barcode_to_particle( ipart );

      if ( part->status() == 1 )
      {
         fDecayer->event.reset();
         Particle py8part(  part->pdg_id(), 93, 0, 0, 0, 0, 0, 0,
                         part->momentum().x(),
                         part->momentum().y(),
                         part->momentum().z(),
                         part->momentum().t(),
                         part->generated_mass() );
         HepMC::GenVertex* ProdVtx = part->production_vertex();
         py8part.vProd( ProdVtx->position().x(), ProdVtx->position().y(), 
                     ProdVtx->position().z(), ProdVtx->position().t() );
         py8part.tau( (fDecayer->particleData).tau0( part->pdg_id() ) );
         fDecayer->event.append( py8part );
         int nentries = fDecayer->event.size();
         if ( !fDecayer->event[nentries-1].mayDecay() ) continue;
         fDecayer->next();
         int nentries1 = fDecayer->event.size();
         if ( nentries1 <= nentries ) continue; //same number of particles, no decays...
            
         part->set_status(2);
            
         Particle& py8daughter = fDecayer->event[nentries]; // the 1st daughter
         HepMC::GenVertex* DecVtx = new HepMC::GenVertex( HepMC::FourVector(py8daughter.xProd(),
                                                          py8daughter.yProd(),
                                                          py8daughter.zProd(),
                                                          py8daughter.tProd()) );

         DecVtx->add_particle_in( part ); // this will cleanup end_vertex if exists, replace with the new one
                                          // I presume (vtx) barcode will be given automatically
            
         HepMC::FourVector pmom( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );
            
         HepMC::GenParticle* daughter =
                             new HepMC::GenParticle( pmom, py8daughter.id(), 1 );
            
         NewBarcode++;
         daughter->suggest_barcode( NewBarcode );
         DecVtx->add_particle_out( daughter );
                    
         for ( ipart=nentries+1; ipart<nentries1; ipart++ )
         {
            py8daughter = fDecayer->event[ipart];
            HepMC::FourVector pmomN( py8daughter.px(), py8daughter.py(), py8daughter.pz(), py8daughter.e() );       
            HepMC::GenParticle* daughterN =
                                new HepMC::GenParticle( pmomN, py8daughter.id(), 1 );
            NewBarcode++;
            daughterN->suggest_barcode( NewBarcode );
            DecVtx->add_particle_out( daughterN );
         }
            
         event().get()->add_vertex( DecVtx );

      }
   } 

   return true;

}
GenRunInfoProduct& gen::Py8GunBase::runInfo ( ) [inline, protected]

Definition at line 49 of file Py8GunBase.h.

References genRunInfo_.

Referenced by Py8GunBase(), and statistics().

{ return genRunInfo_; }
virtual bool gen::Py8GunBase::select ( HepMC::GenEvent *  ) const [inline, virtual]

Definition at line 41 of file Py8GunBase.h.

{ return true;}
void gen::Py8GunBase::setEDMEvent ( edm::Event event) [inline]

Definition at line 39 of file Py8GunBase.h.

References edmEvent_, and event().

{ edmEvent_ = &event; }
void gen::Py8GunBase::statistics ( ) [virtual]

Reimplemented from gen::Py8InterfaceBase.

Definition at line 168 of file Py8GunBase.cc.

References gen::Py8InterfaceBase::fMasterGen, runInfo(), and GenRunInfoProduct::setInternalXSec().

{

   fMasterGen->statistics();

   double xsec = fMasterGen->info.sigmaGen(); // cross section in mb
   xsec *= 1.0e9; // translate to pb (CMS/Gen "convention" as of May 2009)
   runInfo().setInternalXSec(xsec);
   return;

}

Member Data Documentation

Definition at line 64 of file Py8GunBase.h.

Referenced by getEDMEvent(), and setEDMEvent().

double gen::Py8GunBase::fMaxPhi [protected]
double gen::Py8GunBase::fMinPhi [protected]
std::vector<int> gen::Py8GunBase::fPartIDs [protected]
std::auto_ptr<HepMC::GenEvent> gen::Py8GunBase::genEvent_ [private]

Definition at line 61 of file Py8GunBase.h.

Referenced by event(), getGenEvent(), and resetEvent().

Definition at line 62 of file Py8GunBase.h.

Referenced by eventInfo(), getGenEventInfo(), and resetEventInfo().

Definition at line 60 of file Py8GunBase.h.

Referenced by getGenRunInfo(), and runInfo().