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

#include <G4GEMChannel.hh>

Inheritance diagram for G4GEMChannel:
G4VEvaporationChannel 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 G4TritonGEMChannel

Public Member Functions

virtual G4bool BreakUpChain (G4FragmentVector *theResult, G4Fragment *theNucleus)
 
G4FragmentVectorBreakUpFragment (G4Fragment *theNucleus)
 
virtual void Dump () const
 
virtual G4FragmentEmittedFragment (G4Fragment *theNucleus)
 
 G4GEMChannel (G4int theA, G4int theZ, const G4String &aName, G4GEMProbability *aEmissionStrategy)
 
virtual G4double GetEmissionProbability (G4Fragment *theNucleus)
 
virtual G4double GetLifeTime (G4Fragment *theNucleus)
 
virtual void Initialise ()
 
virtual void RDMForced (G4bool)
 
virtual void SetICM (G4bool)
 
void SetLevelDensityParameter (G4VLevelDensityParameter *aLevelDensity)
 
void SetOPTxs (G4int opt)
 
void UseSICB (G4bool use)
 
virtual ~G4GEMChannel ()
 

Protected Attributes

G4int OPTxs
 
G4bool useSICB
 

Private Member Functions

 G4GEMChannel (const G4GEMChannel &right)=delete
 
G4bool operator!= (const G4GEMChannel &right) const =delete
 
const G4GEMChanneloperator= (const G4GEMChannel &right)=delete
 
G4bool operator== (const G4GEMChannel &right) const =delete
 
G4double SampleKineticEnergy (const G4Fragment &fragment)
 

Private Attributes

G4int A
 
G4double CoulombBarrier
 
G4double EmissionProbability
 
G4double EvaporatedMass
 
G4PowfG4pow
 
G4NuclearLevelDatafNucData
 
G4double MaximalKineticEnergy
 
G4bool MyOwnLevelDensity
 
G4int ResidualA
 
G4double ResidualMass
 
G4int ResidualZ
 
G4int secID
 
G4VCoulombBarriertheCoulombBarrierPtr
 
G4GEMProbabilitytheEvaporationProbabilityPtr
 
G4VLevelDensityParametertheLevelDensityPtr
 
G4int Z
 

Detailed Description

Definition at line 47 of file G4GEMChannel.hh.

Constructor & Destructor Documentation

◆ G4GEMChannel() [1/2]

G4GEMChannel::G4GEMChannel ( G4int  theA,
G4int  theZ,
const G4String aName,
G4GEMProbability aEmissionStrategy 
)
explicit

Definition at line 49 of file G4GEMChannel.cc.

50 :
52 A(theA),
53 Z(theZ),
56 theEvaporationProbabilityPtr(aEmissionStrategy),
57 secID(-1)
58{
62 MyOwnLevelDensity = true;
66 ResidualZ = ResidualA = 0;
68 secID = G4PhysicsModelCatalog::GetModelID("model_G4GEMChannel");
69}
G4double MaximalKineticEnergy
G4bool MyOwnLevelDensity
G4VCoulombBarrier * theCoulombBarrierPtr
G4double CoulombBarrier
Definition: G4GEMChannel.hh:98
G4VLevelDensityParameter * theLevelDensityPtr
G4double EvaporatedMass
Definition: G4GEMChannel.hh:96
G4GEMProbability * theEvaporationProbabilityPtr
G4NuclearLevelData * fNucData
G4double EmissionProbability
Definition: G4GEMChannel.hh:99
G4Pow * fG4pow
G4double ResidualMass
Definition: G4GEMChannel.hh:97
void SetCoulomBarrier(const G4VCoulombBarrier *aCoulombBarrierStrategy)
static G4NuclearLevelData * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
static G4int GetModelID(const G4int modelIndex)
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4VEvaporationChannel(const G4String &aName="")
static constexpr double GeV

References A, CoulombBarrier, EvaporatedMass, fG4pow, fNucData, G4Pow::GetInstance(), G4NuclearLevelData::GetInstance(), G4PhysicsModelCatalog::GetModelID(), G4NucleiProperties::GetNuclearMass(), MyOwnLevelDensity, ResidualA, ResidualMass, ResidualZ, secID, G4GEMProbability::SetCoulomBarrier(), theCoulombBarrierPtr, theEvaporationProbabilityPtr, theLevelDensityPtr, and Z.

◆ ~G4GEMChannel()

G4GEMChannel::~G4GEMChannel ( )
virtual

Definition at line 71 of file G4GEMChannel.cc.

72{
75}

References MyOwnLevelDensity, theCoulombBarrierPtr, and theLevelDensityPtr.

◆ G4GEMChannel() [2/2]

G4GEMChannel::G4GEMChannel ( const G4GEMChannel right)
privatedelete

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 G4GEMChannel::Dump ( ) const
virtual

Reimplemented from G4VEvaporationChannel.

Definition at line 254 of file G4GEMChannel.cc.

References G4GEMProbability::Dump(), and theEvaporationProbabilityPtr.

◆ EmittedFragment()

G4Fragment * G4GEMChannel::EmittedFragment ( G4Fragment theNucleus)
virtual

Reimplemented from G4VEvaporationChannel.

Definition at line 130 of file G4GEMChannel.cc.

131{
132 G4Fragment* evFragment = 0;
133 G4double evEnergy = SampleKineticEnergy(*theNucleus) + EvaporatedMass;
134
135 G4ThreeVector momentum = G4RandomDirection()*
136 std::sqrt((evEnergy - EvaporatedMass)*(evEnergy + EvaporatedMass));
137
138 G4LorentzVector EvaporatedMomentum(momentum, evEnergy);
139 G4LorentzVector ResidualMomentum = theNucleus->GetMomentum();
140 EvaporatedMomentum.boost(ResidualMomentum.boostVector());
141
142 evFragment = new G4Fragment(A, Z, EvaporatedMomentum);
143 if ( evFragment != nullptr ) { evFragment->SetCreatorModelID(secID); }
144 ResidualMomentum -= EvaporatedMomentum;
145 theNucleus->SetZandA_asInt(ResidualZ, ResidualA);
146 theNucleus->SetMomentum(ResidualMomentum);
147 theNucleus->SetCreatorModelID(secID);
148
149 return evFragment;
150}
G4ThreeVector G4RandomDirection()
double G4double
Definition: G4Types.hh:83
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
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
G4double SampleKineticEnergy(const G4Fragment &fragment)

References A, CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), EvaporatedMass, G4RandomDirection(), G4Fragment::GetMomentum(), ResidualA, ResidualZ, SampleKineticEnergy(), secID, G4Fragment::SetCreatorModelID(), G4Fragment::SetMomentum(), G4Fragment::SetZandA_asInt(), and Z.

◆ GetEmissionProbability()

G4double G4GEMChannel::GetEmissionProbability ( G4Fragment theNucleus)
virtual

Implements G4VEvaporationChannel.

Definition at line 77 of file G4GEMChannel.cc.

78{
79 G4int anA = fragment->GetA_asInt();
80 G4int aZ = fragment->GetZ_asInt();
81 ResidualA = anA - A;
82 ResidualZ = aZ - Z;
83 /*
84 G4cout << "G4GEMChannel: Z= " << Z << " A= " << A
85 << " FragmentZ= " << aZ << " FragmentA= " << anA
86 << " Zres= " << ResidualZ << " Ares= " << ResidualA
87 << G4endl;
88 */
89 // We only take into account channels which are physically allowed
91
92 // Only channels which are physically allowed are taken into account
93 if (ResidualA >= ResidualZ && ResidualZ >= 0 && ResidualA >= A) {
94
95 //Effective excitation energy
96 G4double ExEnergy = fragment->GetExcitationEnergy()
98 if(ExEnergy > 0.0) {
100 G4double FragmentMass = fragment->GetGroundStateMass();
101 G4double Etot = FragmentMass + ExEnergy;
102 // Coulomb Barrier calculation
105 /*
106 G4cout << "Eexc(MeV)= " << ExEnergy/MeV
107 << " CoulBarrier(MeV)= " << CoulombBarrier/MeV << G4endl;
108 */
110
111 // Maximal Kinetic Energy
113 + EvaporatedMass*EvaporatedMass)/(2.0*Etot)
115
116 //G4cout << "CBarrier(MeV)= " << CoulombBarrier/MeV << G4endl;
117
118 if (MaximalKineticEnergy > 0.0) {
119 // Total emission probability for this channel
122 }
123 }
124 }
125 }
126 //G4cout << "Prob= " << EmissionProbability << G4endl;
127 return EmissionProbability;
128}
int G4int
Definition: G4Types.hh:85
G4PairingCorrection * GetPairingCorrection()
virtual G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U) const =0

References A, CoulombBarrier, EmissionProbability, EvaporatedMass, fNucData, G4Fragment::GetA_asInt(), G4VCoulombBarrier::GetCoulombBarrier(), G4Fragment::GetExcitationEnergy(), G4Fragment::GetGroundStateMass(), G4NucleiProperties::GetNuclearMass(), G4NuclearLevelData::GetPairingCorrection(), G4Fragment::GetZ_asInt(), MaximalKineticEnergy, ResidualA, ResidualMass, ResidualZ, theCoulombBarrierPtr, theEvaporationProbabilityPtr, and Z.

◆ GetLifeTime()

G4double G4VEvaporationChannel::GetLifeTime ( G4Fragment theNucleus)
virtualinherited

Definition at line 50 of file G4VEvaporationChannel.cc.

51{
52 return 0.0;
53}

◆ Initialise()

void G4VEvaporationChannel::Initialise ( )
virtualinherited

◆ operator!=()

G4bool G4GEMChannel::operator!= ( const G4GEMChannel right) const
privatedelete

◆ operator=()

const G4GEMChannel & G4GEMChannel::operator= ( const G4GEMChannel right)
privatedelete

◆ operator==()

G4bool G4GEMChannel::operator== ( const G4GEMChannel right) const
privatedelete

◆ RDMForced()

void G4VEvaporationChannel::RDMForced ( G4bool  )
virtualinherited

Reimplemented in G4PhotonEvaporation.

Definition at line 58 of file G4VEvaporationChannel.cc.

59{}

◆ SampleKineticEnergy()

G4double G4GEMChannel::SampleKineticEnergy ( const G4Fragment fragment)
private

Definition at line 152 of file G4GEMChannel.cc.

154{
155 G4double U = fragment.GetExcitationEnergy();
156
159
160 // ***RESIDUAL***
161 //JMQ (September 2009) the following quantities refer to the RESIDUAL:
163 G4double Ux = (2.5 + 150.0/ResidualA)*MeV;
164 G4double Ex = Ux + delta0;
165 G4double InitialLevelDensity;
166 // ***end RESIDUAL ***
167
168 // ***PARENT***
169 //JMQ (September 2009) the following quantities refer to the PARENT:
170
171 G4double deltaCN = fNucData->GetPairingCorrection(fragment.GetZ_asInt(),
172 fragment.GetA_asInt());
174 fragment.GetZ_asInt(),
175 U-deltaCN);
176 G4double UxCN = (2.5 + 150.0/G4double(fragment.GetA_asInt()))*MeV;
177 G4double ExCN = UxCN + deltaCN;
178 G4double TCN = 1.0/(std::sqrt(aCN/UxCN) - 1.5/UxCN);
179 // ***end PARENT***
180
181 //JMQ quantities calculated for CN in InitialLevelDensity
182 if ( U < ExCN )
183 {
184 G4double E0CN = ExCN - TCN*(G4Log(TCN/MeV) - G4Log(aCN*MeV)/4.0
185 - 1.25*G4Log(UxCN/MeV) + 2.0*std::sqrt(aCN*UxCN));
186 InitialLevelDensity = (pi/12.0)*G4Exp((U-E0CN)/TCN)/TCN;
187 }
188 else
189 {
190 G4double x = U-deltaCN;
191 G4double x1 = std::sqrt(aCN*x);
192 InitialLevelDensity = (pi/12.0)*G4Exp(2*x1)/(x*std::sqrt(x1));
193 }
194
196 //JMQ BIG BUG fixed: hbarc instead of hbar_Planck !!!!
197 // it was also fixed in total probability
198 G4double gg = (2.0*Spin+1.0)*EvaporatedMass/(pi2* hbarc*hbarc);
199 //
200 //JMQ fix on Rb and geometrical cross sections according to Furihata's paper
201 // (JAERI-Data/Code 2001-105, p6)
202 G4double Rb = 0.0;
203 if (A > 4)
204 {
206 G4double Aj = fG4pow->Z13(A);
207 Rb = (1.12*(Aj + Ad) - 0.86*((Aj+Ad)/(Aj*Ad))+2.85)*fermi;
208 }
209 else if (A>1)
210 {
212 G4double Aj = fG4pow->Z13(A);
213 Rb=1.5*(Aj+Ad)*fermi;
214 }
215 else
216 {
218 Rb = 1.5*Ad*fermi;
219 }
220 G4double GeometricalXS = pi*Rb*Rb;
221
222 G4double ConstantFactor = gg*GeometricalXS*Alpha*pi/(InitialLevelDensity*12);
223 //JMQ : this is the assimptotic maximal kinetic energy of the ejectile
224 // (obtained by energy-momentum conservation).
225 // In general smaller than U-theSeparationEnergy
228 G4double Probability;
229
230 for(G4int i=0; i<100; ++i) {
232 G4double edelta = theEnergy-KineticEnergy-delta0;
233 Probability = ConstantFactor*(KineticEnergy + Beta);
234 G4double a =
236 G4double T = 1.0/(std::sqrt(a/Ux) - 1.5/Ux);
237 //JMQ fix in units
238
239 if (theEnergy - KineticEnergy < Ex) {
240 G4double E0 = Ex - T*(G4Log(T) - G4Log(a)*0.25
241 - 1.25*G4Log(Ux) + 2.0*std::sqrt(a*Ux));
242 Probability *= G4Exp((theEnergy-KineticEnergy-E0)/T)/T;
243 } else {
244 G4double e2 = edelta*edelta;
245 Probability *=
246 G4Exp(2*std::sqrt(a*edelta) - 0.25*G4Log(a*edelta*e2*e2));
247 }
248 if(EmissionProbability*G4UniformRand() <= Probability) { break; }
249 }
250
251 return KineticEnergy;
252}
static const G4double e2[44]
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static constexpr double pi2
Definition: G4SIunits.hh:58
static constexpr double fermi
Definition: G4SIunits.hh:83
static constexpr double MeV
Definition: G4SIunits.hh:200
static constexpr double pi
Definition: G4SIunits.hh:55
#define G4UniformRand()
Definition: Randomize.hh:52
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:299
G4int GetZ_asInt() const
Definition: G4Fragment.hh:276
G4int GetA_asInt() const
Definition: G4Fragment.hh:271
G4double CalcAlphaParam(const G4Fragment &) const
G4double CalcBetaParam(const G4Fragment &) const
G4double GetSpin(void) const
G4double Z13(G4int Z) const
Definition: G4Pow.hh:123
virtual G4double LevelDensityParameter(G4int A, G4int Z, G4double U) const =0
float hbarc
Definition: hepunit.py:264

References A, Alpha, G4GEMProbability::CalcAlphaParam(), G4GEMProbability::CalcBetaParam(), CoulombBarrier, e2, EmissionProbability, EvaporatedMass, fermi, fG4pow, fNucData, G4Exp(), G4Log(), G4UniformRand, G4Fragment::GetA_asInt(), G4Fragment::GetExcitationEnergy(), G4NuclearLevelData::GetPairingCorrection(), G4GEMProbability::GetSpin(), G4Fragment::GetZ_asInt(), source.hepunit::hbarc, G4VLevelDensityParameter::LevelDensityParameter(), MaximalKineticEnergy, MeV, pi, pi2, ResidualA, ResidualZ, theEvaporationProbabilityPtr, theLevelDensityPtr, and G4Pow::Z13().

Referenced by EmittedFragment().

◆ SetICM()

void G4VEvaporationChannel::SetICM ( G4bool  )
virtualinherited

Reimplemented in G4PhotonEvaporation.

Definition at line 55 of file G4VEvaporationChannel.cc.

56{}

Referenced by G4NeutronRadCapture::InitialiseModel().

◆ SetLevelDensityParameter()

void G4GEMChannel::SetLevelDensityParameter ( G4VLevelDensityParameter aLevelDensity)
inline

Definition at line 63 of file G4GEMChannel.hh.

64 {
66 theLevelDensityPtr = aLevelDensity;
67 MyOwnLevelDensity = false;
68 }

References MyOwnLevelDensity, and theLevelDensityPtr.

◆ 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

◆ A

G4int G4GEMChannel::A
private

◆ CoulombBarrier

G4double G4GEMChannel::CoulombBarrier
private

Definition at line 98 of file G4GEMChannel.hh.

Referenced by G4GEMChannel(), GetEmissionProbability(), and SampleKineticEnergy().

◆ EmissionProbability

G4double G4GEMChannel::EmissionProbability
private

Definition at line 99 of file G4GEMChannel.hh.

Referenced by GetEmissionProbability(), and SampleKineticEnergy().

◆ EvaporatedMass

G4double G4GEMChannel::EvaporatedMass
private

◆ fG4pow

G4Pow* G4GEMChannel::fG4pow
private

Definition at line 104 of file G4GEMChannel.hh.

Referenced by G4GEMChannel(), and SampleKineticEnergy().

◆ fNucData

G4NuclearLevelData* G4GEMChannel::fNucData
private

Definition at line 116 of file G4GEMChannel.hh.

Referenced by G4GEMChannel(), GetEmissionProbability(), and SampleKineticEnergy().

◆ MaximalKineticEnergy

G4double G4GEMChannel::MaximalKineticEnergy
private

Definition at line 102 of file G4GEMChannel.hh.

Referenced by GetEmissionProbability(), and SampleKineticEnergy().

◆ MyOwnLevelDensity

G4bool G4GEMChannel::MyOwnLevelDensity
private

Definition at line 110 of file G4GEMChannel.hh.

Referenced by G4GEMChannel(), SetLevelDensityParameter(), and ~G4GEMChannel().

◆ OPTxs

G4int G4VEvaporationChannel::OPTxs
protectedinherited

◆ ResidualA

G4int G4GEMChannel::ResidualA
private

◆ ResidualMass

G4double G4GEMChannel::ResidualMass
private

Definition at line 97 of file G4GEMChannel.hh.

Referenced by G4GEMChannel(), and GetEmissionProbability().

◆ ResidualZ

G4int G4GEMChannel::ResidualZ
private

◆ secID

G4int G4GEMChannel::secID
private

Definition at line 119 of file G4GEMChannel.hh.

Referenced by EmittedFragment(), and G4GEMChannel().

◆ theCoulombBarrierPtr

G4VCoulombBarrier* G4GEMChannel::theCoulombBarrierPtr
private

Definition at line 114 of file G4GEMChannel.hh.

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

◆ theEvaporationProbabilityPtr

G4GEMProbability* G4GEMChannel::theEvaporationProbabilityPtr
private

Definition at line 107 of file G4GEMChannel.hh.

Referenced by Dump(), G4GEMChannel(), GetEmissionProbability(), and SampleKineticEnergy().

◆ theLevelDensityPtr

G4VLevelDensityParameter* G4GEMChannel::theLevelDensityPtr
private

◆ useSICB

G4bool G4VEvaporationChannel::useSICB
protectedinherited

Definition at line 94 of file G4VEvaporationChannel.hh.

◆ Z

G4int G4GEMChannel::Z
private

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