Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EvaporationChannel Class Reference

#include <G4EvaporationChannel.hh>

+ Inheritance diagram for G4EvaporationChannel:

Public Member Functions

 G4EvaporationChannel (G4int A, G4int Z, G4EvaporationProbability *)
 
 ~G4EvaporationChannel () override
 
void Initialise () override
 
G4double GetEmissionProbability (G4Fragment *fragment) override
 
G4FragmentEmittedFragment (G4Fragment *theNucleus) override
 
G4double ComputeInverseXSection (G4Fragment *, G4double kinEnergy) override
 
G4double ComputeProbability (G4Fragment *, G4double kinEnergy) override
 
G4int GetZ () const
 
G4int GetA () const
 
G4EvaporationProbabilityGetEvaporationProbability ()
 
 G4EvaporationChannel (const G4EvaporationChannel &right)=delete
 
const G4EvaporationChanneloperator= (const G4EvaporationChannel &right)=delete
 
G4bool operator== (const G4EvaporationChannel &right) const =delete
 
G4bool operator!= (const G4EvaporationChannel &right) const =delete
 
- Public Member Functions inherited from G4VEvaporationChannel
 G4VEvaporationChannel (const G4String &aName="")
 
virtual ~G4VEvaporationChannel ()=default
 
virtual G4double GetLifeTime (G4Fragment *theNucleus)
 
virtual G4bool BreakUpChain (G4FragmentVector *theResult, G4Fragment *theNucleus)
 
G4FragmentVectorBreakUpFragment (G4Fragment *theNucleus)
 
virtual void Dump () const
 
virtual void SetICM (G4bool)
 
virtual void RDMForced (G4bool)
 
void SetOPTxs (G4int)
 
void UseSICB (G4bool)
 
 G4VEvaporationChannel (const G4VEvaporationChannel &right)=delete
 
const G4VEvaporationChanneloperator= (const G4VEvaporationChannel &right)=delete
 
G4bool operator== (const G4VEvaporationChannel &right) const =delete
 
G4bool operator!= (const G4VEvaporationChannel &right) const =delete
 

Additional Inherited Members

- Protected Attributes inherited from G4VEvaporationChannel
G4int OPTxs {3}
 
G4bool useSICB {true}
 

Detailed Description

Definition at line 45 of file G4EvaporationChannel.hh.

Constructor & Destructor Documentation

◆ G4EvaporationChannel() [1/2]

G4EvaporationChannel::G4EvaporationChannel ( G4int A,
G4int Z,
G4EvaporationProbability * aprob )
explicit

Definition at line 53 of file G4EvaporationChannel.cc.

54 :
56 theProbability(aprob),
57 theCoulombBarrier(new G4CoulombBarrier(anA, aZ)),
58 theA(anA), theZ(aZ)
59{
60 secID = G4PhysicsModelCatalog::GetModelID("model_G4EvaporationChannel");
61 evapMass = G4NucleiProperties::GetNuclearMass(theA, theZ);
62 evapMass2 = evapMass*evapMass;
63 theLevelData = G4NuclearLevelData::GetInstance();
64}
static G4NuclearLevelData * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
static G4int GetModelID(const G4int modelIndex)
G4VEvaporationChannel(const G4String &aName="")

◆ ~G4EvaporationChannel()

G4EvaporationChannel::~G4EvaporationChannel ( )
override

Definition at line 66 of file G4EvaporationChannel.cc.

67{
68 delete theCoulombBarrier;
69}

◆ G4EvaporationChannel() [2/2]

G4EvaporationChannel::G4EvaporationChannel ( const G4EvaporationChannel & right)
delete

Member Function Documentation

◆ ComputeInverseXSection()

G4double G4EvaporationChannel::ComputeInverseXSection ( G4Fragment * frag,
G4double kinEnergy )
overridevirtual

Reimplemented from G4VEvaporationChannel.

Definition at line 162 of file G4EvaporationChannel.cc.

164{
165 ComputeProbability(frag, kinEnergy);
166 return theProbability->CrossSection(kinEnergy, bCoulomb);
167}
G4double ComputeProbability(G4Fragment *, G4double kinEnergy) override
G4double CrossSection(G4double K, G4double CB)

◆ ComputeProbability()

G4double G4EvaporationChannel::ComputeProbability ( G4Fragment * frag,
G4double kinEnergy )
overridevirtual

Reimplemented from G4VEvaporationChannel.

Definition at line 169 of file G4EvaporationChannel.cc.

171{
172 G4int fragA = frag->GetA_asInt();
173 G4int fragZ = frag->GetZ_asInt();
174 resA = fragA - theA;
175 resZ = fragZ - theZ;
176 bCoulomb = theCoulombBarrier->GetCoulombBarrier(resA, resZ, 0.0);
177 return theProbability->ComputeProbability(kinEnergy, bCoulomb);
178}
int G4int
Definition G4Types.hh:85
G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U) const override
G4double ComputeProbability(G4double K, G4double CB) override
G4int GetZ_asInt() const
G4int GetA_asInt() const

Referenced by ComputeInverseXSection().

◆ EmittedFragment()

G4Fragment * G4EvaporationChannel::EmittedFragment ( G4Fragment * theNucleus)
overridevirtual

Reimplemented from G4VEvaporationChannel.

Definition at line 140 of file G4EvaporationChannel.cc.

141{
142 G4double ekin = ekinmax;
143 // assumed, that TotalProbability(...) was already called
144 // if value iz zero no possiblity to sample final state
145 if(resA > 4 && theProbability->GetProbability() > 0.0) {
146 ekin = theProbability->SampleEnergy();
147 }
148 ekin = std::max(ekin, 0.0);
149 G4LorentzVector lv0 = theNucleus->GetMomentum();
150 G4LorentzVector lv(std::sqrt(ekin*(ekin + 2.0*evapMass))*G4RandomDirection(),
151 ekin + evapMass);
152 lv.boost(lv0.boostVector());
153
154 G4Fragment* evFragment = new G4Fragment(theA, theZ, lv);
155 evFragment->SetCreatorModelID(secID);
156 lv0 -= lv;
157 theNucleus->SetZAandMomentum(lv0, resZ, resA);
158 theNucleus->SetCreatorModelID(secID);
159 return evFragment;
160}
G4ThreeVector G4RandomDirection()
double G4double
Definition G4Types.hh:83
Hep3Vector boostVector() const
const G4LorentzVector & GetMomentum() const
void SetCreatorModelID(G4int value)
void SetZAandMomentum(const G4LorentzVector &, G4int Z, G4int A, G4int nLambdas=0)

◆ GetA()

G4int G4EvaporationChannel::GetA ( ) const
inline

Definition at line 66 of file G4EvaporationChannel.hh.

66{ return theA; };

◆ GetEmissionProbability()

G4double G4EvaporationChannel::GetEmissionProbability ( G4Fragment * fragment)
overridevirtual

Implements G4VEvaporationChannel.

Definition at line 77 of file G4EvaporationChannel.cc.

78{
79 theProbability->ResetProbability();
80 G4int fragA = fragment->GetA_asInt();
81 G4int fragZ = fragment->GetZ_asInt();
82 resA = fragA - theA;
83 resZ = fragZ - theZ;
84
85 // Only channels which are physically allowed are taken into account
86 if(resA < theA || resA < resZ || resZ < 0 || (resA == theA && resZ < theZ)
87 || ((resA > 1) && (resA == resZ || resZ == 0)))
88 { return 0.0; }
89
90 G4double exEnergy = fragment->GetExcitationEnergy();
91 G4double delta0 = theLevelData->GetPairingCorrection(fragZ,fragA);
92 /*
93 G4cout << "G4EvaporationChannel::Initialize Z= "<<theZ<<" A= "<<theA
94 << " FragZ= " << fragZ << " FragA= " << fragA
95 << " exEnergy= " << exEnergy << " d0= " << delta0 << G4endl;
96 */
97 if(exEnergy < delta0) { return 0.0; }
98
99 G4double fragMass = fragment->GetGroundStateMass();
100 mass = fragMass + exEnergy;
101
102 resMass = G4NucleiProperties::GetNuclearMass(resA, resZ);
103 ekinmax = 0.5*((mass-resMass)*(mass+resMass) + evapMass2)/mass - evapMass;
104
105 G4double elim = 0.0;
106 if(theZ > 0) {
107 bCoulomb = theCoulombBarrier->GetCoulombBarrier(resA, resZ, 0.0);
108
109 // for OPTxs >0 penetration under the barrier is taken into account
110 elim = (0 != OPTxs) ? bCoulomb*0.6 : bCoulomb;
111 }
112 /*
113 G4cout << "exEnergy= " << exEnergy << " Ec= " << bCoulomb
114 << " elim=" << elim << " d0= " << delta0
115 << " Free= " << mass - resMass - evapMass
116 << G4endl;
117 */
118 if(mass <= resMass + evapMass + elim) { return 0.0; }
119
120 G4double ekinmin = 0.0;
121 if(elim > 0.0) {
122 G4double resM = mass - evapMass - elim;
123 ekinmin =
124 std::max(0.5*((mass-resM)*(mass+resM) + evapMass2)/mass - evapMass, 0.0);
125 }
126 /*
127 G4cout << "Emin= " <<ekinmin<<" Emax= "<<ekinmax
128 << " mass= " << mass << " resM= " << resMass
129 << " evapM= " << evapMass << G4endl;
130 */
131 if(ekinmax <= ekinmin) { return 0.0; }
132
133 theProbability->SetDecayKinematics(resZ, resA, resMass, mass);
134 G4double prob = theProbability->TotalProbability(*fragment, ekinmin,
135 ekinmax, bCoulomb,
136 exEnergy - delta0);
137 return prob;
138}
virtual G4double TotalProbability(const G4Fragment &fragment, G4double minKinEnergy, G4double maxKinEnergy, G4double CB, G4double exEnergy)
G4double GetGroundStateMass() const
G4double GetExcitationEnergy() const
G4PairingCorrection * GetPairingCorrection()
void SetDecayKinematics(G4int rZ, G4int rA, G4double rmass, G4double fmass)

◆ GetEvaporationProbability()

G4EvaporationProbability * G4EvaporationChannel::GetEvaporationProbability ( )
inline

Definition at line 68 of file G4EvaporationChannel.hh.

69 { return theProbability; }

◆ GetZ()

G4int G4EvaporationChannel::GetZ ( ) const
inline

Definition at line 64 of file G4EvaporationChannel.hh.

64{ return theZ; };

◆ Initialise()

void G4EvaporationChannel::Initialise ( )
overridevirtual

Reimplemented from G4VEvaporationChannel.

Definition at line 71 of file G4EvaporationChannel.cc.

◆ operator!=()

G4bool G4EvaporationChannel::operator!= ( const G4EvaporationChannel & right) const
delete

◆ operator=()

const G4EvaporationChannel & G4EvaporationChannel::operator= ( const G4EvaporationChannel & right)
delete

◆ operator==()

G4bool G4EvaporationChannel::operator== ( const G4EvaporationChannel & right) const
delete

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