CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ExternalDecayDriver.cc
Go to the documentation of this file.
2 
10 #include "HepMC/GenEvent.h"
12 // LHE Run
15 
16 // LHE Event
19 
20 using namespace gen;
21 using namespace edm;
22 
23 CLHEP::HepRandomEngine* decayRandomEngine = nullptr;
24 
26  : fIsInitialized(false),
27  fTauolaInterface(0),
28  fEvtGenInterface(0),
29  fPhotosInterface(0)
30 {
31  std::vector<std::string> extGenNames =
32  pset.getParameter< std::vector<std::string> >("parameterSets");
33 
34  for (unsigned int ip=0; ip<extGenNames.size(); ++ip ){
35  std::string curSet = extGenNames[ip];
36  if ( curSet == "EvtGen" || curSet == "EvtGenLHC91"){
37  fEvtGenInterface = (EvtGenInterfaceBase*)(EvtGenFactory::get()->create("EvtGenLHC91", pset.getUntrackedParameter< ParameterSet >(curSet)));
41  }
42  else if( curSet == "EvtGen130"){
49  }
50  else if ( curSet == "Tauola" || curSet == "Tauolapp114" ){
51  fTauolaInterface = (TauolaInterfaceBase*)(TauolaFactory::get()->create("Tauolapp114", pset.getUntrackedParameter< ParameterSet >(curSet)));
52  fPhotosInterface = (PhotosInterfaceBase*)(PhotosFactory::get()->create("Photos2155", pset.getUntrackedParameter< ParameterSet >(curSet)));
57  }
58  else if ( curSet == "Photos" || curSet == "Photos2155" ){
59  if ( !fPhotosInterface ){
62  }
63  }
64  else if (curSet == "Photospp355" ){
65  if ( !fPhotosInterface ){
66  fPhotosInterface = (PhotosInterfaceBase*)(PhotosFactory::get()->create("Photospp355", pset.getUntrackedParameter< ParameterSet>(curSet)));
68  }
69  }
70  }
71 }
72 
74 {
75  if ( fEvtGenInterface ) delete fEvtGenInterface;
76  if ( fTauolaInterface ) delete fTauolaInterface;
77  if ( fPhotosInterface ) delete fPhotosInterface;
78 }
79 
80 
81 HepMC::GenEvent* ExternalDecayDriver::decay( HepMC::GenEvent* evt )
82 {
83  if ( !fIsInitialized ) return evt;
84 
85  if ( fEvtGenInterface ){
86  evt = fEvtGenInterface->decay( evt );
87  if ( !evt ) return 0;
88  }
89 
90  if ( fTauolaInterface ){
91  evt = fTauolaInterface->decay( evt );
92  if ( !evt ) return 0;
93  }
94 
95  if ( fPhotosInterface ){
96  evt = fPhotosInterface->apply( evt );
97  if ( !evt ) return 0;
98  }
99 
100  return evt;
101 }
102 
104 {
105 
106  if ( fIsInitialized ) return;
107 
108  if ( fTauolaInterface ){
109  fTauolaInterface->init( es );
110  for ( std::vector<int>::const_iterator i=fTauolaInterface->operatesOnParticles().begin();
111  i!=fTauolaInterface->operatesOnParticles().end(); i++ )
112  fPDGs.push_back( *i );
113  }
114 
115  if ( fEvtGenInterface ){
117  for ( std::vector<int>::const_iterator i=fEvtGenInterface->operatesOnParticles().begin();
118  i!=fEvtGenInterface->operatesOnParticles().end(); i++ )
119  fPDGs.push_back( *i );
120  }
121 
122 
123  if ( fPhotosInterface ){
125  // for tauola++
126  if ( fPhotosInterface ){
127  for ( unsigned int iss=0; iss<fPhotosInterface->specialSettings().size(); iss++ ){
129  }
130  }
131  }
132 
133  fIsInitialized = true;
134 
135  return;
136 }
137 
139 {
141  //if ( fPhotosInterface ) fPhotosInterface->statistics();
142  // similar for EvtGen if needed
143  return;
144 }
145 
146 void ExternalDecayDriver::setRandomEngine(CLHEP::HepRandomEngine* v)
147 {
152 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
virtual void init()=0
virtual void configureOnlyFor(int)=0
TauolaInterfaceBase * fTauolaInterface
CLHEP::HepRandomEngine * decayRandomEngine
virtual void init(const edm::EventSetup &)
PhotosInterfaceBase * fPhotosInterface
static const std::string kTauola
HepMC::GenEvent * decay(HepMC::GenEvent *)
static const std::string kPhotos
double v[5][pyjets_maxn]
virtual void setRandomEngine(CLHEP::HepRandomEngine *decayRandomEngine)=0
void init(const edm::EventSetup &)
static const std::string kPythia6
virtual HepMC::GenEvent * apply(HepMC::GenEvent *evt)
static const std::string kEvtGen
std::vector< std::string > exSharedResources
virtual const std::vector< int > & operatesOnParticles()
std::vector< std::string > fSpecialSettings
static const std::string kFortranInstance
EvtGenInterfaceBase * fEvtGenInterface
ExternalDecayDriver(const edm::ParameterSet &)
virtual HepMC::GenEvent * decay(HepMC::GenEvent *evt)
virtual HepMC::GenEvent * decay(HepMC::GenEvent *evt)
virtual const std::vector< std::string > & specialSettings()
virtual const std::vector< int > & operatesOnParticles()
volatile std::atomic< bool > shutdown_flag false
static const std::string kPythia8
virtual void setRandomEngine(CLHEP::HepRandomEngine *v)=0
void setRandomEngine(CLHEP::HepRandomEngine *)
virtual void setRandomEngine(CLHEP::HepRandomEngine *v)=0
T get(const Candidate &c)
Definition: component.h:55
virtual void avoidTauLeptonicDecays()=0