Geant4-11
Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4EvaporationChannel Class Reference

#include <G4EvaporationChannel.hh>

Inheritance diagram for G4EvaporationChannel:
G4VEvaporationChannel G4AlphaEvaporationChannel G4DeuteronEvaporationChannel G4He3EvaporationChannel G4NeutronEvaporationChannel G4ProtonEvaporationChannel G4TritonEvaporationChannel

Public Member Functions

virtual G4bool BreakUpChain (G4FragmentVector *theResult, G4Fragment *theNucleus)
 
G4FragmentVectorBreakUpFragment (G4Fragment *theNucleus)
 
virtual void Dump () const
 
G4FragmentEmittedFragment (G4Fragment *theNucleus) override
 
 G4EvaporationChannel (G4int A, G4int Z, G4EvaporationProbability *)
 
G4double GetEmissionProbability (G4Fragment *fragment) override
 
virtual G4double GetLifeTime (G4Fragment *theNucleus)
 
void Initialise () override
 
virtual void RDMForced (G4bool)
 
virtual void SetICM (G4bool)
 
void SetOPTxs (G4int opt)
 
void UseSICB (G4bool use)
 
 ~G4EvaporationChannel () override
 

Protected Attributes

G4int OPTxs
 
G4bool useSICB
 

Private Member Functions

 G4EvaporationChannel (const G4EvaporationChannel &right)
 
G4bool operator!= (const G4EvaporationChannel &right) const
 
const G4EvaporationChanneloperator= (const G4EvaporationChannel &right)
 
G4bool operator== (const G4EvaporationChannel &right) const
 

Private Attributes

G4double evapMass
 
G4double evapMass2
 
G4double mass
 
G4int resA
 
G4double resMass
 
G4int resZ
 
G4int secID
 
G4int theA
 
G4CoulombBarriertheCoulombBarrier
 
G4NuclearLevelDatatheLevelData
 
G4EvaporationProbabilitytheProbability
 
G4int theZ
 

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 54 of file G4EvaporationChannel.cc.

55 :
57 theA(anA),
58 theZ(aZ),
59 secID(-1),
60 theProbability(aprob),
62{
63 resA = resZ = 0;
64 secID = G4PhysicsModelCatalog::GetModelID("model_G4EvaporationChannel");
65 mass = resMass = 0.0;
67 //G4cout << "G4EvaporationChannel: Z= " << theZ << " A= " << theA
68 // << " M(GeV)= " << evapMass/GeV << G4endl;
71}
G4NuclearLevelData * theLevelData
G4CoulombBarrier * theCoulombBarrier
G4EvaporationProbability * theProbability
static G4NuclearLevelData * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
static G4int GetModelID(const G4int modelIndex)
G4VEvaporationChannel(const G4String &aName="")

References evapMass, evapMass2, G4NuclearLevelData::GetInstance(), G4PhysicsModelCatalog::GetModelID(), G4NucleiProperties::GetNuclearMass(), mass, resA, resMass, resZ, secID, theA, theLevelData, and theZ.

◆ ~G4EvaporationChannel()

G4EvaporationChannel::~G4EvaporationChannel ( )
override

Definition at line 73 of file G4EvaporationChannel.cc.

74{
75 delete theCoulombBarrier;
76}

References theCoulombBarrier.

◆ G4EvaporationChannel() [2/2]

G4EvaporationChannel::G4EvaporationChannel ( const G4EvaporationChannel right)
private

Member Function Documentation

◆ BreakUpChain()

G4bool G4VEvaporationChannel::BreakUpChain ( G4FragmentVector theResult,
G4Fragment theNucleus 
)
virtualinherited

Reimplemented in G4UnstableFragmentBreakUp, and G4PhotonEvaporation.

Definition at line 66 of file G4VEvaporationChannel.cc.

67{
68 return false;
69}

Referenced by G4VEvaporationChannel::BreakUpFragment().

◆ BreakUpFragment()

G4FragmentVector * G4VEvaporationChannel::BreakUpFragment ( G4Fragment theNucleus)
inlineinherited

Definition at line 106 of file G4VEvaporationChannel.hh.

107{
108 G4FragmentVector* results = new G4FragmentVector();
109 BreakUpChain(results, theNucleus);
110 return results;
111}
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:64
virtual G4bool BreakUpChain(G4FragmentVector *theResult, G4Fragment *theNucleus)

References G4VEvaporationChannel::BreakUpChain().

Referenced by G4NeutronRadCapture::ApplyYourself().

◆ Dump()

void G4VEvaporationChannel::Dump ( ) const
virtualinherited

Reimplemented in G4GEMChannel, and G4GEMChannelVI.

Definition at line 71 of file G4VEvaporationChannel.cc.

72{}

◆ EmittedFragment()

G4Fragment * G4EvaporationChannel::EmittedFragment ( G4Fragment theNucleus)
overridevirtual

Reimplemented from G4VEvaporationChannel.

Definition at line 156 of file G4EvaporationChannel.cc.

157{
158 G4double ekin;
159 // assumed, that TotalProbability(...) was already called
160 // if value iz zero no possiblity to sample final state
161 if(resA <= 4 || theProbability->GetProbability() == 0.0) {
162 ekin = 0.5*(mass*mass - resMass*resMass + evapMass2)/mass - evapMass;
163 } else {
165 }
166 ekin = std::max(ekin, 0.0);
167 G4LorentzVector lv0 = theNucleus->GetMomentum();
168 G4LorentzVector lv(std::sqrt(ekin*(ekin + 2.0*evapMass))*G4RandomDirection(),
169 ekin + evapMass);
170 lv.boost(lv0.boostVector());
171
172 G4Fragment* evFragment = new G4Fragment(theA, theZ, lv);
173 if(evFragment != nullptr) { evFragment->SetCreatorModelID(secID); }
174 lv0 -= lv;
175 theNucleus->SetZandA_asInt(resZ, resA);
176 theNucleus->SetMomentum(lv0);
177 theNucleus->SetCreatorModelID(secID);
178
179 //G4cout << "Residual: Z= " << resZ << " A= " << resA << " Eex= "
180 // << theNucleus->GetExcitationEnergy() << G4endl;
181 return evFragment;
182}
G4ThreeVector G4RandomDirection()
double G4double
Definition: G4Types.hh:83
Hep3Vector boostVector() const
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:323
void SetCreatorModelID(G4int value)
Definition: G4Fragment.hh:428
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:328
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:281
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), evapMass, evapMass2, G4RandomDirection(), G4Fragment::GetMomentum(), mass, G4INCL::Math::max(), resA, resMass, resZ, G4VEmissionProbability::SampleEnergy(), secID, G4Fragment::SetCreatorModelID(), G4Fragment::SetMomentum(), G4Fragment::SetZandA_asInt(), theA, theProbability, and theZ.

◆ GetEmissionProbability()

G4double G4EvaporationChannel::GetEmissionProbability ( G4Fragment fragment)
overridevirtual

Implements G4VEvaporationChannel.

Definition at line 84 of file G4EvaporationChannel.cc.

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

References ekinmin, evapMass, evapMass2, G4Fragment::GetA_asInt(), G4CoulombBarrier::GetCoulombBarrier(), G4Fragment::GetExcitationEnergy(), G4Fragment::GetGroundStateMass(), G4NucleiProperties::GetNuclearMass(), G4NuclearLevelData::GetPairingCorrection(), G4Fragment::GetZ_asInt(), mass, G4INCL::Math::max(), CLHEP::MeV, G4VEvaporationChannel::OPTxs, resA, G4VEmissionProbability::ResetProbability(), resMass, resZ, G4VEmissionProbability::SetDecayKinematics(), theA, theCoulombBarrier, theLevelData, theProbability, theZ, and G4EvaporationProbability::TotalProbability().

◆ GetLifeTime()

G4double G4VEvaporationChannel::GetLifeTime ( G4Fragment theNucleus)
virtualinherited

Definition at line 50 of file G4VEvaporationChannel.cc.

51{
52 return 0.0;
53}

◆ Initialise()

void G4EvaporationChannel::Initialise ( )
overridevirtual

◆ operator!=()

G4bool G4EvaporationChannel::operator!= ( const G4EvaporationChannel right) const
private

◆ operator=()

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

◆ operator==()

G4bool G4EvaporationChannel::operator== ( const G4EvaporationChannel right) const
private

◆ RDMForced()

void G4VEvaporationChannel::RDMForced ( G4bool  )
virtualinherited

Reimplemented in G4PhotonEvaporation.

Definition at line 58 of file G4VEvaporationChannel.cc.

59{}

◆ SetICM()

void G4VEvaporationChannel::SetICM ( G4bool  )
virtualinherited

Reimplemented in G4PhotonEvaporation.

Definition at line 55 of file G4VEvaporationChannel.cc.

56{}

Referenced by G4NeutronRadCapture::InitialiseModel().

◆ SetOPTxs()

void G4VEvaporationChannel::SetOPTxs ( G4int  opt)
inlineinherited

Definition at line 113 of file G4VEvaporationChannel.hh.

114{}

◆ UseSICB()

void G4VEvaporationChannel::UseSICB ( G4bool  use)
inlineinherited

Definition at line 116 of file G4VEvaporationChannel.hh.

117{}

Field Documentation

◆ evapMass

G4double G4EvaporationChannel::evapMass
private

◆ evapMass2

G4double G4EvaporationChannel::evapMass2
private

◆ mass

G4double G4EvaporationChannel::mass
private

◆ OPTxs

G4int G4VEvaporationChannel::OPTxs
protectedinherited

Definition at line 93 of file G4VEvaporationChannel.hh.

Referenced by GetEmissionProbability().

◆ resA

G4int G4EvaporationChannel::resA
private

◆ resMass

G4double G4EvaporationChannel::resMass
private

◆ resZ

G4int G4EvaporationChannel::resZ
private

◆ secID

G4int G4EvaporationChannel::secID
private

Definition at line 74 of file G4EvaporationChannel.hh.

Referenced by EmittedFragment(), and G4EvaporationChannel().

◆ theA

G4int G4EvaporationChannel::theA
private

◆ theCoulombBarrier

G4CoulombBarrier* G4EvaporationChannel::theCoulombBarrier
private

Definition at line 86 of file G4EvaporationChannel.hh.

Referenced by GetEmissionProbability(), and ~G4EvaporationChannel().

◆ theLevelData

G4NuclearLevelData* G4EvaporationChannel::theLevelData
private

Definition at line 88 of file G4EvaporationChannel.hh.

Referenced by G4EvaporationChannel(), and GetEmissionProbability().

◆ theProbability

G4EvaporationProbability* G4EvaporationChannel::theProbability
private

Definition at line 83 of file G4EvaporationChannel.hh.

Referenced by EmittedFragment(), GetEmissionProbability(), and Initialise().

◆ theZ

G4int G4EvaporationChannel::theZ
private

◆ useSICB

G4bool G4VEvaporationChannel::useSICB
protectedinherited

Definition at line 94 of file G4VEvaporationChannel.hh.


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