#include <G4GEMChannel.hh>
Inherits G4VEvaporationChannel.
Inherited by G4AlphaGEMChannel, G4B10GEMChannel, G4B11GEMChannel, G4B12GEMChannel, G4B13GEMChannel, G4B8GEMChannel, G4Be10GEMChannel, G4Be11GEMChannel, G4Be12GEMChannel, G4Be7GEMChannel, G4Be9GEMChannel, G4C10GEMChannel, G4C11GEMChannel, G4C12GEMChannel, G4C13GEMChannel, G4C14GEMChannel, G4C15GEMChannel, G4C16GEMChannel, G4DeuteronGEMChannel, G4F17GEMChannel, G4F18GEMChannel, G4F19GEMChannel, G4F20GEMChannel, G4F21GEMChannel, G4He3GEMChannel, G4He6GEMChannel, G4He8GEMChannel, G4Li6GEMChannel, G4Li7GEMChannel, G4Li8GEMChannel, G4Li9GEMChannel, G4Mg22GEMChannel, G4Mg23GEMChannel, G4Mg24GEMChannel, G4Mg25GEMChannel, G4Mg26GEMChannel, G4Mg27GEMChannel, G4Mg28GEMChannel, G4N12GEMChannel, G4N13GEMChannel, G4N14GEMChannel, G4N15GEMChannel, G4N16GEMChannel, G4N17GEMChannel, G4Na21GEMChannel, G4Na22GEMChannel, G4Na23GEMChannel, G4Na24GEMChannel, G4Na25GEMChannel, G4Ne18GEMChannel, G4Ne19GEMChannel, G4Ne20GEMChannel, G4Ne21GEMChannel, G4Ne22GEMChannel, G4Ne23GEMChannel, G4Ne24GEMChannel, G4NeutronGEMChannel, G4O14GEMChannel, G4O15GEMChannel, G4O16GEMChannel, G4O17GEMChannel, G4O18GEMChannel, G4O19GEMChannel, G4O20GEMChannel, G4ProtonGEMChannel, and G4TritonGEMChannel.
|
| G4GEMChannel (const G4int theA, const G4int theZ, const G4String &aName, G4GEMProbability *aEmissionStrategy, G4VCoulombBarrier *aCoulombBarrier) |
|
virtual | ~G4GEMChannel () |
|
virtual G4double | GetEmissionProbability (G4Fragment *theNucleus) |
|
virtual G4Fragment * | EmittedFragment (G4Fragment *theNucleus) |
|
virtual G4FragmentVector * | BreakUp (const G4Fragment &theNucleus) |
|
virtual void | Dump () const |
|
void | SetLevelDensityParameter (G4VLevelDensityParameter *aLevelDensity) |
|
| G4VEvaporationChannel (const G4String &aName="") |
|
virtual | ~G4VEvaporationChannel () |
|
virtual void | Initialise () |
|
virtual G4double | GetLifeTime (G4Fragment *theNucleus) |
|
virtual G4FragmentVector * | BreakUpFragment (G4Fragment *theNucleus) |
|
virtual G4bool | BreakUpChain (G4FragmentVector *theResult, G4Fragment *theNucleus) |
|
virtual void | SetICM (G4bool) |
|
virtual void | RDMForced (G4bool) |
|
virtual G4double | GetFinalLevelEnergy (G4int Z, G4int A, G4double energy) |
|
virtual G4double | GetUpperLevelEnergy (G4int Z, G4int A) |
|
G4double | GetMaxLevelEnergy (G4int Z, G4int A) |
|
G4double | GetNearestLevelEnergy (G4int Z, G4int A, G4double energy) |
|
void | SetPhotonEvaporation (G4VEvaporationChannel *p) |
|
void | SetOPTxs (G4int opt) |
|
void | UseSICB (G4bool use) |
|
Definition at line 49 of file G4GEMChannel.hh.
◆ G4GEMChannel() [1/2]
Definition at line 46 of file G4GEMChannel.cc.
static G4Pow * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4VEvaporationChannel(const G4String &aName="")
G4VCoulombBarrier * theCoulombBarrierPtr
G4PairingCorrection * pairingCorrection
G4double EmissionProbability
static G4PairingCorrection * GetInstance()
G4VLevelDensityParameter * theLevelDensityPtr
G4GEMProbability * theEvaporationProbabilityPtr
G4double MaximalKineticEnergy
◆ ~G4GEMChannel()
G4GEMChannel::~G4GEMChannel |
( |
| ) |
|
|
virtual |
Definition at line 66 of file G4GEMChannel.cc.
G4VLevelDensityParameter * theLevelDensityPtr
◆ G4GEMChannel() [2/2]
◆ BreakUp()
Implements G4VEvaporationChannel.
Definition at line 144 of file G4GEMChannel.cc.
149 if(frag1) { theResult->push_back(frag1); }
150 theResult->push_back(frag0);
virtual G4Fragment * EmittedFragment(G4Fragment *theNucleus)
std::vector< G4Fragment * > G4FragmentVector
◆ Dump()
void G4GEMChannel::Dump |
( |
| ) |
const |
|
virtual |
◆ EmittedFragment()
Reimplemented from G4VEvaporationChannel.
Definition at line 124 of file G4GEMChannel.cc.
137 ResidualMomentum -= EvaporatedMomentum;
G4ThreeVector IsotropicVector(G4double Magnitude=1.0)
G4double SampleKineticEnergy(const G4Fragment &fragment)
HepLorentzVector & boost(double, double, double)
void SetMomentum(const G4LorentzVector &value)
const G4LorentzVector & GetMomentum() const
Hep3Vector boostVector() const
void SetZandA_asInt(G4int Znew, G4int Anew)
◆ GetEmissionProbability()
Implements G4VEvaporationChannel.
Definition at line 71 of file G4GEMChannel.cc.
73 G4int anA = fragment->GetA_asInt();
74 G4int aZ = fragment->GetZ_asInt();
90 G4double ExEnergy = fragment->GetExcitationEnergy()
94 G4double FragmentMass = fragment->GetGroundStateMass();
95 G4double Etot = FragmentMass + ExEnergy;
G4double GetPairingCorrection(G4int A, G4int Z) const
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4VCoulombBarrier * theCoulombBarrierPtr
G4PairingCorrection * pairingCorrection
G4double EmissionProbability
virtual G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U) const =0
G4GEMProbability * theEvaporationProbabilityPtr
G4double MaximalKineticEnergy
◆ IsotropicVector()
Definition at line 257 of file G4GEMChannel.cc.
262 G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta);
265 Magnitude*std::sin(Phi)*SinTheta,
static const double twopi
◆ operator!=()
◆ operator=()
◆ operator==()
◆ SampleKineticEnergy()
Definition at line 154 of file G4GEMChannel.cc.
180 G4double TCN = 1.0/(std::sqrt(aCN/UxCN) - 1.5/UxCN);
187 - 1.25*
G4Log(UxCN/
MeV) + 2.0*std::sqrt(aCN*UxCN));
188 InitialLevelDensity = (
pi/12.0)*
G4Exp((U-E0CN)/TCN)/TCN;
194 InitialLevelDensity = (
pi/12.0)*
G4Exp(2*x1)/(x*std::sqrt(x1));
210 Rb = (1.12*(Aj + Ad) - 0.86*((Aj+Ad)/(Aj*Ad))+2.85)*
fermi;
216 Rb=1.5*(Aj+Ad)*
fermi;
225 G4double ConstantFactor = gg*GeometricalXS*Alpha*
pi/(InitialLevelDensity*12);
233 for(
G4int i=0; i<100; ++i) {
235 G4double edelta = theEnergy-KineticEnergy-delta0;
236 Probability = ConstantFactor*(KineticEnergy + Beta);
239 G4double T = 1.0/(std::sqrt(a/Ux) - 1.5/Ux);
242 if (theEnergy - KineticEnergy < Ex) {
244 - 1.25*
G4Log(Ux) + 2.0*std::sqrt(a*Ux));
245 Probability *=
G4Exp((theEnergy-KineticEnergy-E0)/T)/T;
249 G4Exp(2*std::sqrt(a*edelta) - 0.25*
G4Log(a*edelta*e2*e2));
254 return KineticEnergy;
G4double CalcAlphaParam(const G4Fragment &) const
G4double GetExcitationEnergy() const
G4double GetPairingCorrection(G4int A, G4int Z) const
G4double CalcBetaParam(const G4Fragment &) const
G4PairingCorrection * pairingCorrection
G4double EmissionProbability
G4double Z13(G4int Z) const
G4double G4Log(G4double x)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
G4VLevelDensityParameter * theLevelDensityPtr
virtual G4double LevelDensityParameter(G4int A, G4int Z, G4double U) const =0
G4GEMProbability * theEvaporationProbabilityPtr
G4double MaximalKineticEnergy
G4double GetSpin(void) const
static const double fermi
◆ SetLevelDensityParameter()
Definition at line 68 of file G4GEMChannel.hh.
G4VLevelDensityParameter * theLevelDensityPtr
◆ CoulombBarrier
◆ EmissionProbability
G4double G4GEMChannel::EmissionProbability |
|
private |
◆ EvaporatedMass
◆ fG4pow
G4Pow* G4GEMChannel::fG4pow |
|
private |
◆ MaximalKineticEnergy
G4double G4GEMChannel::MaximalKineticEnergy |
|
private |
◆ MyOwnLevelDensity
G4bool G4GEMChannel::MyOwnLevelDensity |
|
private |
◆ pairingCorrection
◆ ResidualA
G4int G4GEMChannel::ResidualA |
|
private |
◆ ResidualMass
◆ ResidualZ
G4int G4GEMChannel::ResidualZ |
|
private |
◆ theCoulombBarrierPtr
◆ theEvaporationProbabilityPtr
◆ theLevelDensityPtr
The documentation for this class was generated from the following files: