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 
13 using namespace gen;
14 using namespace edm;
15 
16 CLHEP::HepRandomEngine* decayRandomEngine = nullptr;
17 
19  : fIsInitialized(false),
20  fTauolaInterface(0),
21  fEvtGenInterface(0),
22  fPhotosInterface(0)
23 {
24  std::vector<std::string> extGenNames =
25  pset.getParameter< std::vector<std::string> >("parameterSets");
26 
27  for (unsigned int ip=0; ip<extGenNames.size(); ++ip ){
28  std::string curSet = extGenNames[ip];
29  if ( curSet == "EvtGen" || curSet == "EvtGenLHC91"){
30  fEvtGenInterface = (EvtGenInterfaceBase*)(EvtGenFactory::get()->create("EvtGenLHC91", pset.getUntrackedParameter< ParameterSet >(curSet)));
34  }
35  else if ( curSet == "Tauola" || curSet == "Tauolapp114" ){
36  fTauolaInterface = (TauolaInterfaceBase*)(TauolaFactory::get()->create("Tauolapp114", pset.getUntrackedParameter< ParameterSet >(curSet)));
37  fPhotosInterface = (PhotosInterfaceBase*)(PhotosFactory::get()->create("Photos2155", pset.getUntrackedParameter< ParameterSet >(curSet)));
42  }
43  else if ( curSet == "Photos" || curSet == "Photos2155" ){
44  if ( !fPhotosInterface ){
47  }
48  }
49  }
50 }
51 
53 {
54  if ( fEvtGenInterface ) delete fEvtGenInterface;
55  if ( fTauolaInterface ) delete fTauolaInterface;
56  if ( fPhotosInterface ) delete fPhotosInterface;
57 }
58 
59 HepMC::GenEvent* ExternalDecayDriver::decay( HepMC::GenEvent* evt )
60 {
61 
62  if ( !fIsInitialized ) return evt;
63 
64  if ( fEvtGenInterface )
65  {
66  evt = fEvtGenInterface->decay( evt );
67  if ( !evt ) return 0;
68  }
69 
70  if ( fTauolaInterface )
71  {
72  evt = fTauolaInterface->decay( evt );
73  if ( !evt ) return 0;
74  }
75 
76  if ( fPhotosInterface )
77  {
78  evt = fPhotosInterface->apply( evt );
79  if ( !evt ) return 0;
80  }
81 
82  return evt;
83 }
84 
86 {
87 
88  if ( fIsInitialized ) return;
89 
90  if ( fTauolaInterface )
91  {
92  fTauolaInterface->init( es );
93  for ( std::vector<int>::const_iterator i=fTauolaInterface->operatesOnParticles().begin();
94  i!=fTauolaInterface->operatesOnParticles().end(); i++ )
95  fPDGs.push_back( *i );
96  }
97 
98  if ( fEvtGenInterface )
99  {
101  for ( std::vector<int>::const_iterator i=fEvtGenInterface->operatesOnParticles().begin();
102  i!=fEvtGenInterface->operatesOnParticles().end(); i++ )
103  fPDGs.push_back( *i );
104  }
105 
106 
107  if ( fPhotosInterface )
108  {
110 // for tauola++
111  if ( fPhotosInterface )
112  {
113  for ( unsigned int iss=0; iss<fPhotosInterface->specialSettings().size(); iss++ )
114  {
116  }
117  }
118  }
119 
120 // this is specific to old tauola27 only, because it calls up photos automatically
121 //
122 //
123 // if ( fTauolaInterface )
124 // {
125 // // override !
126 // fSpecialSettings.clear();
127 // fSpecialSettings.push_back( "QED-brem-off:15" );
128 // }
129 
130  fIsInitialized = true;
131 
132  return;
133 }
134 
136 {
138  // similar for EvtGen and/or Photos, if needs be
139  return;
140 }
141 
142 void ExternalDecayDriver::setRandomEngine(CLHEP::HepRandomEngine* v)
143 {
148 }
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
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