Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4INCL::PhaseSpaceRauboldLynch Class Reference

Generate momenta using the RauboldLynch method. More...

#include <G4INCLPhaseSpaceRauboldLynch.hh>

+ Inheritance diagram for G4INCL::PhaseSpaceRauboldLynch:

Public Member Functions

 PhaseSpaceRauboldLynch ()
 
virtual ~PhaseSpaceRauboldLynch ()
 
 PhaseSpaceRauboldLynch (PhaseSpaceRauboldLynch const &other)
 Dummy copy constructor to silence Coverity warning.
 
PhaseSpaceRauboldLynchoperator= (PhaseSpaceRauboldLynch const &rhs)
 Dummy assignment operator to silence Coverity warning.
 
void generate (const G4double sqrtS, ParticleList &particles)
 Generate momenta according to a uniform, Lorentz-invariant phase-space model.
 
G4double getMaxGeneratedWeight () const
 Return the largest generated weight.
 
- Public Member Functions inherited from G4INCL::IPhaseSpaceGenerator
 IPhaseSpaceGenerator ()
 
virtual ~IPhaseSpaceGenerator ()
 

Detailed Description

Generate momenta using the RauboldLynch method.

Definition at line 49 of file G4INCLPhaseSpaceRauboldLynch.hh.

Constructor & Destructor Documentation

◆ PhaseSpaceRauboldLynch() [1/2]

G4INCL::PhaseSpaceRauboldLynch::PhaseSpaceRauboldLynch ( )

Definition at line 180 of file G4INCLPhaseSpaceRauboldLynch.cc.

180 :
181 nParticles(0),
182 sqrtS(0.),
183 availableEnergy(0.),
184 maxGeneratedWeight(0.)
185 {
186 std::vector<G4double> wMaxMasslessXV(wMaxMasslessX, wMaxMasslessX + wMaxNE);
187 std::vector<G4double> wMaxMasslessYV(wMaxMasslessY, wMaxMasslessY + wMaxNE);
188 wMaxMassless = new InterpolationTable(wMaxMasslessXV, wMaxMasslessYV);
189 std::vector<G4double> wMaxCorrectionXV(wMaxCorrectionX, wMaxCorrectionX + wMaxNE);
190 std::vector<G4double> wMaxCorrectionYV(wMaxCorrectionY, wMaxCorrectionY + wMaxNE);
191 wMaxCorrection = new InterpolationTable(wMaxCorrectionXV, wMaxCorrectionYV);
192
193 // initialize the precalculated coefficients
194 prelog[0] = 0.;
195 for(size_t i=1; i<wMaxNP; ++i) {
196 prelog[i] = -std::log(G4double(i));
197 }
198 }
double G4double
Definition G4Types.hh:83

◆ ~PhaseSpaceRauboldLynch()

G4INCL::PhaseSpaceRauboldLynch::~PhaseSpaceRauboldLynch ( )
virtual

Definition at line 200 of file G4INCLPhaseSpaceRauboldLynch.cc.

200 {
201 delete wMaxMassless;
202 delete wMaxCorrection;
203 }

◆ PhaseSpaceRauboldLynch() [2/2]

G4INCL::PhaseSpaceRauboldLynch::PhaseSpaceRauboldLynch ( PhaseSpaceRauboldLynch const & other)

Dummy copy constructor to silence Coverity warning.

Member Function Documentation

◆ generate()

void G4INCL::PhaseSpaceRauboldLynch::generate ( const G4double sqrtS,
ParticleList & particles )
virtual

Generate momenta according to a uniform, Lorentz-invariant phase-space model.

This function will assign momenta to the particles in the list that is passed as an argument. The event is generated in the CM frame.

Parameters
sqrtStotal centre-of-mass energy of the system
particleslist of particles

Implements G4INCL::IPhaseSpaceGenerator.

Definition at line 205 of file G4INCLPhaseSpaceRauboldLynch.cc.

205 {
206 maxGeneratedWeight = 0.;
207
208 sqrtS = sqs;
209
210 // initialize the structures containing the particle masses
211 initialize(particles);
212
213 // initialize the maximum weight
214 const G4double weightMax = computeMaximumWeightParam();
215
216 const G4int maxIter = 500;
217 G4int iter = 0;
218 G4double weight, r;
219 do {
220 weight = computeWeight();
221 maxGeneratedWeight = std::max(weight, maxGeneratedWeight);
222// assert(weight<=weightMax);
223 r = Random::shoot();
224 } while(++iter<maxIter && r*weightMax>weight); /* Loop checking, 10.07.2015, D.Mancusi */
225
226#ifndef INCLXX_IN_GEANT4_MODE
227 if(iter>=maxIter) {
228 INCL_WARN("Number of tries exceeded in PhaseSpaceRauboldLynch::generate()\n"
229 << "nParticles=" << nParticles << ", sqrtS=" << sqrtS << '\n');
230 }
231#endif
232
233 generateEvent(particles);
234 }
#define INCL_WARN(x)
int G4int
Definition G4Types.hh:85
G4double shoot()

◆ getMaxGeneratedWeight()

G4double G4INCL::PhaseSpaceRauboldLynch::getMaxGeneratedWeight ( ) const

Return the largest generated weight.

Definition at line 345 of file G4INCLPhaseSpaceRauboldLynch.cc.

345 {
346 return maxGeneratedWeight;
347 }

◆ operator=()

PhaseSpaceRauboldLynch & G4INCL::PhaseSpaceRauboldLynch::operator= ( PhaseSpaceRauboldLynch const & rhs)

Dummy assignment operator to silence Coverity warning.


The documentation for this class was generated from the following files: