Geant4  10.02.p03
G4PreCompoundModel Class Reference

#include <G4PreCompoundModel.hh>

Inheritance diagram for G4PreCompoundModel:
Collaboration diagram for G4PreCompoundModel:

Public Member Functions

 G4PreCompoundModel (G4ExcitationHandler *ptr=0)
 
virtual ~G4PreCompoundModel ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &thePrimary, G4Nucleus &theNucleus)
 
virtual G4ReactionProductVectorDeExcite (G4Fragment &aFragment)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void DeExciteModelDescription (std::ostream &outFile) const
 
void UseHETCEmission ()
 
void UseDefaultEmission ()
 
void UseGNASHTransition ()
 
void UseDefaultTransition ()
 
void SetOPTxs (G4int opt)
 
void UseSICB ()
 
void UseNGB ()
 
void UseSCO ()
 
void UseCEMtr ()
 
- Public Member Functions inherited from G4VPreCompoundModel
 G4VPreCompoundModel (G4ExcitationHandler *ptr=0, const G4String &modelName="PrecompoundModel")
 
virtual ~G4VPreCompoundModel ()
 
void SetExcitationHandler (G4ExcitationHandler *ptr)
 
G4ExcitationHandlerGetExcitationHandler () const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &, G4Nucleus &)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
const G4HadronicInteractionGetMyPointer () const
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
G4bool operator== (const G4HadronicInteraction &right) const
 
G4bool operator!= (const G4HadronicInteraction &right) const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 

Private Member Functions

void PerformEquilibriumEmission (const G4Fragment &aFragment, G4ReactionProductVector *theResult) const
 
 G4PreCompoundModel (const G4PreCompoundModel &)
 
const G4PreCompoundModeloperator= (const G4PreCompoundModel &right)
 
G4bool operator== (const G4PreCompoundModel &right) const
 
G4bool operator!= (const G4PreCompoundModel &right) const
 

Private Attributes

G4PreCompoundEmissiontheEmission
 
G4VPreCompoundTransitionstheTransition
 
const G4ParticleDefinitionproton
 
const G4ParticleDefinitionneutron
 
G4double fLevelDensity
 
G4bool useHETCEmission
 
G4bool useGNASHTransition
 
G4int OPTxs
 
G4bool useSICB
 
G4bool useNGB
 
G4bool useSCO
 
G4bool useCEMtr
 
G4int maxZ
 
G4int maxA
 
G4HadFinalState theResult
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 64 of file G4PreCompoundModel.hh.

Constructor & Destructor Documentation

◆ G4PreCompoundModel() [1/2]

G4PreCompoundModel::G4PreCompoundModel ( G4ExcitationHandler ptr = 0)

Definition at line 70 of file G4PreCompoundModel.cc.

71  : G4VPreCompoundModel(ptr,"PRECO"), useHETCEmission(false),
72  useGNASHTransition(false), OPTxs(3), useSICB(false),
73  useNGB(false), useSCO(false), useCEMtr(true), maxZ(3), maxA(5)
74 {
75  if(!ptr) { SetExcitationHandler(new G4ExcitationHandler()); }
77 
78  // 12/pi2 factor is used in real computation
80 
83  else { theEmission->SetDefaultModel(); }
86 
88  else { theTransition = new G4PreCompoundTransitions(); }
91 
94 }
G4PreCompoundEmission * theEmission
const G4ParticleDefinition * neutron
void SetExcitationHandler(G4ExcitationHandler *ptr)
G4VPreCompoundTransitions * theTransition
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static const double pi2
Definition: SystemOfUnits.h:56
const G4ParticleDefinition * proton
G4VPreCompoundModel(G4ExcitationHandler *ptr=0, const G4String &modelName="PrecompoundModel")
Here is the call graph for this function:

◆ ~G4PreCompoundModel()

G4PreCompoundModel::~G4PreCompoundModel ( )
virtual

Definition at line 98 of file G4PreCompoundModel.cc.

99 {
100  delete theEmission;
101  delete theTransition;
102  delete GetExcitationHandler();
103 }
G4PreCompoundEmission * theEmission
G4VPreCompoundTransitions * theTransition
G4ExcitationHandler * GetExcitationHandler() const
Here is the call graph for this function:

◆ G4PreCompoundModel() [2/2]

G4PreCompoundModel::G4PreCompoundModel ( const G4PreCompoundModel )
private

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4PreCompoundModel::ApplyYourself ( const G4HadProjectile thePrimary,
G4Nucleus theNucleus 
)
virtual

Implements G4VPreCompoundModel.

Definition at line 108 of file G4PreCompoundModel.cc.

110 {
111  const G4ParticleDefinition* primary = thePrimary.GetDefinition();
112  if(primary != neutron && primary != proton) {
114  ed << "G4PreCompoundModel is used for ";
115  if(primary) { ed << primary->GetParticleName(); }
116  G4Exception("G4PreCompoundModel::ApplyYourself()","had0033",FatalException,
117  ed,"");
118  return 0;
119  }
120 
121  G4int Zp = 0;
122  G4int Ap = 1;
123  if(primary == proton) { Zp = 1; }
124 
125  G4int A = theNucleus.GetA_asInt();
126  G4int Z = theNucleus.GetZ_asInt();
127 
128  //G4cout << "### G4PreCompoundModel::ApplyYourself: A= " << A << " Z= " << Z
129  // << " Ap= " << Ap << " Zp= " << Zp << G4endl;
130  // 4-Momentum
131  G4LorentzVector p = thePrimary.Get4Momentum();
133  p += G4LorentzVector(0.0,0.0,0.0,mass);
134  //G4cout << "Primary 4-mom " << p << " mass= " << mass << G4endl;
135 
136  // prepare fragment
137  G4Fragment anInitialState(A + Ap, Z + Zp, p);
138  anInitialState.SetNumberOfExcitedParticle(2, 1);
139  anInitialState.SetNumberOfHoles(1,0);
140  anInitialState.SetCreationTime(thePrimary.GetGlobalTime());
141 
142  // call excitation handler
143  G4ReactionProductVector * result = DeExcite(anInitialState);
144 
145  // fill particle change
146  theResult.Clear();
148  for(G4ReactionProductVector::iterator i= result->begin();
149  i != result->end(); ++i)
150  {
151  G4DynamicParticle * aNew =
152  new G4DynamicParticle((*i)->GetDefinition(),
153  (*i)->GetTotalEnergy(),
154  (*i)->GetMomentum());
155  delete (*i);
156  theResult.AddSecondary(aNew);
157  }
158  delete result;
159 
160  //return the filled particle change
161  return &theResult;
162 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4LorentzVector & Get4Momentum() const
const G4ParticleDefinition * neutron
G4int GetA_asInt() const
Definition: G4Nucleus.hh:109
int G4int
Definition: G4Types.hh:78
G4HadFinalState theResult
void SetStatusChange(G4HadFinalStateStatus aS)
std::vector< G4ReactionProduct * > G4ReactionProductVector
const G4String & GetParticleName() const
double A(double temperature)
Float_t Z
G4double GetGlobalTime() const
virtual G4ReactionProductVector * DeExcite(G4Fragment &aFragment)
G4int GetZ_asInt() const
Definition: G4Nucleus.hh:115
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
const G4ParticleDefinition * GetDefinition() const
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
const G4ParticleDefinition * proton
double G4double
Definition: G4Types.hh:76
CLHEP::HepLorentzVector G4LorentzVector
Here is the call graph for this function:

◆ DeExcite()

G4ReactionProductVector * G4PreCompoundModel::DeExcite ( G4Fragment aFragment)
virtual

Implements G4VPreCompoundModel.

Definition at line 166 of file G4PreCompoundModel.cc.

167 {
169  G4double Eex = aFragment.GetExcitationEnergy();
170  G4int Z = aFragment.GetZ_asInt();
171  G4int A = aFragment.GetA_asInt();
172 
173  //G4cout << "### G4PreCompoundModel::DeExcite" << G4endl;
174  //G4cout << aFragment << G4endl;
175 
176  // Perform Equilibrium Emission
177  if ((Z < maxZ && A < maxA) || Eex < MeV /*|| Eex > 3.*MeV*A*/) {
178  PerformEquilibriumEmission(aFragment, Result);
179  return Result;
180  }
181 
182  // main loop
183  G4int count = 0;
184  const G4int countmax = 1000;
185  for (;;) {
186  //G4cout << "### PreCompound loop over fragment" << G4endl;
187  //G4cout << aFragment << G4endl;
188  G4int EquilibriumExcitonNumber =
189  G4lrint(std::sqrt(aFragment.GetExcitationEnergy()
190  *aFragment.GetA_asInt()*fLevelDensity));
191  //
192  // G4cout<<"Neq="<<EquilibriumExcitonNumber<<G4endl;
193  //
194  // J. M. Quesada (Jan. 08) equilibrium hole number could be used as preeq.
195  // evap. delimiter (IAEA report)
196 
197  // Loop for transitions, it is performed while there are
198  // preequilibrium transitions.
199  G4bool ThereIsTransition = false;
200 
201  // G4cout<<"----------------------------------------"<<G4endl;
202  // G4double NP=aFragment.GetNumberOfParticles();
203  // G4double NH=aFragment.GetNumberOfHoles();
204  // G4double NE=aFragment.GetNumberOfExcitons();
205  // G4cout<<" Ex. Energy="<<aFragment.GetExcitationEnergy()<<G4endl;
206  // G4cout<<"N. excitons="<<NE<<" N. Part="<<NP<<"N. Holes ="<<NH<<G4endl;
207  do {
208  ++count;
209  //G4cout<<"transition number .."<<count
210  // <<" n ="<<aFragment.GetNumberOfExcitons()<<G4endl;
211  G4bool go_ahead = false;
212  // soft cutoff criterium as an "ad-hoc" solution to force
213  // increase in evaporation
214  G4int test = aFragment.GetNumberOfExcitons();
215  if (test <= EquilibriumExcitonNumber) { go_ahead=true; }
216 
217  //J. M. Quesada (Apr. 08): soft-cutoff switched off by default
218  if (useSCO && go_ahead)
219  {
220  G4double x = G4double(test)/G4double(EquilibriumExcitonNumber) - 1;
221  if( G4UniformRand() < 1.0 - G4Exp(-x*x/0.32) ) { go_ahead = false; }
222  }
223 
224  // JMQ: WARNING: CalculateProbability MUST be called prior to Get!!
225  // (O values would be returned otherwise)
226  G4double TotalTransitionProbability =
231  //G4cout<<"#0 P1="<<P1<<" P2="<<P2<<" P3="<<P3<<G4endl;
232 
233  //J.M. Quesada (May 2008) Physical criterium (lamdas) PREVAILS over
234  // approximation (critical exciton number)
235  //V.Ivanchenko (May 2011) added check on number of nucleons
236  // to send a fragment to FermiBreakUp
237  if(!go_ahead || P1 <= P2+P3 ||
238  (aFragment.GetZ_asInt() < maxZ && aFragment.GetA_asInt() < maxA) )
239  {
240  //G4cout<<"#4 EquilibriumEmission"<<G4endl;
241  PerformEquilibriumEmission(aFragment,Result);
242  return Result;
243  }
244  else
245  {
246  //
247  // Check if number of excitons is greater than 0
248  // else perform equilibrium emission
249  if (aFragment.GetNumberOfExcitons() <= 0)
250  {
251  PerformEquilibriumEmission(aFragment,Result);
252  return Result;
253  }
254 
255  G4double TotalEmissionProbability =
256  theEmission->GetTotalProbability(aFragment);
257  //
258  //G4cout<<"#1 TotalEmissionProbability="<<TotalEmissionProbability
259  // <<" Nex= " <<aFragment.GetNumberOfExcitons()<<G4endl;
260  //J.M.Quesada (May 08) this has already been done in order to decide
261  // what to do (preeq-eq)
262  // Sum of all probabilities
263  G4double TotalProbability = TotalEmissionProbability
264  + TotalTransitionProbability;
265 
266  // Select subprocess
267  if (TotalProbability*G4UniformRand() > TotalEmissionProbability)
268  {
269  //G4cout<<"#2 Transition"<<G4endl;
270  // It will be transition to state with a new number of excitons
271  ThereIsTransition = true;
272  // Perform the transition
273  theTransition->PerformTransition(aFragment);
274  }
275  else
276  {
277  //G4cout<<"#3 Emission"<<G4endl;
278  // It will be fragment emission
279  ThereIsTransition = false;
280  Result->push_back(theEmission->PerformEmission(aFragment));
281  }
282  }
283  // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
284  } while (ThereIsTransition); // end of do loop
285 
286  // stop if too many iterations
287  if(count >= countmax) {
289  ed << "G4PreCompoundModel loop over " << countmax << " iterations; "
290  << "current G4Fragment: \n" << aFragment;
291  G4Exception("G4PreCompoundModel::DeExcite()","had0034",JustWarning,
292  ed,"");
293  PerformEquilibriumEmission(aFragment, Result);
294  return Result;
295  }
296  } // end of for (;;) loop
297  return Result;
298 }
G4int GetNumberOfExcitons() const
Definition: G4Fragment.hh:322
virtual void PerformTransition(G4Fragment &aFragment)=0
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:273
static const double MeV
Definition: G4SIunits.hh:211
static const G4double * P1[nN]
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4PreCompoundEmission * theEmission
G4int GetA_asInt() const
Definition: G4Fragment.hh:256
int G4int
Definition: G4Types.hh:78
std::vector< G4ReactionProduct * > G4ReactionProductVector
#define G4UniformRand()
Definition: Randomize.hh:97
double A(double temperature)
G4int GetZ_asInt() const
Definition: G4Fragment.hh:261
G4VPreCompoundTransitions * theTransition
Float_t Z
bool G4bool
Definition: G4Types.hh:79
void PerformEquilibriumEmission(const G4Fragment &aFragment, G4ReactionProductVector *theResult) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
virtual G4double CalculateProbability(const G4Fragment &aFragment)=0
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4ReactionProduct * PerformEmission(G4Fragment &aFragment)
int G4lrint(double ad)
Definition: templates.hh:163
static const G4double * P2[nN]
double G4double
Definition: G4Types.hh:76
G4double GetTotalProbability(const G4Fragment &aFragment)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DeExciteModelDescription()

void G4PreCompoundModel::DeExciteModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VPreCompoundModel.

Definition at line 386 of file G4PreCompoundModel.cc.

387 {
388  outFile << "description of precompound model as used with DeExcite()"
389  << "\n";
390 }

◆ ModelDescription()

void G4PreCompoundModel::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 363 of file G4PreCompoundModel.cc.

364 {
365  outFile << "The GEANT4 precompound model is considered as an extension of the\n"
366  << "hadron kinetic model. It gives a possibility to extend the low energy range\n"
367  << "of the hadron kinetic model for nucleon-nucleus inelastic collision and it \n"
368  << "provides a ”smooth” transition from kinetic stage of reaction described by the\n"
369  << "hadron kinetic model to the equilibrium stage of reaction described by the\n"
370  << "equilibrium deexcitation models.\n"
371  << "The initial information for calculation of pre-compound nuclear stage\n"
372  << "consists of the atomic mass number A, charge Z of residual nucleus, its\n"
373  << "four momentum P0 , excitation energy U and number of excitons n, which equals\n"
374  << "the sum of the number of particles p (from them p_Z are charged) and the number of\n"
375  << "holes h.\n"
376  << "At the preequilibrium stage of reaction, we follow the exciton model approach in ref. [1],\n"
377  << "taking into account the competition among all possible nuclear transitions\n"
378  << "with ∆n = +2, −2, 0 (which are defined by their associated transition probabilities) and\n"
379  << "the emission of neutrons, protons, deutrons, thritium and helium nuclei (also defined by\n"
380  << "their associated emission probabilities according to exciton model)\n"
381  << "\n"
382  << "[1] K.K. Gudima, S.G. Mashnik, V.D. Toneev, Nucl. Phys. A401 329 (1983)\n"
383  << "\n";
384 }

◆ operator!=()

G4bool G4PreCompoundModel::operator!= ( const G4PreCompoundModel right) const
private

◆ operator=()

const G4PreCompoundModel& G4PreCompoundModel::operator= ( const G4PreCompoundModel right)
private

◆ operator==()

G4bool G4PreCompoundModel::operator== ( const G4PreCompoundModel right) const
private

◆ PerformEquilibriumEmission()

void G4PreCompoundModel::PerformEquilibriumEmission ( const G4Fragment aFragment,
G4ReactionProductVector theResult 
) const
inlineprivate

Definition at line 137 of file G4PreCompoundModel.hh.

139 {
140  G4ReactionProductVector* theEquilibriumResult =
141  GetExcitationHandler()->BreakItUp(aFragment);
142  Result->insert(Result->end(),theEquilibriumResult->begin(), theEquilibriumResult->end());
143  delete theEquilibriumResult;
144 }
G4ReactionProductVector * BreakItUp(const G4Fragment &theInitialState)
std::vector< G4ReactionProduct * > G4ReactionProductVector
G4ExcitationHandler * GetExcitationHandler() const
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetOPTxs()

void G4PreCompoundModel::SetOPTxs ( G4int  opt)

Definition at line 332 of file G4PreCompoundModel.cc.

333 {
334  OPTxs = opt;
336 }
G4PreCompoundEmission * theEmission
Here is the call graph for this function:

◆ UseCEMtr()

void G4PreCompoundModel::UseCEMtr ( )

Definition at line 354 of file G4PreCompoundModel.cc.

355 {
356  useCEMtr = true;
357 }

◆ UseDefaultEmission()

void G4PreCompoundModel::UseDefaultEmission ( )

Definition at line 310 of file G4PreCompoundModel.cc.

311 {
312  useHETCEmission = false;
314 }
G4PreCompoundEmission * theEmission
Here is the call graph for this function:

◆ UseDefaultTransition()

void G4PreCompoundModel::UseDefaultTransition ( )

Definition at line 324 of file G4PreCompoundModel.cc.

Here is the call graph for this function:

◆ UseGNASHTransition()

void G4PreCompoundModel::UseGNASHTransition ( )

Definition at line 316 of file G4PreCompoundModel.cc.

Here is the call graph for this function:

◆ UseHETCEmission()

void G4PreCompoundModel::UseHETCEmission ( )

Definition at line 304 of file G4PreCompoundModel.cc.

305 {
306  useHETCEmission = true;
308 }
G4PreCompoundEmission * theEmission
Here is the call graph for this function:

◆ UseNGB()

void G4PreCompoundModel::UseNGB ( )

Definition at line 344 of file G4PreCompoundModel.cc.

345 {
346  useNGB = true;
347 }

◆ UseSCO()

void G4PreCompoundModel::UseSCO ( )

Definition at line 349 of file G4PreCompoundModel.cc.

350 {
351  useSCO = true;
352 }

◆ UseSICB()

void G4PreCompoundModel::UseSICB ( )

Definition at line 338 of file G4PreCompoundModel.cc.

339 {
340  useSICB = true;
342 }
G4PreCompoundEmission * theEmission
Here is the call graph for this function:

Member Data Documentation

◆ fLevelDensity

G4double G4PreCompoundModel::fLevelDensity
private

Definition at line 115 of file G4PreCompoundModel.hh.

◆ maxA

G4int G4PreCompoundModel::maxA
private

Definition at line 130 of file G4PreCompoundModel.hh.

◆ maxZ

G4int G4PreCompoundModel::maxZ
private

Definition at line 129 of file G4PreCompoundModel.hh.

◆ neutron

const G4ParticleDefinition* G4PreCompoundModel::neutron
private

Definition at line 113 of file G4PreCompoundModel.hh.

◆ OPTxs

G4int G4PreCompoundModel::OPTxs
private

Definition at line 121 of file G4PreCompoundModel.hh.

◆ proton

const G4ParticleDefinition* G4PreCompoundModel::proton
private

Definition at line 112 of file G4PreCompoundModel.hh.

◆ theEmission

G4PreCompoundEmission* G4PreCompoundModel::theEmission
private

Definition at line 109 of file G4PreCompoundModel.hh.

◆ theResult

G4HadFinalState G4PreCompoundModel::theResult
private

Definition at line 132 of file G4PreCompoundModel.hh.

◆ theTransition

G4VPreCompoundTransitions* G4PreCompoundModel::theTransition
private

Definition at line 110 of file G4PreCompoundModel.hh.

◆ useCEMtr

G4bool G4PreCompoundModel::useCEMtr
private

Definition at line 127 of file G4PreCompoundModel.hh.

◆ useGNASHTransition

G4bool G4PreCompoundModel::useGNASHTransition
private

Definition at line 118 of file G4PreCompoundModel.hh.

◆ useHETCEmission

G4bool G4PreCompoundModel::useHETCEmission
private

Definition at line 117 of file G4PreCompoundModel.hh.

◆ useNGB

G4bool G4PreCompoundModel::useNGB
private

Definition at line 125 of file G4PreCompoundModel.hh.

◆ useSCO

G4bool G4PreCompoundModel::useSCO
private

Definition at line 126 of file G4PreCompoundModel.hh.

◆ useSICB

G4bool G4PreCompoundModel::useSICB
private

Definition at line 124 of file G4PreCompoundModel.hh.


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