34#define INCLXX_IN_GEANT4_MODE 1
47 const G4double PhaseSpaceRauboldLynch::wMaxMasslessX[PhaseSpaceRauboldLynch::wMaxNE] = {
80 const G4double PhaseSpaceRauboldLynch::wMaxMasslessY[PhaseSpaceRauboldLynch::wMaxNE] = {
113 const G4double PhaseSpaceRauboldLynch::wMaxCorrectionX[PhaseSpaceRauboldLynch::wMaxNE] = {
145 const G4double PhaseSpaceRauboldLynch::wMaxCorrectionY[PhaseSpaceRauboldLynch::wMaxNE] = {
178 const G4double PhaseSpaceRauboldLynch::wMaxInterpolationMargin = std::log(1.5);
184 maxGeneratedWeight(0.)
186 std::vector<G4double> wMaxMasslessXV(wMaxMasslessX, wMaxMasslessX + wMaxNE);
187 std::vector<G4double> wMaxMasslessYV(wMaxMasslessY, wMaxMasslessY + wMaxNE);
189 std::vector<G4double> wMaxCorrectionXV(wMaxCorrectionX, wMaxCorrectionX + wMaxNE);
190 std::vector<G4double> wMaxCorrectionYV(wMaxCorrectionY, wMaxCorrectionY + wMaxNE);
195 for(
size_t i=1; i<wMaxNP; ++i) {
202 delete wMaxCorrection;
206 maxGeneratedWeight = 0.;
211 initialize(particles);
214 const G4double weightMax = computeMaximumWeightParam();
216 const G4int maxIter = 500;
220 weight = computeWeight();
221 maxGeneratedWeight = std::max(weight, maxGeneratedWeight);
224 }
while(++iter<maxIter && r*weightMax>weight);
226#ifndef INCLXX_IN_GEANT4_MODE
228 INCL_WARN(
"Number of tries exceeded in PhaseSpaceRauboldLynch::generate()\n"
229 <<
"nParticles=" << nParticles <<
", sqrtS=" << sqrtS <<
'\n');
233 generateEvent(particles);
236 void PhaseSpaceRauboldLynch::initialize(
ParticleList &particles) {
237 nParticles = particles.size();
241 masses.resize(nParticles);
242 sumMasses.resize(nParticles);
243 std::transform(particles.begin(), particles.end(), masses.begin(), std::mem_fn(&
Particle::getMass));
244 std::partial_sum(masses.begin(), masses.end(), sumMasses.begin());
247 availableEnergy = sqrtS-sumMasses[nParticles-1];
249 if(availableEnergy<0.)
250 availableEnergy = 0.;
252 rnd.resize(nParticles);
253 invariantMasses.resize(nParticles);
254 momentaCM.resize(nParticles-1);
257 G4double PhaseSpaceRauboldLynch::computeMaximumWeightNaive() {
262 for(
size_t i=1; i<nParticles; i++) {
263 eMMin += masses[i-1];
270 G4double PhaseSpaceRauboldLynch::computeMaximumWeightParam() {
271#ifndef INCLXX_IN_GEANT4_MODE
272 if(nParticles>=wMaxNP) {
273 INCL_WARN(
"The requested number of particles (" << nParticles <<
") requires extrapolation the tables in PhaseSpaceRauboldLynch. YMMV." <<
'\n');
276 if(availableEnergy>wMaxMasslessX[wMaxNE-1] || availableEnergy<wMaxMasslessX[0] ) {
277 INCL_WARN(
"The requested available energy (" << availableEnergy <<
" MeV) requires extrapolation the tables in PhaseSpaceRauboldLynch. YMMV." <<
'\n');
282 const G4double logMassless = ((*wMaxMassless)(availableEnergy)+prelog[nParticles])*(nParticles-1);
283 const G4double reducedSqrtS = availableEnergy/sumMasses[nParticles-1];
284 const G4double correction = (*wMaxCorrection)(reducedSqrtS);
285 const G4double wMax = std::exp(correction*
G4double(nParticles-1) + logMassless + wMaxInterpolationMargin);
289 return computeMaximumWeightNaive();
292 G4double PhaseSpaceRauboldLynch::computeWeight() {
295 for(
size_t i=1; i<nParticles-1; ++i)
297 rnd[nParticles-1] = 1.;
298 std::sort(rnd.begin()+1, rnd.begin()+nParticles-1);
301 for(
size_t i=0; i<nParticles; ++i)
302 invariantMasses[i] = rnd[i]*availableEnergy + sumMasses[i];
306 momentaCM[0] = weight;
307 for(
size_t i=1; i<nParticles-1; ++i) {
310 if(invariantMasses[i+1]-invariantMasses[i]-masses[i+1] < 0.) momentumCM = 0.;
312 momentaCM[i] = momentumCM;
313 weight *= momentumCM;
319 void PhaseSpaceRauboldLynch::generateEvent(ParticleList &particles) {
320 Particle *p = particles[0];
322 p->setMomentum(momentum);
323 p->adjustEnergyFromMomentum();
327 for(
size_t i=1; i<nParticles; ++i) {
329 p->setMomentum(-momentum);
330 p->adjustEnergyFromMomentum();
337 const G4double iM = invariantMasses[i];
338 const G4double recoilE = std::sqrt(momentum.mag2() + iM*iM);
339 boostV = -momentum/recoilE;
340 for(
size_t j=0; j<=i; ++j)
341 particles[j]->boost(boostV);
346 return maxGeneratedWeight;
Class for interpolating the of a 1-dimensional function.
G4double getMass() const
Get the cached particle mass.
G4double getMaxGeneratedWeight() const
Return the largest generated weight.
virtual ~PhaseSpaceRauboldLynch()
void generate(const G4double sqrtS, ParticleList &particles)
Generate momenta according to a uniform, Lorentz-invariant phase-space model.
G4double momentumInCM(Particle const *const p1, Particle const *const p2)
gives the momentum in the CM frame of two particles.
ThreeVector normVector(G4double norm=1.)