Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | Friends
G4VPreCompoundFragment Class Referenceabstract

#include <G4VPreCompoundFragment.hh>

Inheritance diagram for G4VPreCompoundFragment:
G4HETCFragment G4PreCompoundFragment G4HETCChargedFragment G4HETCNeutron G4PreCompoundIon G4PreCompoundNucleon G4HETCAlpha G4HETCDeuteron G4HETCHe3 G4HETCProton G4HETCTriton G4PreCompoundAlpha G4PreCompoundDeuteron G4PreCompoundHe3 G4PreCompoundTriton G4PreCompoundNeutron G4PreCompoundProton

Public Member Functions

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

Protected Member Functions

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

Protected Attributes

G4NuclearLevelDatafNucData
 
G4Powg4calc
 
G4int OPTxs
 
G4int theA
 
G4double theBindingEnergy
 
G4double theCoulombBarrier
 
G4double theEmissionProbability
 
G4int theFragA
 
G4int theFragZ
 
G4double theMass
 
G4double theMaxKinEnergy
 
G4double theMinKinEnergy
 
G4DeexPrecoParameterstheParameters
 
G4double theReducedMass
 
G4int theResA
 
G4double theResA13
 
G4double theResMass
 
G4int theResZ
 
G4int theZ
 
G4bool useSICB
 

Private Attributes

const G4ParticleDefinitionparticle
 
G4VCoulombBarriertheCoulombBarrierPtr
 
G4LorentzVector theMomentum
 

Friends

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

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 39 of file G4VPreCompoundFragment.cc.

41 : particle(part), theCoulombBarrierPtr(aCoulombBarrier),
42 theMomentum(0.,0.,0.,0.),
47 theReducedMass(0.0),
49 OPTxs(3),useSICB(true)
50{
55 theResA13 = 0.0;
56}
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
G4double GetPDGCharge() const
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4NuclearLevelData * fNucData
G4VCoulombBarrier * theCoulombBarrierPtr
G4DeexPrecoParameters * theParameters
const G4ParticleDefinition * particle
int G4lrint(double ad)
Definition: templates.hh:134

References fNucData, g4calc, G4Pow::GetInstance(), G4NuclearLevelData::GetInstance(), G4NuclearLevelData::GetParameters(), G4ParticleDefinition::GetPDGMass(), particle, theMass, theParameters, and theResA13.

◆ ~G4VPreCompoundFragment()

G4VPreCompoundFragment::~G4VPreCompoundFragment ( )
virtual

Definition at line 58 of file G4VPreCompoundFragment.cc.

59{}

◆ G4VPreCompoundFragment() [2/2]

G4VPreCompoundFragment::G4VPreCompoundFragment ( const G4VPreCompoundFragment right)
delete

Member Function Documentation

◆ CalcEmissionProbability()

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

Implemented in G4HETCFragment, and G4PreCompoundFragment.

◆ GetA()

G4int G4VPreCompoundFragment::GetA ( ) const
inline

◆ GetAlpha()

virtual G4double G4VPreCompoundFragment::GetAlpha ( ) const
protectedpure virtual

◆ GetBeta()

virtual G4double G4VPreCompoundFragment::GetBeta ( ) const
protectedpure virtual

◆ GetBindingEnergy()

G4double G4VPreCompoundFragment::GetBindingEnergy ( ) const
inline

◆ GetEmissionProbability()

G4double G4VPreCompoundFragment::GetEmissionProbability ( ) const
inline

◆ GetEnergyThreshold()

G4double G4VPreCompoundFragment::GetEnergyThreshold ( ) const
inline

◆ GetMomentum()

const G4LorentzVector & G4VPreCompoundFragment::GetMomentum ( ) const
inline

◆ GetNuclearMass()

G4double G4VPreCompoundFragment::GetNuclearMass ( ) const
inline

◆ GetReactionProduct()

G4ReactionProduct * G4VPreCompoundFragment::GetReactionProduct ( ) const
inline

◆ GetRestA()

G4int G4VPreCompoundFragment::GetRestA ( ) const
inline

◆ GetRestNuclearMass()

G4double G4VPreCompoundFragment::GetRestNuclearMass ( ) const
inline

◆ GetRestZ()

G4int G4VPreCompoundFragment::GetRestZ ( ) const
inline

◆ GetZ()

G4int G4VPreCompoundFragment::GetZ ( ) const
inline

◆ Initialize()

void G4VPreCompoundFragment::Initialize ( const G4Fragment aFragment)

Definition at line 79 of file G4VPreCompoundFragment.cc.

80{
81 theFragA = aFragment.GetA_asInt();
82 theFragZ = aFragment.GetZ_asInt();
85
87 if ((theResA < theResZ) || (theResA < theA) || (theResZ < theZ)) {
88 return;
89 }
90
93 GetCoulombBarrier(theResA,theResZ,aFragment.GetExcitationEnergy());
94
96
97 // Calculate masses
100
101 // Compute Binding Energies for fragments
102 // needed to separate a fragment from the nucleus
104
105 // Compute Maximal Kinetic Energy which can be carried by fragments
106 // after separation - the true assimptotic value
107 G4double Ecm = aFragment.GetMomentum().m();
108 G4double twoEcm = Ecm + Ecm;
110 /twoEcm - theMass,0.0);
111 theMinKinEnergy = (elim == 0.0) ? 0.0 :
112 std::max(((theMass+elim)*(twoEcm-theMass-elim) +
113 theMass*theMass)/twoEcm - theMass,0.0);
114}
double G4double
Definition: G4Types.hh:83
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:304
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:299
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:323
G4int GetZ_asInt() const
Definition: G4Fragment.hh:276
G4int GetA_asInt() const
Definition: G4Fragment.hh:271
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

References g4calc, G4Fragment::GetA_asInt(), G4Fragment::GetExcitationEnergy(), G4Fragment::GetGroundStateMass(), G4Fragment::GetMomentum(), G4NucleiProperties::GetNuclearMass(), G4Fragment::GetZ_asInt(), CLHEP::HepLorentzVector::m(), G4INCL::Math::max(), OPTxs, theA, theBindingEnergy, theCoulombBarrier, theCoulombBarrierPtr, theFragA, theFragZ, theMass, theMaxKinEnergy, theMinKinEnergy, theReducedMass, theResA, theResA13, theResMass, theResZ, theZ, and G4Pow::Z13().

◆ IsItPossible()

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

◆ 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 aFragment)
pure virtual

◆ SetMomentum()

void G4VPreCompoundFragment::SetMomentum ( const G4LorentzVector value)
inline

◆ SetOPTxs()

void G4VPreCompoundFragment::SetOPTxs ( G4int  )
inline

◆ UseSICB()

void G4VPreCompoundFragment::UseSICB ( G4bool  )
inline

Friends And Related Function Documentation

◆ operator<< [1/2]

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

Definition at line 61 of file G4VPreCompoundFragment.cc.

63{
64 out << &theFragment;
65 return out;
66}

◆ operator<< [2/2]

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

Definition at line 68 of file G4VPreCompoundFragment.cc.

70{
71 out
72 << "PreCompoundModel Emitted Fragment: Z= " << theFragment->GetZ()
73 << " A= " << theFragment->GetA()
74 << " Mass(GeV)= " << theFragment->GetNuclearMass()/CLHEP::GeV;
75 return out;
76}
G4double GetNuclearMass() const
static constexpr double GeV

Field Documentation

◆ fNucData

G4NuclearLevelData* G4VPreCompoundFragment::fNucData
protected

◆ g4calc

G4Pow* G4VPreCompoundFragment::g4calc
protected

◆ OPTxs

G4int G4VPreCompoundFragment::OPTxs
protected

Definition at line 161 of file G4VPreCompoundFragment.hh.

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

◆ particle

const G4ParticleDefinition* G4VPreCompoundFragment::particle
private

◆ theA

G4int G4VPreCompoundFragment::theA
protected

◆ theBindingEnergy

G4double G4VPreCompoundFragment::theBindingEnergy
protected

◆ theCoulombBarrier

G4double G4VPreCompoundFragment::theCoulombBarrier
protected

◆ theCoulombBarrierPtr

G4VCoulombBarrier* G4VPreCompoundFragment::theCoulombBarrierPtr
private

Definition at line 132 of file G4VPreCompoundFragment.hh.

Referenced by Initialize().

◆ theEmissionProbability

G4double G4VPreCompoundFragment::theEmissionProbability
protected

◆ theFragA

G4int G4VPreCompoundFragment::theFragA
protected

◆ theFragZ

G4int G4VPreCompoundFragment::theFragZ
protected

◆ theMass

G4double G4VPreCompoundFragment::theMass
protected

Definition at line 155 of file G4VPreCompoundFragment.hh.

Referenced by G4VPreCompoundFragment(), and Initialize().

◆ theMaxKinEnergy

G4double G4VPreCompoundFragment::theMaxKinEnergy
protected

◆ theMinKinEnergy

G4double G4VPreCompoundFragment::theMinKinEnergy
protected

◆ theMomentum

G4LorentzVector G4VPreCompoundFragment::theMomentum
private

Definition at line 134 of file G4VPreCompoundFragment.hh.

◆ theParameters

G4DeexPrecoParameters* G4VPreCompoundFragment::theParameters
protected

◆ theReducedMass

G4double G4VPreCompoundFragment::theReducedMass
protected

◆ theResA

G4int G4VPreCompoundFragment::theResA
protected

◆ theResA13

G4double G4VPreCompoundFragment::theResA13
protected

◆ theResMass

G4double G4VPreCompoundFragment::theResMass
protected

Definition at line 153 of file G4VPreCompoundFragment.hh.

Referenced by Initialize().

◆ theResZ

G4int G4VPreCompoundFragment::theResZ
protected

◆ theZ

G4int G4VPreCompoundFragment::theZ
protected

◆ useSICB

G4bool G4VPreCompoundFragment::useSICB
protected

Definition at line 163 of file G4VPreCompoundFragment.hh.


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