Geant4-11
Public Member Functions | Protected Member Functions | Private Attributes
G4HadPhaseSpaceGenbod Class Reference

#include <G4HadPhaseSpaceGenbod.hh>

Inheritance diagram for G4HadPhaseSpaceGenbod:
G4VHadPhaseSpaceAlgorithm G4VHadDecayAlgorithm

Public Member Functions

 G4HadPhaseSpaceGenbod (G4int verbose=0)
 
void Generate (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
const G4StringGetName () const
 
G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int verbose)
 
virtual ~G4HadPhaseSpaceGenbod ()
 

Protected Member Functions

G4bool AcceptEvent () const
 
void AccumulateFinalState (size_t i, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
G4double ComputeWeight () const
 
void ComputeWeightScale (const std::vector< G4double > &masses)
 
void FillEnergySteps (G4double initialMass, const std::vector< G4double > &masses)
 
void FillRandomBuffer ()
 
void GenerateMomenta (const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
virtual void GenerateMultiBody (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
virtual void GenerateTwoBody (G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
 
void Initialize (G4double initialMass, const std::vector< G4double > &masses)
 
virtual G4bool IsDecayAllowed (G4double initialMass, const std::vector< G4double > &masses) const
 
void PrintVector (const std::vector< G4double > &v, const G4String &name, std::ostream &os) const
 
G4double TwoBodyMomentum (G4double M0, G4double M1, G4double M2) const
 
G4double UniformPhi () const
 
G4double UniformTheta () const
 
G4ThreeVector UniformVector (G4double mag=1.) const
 

Private Attributes

G4double massExcess
 
std::vector< G4doublemeff
 
std::vector< G4doublemsq
 
std::vector< G4doublemsum
 
G4String name
 
size_t nFinal
 
G4int nTrials
 
std::vector< G4doublepd
 
std::vector< G4doublerndm
 
G4double totalMass
 
G4int verboseLevel
 
G4double weightMax
 

Detailed Description

Definition at line 37 of file G4HadPhaseSpaceGenbod.hh.

Constructor & Destructor Documentation

◆ G4HadPhaseSpaceGenbod()

G4HadPhaseSpaceGenbod::G4HadPhaseSpaceGenbod ( G4int  verbose = 0)

Definition at line 51 of file G4HadPhaseSpaceGenbod.cc.

52 : G4VHadPhaseSpaceAlgorithm("G4HadPhaseSpaceGenbod",verbose),
53 nFinal(0), totalMass(0.), massExcess(0.), weightMax(0.), nTrials(0) {;}
G4VHadPhaseSpaceAlgorithm(const G4String &algName, G4int verbose=0)

◆ ~G4HadPhaseSpaceGenbod()

virtual G4HadPhaseSpaceGenbod::~G4HadPhaseSpaceGenbod ( )
inlinevirtual

Definition at line 40 of file G4HadPhaseSpaceGenbod.hh.

40{;}

Member Function Documentation

◆ AcceptEvent()

G4bool G4HadPhaseSpaceGenbod::AcceptEvent ( ) const
protected

Definition at line 170 of file G4HadPhaseSpaceGenbod.cc.

170 {
171 if (GetVerboseLevel()>1)
172 G4cout << GetName() << "::AcceptEvent? " << nTrials << G4endl;
173
174 return (G4UniformRand() <= ComputeWeight());
175}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
const G4String & GetName() const

References ComputeWeight(), G4cout, G4endl, G4UniformRand, G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::GetVerboseLevel(), and nTrials.

Referenced by GenerateMultiBody().

◆ AccumulateFinalState()

void G4HadPhaseSpaceGenbod::AccumulateFinalState ( size_t  i,
const std::vector< G4double > &  masses,
std::vector< G4LorentzVector > &  finalState 
)
protected

Definition at line 196 of file G4HadPhaseSpaceGenbod.cc.

199 {
200 if (GetVerboseLevel()>2)
201 G4cout << GetName() << "::AccumulateFinalState " << i << G4endl;
202
203 if (i==0) { // First final state particle left alone
204 finalState[i].setVectM(G4ThreeVector(0.,pd[i],0.),masses[i]);
205 return;
206 }
207
208 finalState[i].setVectM(G4ThreeVector(0.,-pd[i-1],0.),masses[i]);
209 G4double phi = G4UniformRand() * twopi;
210 G4double theta = std::acos(2.*G4UniformRand() - 1.);
211
212 if (GetVerboseLevel() > 2) {
213 G4cout << " initialized Py " << -pd[i-1] << " phi " << phi
214 << " theta " << theta << G4endl;
215 }
216
217 G4double esys=0.,beta=0.,gamma=1.;
218 if (i < nFinal-1) { // Do not boost final particle
219 esys = std::sqrt(pd[i]*pd[i]+meff[i]*meff[i]);
220 beta = pd[i] / esys;
221 gamma = esys / meff[i];
222
223 if (GetVerboseLevel()>2)
224 G4cout << " esys " << esys << " beta " << beta << " gamma " << gamma
225 << G4endl;
226 }
227
228 for (size_t j=0; j<=i; j++) { // Accumulate rotations
229 finalState[j].rotateZ(theta).rotateY(phi);
230 finalState[j].setY(gamma*(finalState[j].y() + beta*finalState[j].e()));
231 if (GetVerboseLevel()>2) G4cout << " j " << j << " " << finalState[j] << G4endl;
232 }
233}
static constexpr double twopi
Definition: G4SIunits.hh:56
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
std::vector< G4double > pd
std::vector< G4double > meff

References anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta, G4cout, G4endl, G4UniformRand, G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::GetVerboseLevel(), meff, nFinal, pd, and twopi.

Referenced by GenerateMomenta().

◆ ComputeWeight()

G4double G4HadPhaseSpaceGenbod::ComputeWeight ( ) const
protected

Definition at line 163 of file G4HadPhaseSpaceGenbod.cc.

163 {
164 if (GetVerboseLevel()>1) G4cout << GetName() << "::ComputeWeight" << G4endl;
165
166 return (std::accumulate(pd.begin(), pd.end(), 1./weightMax,
167 std::multiplies<G4double>()));
168}

References G4cout, G4endl, G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::GetVerboseLevel(), pd, and weightMax.

Referenced by AcceptEvent().

◆ ComputeWeightScale()

void G4HadPhaseSpaceGenbod::ComputeWeightScale ( const std::vector< G4double > &  masses)
protected

Definition at line 148 of file G4HadPhaseSpaceGenbod.cc.

148 {
149 if (GetVerboseLevel()>1)
150 G4cout << GetName() << "::ComputeWeightScale" << G4endl;
151
152 weightMax = 1.;
153 for (size_t i=1; i<nFinal; i++) {
154 weightMax *= TwoBodyMomentum(massExcess+msum[i], msum[i-1], masses[i]);
155 }
156
157 if (GetVerboseLevel()>2) G4cout << " weightMax = " << weightMax << G4endl;
158}
std::vector< G4double > msum
G4double TwoBodyMomentum(G4double M0, G4double M1, G4double M2) const

References G4cout, G4endl, G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::GetVerboseLevel(), massExcess, msum, nFinal, G4VHadDecayAlgorithm::TwoBodyMomentum(), and weightMax.

Referenced by Initialize().

◆ FillEnergySteps()

void G4HadPhaseSpaceGenbod::FillEnergySteps ( G4double  initialMass,
const std::vector< G4double > &  masses 
)
protected

Definition at line 123 of file G4HadPhaseSpaceGenbod.cc.

124 {
125 if (GetVerboseLevel()>1) G4cout << GetName() << "::FillEnergySteps" << G4endl;
126
127 meff.clear();
128 pd.clear();
129
130 meff.push_back(masses[0]);
131 for (size_t i=1; i<nFinal-1; i++) {
132 meff.push_back(rndm[i-1]*massExcess + msum[i]);
133 pd.push_back(TwoBodyMomentum(meff[i], meff[i-1], masses[i]));
134 }
135 meff.push_back(initialMass);
136 pd.push_back(TwoBodyMomentum(meff[nFinal-1], meff[nFinal-2], masses[nFinal-1]));
137
138 if (GetVerboseLevel()>2) {
139 PrintVector(meff,"meff",G4cout);
140 PrintVector(pd,"pd",G4cout);
141 }
142}
std::vector< G4double > rndm
void PrintVector(const std::vector< G4double > &v, const G4String &name, std::ostream &os) const

References G4cout, G4endl, G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::GetVerboseLevel(), massExcess, meff, msum, nFinal, pd, G4VHadDecayAlgorithm::PrintVector(), rndm, and G4VHadDecayAlgorithm::TwoBodyMomentum().

Referenced by GenerateMultiBody().

◆ FillRandomBuffer()

void G4HadPhaseSpaceGenbod::FillRandomBuffer ( )
protected

Definition at line 110 of file G4HadPhaseSpaceGenbod.cc.

110 {
111 if (GetVerboseLevel()>1) G4cout << GetName() << "::FillRandomBuffer" << G4endl;
112
113 rndm.resize(nFinal-2,0.); // Final states generated in sorted order
114 std::generate(rndm.begin(), rndm.end(), uniformRand);
115 std::sort(rndm.begin(), rndm.end());
116 if (GetVerboseLevel()>2) PrintVector(rndm, "rndm", G4cout);
117}
void generate(const G4double sqrtS, ParticleList &particles)
Generate an event in the CM system.

References G4cout, G4endl, G4INCL::PhaseSpaceGenerator::generate(), G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::GetVerboseLevel(), nFinal, G4VHadDecayAlgorithm::PrintVector(), rndm, and anonymous_namespace{G4HadPhaseSpaceGenbod.cc}::uniformRand().

Referenced by GenerateMultiBody().

◆ Generate()

void G4VHadDecayAlgorithm::Generate ( G4double  initialMass,
const std::vector< G4double > &  masses,
std::vector< G4LorentzVector > &  finalState 
)
inherited

Definition at line 48 of file G4VHadDecayAlgorithm.cc.

50 {
51 if (verboseLevel) G4cout << GetName() << "::Generate" << G4endl;
52
53 // Initialization and sanity check
54 finalState.clear();
55 if (!IsDecayAllowed(initialMass, masses)) return;
56
57 // Allow different procedures for two-body or N-body distributions
58 if (masses.size() == 2U)
59 GenerateTwoBody(initialMass, masses, finalState);
60 else
61 GenerateMultiBody(initialMass, masses, finalState);
62}
virtual void GenerateTwoBody(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)=0
virtual G4bool IsDecayAllowed(G4double initialMass, const std::vector< G4double > &masses) const
virtual void GenerateMultiBody(G4double initialMass, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)=0

References G4cout, G4endl, G4VHadDecayAlgorithm::GenerateMultiBody(), G4VHadDecayAlgorithm::GenerateTwoBody(), G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::IsDecayAllowed(), and G4VHadDecayAlgorithm::verboseLevel.

Referenced by G4HadDecayGenerator::Generate().

◆ GenerateMomenta()

void G4HadPhaseSpaceGenbod::GenerateMomenta ( const std::vector< G4double > &  masses,
std::vector< G4LorentzVector > &  finalState 
)
protected

Definition at line 180 of file G4HadPhaseSpaceGenbod.cc.

182 {
183 if (GetVerboseLevel()>1) G4cout << GetName() << "::GenerateMomenta" << G4endl;
184
185 finalState.resize(nFinal); // Preallocate vectors for convenience below
186
187 for (size_t i=0; i<nFinal; i++) {
188 AccumulateFinalState(i, masses, finalState);
189 if (GetVerboseLevel()>2)
190 G4cout << " finalState[" << i << "] " << finalState[i] << G4endl;
191 }
192}
void AccumulateFinalState(size_t i, const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)

References AccumulateFinalState(), G4cout, G4endl, G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::GetVerboseLevel(), and nFinal.

Referenced by GenerateMultiBody().

◆ GenerateMultiBody()

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

Implements G4VHadDecayAlgorithm.

Definition at line 58 of file G4HadPhaseSpaceGenbod.cc.

61 {
62 if (GetVerboseLevel()) G4cout << GetName() << "::GenerateMultiBody" << G4endl;
63
64 finalState.clear();
65
66 Initialize(initialMass, masses);
67
68 const G4int maxNumberOfLoops = 10000;
69 nTrials = 0;
70 do { // Apply accept/reject to get distribution
71 ++nTrials;
73 FillEnergySteps(initialMass, masses);
74 } while ( (!AcceptEvent()) && nTrials < maxNumberOfLoops ); /* Loop checking, 02.11.2015, A.Ribon */
75 if ( nTrials >= maxNumberOfLoops ) {
77 ed << " Failed sampling after maxNumberOfLoops attempts : forced exit" << G4endl;
78 G4Exception( " G4HadPhaseSpaceGenbod::GenerateMultiBody ", "HAD_GENBOD_001", FatalException, ed );
79 }
80 GenerateMomenta(masses, finalState);
81}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
int G4int
Definition: G4Types.hh:85
void GenerateMomenta(const std::vector< G4double > &masses, std::vector< G4LorentzVector > &finalState)
void Initialize(G4double initialMass, const std::vector< G4double > &masses)
void FillEnergySteps(G4double initialMass, const std::vector< G4double > &masses)

References AcceptEvent(), FatalException, FillEnergySteps(), FillRandomBuffer(), G4cout, G4endl, G4Exception(), GenerateMomenta(), G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::GetVerboseLevel(), Initialize(), and nTrials.

◆ GenerateTwoBody()

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

Implements G4VHadDecayAlgorithm.

Definition at line 50 of file G4VHadPhaseSpaceAlgorithm.cc.

53 {
54 if (GetVerboseLevel()>1)
55 G4cout << " >>> G4HadDecayGenerator::FillTwoBody" << G4endl;
56
57 // Initialization and sanity check
58 finalState.clear();
59 if (masses.size() != 2U) return; // Should not have been called
60
61 // Momentum of final state (energy balance has already been checked)
62 G4double p = TwoBodyMomentum(initialMass,masses[0],masses[1]);
63 if (GetVerboseLevel()>2) G4cout << " finalState momentum = " << p << G4endl;
64
65 finalState.resize(2); // Allows filling by index
66 finalState[0].setVectM(UniformVector(p), masses[0]);
67 finalState[1].setVectM(-finalState[0].vect(), masses[1]);
68}
G4ThreeVector UniformVector(G4double mag=1.) const

References G4cout, G4endl, G4VHadDecayAlgorithm::GetVerboseLevel(), G4VHadDecayAlgorithm::TwoBodyMomentum(), and G4VHadPhaseSpaceAlgorithm::UniformVector().

◆ GetName()

const G4String & G4VHadDecayAlgorithm::GetName ( ) const
inlineinherited

◆ GetVerboseLevel()

G4int G4VHadDecayAlgorithm::GetVerboseLevel ( ) const
inlineinherited

◆ Initialize()

void G4HadPhaseSpaceGenbod::Initialize ( G4double  initialMass,
const std::vector< G4double > &  masses 
)
protected

Definition at line 83 of file G4HadPhaseSpaceGenbod.cc.

84 {
85 if (GetVerboseLevel()>1) G4cout << GetName() << "::Initialize" << G4endl;
86
87 nFinal = masses.size();
88 msum.resize(nFinal, 0.); // Initialize buffers for filling
89 msq.resize(nFinal, 0.);
90
91 std::partial_sum(masses.begin(), masses.end(), msum.begin());
92 std::transform(masses.begin(), masses.end(), masses.begin(), msq.begin(),
93 std::multiplies<G4double>());
94 totalMass = msum.back();
95 massExcess = initialMass - totalMass;
96
97 if (GetVerboseLevel()>2) {
98 PrintVector(msum, "msum", G4cout);
99 PrintVector(msq, "msq", G4cout);
100 G4cout << " totalMass " << totalMass << " massExcess " << massExcess
101 << G4endl;
102 }
103
104 ComputeWeightScale(masses);
105}
std::vector< G4double > msq
void ComputeWeightScale(const std::vector< G4double > &masses)
G4bool transform(G4String &input, const G4String &type)

References ComputeWeightScale(), G4cout, G4endl, G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::GetVerboseLevel(), massExcess, msq, msum, nFinal, G4VHadDecayAlgorithm::PrintVector(), totalMass, and G4coutFormatters::anonymous_namespace{G4coutFormatters.cc}::transform().

Referenced by GenerateMultiBody().

◆ IsDecayAllowed()

G4bool G4VHadDecayAlgorithm::IsDecayAllowed ( G4double  initialMass,
const std::vector< G4double > &  masses 
) const
protectedvirtualinherited

Definition at line 67 of file G4VHadDecayAlgorithm.cc.

69 {
70 G4bool okay =
71 (initialMass > 0. && masses.size() >= 2 &&
72 initialMass >= std::accumulate(masses.begin(),masses.end(),0.));
73
74 if (verboseLevel) {
75 G4cout << GetName() << "::IsDecayAllowed? initialMass " << initialMass
76 << " " << masses.size() << " masses sum "
77 << std::accumulate(masses.begin(),masses.end(),0.) << G4endl;
78
79 if (verboseLevel>1) PrintVector(masses," ",G4cout);
80
81 G4cout << " Returning " << okay << G4endl;
82 }
83
84 return okay;
85}
bool G4bool
Definition: G4Types.hh:86

References G4cout, G4endl, G4VHadDecayAlgorithm::GetName(), G4VHadDecayAlgorithm::PrintVector(), and G4VHadDecayAlgorithm::verboseLevel.

Referenced by G4VHadDecayAlgorithm::Generate().

◆ PrintVector()

void G4VHadDecayAlgorithm::PrintVector ( const std::vector< G4double > &  v,
const G4String name,
std::ostream &  os 
) const
protectedinherited

Definition at line 121 of file G4VHadDecayAlgorithm.cc.

123 {
124 os << " " << vname << "(" << v.size() << ") ";
125 std::copy(v.begin(), v.end(), std::ostream_iterator<G4double>(os, " "));
126 os << std::endl;
127}
void copy(G4double dst[], const G4double src[], size_t size=G4FieldTrack::ncompSVEC)
Definition: G4FieldUtils.cc:98

References field_utils::copy().

Referenced by FillEnergySteps(), FillRandomBuffer(), G4HadPhaseSpaceNBodyAsai::GenerateMultiBody(), Initialize(), and G4VHadDecayAlgorithm::IsDecayAllowed().

◆ SetVerboseLevel()

virtual void G4VHadDecayAlgorithm::SetVerboseLevel ( G4int  verbose)
inlinevirtualinherited

◆ TwoBodyMomentum()

G4double G4VHadDecayAlgorithm::TwoBodyMomentum ( G4double  M0,
G4double  M1,
G4double  M2 
) const
protectedinherited

Definition at line 90 of file G4VHadDecayAlgorithm.cc.

91 {
92 G4double PSQ = (M0+M1+M2)*(M0+M1-M2)*(M0-M1+M2)*(M0-M1-M2);
93 if (PSQ < 0.) {
94 G4cout << GetName() << ": problem of decay of M(GeV) " << M0/GeV
95 << " to M1(GeV) " << M1/GeV << " and M2(GeV) " << M2/GeV
96 << " PSQ(MeV) " << PSQ/MeV << " < 0" << G4endl;
97 // exception only if the problem is numerically significant
98 if (PSQ < -CLHEP::eV) {
99 throw G4HadronicException(__FILE__, __LINE__,"Error in decay kinematics");
100 }
101
102 PSQ = 0.;
103 }
104
105 return std::sqrt(PSQ)/(2.*M0);
106}
static constexpr double GeV
Definition: G4SIunits.hh:203
static constexpr double MeV
Definition: G4SIunits.hh:200
static constexpr double eV

References CLHEP::eV, G4cout, G4endl, G4VHadDecayAlgorithm::GetName(), GeV, and MeV.

Referenced by ComputeWeightScale(), FillEnergySteps(), G4CascadeFinalStateAlgorithm::FillUsingKopylov(), G4HadPhaseSpaceKopylov::GenerateMultiBody(), G4HadPhaseSpaceNBodyAsai::GenerateMultiBody(), G4CascadeFinalStateAlgorithm::GenerateTwoBody(), and G4VHadPhaseSpaceAlgorithm::GenerateTwoBody().

◆ UniformPhi()

G4double G4VHadDecayAlgorithm::UniformPhi ( ) const
protectedinherited

◆ UniformTheta()

G4double G4VHadDecayAlgorithm::UniformTheta ( ) const
protectedinherited

Definition at line 110 of file G4VHadDecayAlgorithm.cc.

110 {
111 return std::acos(2.0*G4UniformRand() - 1.0);
112}

References G4UniformRand.

Referenced by G4CascadeFinalStateAlgorithm::FillUsingKopylov(), and G4VHadPhaseSpaceAlgorithm::UniformVector().

◆ UniformVector()

G4ThreeVector G4VHadPhaseSpaceAlgorithm::UniformVector ( G4double  mag = 1.) const
protectedinherited

Definition at line 73 of file G4VHadPhaseSpaceAlgorithm.cc.

73 {
74 // FIXME: Should this be made a static thread-local buffer?
77 return v;
78}
void setRThetaPhi(double r, double theta, double phi)
G4double UniformTheta() const

References CLHEP::Hep3Vector::setRThetaPhi(), G4VHadDecayAlgorithm::UniformPhi(), and G4VHadDecayAlgorithm::UniformTheta().

Referenced by G4HadPhaseSpaceKopylov::GenerateMultiBody(), G4HadPhaseSpaceNBodyAsai::GenerateMultiBody(), and G4VHadPhaseSpaceAlgorithm::GenerateTwoBody().

Field Documentation

◆ massExcess

G4double G4HadPhaseSpaceGenbod::massExcess
private

Definition at line 71 of file G4HadPhaseSpaceGenbod.hh.

Referenced by ComputeWeightScale(), FillEnergySteps(), and Initialize().

◆ meff

std::vector<G4double> G4HadPhaseSpaceGenbod::meff
private

Definition at line 78 of file G4HadPhaseSpaceGenbod.hh.

Referenced by AccumulateFinalState(), and FillEnergySteps().

◆ msq

std::vector<G4double> G4HadPhaseSpaceGenbod::msq
private

Definition at line 76 of file G4HadPhaseSpaceGenbod.hh.

Referenced by Initialize().

◆ msum

std::vector<G4double> G4HadPhaseSpaceGenbod::msum
private

Definition at line 75 of file G4HadPhaseSpaceGenbod.hh.

Referenced by ComputeWeightScale(), FillEnergySteps(), and Initialize().

◆ name

G4String G4VHadDecayAlgorithm::name
privateinherited

◆ nFinal

size_t G4HadPhaseSpaceGenbod::nFinal
private

◆ nTrials

G4int G4HadPhaseSpaceGenbod::nTrials
private

Definition at line 73 of file G4HadPhaseSpaceGenbod.hh.

Referenced by AcceptEvent(), and GenerateMultiBody().

◆ pd

std::vector<G4double> G4HadPhaseSpaceGenbod::pd
private

Definition at line 79 of file G4HadPhaseSpaceGenbod.hh.

Referenced by AccumulateFinalState(), ComputeWeight(), and FillEnergySteps().

◆ rndm

std::vector<G4double> G4HadPhaseSpaceGenbod::rndm
private

Definition at line 77 of file G4HadPhaseSpaceGenbod.hh.

Referenced by FillEnergySteps(), and FillRandomBuffer().

◆ totalMass

G4double G4HadPhaseSpaceGenbod::totalMass
private

Definition at line 70 of file G4HadPhaseSpaceGenbod.hh.

Referenced by Initialize().

◆ verboseLevel

G4int G4VHadDecayAlgorithm::verboseLevel
privateinherited

◆ weightMax

G4double G4HadPhaseSpaceGenbod::weightMax
private

Definition at line 72 of file G4HadPhaseSpaceGenbod.hh.

Referenced by ComputeWeight(), and ComputeWeightScale().


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