Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CascadeFinalStateAlgorithm Class Reference

#include <G4CascadeFinalStateAlgorithm.hh>

+ Inheritance diagram for G4CascadeFinalStateAlgorithm:

Public Member Functions

 G4CascadeFinalStateAlgorithm ()
 
virtual ~G4CascadeFinalStateAlgorithm ()
 
virtual void SetVerboseLevel (G4int verbose)
 
void Configure (G4InuclElementaryParticle *bullet, G4InuclElementaryParticle *target, const std::vector< G4int > &particle_kinds)
 
- Public Member Functions inherited from G4VHadDecayAlgorithm
 G4VHadDecayAlgorithm (const G4String &algName, G4int verbose=0)
 
virtual ~G4VHadDecayAlgorithm ()
 
void Generate (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
G4int GetVerboseLevel () const
 
const G4StringGetName () const
 

Protected Member Functions

virtual void GenerateTwoBody (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
virtual void GenerateMultiBody (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
void SaveKinematics (G4InuclElementaryParticle *bullet, G4InuclElementaryParticle *target)
 
void ChooseGenerators (G4int is, G4int fs)
 
void FillMagnitudes (G4double initialMass, const std::vector< G4double > &masses)
 
G4bool satisfyTriangle (const std::vector< G4double > &pmod) const
 
void FillDirections (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
void FillDirThreeBody (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
void FillDirManyBody (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
G4double GenerateCosTheta (G4int ptype, G4double pmod) const
 
void FillUsingKopylov (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
G4double BetaKopylov (G4int K) const
 
- Protected Member Functions inherited from G4VHadDecayAlgorithm
virtual G4bool IsDecayAllowed (G4double initialMass, const std::vector< G4double > &masses) const
 
G4double TwoBodyMomentum (G4double M0, G4double M1, G4double M2) const
 
G4double UniformTheta () const
 
G4double UniformPhi () const
 
void PrintVector (const std::vector< G4double > &v, const G4String &name, std::ostream &os) const
 

Detailed Description

Definition at line 48 of file G4CascadeFinalStateAlgorithm.hh.

Constructor & Destructor Documentation

◆ G4CascadeFinalStateAlgorithm()

G4CascadeFinalStateAlgorithm::G4CascadeFinalStateAlgorithm ( )

Definition at line 74 of file G4CascadeFinalStateAlgorithm.cc.

75 : G4VHadDecayAlgorithm("G4CascadeFinalStateAlgorithm"),
76 momDist(0), angDist(0), multiplicity(0), bullet_ekin(0.) {;}
G4VHadDecayAlgorithm(const G4String &algName, G4int verbose=0)

◆ ~G4CascadeFinalStateAlgorithm()

G4CascadeFinalStateAlgorithm::~G4CascadeFinalStateAlgorithm ( )
virtual

Definition at line 78 of file G4CascadeFinalStateAlgorithm.cc.

78{;}

Member Function Documentation

◆ BetaKopylov()

G4double G4CascadeFinalStateAlgorithm::BetaKopylov ( G4int K) const
protected

Definition at line 481 of file G4CascadeFinalStateAlgorithm.cc.

481 {
482 G4Pow* g4pow = G4Pow::GetInstance();
483
484 G4int N = 3*K - 5;
485 G4double xN = G4double(N);
486 G4double Fmax = std::sqrt(g4pow->powN(xN/(xN+1.),N)/(xN+1.));
487
488 G4double F, chi;
489 do { /* Loop checking 08.06.2015 MHK */
490 chi = G4UniformRand();
491 F = std::sqrt(g4pow->powN(chi,N)*(1.-chi));
492 } while ( Fmax*G4UniformRand() > F);
493 return chi;
494}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double powN(G4double x, G4int n) const
Definition G4Pow.cc:162
#define N
Definition crc32.c:57

Referenced by FillUsingKopylov().

◆ ChooseGenerators()

void G4CascadeFinalStateAlgorithm::ChooseGenerators ( G4int is,
G4int fs )
protected

Definition at line 135 of file G4CascadeFinalStateAlgorithm.cc.

135 {
136 if (GetVerboseLevel()>1)
137 G4cout << " >>> " << GetName() << "::ChooseGenerators"
138 << " is " << is << " fs " << fs << G4endl;
139
140 // Get generators for momentum and angle
141 if (G4CascadeParameters::usePhaseSpace()) momDist = 0;
142 else momDist = G4MultiBodyMomentumDist::GetDist(is, multiplicity);
143
144 if (fs > 0 && multiplicity == 2) {
145 G4int kw = (fs==is) ? 1 : 2;
146 angDist = G4TwoBodyAngularDist::GetDist(is, fs, kw);
147 } else if (multiplicity == 3) {
148 angDist = G4TwoBodyAngularDist::GetDist(is);
149 } else {
150 angDist = 0;
151 }
152
153 if (GetVerboseLevel()>1) {
154 G4cout << " " << (momDist?momDist->GetName().c_str():"") << " "
155 << (angDist?angDist->GetName().c_str():"") << G4endl;
156 }
157}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4bool usePhaseSpace()
static const G4VMultiBodyMomDst * GetDist(G4int is, G4int mult)
static const G4VTwoBodyAngDst * GetDist(G4int is, G4int fs, G4int kw)
const G4String & GetName() const

Referenced by Configure().

◆ Configure()

void G4CascadeFinalStateAlgorithm::Configure ( G4InuclElementaryParticle * bullet,
G4InuclElementaryParticle * target,
const std::vector< G4int > & particle_kinds )

Definition at line 90 of file G4CascadeFinalStateAlgorithm.cc.

93 {
94 if (GetVerboseLevel()>1)
95 G4cout << " >>> " << GetName() << "::Configure" << G4endl;
96
97 // Identify initial and final state (if two-body) for algorithm selection
98 multiplicity = (G4int)particle_kinds.size();
99 G4int is = bullet->type() * target->type();
100 G4int fs = (multiplicity==2) ? particle_kinds[0]*particle_kinds[1] : 0;
101
102 ChooseGenerators(is, fs);
103
104 // Save kinematics for use with distributions
105 SaveKinematics(bullet, target);
106
107 // Save particle types for use with distributions
108 kinds = particle_kinds;
109}
void SaveKinematics(G4InuclElementaryParticle *bullet, G4InuclElementaryParticle *target)

Referenced by G4CascadeFinalStateGenerator::Configure().

◆ FillDirections()

void G4CascadeFinalStateAlgorithm::FillDirections ( G4double initialMass,
const std::vector< G4double > & masses,
std::vector< G4LorentzVector > & finalState )
protected

Definition at line 311 of file G4CascadeFinalStateAlgorithm.cc.

313 {
314 if (GetVerboseLevel()>1)
315 G4cout << " >>> " << GetName() << "::FillDirections" << G4endl;
316
317 finalState.clear(); // Initialization and sanity check
318 if ((G4int)modules.size() != multiplicity) return;
319
320 // Different order of processing for three vs. N body
321 if (multiplicity == 3)
322 FillDirThreeBody(initialMass, masses, finalState);
323 else
324 FillDirManyBody(initialMass, masses, finalState);
325}
void FillDirThreeBody(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
void FillDirManyBody(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)

Referenced by GenerateMultiBody().

◆ FillDirManyBody()

void G4CascadeFinalStateAlgorithm::FillDirManyBody ( G4double initialMass,
const std::vector< G4double > & masses,
std::vector< G4LorentzVector > & finalState )
protected

Definition at line 360 of file G4CascadeFinalStateAlgorithm.cc.

362 {
363 if (GetVerboseLevel()>1)
364 G4cout << " >>> " << GetName() << "::FillDirManyBody" << G4endl;
365
366 // Fill all but the last two particles randomly
367 G4double costh = 0.;
368
369 finalState.resize(multiplicity);
370
371 for (G4int i=0; i<multiplicity-2; i++) {
372 costh = GenerateCosTheta(kinds[i], modules[i]);
373 finalState[i] = generateWithFixedTheta(costh, modules[i], masses[i]);
374 finalState[i] = toSCM.rotate(finalState[i]); // Align target axis
375 }
376
377 // Total momentum so far, to compute recoil of last two particles
378 G4LorentzVector psum =
379 std::accumulate(finalState.begin(), finalState.end()-2, G4LorentzVector());
380 G4double pmod = psum.rho();
381
382 costh = -0.5 * (pmod*pmod +
383 modules[multiplicity-2]*modules[multiplicity-2] -
384 modules[multiplicity-1]*modules[multiplicity-1])
385 / pmod / modules[multiplicity-2];
386
387 if (GetVerboseLevel() > 2) G4cout << " ct last " << costh << G4endl;
388
389 if (std::fabs(costh) >= maxCosTheta) { // Bad kinematics; abort generation
390 finalState.clear();
391 return;
392 }
393
394 // Report success
395 if (GetVerboseLevel()>2) G4cout << " ok for mult " << multiplicity << G4endl;
396
397 // First particle is at fixed angle to recoil system
398 finalState[multiplicity-2] =
399 generateWithFixedTheta(costh, modules[multiplicity-2],
400 masses[multiplicity-2]);
401 finalState[multiplicity-2] = toSCM.rotate(psum, finalState[multiplicity-2]);
402
403 // Remaining particle is constrained to recoil from entire rest of system
404 finalState[multiplicity-1].set(0.,0.,0.,initialMass);
405 finalState[multiplicity-1] -= psum + finalState[multiplicity-2];
406}
CLHEP::HepLorentzVector G4LorentzVector
G4double GenerateCosTheta(G4int ptype, G4double pmod) const
G4LorentzVector generateWithFixedTheta(G4double ct, G4double p, G4double mass=0.)

Referenced by FillDirections().

◆ FillDirThreeBody()

void G4CascadeFinalStateAlgorithm::FillDirThreeBody ( G4double initialMass,
const std::vector< G4double > & masses,
std::vector< G4LorentzVector > & finalState )
protected

Definition at line 327 of file G4CascadeFinalStateAlgorithm.cc.

329 {
330 if (GetVerboseLevel()>1)
331 G4cout << " >>> " << GetName() << "::FillDirThreeBody" << G4endl;
332
333 finalState.resize(3);
334
335 G4double costh = GenerateCosTheta(kinds[2], modules[2]);
336 finalState[2] = generateWithFixedTheta(costh, modules[2], masses[2]);
337 finalState[2] = toSCM.rotate(finalState[2]); // Align target axis
338
339 // Generate direction of first particle
340 costh = -0.5 * (modules[2]*modules[2] + modules[0]*modules[0] -
341 modules[1]*modules[1]) / modules[2] / modules[0];
342
343 if (std::fabs(costh) >= maxCosTheta) { // Bad kinematics; abort generation
344 finalState.clear();
345 return;
346 }
347
348 // Report success
349 if (GetVerboseLevel()>2) G4cout << " ok for mult 3" << G4endl;
350
351 // First particle is at fixed angle to recoil system
352 finalState[0] = generateWithFixedTheta(costh, modules[0], masses[0]);
353 finalState[0] = toSCM.rotate(finalState[2], finalState[0]);
354
355 // Remaining particle is constrained to recoil from entire rest of system
356 finalState[1].set(0.,0.,0.,initialMass);
357 finalState[1] -= finalState[0] + finalState[2];
358}

Referenced by FillDirections().

◆ FillMagnitudes()

void G4CascadeFinalStateAlgorithm::FillMagnitudes ( G4double initialMass,
const std::vector< G4double > & masses )
protected

Definition at line 226 of file G4CascadeFinalStateAlgorithm.cc.

227 {
228 if (GetVerboseLevel()>1)
229 G4cout << " >>> " << GetName() << "::FillMagnitudes" << G4endl;
230
231 modules.clear(); // Initialization and sanity checks
232 if (!momDist) return;
233
234 modules.resize(multiplicity,0.); // Pre-allocate to avoid resizing
235
236 G4double mass_last = masses.back();
237 G4double pmod = 0.;
238
239 if (GetVerboseLevel() > 3){
240 G4cout << " knd_last " << kinds.back() << " mass_last "
241 << mass_last << G4endl;
242 }
243
244 G4int itry = -1;
245 while (++itry < itry_max) { /* Loop checking 08.06.2015 MHK */
246 if (GetVerboseLevel() > 3) {
247 G4cout << " itry in fillMagnitudes " << itry << G4endl;
248 }
249
250 G4double eleft = initialMass;
251
252 G4int i; // For access outside of loop
253 for (i=0; i < multiplicity-1; i++) {
254 pmod = momDist->GetMomentum(kinds[i], bullet_ekin);
255
256 if (pmod < small) break;
257 eleft -= std::sqrt(pmod*pmod + masses[i]*masses[i]);
258
259 if (GetVerboseLevel() > 3) {
260 G4cout << " kp " << kinds[i] << " pmod " << pmod
261 << " mass2 " << masses[i]*masses[i] << " eleft " << eleft
262 << "\n x1 " << eleft - mass_last << G4endl;
263 }
264
265 if (eleft <= mass_last) break;
266
267 modules[i] = pmod;
268 }
269
270 if (i < multiplicity-1) continue; // Failed to generate full kinematics
271
272 G4double plast = eleft * eleft - masses.back()*masses.back();
273 if (GetVerboseLevel() > 2) G4cout << " plast ** 2 " << plast << G4endl;
274
275 if (plast <= small) continue; // Not enough momentum left over
276
277 plast = std::sqrt(plast); // Final momentum is what's left over
278 modules.back() = plast;
279
280 if (multiplicity > 3 || satisfyTriangle(modules)) break; // Successful
281 } // while (itry < itry_max)
282
283 if (itry >= itry_max) { // Too many attempts
284 if (GetVerboseLevel() > 2)
285 G4cerr << " Unable to generate momenta for multiplicity "
286 << multiplicity << G4endl;
287
288 modules.clear(); // Something went wrong, throw away partial
289 }
290}
G4GLOB_DLL std::ostream G4cerr
G4bool satisfyTriangle(const std::vector< G4double > &pmod) const

Referenced by GenerateMultiBody().

◆ FillUsingKopylov()

void G4CascadeFinalStateAlgorithm::FillUsingKopylov ( G4double initialMass,
const std::vector< G4double > & masses,
std::vector< G4LorentzVector > & finalState )
protected

Definition at line 437 of file G4CascadeFinalStateAlgorithm.cc.

440 {
441 if (GetVerboseLevel()>2)
442 G4cout << " >>> " << GetName() << "::FillUsingKopylov" << G4endl;
443
444 finalState.clear();
445
446 std::size_t N = masses.size();
447 finalState.resize(N);
448
449 G4double mtot = std::accumulate(masses.begin(), masses.end(), 0.0);
450 G4double mu = mtot;
451 G4double Mass = initialMass;
452 G4double T = Mass-mtot;
453 G4double recoilMass = 0.0;
454 G4ThreeVector momV, boostV; // Buffers to reduce memory churn
455 G4LorentzVector recoil(0.0,0.0,0.0,Mass);
456
457 for (std::size_t k=N-1; k>0; --k) {
458 mu -= masses[k];
459 T *= (k>1) ? BetaKopylov((G4int)k) : 0.;
460
461 recoilMass = mu + T;
462
463 boostV = recoil.boostVector(); // Previous system's rest frame
464
465 // Create momentum with a random direction isotropically distributed
466 // FIXME: Should theta distribution should use Bertini fit function?
467 momV.setRThetaPhi(TwoBodyMomentum(Mass,masses[k],recoilMass),
469
470 finalState[k].setVectM(momV,masses[k]);
471 recoil.setVectM(-momV,recoilMass);
472
473 finalState[k].boost(boostV);
474 recoil.boost(boostV);
475 Mass = recoilMass;
476 }
477
478 finalState[0] = recoil;
479}
CLHEP::Hep3Vector G4ThreeVector
void setRThetaPhi(double r, double theta, double phi)
G4double TwoBodyMomentum(G4double M0, G4double M1, G4double M2) const

Referenced by GenerateMultiBody().

◆ GenerateCosTheta()

G4double G4CascadeFinalStateAlgorithm::GenerateCosTheta ( G4int ptype,
G4double pmod ) const
protected

Definition at line 411 of file G4CascadeFinalStateAlgorithm.cc.

412 {
413 if (GetVerboseLevel() > 2) {
414 G4cout << " >>> " << GetName() << "::GenerateCosTheta " << ptype
415 << " " << pmod << G4endl;
416 }
417
418 if (multiplicity == 3) { // Use distribution for three-body
419 return angDist->GetCosTheta(bullet_ekin, ptype);
420 }
421
422 // Throw multi-body distribution
423
424 // Sample costheta directly from exp(-a*pmod*(1 - costheta) )
425 // Previous method sampled from a*sintheta*exp(-a*sintheta),
426 // converted to costheta and (incorrectly) reflected around 180 degrees
427 //
428 G4double p0 = ptype < 3 ? 0.36 : 0.25; // 0.36 for nucleon, 0.25 for all others
429 G4double alf = 3.*pmod/p0;
430 return G4Log(G4UniformRand()*(G4Exp(2.*alf) - 1.) + 1.)/alf - 1.;
431}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
G4double G4Log(G4double x)
Definition G4Log.hh:227

Referenced by FillDirManyBody(), and FillDirThreeBody().

◆ GenerateMultiBody()

void G4CascadeFinalStateAlgorithm::GenerateMultiBody ( G4double initialMass,
const std::vector< G4double > & masses,
std::vector< G4LorentzVector > & finalState )
protectedvirtual

Implements G4VHadDecayAlgorithm.

Definition at line 203 of file G4CascadeFinalStateAlgorithm.cc.

205 {
206 if (GetVerboseLevel()>1)
207 G4cout << " >>> " << GetName() << "::GenerateMultiBody" << G4endl;
208
210 FillUsingKopylov(initialMass, masses, finalState);
211 return;
212 }
213
214 finalState.clear(); // Initialization and sanity checks
215 if (multiplicity < 3) return;
216 if (!momDist) return;
217
218 G4int itry = -1; /* Loop checking 08.06.2015 MHK */
219 while ((G4int)finalState.size() != multiplicity && ++itry < itry_max) {
220 FillMagnitudes(initialMass, masses);
221 FillDirections(initialMass, masses, finalState);
222 }
223}
void FillUsingKopylov(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
void FillDirections(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
void FillMagnitudes(G4double initialMass, const std::vector< G4double > &masses)

◆ GenerateTwoBody()

void G4CascadeFinalStateAlgorithm::GenerateTwoBody ( G4double initialMass,
const std::vector< G4double > & masses,
std::vector< G4LorentzVector > & finalState )
protectedvirtual

Implements G4VHadDecayAlgorithm.

Definition at line 162 of file G4CascadeFinalStateAlgorithm.cc.

164 {
165 if (GetVerboseLevel()>1)
166 G4cout << " >>> " << GetName() << "::GenerateTwoBody" << G4endl;
167
168 finalState.clear(); // Initialization and sanity checks
169
170 if (multiplicity != 2) return;
171
172 // Generate momentum vector in CMS for back-to-back particles
173 G4double pscm = TwoBodyMomentum(initialMass, masses[0], masses[1]);
174
175 G4double costh = angDist ? angDist->GetCosTheta(bullet_ekin, pscm)
176 : (2.*G4UniformRand() - 1.);
177
178 mom.setRThetaPhi(pscm, std::acos(costh), UniformPhi());
179
180 if (GetVerboseLevel()>3) { // Copied from old G4EPCollider
181 G4cout << " Particle kinds = " << kinds[0] << " , " << kinds[1]
182 << "\n pmod " << pscm
183 << "\n before rotation px " << mom.x() << " py " << mom.y()
184 << " pz " << mom.z() << G4endl;
185 }
186
187 finalState.resize(2); // Allows filling by index
188
189 finalState[0].setVectM(mom, masses[0]);
190 finalState[0] = toSCM.rotate(finalState[0]);
191
192 if (GetVerboseLevel()>3) { // Copied from old G4EPCollider
193 G4cout << " after rotation px " << finalState[0].x() << " py "
194 << finalState[0].y() << " pz " << finalState[0].z() << G4endl;
195 }
196
197 finalState[1].setVectM(-finalState[0].vect(), masses[1]);
198}

◆ satisfyTriangle()

G4bool G4CascadeFinalStateAlgorithm::satisfyTriangle ( const std::vector< G4double > & pmod) const
protected

Definition at line 294 of file G4CascadeFinalStateAlgorithm.cc.

295 {
296 if (GetVerboseLevel() > 3)
297 G4cout << " >>> " << GetName() << "::satisfyTriangle" << G4endl;
298
299 return ( (pmod.size() != 3) ||
300 !(pmod[0] < std::fabs(pmod[1] - pmod[2]) ||
301 pmod[0] > pmod[1] + pmod[2] ||
302 pmod[1] < std::fabs(pmod[0] - pmod[2]) ||
303 pmod[1] > pmod[0] + pmod[2] ||
304 pmod[2] < std::fabs(pmod[0] - pmod[1]) ||
305 pmod[2] > pmod[1] + pmod[0])
306 );
307}

Referenced by FillMagnitudes().

◆ SaveKinematics()

void G4CascadeFinalStateAlgorithm::SaveKinematics ( G4InuclElementaryParticle * bullet,
G4InuclElementaryParticle * target )
protected

Definition at line 113 of file G4CascadeFinalStateAlgorithm.cc.

115 {
116 if (GetVerboseLevel()>1)
117 G4cout << " >>> " << GetName() << "::SaveKinematics" << G4endl;
118
119 if (target->nucleon()) { // Which particle originated in nucleus?
120 toSCM.setBullet(bullet);
121 toSCM.setTarget(target);
122 } else {
123 toSCM.setBullet(target);
124 toSCM.setTarget(bullet);
125 }
126
127 toSCM.toTheCenterOfMass();
128
129 bullet_ekin = toSCM.getKinEnergyInTheTRS();
130}

Referenced by Configure().

◆ SetVerboseLevel()

void G4CascadeFinalStateAlgorithm::SetVerboseLevel ( G4int verbose)
virtual

Reimplemented from G4VHadDecayAlgorithm.

Definition at line 80 of file G4CascadeFinalStateAlgorithm.cc.

80 {
84 toSCM.setVerbose(verbose);
85}
static void setVerboseLevel(G4int vb=0)
static void setVerboseLevel(G4int vb=0)
virtual void SetVerboseLevel(G4int verbose)

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