Geant4 10.7.0
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
 
- Public Member Functions inherited from G4VEvaporationChannel
 G4VEvaporationChannel (const G4String &aName="")
 
virtual ~G4VEvaporationChannel ()
 
virtual G4double GetEmissionProbability (G4Fragment *theNucleus)=0
 
virtual void Initialise ()
 
virtual G4double GetLifeTime (G4Fragment *theNucleus)
 
virtual G4FragmentEmittedFragment (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 opt)
 
void UseSICB (G4bool use)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEvaporationChannel
G4int OPTxs
 
G4bool useSICB
 

Detailed Description

Definition at line 45 of file G4EvaporationChannel.hh.

Constructor & Destructor Documentation

◆ G4EvaporationChannel()

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

Definition at line 53 of file G4EvaporationChannel.cc.

54 :
56 theA(anA),
57 theZ(aZ),
58 theProbability(aprob),
59 theCoulombBarrier(new G4CoulombBarrier(anA, aZ))
60{
61 resA = resZ = 0;
62 mass = resMass = 0.0;
63 evapMass = G4NucleiProperties::GetNuclearMass(theA, theZ);
64 //G4cout << "G4EvaporationChannel: Z= " << theZ << " A= " << theA
65 // << " M(GeV)= " << evapMass/GeV << G4endl;
66 evapMass2 = evapMass*evapMass;
67 theLevelData = G4NuclearLevelData::GetInstance();
68}
static G4NuclearLevelData * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)

◆ ~G4EvaporationChannel()

G4EvaporationChannel::~G4EvaporationChannel ( )
override

Definition at line 70 of file G4EvaporationChannel.cc.

71{
72 delete theCoulombBarrier;
73}

Member Function Documentation

◆ EmittedFragment()

G4Fragment * G4EvaporationChannel::EmittedFragment ( G4Fragment theNucleus)
overridevirtual

Reimplemented from G4VEvaporationChannel.

Definition at line 153 of file G4EvaporationChannel.cc.

154{
155 G4double ekin;
156 // assumed, that TotalProbability(...) was already called
157 // if value iz zero no possiblity to sample final state
158 if(resA <= 4 || theProbability->GetProbability() == 0.0) {
159 ekin = 0.5*(mass*mass - resMass*resMass + evapMass2)/mass - evapMass;
160 } else {
161 ekin = theProbability->SampleEnergy();
162 }
163 ekin = std::max(ekin, 0.0);
164 G4LorentzVector lv0 = theNucleus->GetMomentum();
165 G4LorentzVector lv(std::sqrt(ekin*(ekin + 2.0*evapMass))*G4RandomDirection(),
166 ekin + evapMass);
167 lv.boost(lv0.boostVector());
168
169 G4Fragment* evFragment = new G4Fragment(theA, theZ, lv);
170 lv0 -= lv;
171 theNucleus->SetZandA_asInt(resZ, resA);
172 theNucleus->SetMomentum(lv0);
173
174 //G4cout << "Residual: Z= " << resZ << " A= " << resA << " Eex= "
175 // << theNucleus->GetExcitationEnergy() << G4endl;
176 return evFragment;
177}
G4ThreeVector G4RandomDirection()
double G4double
Definition: G4Types.hh:83
Hep3Vector boostVector() const
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:299
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:304
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:268

◆ GetEmissionProbability()

G4double G4EvaporationChannel::GetEmissionProbability ( G4Fragment fragment)
overridevirtual

Implements G4VEvaporationChannel.

Definition at line 81 of file G4EvaporationChannel.cc.

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

◆ Initialise()

void G4EvaporationChannel::Initialise ( )
overridevirtual

Reimplemented from G4VEvaporationChannel.

Definition at line 75 of file G4EvaporationChannel.cc.


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