#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 (G4int theA, G4int theZ, const G4String &aName, G4GEMProbability *aEmissionStrategy) |
|
virtual | ~G4GEMChannel () |
|
virtual G4double | GetEmissionProbability (G4Fragment *theNucleus) |
|
virtual G4Fragment * | EmittedFragment (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 G4bool | BreakUpChain (G4FragmentVector *theResult, G4Fragment *theNucleus) |
|
G4FragmentVector * | BreakUpFragment (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 47 of file G4GEMChannel.hh.
Definition at line 48 of file G4GEMChannel.cc.
53 theEvaporationProbabilityPtr(aEmissionStrategy),
54 EmissionProbability(0.0),
60 MyOwnLevelDensity =
true;
62 ResidualMass = CoulombBarrier = 0.0;
64 ResidualZ = ResidualA = 0;
static G4Pow * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4VEvaporationChannel(const G4String &aName="")
void SetCoulomBarrier(const G4VCoulombBarrier *aCoulombBarrierStrategy)
double A(double temperature)
static G4PairingCorrection * GetInstance()
static constexpr double GeV
G4GEMChannel::~G4GEMChannel |
( |
| ) |
|
|
virtual |
Definition at line 68 of file G4GEMChannel.cc.
70 if (MyOwnLevelDensity) {
delete theLevelDensityPtr; }
71 delete theCoulombBarrierPtr;
void G4GEMChannel::Dump |
( |
| ) |
const |
|
virtual |
Reimplemented from G4VEvaporationChannel.
Definition at line 127 of file G4GEMChannel.cc.
130 G4double evEnergy = SampleKineticEnergy(*theNucleus) + EvaporatedMass;
133 (std::sqrt((evEnergy - EvaporatedMass)*(evEnergy + EvaporatedMass))));
140 ResidualMomentum -= EvaporatedMomentum;
Hep3Vector boostVector() const
double A(double temperature)
const G4LorentzVector & GetMomentum() const
HepLorentzVector & boost(double, double, double)
void SetMomentum(const G4LorentzVector &value)
void SetZandA_asInt(G4int Znew, G4int Anew)
Implements G4VEvaporationChannel.
Definition at line 74 of file G4GEMChannel.cc.
76 G4int anA = fragment->GetA_asInt();
77 G4int aZ = fragment->GetZ_asInt();
87 EmissionProbability = 0.0;
90 if (ResidualA >= ResidualZ && ResidualZ > 0 && ResidualA >= A) {
93 G4double ExEnergy = fragment->GetExcitationEnergy()
97 G4double FragmentMass = fragment->GetGroundStateMass();
98 G4double Etot = FragmentMass + ExEnergy;
106 if(Etot > ResidualMass + EvaporatedMass + CoulombBarrier) {
109 MaximalKineticEnergy = ((Etot-ResidualMass)*(Etot+ResidualMass)
110 + EvaporatedMass*EvaporatedMass)/(2.0*Etot)
111 - EvaporatedMass - CoulombBarrier;
115 if (MaximalKineticEnergy > 0.0) {
117 EmissionProbability = theEvaporationProbabilityPtr->
118 EmissionProbability(*fragment, MaximalKineticEnergy);
124 return EmissionProbability;
static G4double GetNuclearMass(const G4double A, const G4double Z)
double A(double temperature)
G4double GetPairingCorrection(G4int A, G4int Z) const
virtual G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U) const =0
Definition at line 63 of file G4GEMChannel.hh.
65 if (MyOwnLevelDensity) {
delete theLevelDensityPtr; }
66 theLevelDensityPtr = aLevelDensity;
67 MyOwnLevelDensity =
false;
The documentation for this class was generated from the following files:
- source/geant4.10.03.p03/source/processes/hadronic/models/de_excitation/gem_evaporation/include/G4GEMChannel.hh
- source/geant4.10.03.p03/source/processes/hadronic/models/de_excitation/gem_evaporation/src/G4GEMChannel.cc