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

#include <G4VPreCompoundFragment.hh>

+ Inheritance diagram for G4VPreCompoundFragment:

Public Member Functions

 G4VPreCompoundFragment (const G4ParticleDefinition *, G4VCoulombBarrier *)
 
virtual ~G4VPreCompoundFragment ()
 
void Initialize (const G4Fragment &aFragment)
 
virtual G4double CalcEmissionProbability (const G4Fragment &)=0
 
virtual G4double SampleKineticEnergy (const G4Fragment &)=0
 
G4bool IsItPossible (const G4Fragment &aFragment) const
 
G4ReactionProductGetReactionProduct () const
 
G4int GetA () const
 
G4int GetZ () const
 
G4int GetRestA () const
 
G4int GetRestZ () const
 
G4double GetBindingEnergy () const
 
G4double GetEnergyThreshold () const
 
G4double GetEmissionProbability () const
 
G4double GetNuclearMass () const
 
G4double GetRestNuclearMass () const
 
const G4LorentzVectorGetMomentum () const
 
void SetMomentum (const G4LorentzVector &lv)
 
void SetOPTxs (G4int opt)
 
void UseSICB (G4bool use)
 
 G4VPreCompoundFragment (const G4VPreCompoundFragment &right)=delete
 
const G4VPreCompoundFragmentoperator= (const G4VPreCompoundFragment &right)=delete
 
G4bool operator== (const G4VPreCompoundFragment &right) const =delete
 
G4bool operator!= (const G4VPreCompoundFragment &right) const =delete
 

Protected Member Functions

virtual G4double GetAlpha () const =0
 
virtual G4double GetBeta () const
 

Protected Attributes

G4NuclearLevelDatafNucData
 
G4DeexPrecoParameterstheParameters
 
G4Powg4calc
 
G4int theA
 
G4int theZ
 
G4int theResA {0}
 
G4int theResZ {0}
 
G4int theFragA {0}
 
G4int theFragZ {0}
 
G4double theResA13 {0.0}
 
G4double theBindingEnergy {0.0}
 
G4double theMinKinEnergy {0.0}
 
G4double theMaxKinEnergy {0.0}
 
G4double theResMass {0.0}
 
G4double theReducedMass {0.0}
 
G4double theMass
 
G4double theEmissionProbability {0.0}
 
G4double theCoulombBarrier {0.0}
 
G4int OPTxs {3}
 
G4bool useSICB {true}
 

Friends

std::ostream & operator<< (std::ostream &out, const G4VPreCompoundFragment *theFragment)
 
std::ostream & operator<< (std::ostream &out, const G4VPreCompoundFragment &theFragment)
 

Detailed Description

Definition at line 55 of file G4VPreCompoundFragment.hh.

Constructor & Destructor Documentation

◆ G4VPreCompoundFragment() [1/2]

G4VPreCompoundFragment::G4VPreCompoundFragment ( const G4ParticleDefinition * part,
G4VCoulombBarrier * aCoulombBarrier )
explicit

Definition at line 40 of file G4VPreCompoundFragment.cc.

42 : theA(part->GetBaryonNumber()),
43 theZ(G4lrint(part->GetPDGCharge()/CLHEP::eplus)),
44 particle(part),
45 theCoulombBarrierPtr(aCoulombBarrier)
46{
47 theMass = particle->GetPDGMass();
51}
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4DeexPrecoParameters * theParameters
int G4lrint(double ad)
Definition templates.hh:134

◆ ~G4VPreCompoundFragment()

G4VPreCompoundFragment::~G4VPreCompoundFragment ( )
virtual

Definition at line 53 of file G4VPreCompoundFragment.cc.

54{
55 delete theCoulombBarrierPtr;
56}

◆ G4VPreCompoundFragment() [2/2]

G4VPreCompoundFragment::G4VPreCompoundFragment ( const G4VPreCompoundFragment & right)
delete

Member Function Documentation

◆ CalcEmissionProbability()

virtual G4double G4VPreCompoundFragment::CalcEmissionProbability ( const G4Fragment & )
pure virtual

Implemented in G4HETCFragment, and G4PreCompoundFragment.

◆ GetA()

G4int G4VPreCompoundFragment::GetA ( ) const
inline

◆ GetAlpha()

◆ GetBeta()

virtual G4double G4VPreCompoundFragment::GetBeta ( ) const
inlineprotectedvirtual

Reimplemented in G4HETCNeutron, and G4PreCompoundNeutron.

Definition at line 130 of file G4VPreCompoundFragment.hh.

◆ GetBindingEnergy()

G4double G4VPreCompoundFragment::GetBindingEnergy ( ) const
inline

Definition at line 98 of file G4VPreCompoundFragment.hh.

◆ GetEmissionProbability()

G4double G4VPreCompoundFragment::GetEmissionProbability ( ) const
inline

Definition at line 105 of file G4VPreCompoundFragment.hh.

◆ GetEnergyThreshold()

G4double G4VPreCompoundFragment::GetEnergyThreshold ( ) const
inline

◆ GetMomentum()

const G4LorentzVector & G4VPreCompoundFragment::GetMomentum ( ) const
inline

Definition at line 111 of file G4VPreCompoundFragment.hh.

111{ return theMomentum; }

Referenced by GetReactionProduct().

◆ GetNuclearMass()

G4double G4VPreCompoundFragment::GetNuclearMass ( ) const
inline

Definition at line 107 of file G4VPreCompoundFragment.hh.

107{ return theMass; }

Referenced by G4PreCompoundEmission::PerformEmission().

◆ GetReactionProduct()

G4ReactionProduct * G4VPreCompoundFragment::GetReactionProduct ( ) const
inline

Definition at line 174 of file G4VPreCompoundFragment.hh.

175{
176 G4ReactionProduct* theReactionProduct = new G4ReactionProduct(particle);
177 theReactionProduct->SetMomentum(GetMomentum().vect());
178 theReactionProduct->SetTotalEnergy(GetMomentum().e());
179 return theReactionProduct;
180}
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
const G4LorentzVector & GetMomentum() const

Referenced by G4PreCompoundEmission::PerformEmission().

◆ GetRestA()

G4int G4VPreCompoundFragment::GetRestA ( ) const
inline

◆ GetRestNuclearMass()

G4double G4VPreCompoundFragment::GetRestNuclearMass ( ) const
inline

Definition at line 109 of file G4VPreCompoundFragment.hh.

◆ GetRestZ()

G4int G4VPreCompoundFragment::GetRestZ ( ) const
inline

◆ GetZ()

G4int G4VPreCompoundFragment::GetZ ( ) const
inline

Definition at line 92 of file G4VPreCompoundFragment.hh.

92{ return theZ; }

Referenced by G4PreCompoundEmission::PerformEmission().

◆ Initialize()

void G4VPreCompoundFragment::Initialize ( const G4Fragment & aFragment)

Definition at line 76 of file G4VPreCompoundFragment.cc.

77{
78 theFragA = aFragment.GetA_asInt();
79 theFragZ = aFragment.GetZ_asInt();
82
84 if ((theResA < theResZ) || (theResA < theA) || (theResZ < theZ)) {
85 return;
86 }
87
89
90 if (0 < theZ) {
91 theCoulombBarrier = theCoulombBarrierPtr->
92 GetCoulombBarrier(theResA, theResZ, aFragment.GetExcitationEnergy());
93 }
95
96 // Calculate masses
99
100 // Compute Binding Energies for fragments
101 // needed to separate a fragment from the nucleus
103
104 // Compute Maximal Kinetic Energy which can be carried by fragments
105 // after separation - the true assimptotic value
106 G4double Ecm = aFragment.GetMomentum().m();
107 G4double twoEcm = Ecm + Ecm;
108 theMaxKinEnergy = std::max(((Ecm-theResMass)*(Ecm+theResMass) +
109 theMass*theMass)/twoEcm - theMass, 0.0);
110 theMinKinEnergy = (elim == 0.0) ? 0.0 :
111 std::max(((theMass+elim)*(twoEcm-theMass-elim) +
112 theMass*theMass)/twoEcm - theMass, 0.0);
113}
double G4double
Definition G4Types.hh:83
G4double GetGroundStateMass() const
G4double GetExcitationEnergy() const
const G4LorentzVector & GetMomentum() const
G4int GetZ_asInt() const
G4int GetA_asInt() const
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4double Z13(G4int Z) const
Definition G4Pow.hh:123
T max(const T t1, const T t2)
brief Return the largest of the two arguments

◆ IsItPossible()

G4bool G4VPreCompoundFragment::IsItPossible ( const G4Fragment & aFragment) const
inline

Definition at line 167 of file G4VPreCompoundFragment.hh.

168{
169 G4int pplus = aFragment.GetNumberOfCharged();
170 G4int pneut = aFragment.GetNumberOfParticles()-pplus;
171 return (pneut >= theA - theZ && pplus >= theZ && theMaxKinEnergy > 0.0);
172}
int G4int
Definition G4Types.hh:85
G4int GetNumberOfParticles() const
G4int GetNumberOfCharged() const

◆ operator!=()

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

◆ operator=()

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

◆ operator==()

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

◆ SampleKineticEnergy()

virtual G4double G4VPreCompoundFragment::SampleKineticEnergy ( const G4Fragment & )
pure virtual

◆ SetMomentum()

void G4VPreCompoundFragment::SetMomentum ( const G4LorentzVector & lv)
inline

Definition at line 113 of file G4VPreCompoundFragment.hh.

113{ theMomentum = lv; }

Referenced by G4PreCompoundEmission::PerformEmission().

◆ SetOPTxs()

void G4VPreCompoundFragment::SetOPTxs ( G4int opt)
inline

Definition at line 116 of file G4VPreCompoundFragment.hh.

116{ OPTxs = opt; }

◆ UseSICB()

void G4VPreCompoundFragment::UseSICB ( G4bool use)
inline

Definition at line 118 of file G4VPreCompoundFragment.hh.

Friends And Related Symbol Documentation

◆ operator<< [1/2]

std::ostream & operator<< ( std::ostream & out,
const G4VPreCompoundFragment & theFragment )
friend

Definition at line 58 of file G4VPreCompoundFragment.cc.

60{
61 out << &theFragment;
62 return out;
63}

◆ operator<< [2/2]

std::ostream & operator<< ( std::ostream & out,
const G4VPreCompoundFragment * theFragment )
friend

Definition at line 65 of file G4VPreCompoundFragment.cc.

67{
68 out
69 << "PreCompoundModel Emitted Fragment: Z= " << theFragment->GetZ()
70 << " A= " << theFragment->GetA()
71 << " Mass(GeV)= " << theFragment->GetNuclearMass()/CLHEP::GeV;
72 return out;
73}

Member Data Documentation

◆ fNucData

◆ g4calc

◆ OPTxs

G4int G4VPreCompoundFragment::OPTxs {3}
protected

Definition at line 155 of file G4VPreCompoundFragment.hh.

155{3};

Referenced by G4PreCompoundFragment::CrossSection(), Initialize(), and SetOPTxs().

◆ theA

◆ theBindingEnergy

G4double G4VPreCompoundFragment::theBindingEnergy {0.0}
protected

◆ theCoulombBarrier

◆ theEmissionProbability

G4double G4VPreCompoundFragment::theEmissionProbability {0.0}
protected

◆ theFragA

◆ theFragZ

◆ theMass

G4double G4VPreCompoundFragment::theMass
protected

Definition at line 149 of file G4VPreCompoundFragment.hh.

Referenced by G4VPreCompoundFragment(), GetNuclearMass(), and Initialize().

◆ theMaxKinEnergy

◆ theMinKinEnergy

G4double G4VPreCompoundFragment::theMinKinEnergy {0.0}
protected

◆ theParameters

G4DeexPrecoParameters* G4VPreCompoundFragment::theParameters
protected

◆ theReducedMass

G4double G4VPreCompoundFragment::theReducedMass {0.0}
protected

◆ theResA

◆ theResA13

◆ theResMass

G4double G4VPreCompoundFragment::theResMass {0.0}
protected

Definition at line 147 of file G4VPreCompoundFragment.hh.

147{0.0};

Referenced by GetRestNuclearMass(), and Initialize().

◆ theResZ

◆ theZ

◆ useSICB

G4bool G4VPreCompoundFragment::useSICB {true}
protected

Definition at line 157 of file G4VPreCompoundFragment.hh.

157{true};

Referenced by UseSICB().


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