Geant4 11.3.0
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 ()
 
G4bool Initialize (const G4Fragment &aFragment)
 
virtual G4double CalcEmissionProbability (const G4Fragment &)=0
 
virtual G4double SampleKineticEnergy (const G4Fragment &)=0
 
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)
 
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
 
G4InterfaceToXSfXSection {nullptr}
 
G4int theA
 
G4int theZ
 
G4int theResA {0}
 
G4int theResZ {0}
 
G4int theFragA {0}
 
G4int theFragZ {0}
 
G4int OPTxs
 
G4int index {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}
 
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 56 of file G4VPreCompoundFragment.hh.

Constructor & Destructor Documentation

◆ G4VPreCompoundFragment() [1/2]

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

Definition at line 41 of file G4VPreCompoundFragment.cc.

43 : theA(part->GetBaryonNumber()),
44 theZ(G4lrint(part->GetPDGCharge()/CLHEP::eplus)),
45 particle(part),
46 theCoulombBarrierPtr(aCoulombBarrier)
47{
48 theMass = particle->GetPDGMass();
50 theParameters = fNucData->GetParameters();
51 OPTxs = theParameters->GetDeexModelType();
53
54 if (1 == theZ && 1 == theA) { index = 1; }
55 else if (1 == theZ && 2 == theA) { index = 2; }
56 else if (1 == theZ && 3 == theA) { index = 3; }
57 else if (2 == theZ && 3 == theA) { index = 4; }
58 else if (2 == theZ && 4 == theA) { index = 5; }
59
60 if (OPTxs == 1) {
61 fXSection = new G4InterfaceToXS(particle, index);
62 }
63}
static G4NuclearLevelData * GetInstance()
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4DeexPrecoParameters * theParameters
int G4lrint(double ad)
Definition templates.hh:134

Referenced by G4HETCFragment::G4HETCFragment(), G4PreCompoundFragment::G4PreCompoundFragment(), G4VPreCompoundFragment(), operator!=(), operator<<, operator<<, operator=(), and operator==().

◆ ~G4VPreCompoundFragment()

G4VPreCompoundFragment::~G4VPreCompoundFragment ( )
virtual

Definition at line 65 of file G4VPreCompoundFragment.cc.

66{
67 delete theCoulombBarrierPtr;
68 delete fXSection;
69}

◆ 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 129 of file G4VPreCompoundFragment.hh.

◆ GetBindingEnergy()

G4double G4VPreCompoundFragment::GetBindingEnergy ( ) const
inline

Definition at line 97 of file G4VPreCompoundFragment.hh.

◆ GetEmissionProbability()

G4double G4VPreCompoundFragment::GetEmissionProbability ( ) const
inline

Definition at line 104 of file G4VPreCompoundFragment.hh.

◆ GetEnergyThreshold()

G4double G4VPreCompoundFragment::GetEnergyThreshold ( ) const
inline

◆ GetMomentum()

const G4LorentzVector & G4VPreCompoundFragment::GetMomentum ( ) const
inline

Definition at line 110 of file G4VPreCompoundFragment.hh.

110{ return theMomentum; }

Referenced by GetReactionProduct().

◆ GetNuclearMass()

G4double G4VPreCompoundFragment::GetNuclearMass ( ) const
inline

Definition at line 106 of file G4VPreCompoundFragment.hh.

106{ return theMass; }

Referenced by operator<<, and G4PreCompoundEmission::PerformEmission().

◆ GetReactionProduct()

G4ReactionProduct * G4VPreCompoundFragment::GetReactionProduct ( ) const
inline

Definition at line 167 of file G4VPreCompoundFragment.hh.

168{
169 G4ReactionProduct* theReactionProduct = new G4ReactionProduct(particle);
170 theReactionProduct->SetMomentum(GetMomentum().vect());
171 theReactionProduct->SetTotalEnergy(GetMomentum().e());
172 return theReactionProduct;
173}
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 108 of file G4VPreCompoundFragment.hh.

◆ GetRestZ()

G4int G4VPreCompoundFragment::GetRestZ ( ) const
inline

◆ GetZ()

G4int G4VPreCompoundFragment::GetZ ( ) const
inline

Definition at line 91 of file G4VPreCompoundFragment.hh.

91{ return theZ; }

Referenced by operator<<, and G4PreCompoundEmission::PerformEmission().

◆ Initialize()

G4bool G4VPreCompoundFragment::Initialize ( const G4Fragment & aFragment)

Definition at line 89 of file G4VPreCompoundFragment.cc.

90{
91 theFragA = aFragment.GetA_asInt();
92 theFragZ = aFragment.GetZ_asInt();
95
97 if ((theResA < theResZ) || (theResA < theA) || (theResZ < theZ)
98 || (theResA == theA && theResZ < theZ)
99 || ((theResA > 1) && (theResA == theResZ || theResZ == 0))) {
100 return false;
101 }
103 G4double Ecm = aFragment.GetMomentum().m();
104 if (Ecm <= theResMass + theMass) { return 0.0; }
105
106 theResA13 = g4calc->Z13(theResA);
107
108 G4double elim = 0.0;
109 if (0 < theZ) {
110 theCoulombBarrier = theCoulombBarrierPtr->
111 GetCoulombBarrier(theResA, theResZ, aFragment.GetExcitationEnergy());
112 elim = (0 < OPTxs) ? theCoulombBarrier*0.5 : theCoulombBarrier;
113 }
114
115 // Compute Maximal Kinetic Energy which can be carried by fragments
116 // after separation - the true assimptotic value
118 0.5*((Ecm - theResMass)*(Ecm + theResMass) + theMass*theMass)/Ecm - theMass;
119 G4double resM = Ecm - theMass - elim;
120 if (resM < theResMass) { return false; }
122 0.5*((Ecm - resM)*(Ecm + resM) + theMass*theMass)/Ecm - theMass;
123
124 if (theMinKinEnergy >= theMaxKinEnergy) { return false; }
125 // Calculate masses
127
128 // Compute Binding Energies for fragments
129 // needed to separate a fragment from the nucleus
131 return true;
132}
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)

Referenced by G4PreCompoundFragment::CalcEmissionProbability().

◆ 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 112 of file G4VPreCompoundFragment.hh.

112{ theMomentum = lv; }

Referenced by G4PreCompoundEmission::PerformEmission().

◆ SetOPTxs()

void G4VPreCompoundFragment::SetOPTxs ( G4int )
inline

Definition at line 115 of file G4VPreCompoundFragment.hh.

115{}

◆ UseSICB()

void G4VPreCompoundFragment::UseSICB ( G4bool use)
inline

Definition at line 117 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 71 of file G4VPreCompoundFragment.cc.

73{
74 out << &theFragment;
75 return out;
76}

◆ operator<< [2/2]

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

Definition at line 78 of file G4VPreCompoundFragment.cc.

80{
81 out
82 << "PreCompoundModel Emitted Fragment: Z= " << theFragment->GetZ()
83 << " A= " << theFragment->GetA()
84 << " Mass(GeV)= " << theFragment->GetNuclearMass()/CLHEP::GeV;
85 return out;
86}

Member Data Documentation

◆ fNucData

◆ fXSection

G4InterfaceToXS* G4VPreCompoundFragment::fXSection {nullptr}
protected

◆ g4calc

◆ index

G4int G4VPreCompoundFragment::index {0}
protected

Definition at line 144 of file G4VPreCompoundFragment.hh.

144{0};

Referenced by G4PreCompoundFragment::CrossSection(), and G4VPreCompoundFragment().

◆ OPTxs

G4int G4VPreCompoundFragment::OPTxs
protected

◆ 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 152 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 150 of file G4VPreCompoundFragment.hh.

150{0.0};

Referenced by GetRestNuclearMass(), and Initialize().

◆ theResZ

◆ theZ

◆ useSICB

G4bool G4VPreCompoundFragment::useSICB {true}
protected

Definition at line 158 of file G4VPreCompoundFragment.hh.

158{true};

Referenced by UseSICB().


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