Geant4  10.02.p03
G4PreCompoundFragment Class Referenceabstract

#include <G4PreCompoundFragment.hh>

Inheritance diagram for G4PreCompoundFragment:
Collaboration diagram for G4PreCompoundFragment:

Public Member Functions

 G4PreCompoundFragment (const G4ParticleDefinition *, G4VCoulombBarrier *aCoulombBarrier)
 
virtual ~G4PreCompoundFragment ()
 
G4double CalcEmissionProbability (const G4Fragment &aFragment)
 
G4double SampleKineticEnergy (const G4Fragment &aFragment)
 
- Public Member Functions inherited from G4VPreCompoundFragment
 G4VPreCompoundFragment (const G4ParticleDefinition *, G4VCoulombBarrier *aCoulombBarrier)
 
virtual ~G4VPreCompoundFragment ()
 
void Initialize (const G4Fragment &aFragment)
 
G4bool IsItPossible (const G4Fragment &aFragment) const
 
G4ReactionProductGetReactionProduct () const
 
G4int GetA () const
 
G4int GetZ () const
 
G4int GetRestA () const
 
G4int GetRestZ () const
 
G4double GetBindingEnergy () const
 
G4double GetEnergyThreshold () const
 
G4double GetEmissionProbability () const
 
G4double GetNuclearMass () const
 
G4double GetRestNuclearMass () const
 
const G4LorentzVectorGetMomentum () const
 
void SetMomentum (const G4LorentzVector &value)
 
void SetOPTxs (G4int)
 
void UseSICB (G4bool)
 

Protected Member Functions

virtual G4double GetAlpha () const =0
 
virtual G4double GetBeta () const =0
 
G4double CrossSection (G4double ekin) const
 
virtual G4double ProbabilityDistributionFunction (G4double K, const G4Fragment &aFragment)=0
 

Private Member Functions

G4double IntegrateEmissionProbability (G4double Low, G4double Up, const G4Fragment &aFragment)
 
G4double GetOpt0 (G4double ekin) const
 
 G4PreCompoundFragment ()
 
 G4PreCompoundFragment (const G4PreCompoundFragment &right)
 
const G4PreCompoundFragmentoperator= (const G4PreCompoundFragment &right)
 
G4int operator== (const G4PreCompoundFragment &right) const
 
G4int operator!= (const G4PreCompoundFragment &right) const
 

Private Attributes

G4int index
 
G4double muu
 
G4double probability [12]
 
G4double probmax
 

Additional Inherited Members

- Protected Attributes inherited from G4VPreCompoundFragment
G4PreCompoundParameterstheParameters
 
G4Powg4pow
 
G4int theA
 
G4int theZ
 
G4int theResA
 
G4int theResZ
 
G4int theFragA
 
G4int theFragZ
 
G4double theResA13
 
G4double theBindingEnergy
 
G4double theMaxKinEnergy
 
G4double theResMass
 
G4double theReducedMass
 
G4double theMass
 
G4double theEmissionProbability
 
G4double theCoulombBarrier
 
G4int OPTxs
 
G4bool useSICB
 

Detailed Description

Definition at line 44 of file G4PreCompoundFragment.hh.

Constructor & Destructor Documentation

◆ G4PreCompoundFragment() [1/3]

G4PreCompoundFragment::G4PreCompoundFragment ( const G4ParticleDefinition part,
G4VCoulombBarrier aCoulombBarrier 
)

Definition at line 70 of file G4PreCompoundFragment.cc.

72  : G4VPreCompoundFragment(part,aCoulombBarrier)
73 {
74  muu = probmax = 0.0;
75  index = 0;
76  if(1 == theZ) { index = theA; }
77  else { index = theA + 1; }
78  for(G4int i=0; i<12; ++i) { probability[i] = 0.0; }
79 }
int G4int
Definition: G4Types.hh:78

◆ ~G4PreCompoundFragment()

G4PreCompoundFragment::~G4PreCompoundFragment ( )
virtual

Definition at line 81 of file G4PreCompoundFragment.cc.

82 {}
Here is the call graph for this function:

◆ G4PreCompoundFragment() [2/3]

G4PreCompoundFragment::G4PreCompoundFragment ( )
private

◆ G4PreCompoundFragment() [3/3]

G4PreCompoundFragment::G4PreCompoundFragment ( const G4PreCompoundFragment right)
private

Member Function Documentation

◆ CalcEmissionProbability()

G4double G4PreCompoundFragment::CalcEmissionProbability ( const G4Fragment aFragment)
virtual

Implements G4VPreCompoundFragment.

Definition at line 85 of file G4PreCompoundFragment.cc.

86 {
87  //G4cout << theCoulombBarrier << " " << GetMaximalKineticEnergy() << G4endl;
88  // If theCoulombBarrier effect is included in the emission probabilities
89  // Coulomb barrier is the lower limit of integration over kinetic energy
90 
91  G4double LowerLimit = 0.0;
92  if(OPTxs==0 || useSICB) { LowerLimit = theCoulombBarrier; }
93 
94  if (theMaxKinEnergy <= LowerLimit)
95  {
97  return 0.0;
98  }
99 
100  // compute power once
101  if(OPTxs <= 2) {
103  } else {
105  }
106 
108  IntegrateEmissionProbability(LowerLimit,theMaxKinEnergy,aFragment);
109  /*
110  G4cout << "## G4PreCompoundFragment::CalcEmisProb "
111  << "Z= " << aFragment.GetZ_asInt()
112  << " A= " << aFragment.GetA_asInt()
113  << " Elow= " << LowerLimit/MeV
114  << " Eup= " << UpperLimit/MeV
115  << " prob= " << theEmissionProbability
116  << G4endl;
117  */
118  return theEmissionProbability;
119 }
static G4double ComputePowerParameter(G4int resA, G4int idx)
static G4double ComputePowerParameter(G4int resA, G4int idx)
G4double IntegrateEmissionProbability(G4double Low, G4double Up, const G4Fragment &aFragment)
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ CrossSection()

G4double G4PreCompoundFragment::CrossSection ( G4double  ekin) const
protected

Definition at line 142 of file G4PreCompoundFragment.cc.

143 {
144  G4double res;
145  if(OPTxs == 0) {
146  res = GetOpt0(ekin);
147  } else if(OPTxs <= 2) {
149  index, theZ,
150  theResZ, theResA);
151  } else {
153  index, theZ, theA,
154  theResZ, theResA);
155  }
156  return res;
157 }
static G4double ComputeCrossSection(G4double K, G4double resA13, G4double amu1, G4int idx, G4int Z, G4int resZ, G4int resA)
G4double GetOpt0(G4double ekin) const
double G4double
Definition: G4Types.hh:76
static G4double ComputeCrossSection(G4double K, G4double resA13, G4double amu1, G4int idx, G4int Z, G4int A, G4int resZ, G4int resA)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetAlpha()

virtual G4double G4PreCompoundFragment::GetAlpha ( ) const
protectedpure virtual

Implements G4VPreCompoundFragment.

Implemented in G4PreCompoundAlpha, G4PreCompoundDeuteron, G4PreCompoundHe3, G4PreCompoundTriton, G4PreCompoundNeutron, and G4PreCompoundProton.

Here is the caller graph for this function:

◆ GetBeta()

virtual G4double G4PreCompoundFragment::GetBeta ( ) const
protectedpure virtual

Implements G4VPreCompoundFragment.

Implemented in G4PreCompoundIon, G4PreCompoundNeutron, and G4PreCompoundProton.

Here is the caller graph for this function:

◆ GetOpt0()

G4double G4PreCompoundFragment::GetOpt0 ( G4double  ekin) const
private

Definition at line 160 of file G4PreCompoundFragment.cc.

161 {
163  // cross section is now given in mb (r0 is in mm) for the sake of consistency
164  //with the rest of the options
165  return 1.e+25*CLHEP::pi*r0*r0*theResA13*GetAlpha()*(1.+GetBeta()/ekin);
166 }
virtual G4double GetAlpha() const =0
virtual G4double GetBeta() const =0
static const double pi
Definition: SystemOfUnits.h:53
G4PreCompoundParameters * theParameters
double G4double
Definition: G4Types.hh:76
const G4double r0
Here is the call graph for this function:
Here is the caller graph for this function:

◆ IntegrateEmissionProbability()

G4double G4PreCompoundFragment::IntegrateEmissionProbability ( G4double  Low,
G4double  Up,
const G4Fragment aFragment 
)
private

Definition at line 122 of file G4PreCompoundFragment.cc.

124 {
125  G4double sum = 0.0;
126  G4double del = (up - low)*0.5;
127  G4double avr = (up + low)*0.5;
128 
129  G4double e, y;
130  probmax = 0.0;
131 
132  for (G4int i=0; i<NPOINTSGL; ++i) {
133  e = del*x[i] + avr;
134  y = ProbabilityDistributionFunction(e, aFragment);
135  probability[i] = y;
136  probmax = std::max(probmax, y);
137  sum += w[i]*y;
138  }
139  return sum*del;
140 }
const G4double w[NPOINTSGL]
int G4int
Definition: G4Types.hh:78
Double_t y
const G4int NPOINTSGL
virtual G4double ProbabilityDistributionFunction(G4double K, const G4Fragment &aFragment)=0
const G4double x[NPOINTSGL]
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ operator!=()

G4int G4PreCompoundFragment::operator!= ( const G4PreCompoundFragment right) const
private

◆ operator=()

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

◆ operator==()

G4int G4PreCompoundFragment::operator== ( const G4PreCompoundFragment right) const
private

◆ ProbabilityDistributionFunction()

virtual G4double G4PreCompoundFragment::ProbabilityDistributionFunction ( G4double  K,
const G4Fragment aFragment 
)
protectedpure virtual

Implemented in G4PreCompoundIon, and G4PreCompoundNucleon.

Here is the caller graph for this function:

◆ SampleKineticEnergy()

G4double G4PreCompoundFragment::SampleKineticEnergy ( const G4Fragment aFragment)
virtual

Implements G4VPreCompoundFragment.

Definition at line 168 of file G4PreCompoundFragment.cc.

169 {
170  //let's keep this way for consistency with CalcEmissionProbability method
171  G4double limit = 0.0;
172  if(useSICB) { limit = theCoulombBarrier; }
173 
174  if(theMaxKinEnergy <= limit) { return 0.0; }
175 
176  G4double delta = theMaxKinEnergy - limit;
177  static const G4double toler = 1.25;
178  probmax *= toler;
179  G4double prob, T(0.0);
180  CLHEP::HepRandomEngine* rndm = G4Random::getTheEngine();
181  G4int i;
182  for(i=0; i<100; ++i) {
183  T = limit + delta*rndm->flat();
184  prob = ProbabilityDistributionFunction(T, fragment);
185  /*
186  if(prob > probmax) {
187  G4cout << "G4PreCompoundFragment WARNING: prob= " << prob
188  << " probmax= " << probmax << G4endl;
189  G4cout << "i= " << i << " Z= " << theZ << " A= " << theA
190  << " resZ= " << theResZ << " resA= " << theResA << "\n"
191  << " T= " << T << " Tmax= " << theMaxKinEnergy
192  << " Tmin= " << limit
193  << G4endl;
194  for(G4int i=0; i<N; ++i) { G4cout << " " << probability[i]; }
195  G4cout << G4endl;
196  }
197  */
198  // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
199  if(probmax*rndm->flat() <= prob) { break; }
200  }
201  //G4cout << "G4PreCompoundFragment: i= " << i << " T(MeV)= " << T
202  //<< " Emin(MeV)= " << limit << " Emax= " << theMaxKinEnergy << G4endl;
203  return T;
204 }
virtual double flat()=0
int G4int
Definition: G4Types.hh:78
virtual G4double ProbabilityDistributionFunction(G4double K, const G4Fragment &aFragment)=0
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

Member Data Documentation

◆ index

G4int G4PreCompoundFragment::index
private

Definition at line 92 of file G4PreCompoundFragment.hh.

◆ muu

G4double G4PreCompoundFragment::muu
private

Definition at line 94 of file G4PreCompoundFragment.hh.

◆ probability

G4double G4PreCompoundFragment::probability[12]
private

Definition at line 95 of file G4PreCompoundFragment.hh.

◆ probmax

G4double G4PreCompoundFragment::probmax
private

Definition at line 96 of file G4PreCompoundFragment.hh.


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