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 
25 
27  : fIsInitialized(false),
28  fTauolaInterface(0),
29  fEvtGenInterface(0),
30  fPhotosInterface(0)
31 {
32  std::vector<std::string> extGenNames =
33  pset.getParameter< std::vector<std::string> >("parameterSets");
34 
35  for (unsigned int ip=0; ip<extGenNames.size(); ++ip ) {
36  std::string curSet = extGenNames[ip];
37  if ( curSet == "EvtGen") {
42  }
43  else if( curSet == "EvtGen1" || curSet == "EvtGen130" ) {
50  }
51  else if ( curSet == "Tauola" || curSet == "Tauolapp" || curSet == "Tauolapp114" ) {
52  fTauolaInterface = (TauolaInterfaceBase*)(TauolaFactory::get()->create("Tauolapp114", pset.getUntrackedParameter< ParameterSet >(curSet)));
53  fPhotosInterface = (PhotosInterfaceBase*)(PhotosFactory::get()->create("Photos2155", pset.getUntrackedParameter< ParameterSet >(curSet)));
58  }
59  else if ( curSet == "Photos" || curSet == "Photos2155" ) {
60  if ( !fPhotosInterface ) {
63  }
64  }
65  else if (curSet == "Photospp" || curSet == "Photospp356" ) {
66  if ( !fPhotosInterface ) {
67  fPhotosInterface = (PhotosInterfaceBase*)(PhotosFactory::get()->create("Photospp356", pset.getUntrackedParameter< ParameterSet>(curSet)));
69  }
70  }
71  }
72 }
73 
74 
76 {
77  if ( fEvtGenInterface ) delete fEvtGenInterface;
78  if ( fTauolaInterface ) delete fTauolaInterface;
79  if ( fPhotosInterface ) delete fPhotosInterface;
80 }
81 
82 
83 HepMC::GenEvent* ExternalDecayDriver::decay(HepMC::GenEvent* evt, lhef::LHEEvent *lheEvent){
85  return decay(evt);
86 }
87 
88 
89 HepMC::GenEvent* ExternalDecayDriver::decay( HepMC::GenEvent* evt )
90 {
91  if ( !fIsInitialized ) return evt;
92 
93  if ( fEvtGenInterface ) {
94  evt = fEvtGenInterface->decay( evt );
95  if ( !evt ) return 0;
96  }
97 
98  if ( fTauolaInterface ) {
99  evt = fTauolaInterface->decay( evt );
100  if ( !evt ) return 0;
101  }
102 
103  if ( fPhotosInterface ) {
104  evt = fPhotosInterface->apply( evt );
105  if ( !evt ) return 0;
106  }
107 
108  return evt;
109 }
110 
111 
113 {
114  if ( fIsInitialized ) return;
115 
116  if ( fTauolaInterface ) {
117  fTauolaInterface->init( es );
118  for ( std::vector<int>::const_iterator i=fTauolaInterface->operatesOnParticles().begin();
119  i!=fTauolaInterface->operatesOnParticles().end(); i++ )
120  fPDGs.push_back( *i );
121  }
122 
123  if ( fEvtGenInterface ) {
125  for ( std::vector<int>::const_iterator i=fEvtGenInterface->operatesOnParticles().begin();
126  i!=fEvtGenInterface->operatesOnParticles().end(); i++ )
127  fPDGs.push_back( *i );
128  for ( unsigned int iss=0; iss<fEvtGenInterface->specialSettings().size(); iss++ ) {
130  }
131 
132  }
133 
134  if ( fPhotosInterface ) {
136  // for tauola++
137  if ( fPhotosInterface ) {
138  for ( unsigned int iss=0; iss<fPhotosInterface->specialSettings().size(); iss++ ){
140  }
141  }
142  }
143 
144  fIsInitialized = true;
145 
146  return;
147 }
148 
149 
151 {
154  // similar for EvtGen if needed
155  return;
156 }
157 
158 
159 void ExternalDecayDriver::setRandomEngine(CLHEP::HepRandomEngine* v)
160 {
165 }
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
static const std::string kPhotos
double v[5][pyjets_maxn]
virtual void setRandomEngine(CLHEP::HepRandomEngine *decayRandomEngine)=0
void init(const edm::EventSetup &)
def gen
run2 Cosmic #### Run 256259 @ 0T 2015C### Run 272133 @ 3.8T 2016B###
static const std::string kPythia6
virtual HepMC::GenEvent * apply(HepMC::GenEvent *evt)
static const std::string kEvtGen
virtual const std::vector< std::string > & specialSettings()
std::vector< std::string > exSharedResources
virtual const std::vector< int > & operatesOnParticles()
std::vector< std::string > fSpecialSettings
static const std::string kFortranInstance
EvtGenInterfaceBase * fEvtGenInterface
HepMC::GenEvent * decay(HepMC::GenEvent *evt)
ExternalDecayDriver(const edm::ParameterSet &)
virtual void SetLHE(lhef::LHEEvent *l)
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