22 #include "EvtGen_i/EvtGenExternal/EvtPhotosEngine.hh"
24 #include "EvtGenBase/EvtPDL.hh"
25 #include "EvtGenBase/EvtPhotonParticle.hh"
26 #include "EvtGenBase/EvtRandom.hh"
27 #include "EvtGenBase/EvtReport.hh"
28 #include "EvtGenBase/EvtVector4R.hh"
39 EvtPhotosEngine::EvtPhotosEngine( std::string photonType,
bool useEvtGenRandom )
41 m_photonType = photonType;
42 m_gammaId = EvtId( -1, -1 );
46 EvtGenReport( EVTGEN_INFO,
"EvtGen" ) <<
"Setting up PHOTOS." << endl;
48 if ( useEvtGenRandom ==
true ) {
49 EvtGenReport( EVTGEN_INFO,
"EvtGen" )
50 <<
"Using EvtGen random number engine also for Photos++" << endl;
52 Photospp::Photos::setRandomGenerator( EvtRandom::Flat );
58 Photospp::Photos::maxWtInterference(
60 Photospp::Photos::setInterference(
true );
61 Photospp::Photos::setExponentiation(
true );
65 Photospp::Photos::setInfraredCutOff( 1.0
e-7 );
67 m_initialised =
false;
70 void EvtPhotosEngine::initialise()
72 if ( m_initialised ==
false ) {
73 m_gammaId = EvtPDL::getId( m_photonType );
75 if ( m_gammaId == EvtId( -1, -1 ) ) {
76 EvtGenReport( EVTGEN_INFO,
"EvtGen" )
77 <<
"Error in EvtPhotosEngine. Do not recognise the photon type "
78 << m_photonType <<
". Setting this to \"gamma\". " << endl;
79 m_gammaId = EvtPDL::getId(
"gamma" );
82 m_gammaPDG = EvtPDL::getStdHep( m_gammaId );
83 m_mPhoton = EvtPDL::getMeanMass( m_gammaId );
89 bool EvtPhotosEngine::doDecay( EvtParticle* theMother )
91 if ( m_initialised ==
false ) {
95 if ( theMother == 0 ) {
109 int nDaug( theMother->getNDaug() );
110 if ( nDaug == 0 || nDaug >= 10 ) {
119 theEvent->add_vertex( theVertex );
122 GenParticlePtr hepMCMother = this->createGenParticle( theMother,
true );
123 theVertex->add_particle_in( hepMCMother );
127 int iDaug( 0 ), nGamma( 0 );
128 for ( iDaug = 0; iDaug < nDaug; iDaug++ ) {
129 EvtParticle* theDaughter = theMother->getDaug( iDaug );
130 GenParticlePtr hepMCDaughter = this->createGenParticle( theDaughter,
132 theVertex->add_particle_out( hepMCDaughter );
135 int daugId = theDaughter->getPDGId();
136 if ( daugId == m_gammaPDG ) {
145 Photospp::PhotosHepMC3Event photosEvent( theEvent.get() );
147 Photospp::PhotosHepMCEvent photosEvent( theEvent.get() );
151 photosEvent.process();
154 int nPhotons = this->getNumberOfPhotons( theVertex );
157 int nDiffPhotons = nPhotons - nGamma;
160 if ( nDiffPhotons > 0 ) {
166 for (
auto outParticle : theVertex->particles_out() ) {
168 HepMC::GenVertex::particles_out_const_iterator outIter;
169 for ( outIter = theVertex->particles_out_const_begin();
170 outIter != theVertex->particles_out_const_end(); ++outIter ) {
176 double px( 0.0 ),
py( 0.0 ),
pz( 0.0 );
179 if ( outParticle != 0 ) {
180 FourVector HepMCP4 = outParticle->momentum();
184 pdgId = outParticle->pdg_id();
190 if ( iLoop < nDaug ) {
192 EvtParticle* daugParticle = theMother->getDaug( iLoop );
193 if ( daugParticle != 0 ) {
197 double mass = daugParticle->mass();
202 daugParticle->setP4WithFSR( newP4 );
205 }
else if ( pdgId == m_gammaPDG ) {
212 EvtPhotonParticle*
gamma =
new EvtPhotonParticle();
213 gamma->init( m_gammaId, newP4 );
215 gamma->setFSRP4toZero();
217 gamma->addDaug( theMother );
219 gamma->setAttribute(
"FSR", 1 );
220 gamma->setAttribute(
"ISR", 0 );
234 GenParticlePtr EvtPhotosEngine::createGenParticle( EvtParticle* theParticle,
238 if ( theParticle == 0 ) {
243 EvtVector4R p4( 0.0, 0.0, 0.0, 0.0 );
245 if ( incoming ==
true ) {
246 p4 = theParticle->getP4Restframe();
248 p4 = theParticle->getP4();
252 double E = p4.get( 0 );
253 double px = p4.get( 1 );
254 double py = p4.get( 2 );
255 double pz = p4.get( 3 );
257 FourVector hepMC_p4(
px,
py,
pz,
E );
259 int PDGInt = EvtPDL::getStdHep( theParticle->getId() );
263 int status = Photospp::PhotosParticle::HISTORY;
264 if ( incoming ==
false ) {
265 status = Photospp::PhotosParticle::STABLE;
273 int EvtPhotosEngine::getNumberOfPhotons(
const GenVertexPtr theVertex )
const
285 for (
auto outParticle : theVertex->particles_out() ) {
287 HepMC::GenVertex::particles_out_const_iterator outIter;
288 for ( outIter = theVertex->particles_out_const_begin();
289 outIter != theVertex->particles_out_const_end(); ++outIter ) {
296 if ( outParticle != 0 ) {
297 pdgId = outParticle->pdg_id();
301 if ( pdgId == m_gammaPDG ) {