ATLAS Offline Software
Loading...
Searching...
No Matches
EvtPhotosEngine.cxx
Go to the documentation of this file.
1
2/***********************************************************************
3* Copyright 1998-2024 CERN for the benefit of the EvtGen authors *
4* *
5* This file is part of EvtGen. *
6* *
7* EvtGen is free software: you can redistribute it and/or modify *
8* it under the terms of the GNU General Public License as published by *
9* the Free Software Foundation, either version 3 of the License, or *
10* (at your option) any later version. *
11* *
12* EvtGen is distributed in the hope that it will be useful, *
13* but WITHOUT ANY WARRANTY; without even the implied warranty of *
14* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15* GNU General Public License for more details. *
16* *
17* You should have received a copy of the GNU General Public License *
18* along with EvtGen. If not, see <https://www.gnu.org/licenses/>. *
19***********************************************************************/
20
21
22#include "EvtGen_i/EvtGenExternal/EvtPhotosEngine.hh"
23
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"
29
32
33#include <iostream>
34#include <sstream>
35#include <vector>
36
37using std::endl;
38
39EvtPhotosEngine::EvtPhotosEngine( const std::string& photonType, bool useEvtGenRandom ) :
40 m_photonType(photonType),
41 m_gammaId(EvtId( -1, -1 )),
42 m_gammaPDG(22), // default photon pdg integer
43 m_mPhoton(0.0)
44{
45
46 EvtGenReport( EVTGEN_INFO, "EvtGen" ) << "Setting up PHOTOS." << endl;
47
48 if ( useEvtGenRandom == true ) {
49 EvtGenReport( EVTGEN_INFO, "EvtGen" )
50 << "Using EvtGen random number engine also for Photos++" << endl;
51
52 Photospp::Photos::setRandomGenerator( EvtRandom::Flat );
53 }
54
55 Photospp::Photos::initialize();
56
57 // Increase the maximum possible value of the interference weight
58 Photospp::Photos::maxWtInterference(
59 64.0 ); // 2^n, where n = number of charges (+,-)
60 Photospp::Photos::setInterference( true );
61 Photospp::Photos::setExponentiation( true ); // Sets the infrared cutoff at 1e-7
62 // Reset the minimum photon energy, if required, in units of half of the decaying particle mass.
63 // This must be done after exponentiation! Keep the cut at 1e-7, i.e. 0.1 keV at the 1 GeV scale,
64 // which is appropriate for B decays
65 Photospp::Photos::setInfraredCutOff( 1.0e-7 );
66
67 m_initialised = false;
68}
69
70void EvtPhotosEngine::initialise()
71{
72 if ( m_initialised == false ) {
73 m_gammaId = EvtPDL::getId( m_photonType );
74
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" );
80 }
81
82 m_gammaPDG = EvtPDL::getStdHep( m_gammaId );
83 m_mPhoton = EvtPDL::getMeanMass( m_gammaId );
84
85 m_initialised = true;
86 }
87}
88
89bool EvtPhotosEngine::doDecay( EvtParticle* theMother )
90{
91 if ( m_initialised == false ) {
92 this->initialise();
93 }
94
95 if ( theMother == 0 ) {
96 return false;
97 }
98
99 // Create a dummy HepMC GenEvent containing a single vertex, with the mother
100 // assigned as the incoming particle and its daughters as outgoing particles.
101 // We then pass this event to Photos for processing.
102 // It will return a modified version of the event, updating the momentum of
103 // the original particles and will contain any new photon particles.
104 // We add these extra photons to the mother particle daughter list.
105
106 // Skip running Photos if the particle has no daughters, since we can't add FSR.
107 // Also skip Photos if the particle has too many daughters (>= 10) to avoid a problem
108 // with a hard coded upper limit in the PHOENE subroutine.
109 int nDaug( theMother->getNDaug() );
110 if ( nDaug == 0 || nDaug >= 10 ) {
111 return false;
112 }
113
114 // Create the dummy event.
115 auto theEvent = std::make_unique<GenEvent>( Units::GEV, Units::MM );
116
117 // Create the decay "vertex".
118 GenVertexPtr theVertex = newGenVertexPtr();
119 theEvent->add_vertex( theVertex );
120
121 // Add the mother particle as the incoming particle to the vertex.
122 GenParticlePtr hepMCMother = this->createGenParticle( theMother, true );
123 theVertex->add_particle_in( std::move(hepMCMother) );
124
125 // Find all daughter particles and assign them as outgoing particles to the vertex.
126 // Keep track of the number of photons already in the decay (e.g. we may have B -> K* gamma)
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,
131 false );
132 theVertex->add_particle_out( std::move(hepMCDaughter) );
133
134 if ( theDaughter ) {
135 int daugId = theDaughter->getPDGId();
136 if ( daugId == m_gammaPDG ) {
137 nGamma++;
138 }
139 }
140 }
141
142 // Now pass the event to Photos for processing
143 // Create a Photos event object
144#ifdef HEPMC3
145 Photospp::PhotosHepMC3Event photosEvent( theEvent.get() );
146#else
147 Photospp::PhotosHepMCEvent photosEvent( theEvent.get() );
148#endif
149
150 // Run the Photos algorithm
151 photosEvent.process();
152
153 // Find the number of (outgoing) photons in the event
154 int nPhotons = this->getNumberOfPhotons( theVertex );
155
156 // See if Photos has created additional photons. If not, do nothing extra
157 int nDiffPhotons = nPhotons - nGamma;
158 int iLoop( 0 );
159
160 if ( nDiffPhotons > 0 ) {
161 // We have extra particles from Photos; these would have been appended
162 // to the outgoing particle list
163
164 // Get the iterator of outgoing particles for this vertex
165#ifdef HEPMC3
166 for ( auto outParticle : theVertex->particles_out() ) {
167#else
168 HepMC::GenVertex::particles_out_const_iterator outIter;
169 for ( outIter = theVertex->particles_out_const_begin();
170 outIter != theVertex->particles_out_const_end(); ++outIter ) {
171 // Get the next HepMC GenParticle
172 HepMC::GenParticle* outParticle = *outIter;
173#endif
174
175 // Get the three-momentum Photos result for this particle, and the PDG id
176 double px( 0.0 ), py( 0.0 ), pz( 0.0 );
177 int pdgId( 0 );
178
179 if ( outParticle != 0 ) {
180 FourVector HepMCP4 = outParticle->momentum();
181 px = HepMCP4.px();
182 py = HepMCP4.py();
183 pz = HepMCP4.pz();
184 pdgId = outParticle->pdg_id();
185 }
186
187 // Create an empty 4-momentum vector for the new/modified daughters
188 EvtVector4R newP4;
189
190 if ( iLoop < nDaug ) {
191 // Original daughters
192 EvtParticle* daugParticle = theMother->getDaug( iLoop );
193 if ( daugParticle != 0 ) {
194 // Keep the original particle mass, but set the three-momentum
195 // according to what Photos has modified. However, this will
196 // violate energy conservation (from what Photos has provided).
197 double mass = daugParticle->mass();
198 double energy = sqrt( mass * mass + px * px + py * py +
199 pz * pz );
200 newP4.set( energy, px, py, pz );
201 // Set the new four-momentum (FSR applied)
202 daugParticle->setP4WithFSR( newP4 );
203 }
204
205 } else if ( pdgId == m_gammaPDG ) {
206 // Extra photon particle. Setup the four-momentum object
207 double energy = sqrt( m_mPhoton * m_mPhoton + px * px + py * py +
208 pz * pz );
209 newP4.set( energy, px, py, pz );
210
211 // Create a new photon particle and add it to the list of daughters
212 EvtPhotonParticle* gamma = new EvtPhotonParticle();
213 gamma->init( m_gammaId, newP4 );
214 // Set the pre-FSR photon momentum to zero
215 gamma->setFSRP4toZero();
216 // Let the mother know about this new photon
217 gamma->addDaug( theMother );
218 // Set its particle attribute to specify it is a FSR photon
219 gamma->setAttribute( "FSR", 1 ); // it is a FSR photon
220 gamma->setAttribute( "ISR", 0 ); // it is not an ISR photon
221 }
222
223 // Increment the loop counter for detecting additional photon particles
224 iLoop++;
225 }
226 }
227
228 // Cleanup
229 theEvent->clear();
230
231 return true;
232}
233
234GenParticlePtr EvtPhotosEngine::createGenParticle( EvtParticle* theParticle,
235 bool incoming )
236{
237 // Method to create an HepMC::GenParticle version of the given EvtParticle.
238 if ( theParticle == 0 ) {
239 return 0;
240 }
241
242 // Get the 4-momentum (E, px, py, pz) for the EvtParticle
243 EvtVector4R p4( 0.0, 0.0, 0.0, 0.0 );
244
245 if ( incoming == true ) {
246 p4 = theParticle->getP4Restframe();
247 } else {
248 p4 = theParticle->getP4();
249 }
250
251 // Convert this to the HepMC 4-momentum
252 double E = p4.get( 0 );
253 double px = p4.get( 1 );
254 double py = p4.get( 2 );
255 double pz = p4.get( 3 );
256
257 FourVector hepMC_p4( px, py, pz, E );
258
259 int PDGInt = EvtPDL::getStdHep( theParticle->getId() );
260
261 // Set the status flag for the particle. This is required, otherwise Photos++
262 // will crash from out-of-bounds array index problems.
263 int status = Photospp::PhotosParticle::HISTORY;
264 if ( incoming == false ) {
265 status = Photospp::PhotosParticle::STABLE;
266 }
267
268 GenParticlePtr genParticle = newGenParticlePtr( hepMC_p4, PDGInt, status );
269
270 return genParticle;
271}
272
273int EvtPhotosEngine::getNumberOfPhotons( const GenVertexPtr theVertex ) const
274{
275 // Find the number of photons from the outgoing particle list
276
277 if ( !theVertex ) {
278 return 0;
279 }
280
281 int nPhotons( 0 );
282
283 // Get the iterator of outgoing particles for this vertex
284#ifdef HEPMC3
285 for ( auto outParticle : theVertex->particles_out() ) {
286#else
287 HepMC::GenVertex::particles_out_const_iterator outIter;
288 for ( outIter = theVertex->particles_out_const_begin();
289 outIter != theVertex->particles_out_const_end(); ++outIter ) {
290 // Get the next HepMC GenParticle
291 HepMC::GenParticle* outParticle = *outIter;
292#endif
293
294 // Get the PDG id
295 int pdgId( 0 );
296 if ( outParticle != 0 ) {
297 pdgId = outParticle->pdg_id();
298 }
299
300 // Keep track of how many photons there are
301 if ( pdgId == m_gammaPDG ) {
302 nPhotons++;
303 }
304 }
305
306 return nPhotons;
307}
308
HepMC::GenVertex * GenVertexPtr
Definition GenVertex.h:59
GenVertexPtr newGenVertexPtr(const HepMC::FourVector &pos=HepMC::FourVector(0.0, 0.0, 0.0, 0.0), const int i=0)
Definition GenVertex.h:64
GenParticlePtr newGenParticlePtr(const HepMC::FourVector &mom=HepMC::FourVector(0.0, 0.0, 0.0, 0.0), int pid=0, int status=0)
Definition GenParticle.h:39
GenParticle * GenParticlePtr
Definition GenParticle.h:37
status
Definition merge.py:16