Geant4  10.02.p03
G4VPreCompoundFragment Class Referenceabstract

#include <G4VPreCompoundFragment.hh>

Inheritance diagram for G4VPreCompoundFragment:
Collaboration diagram for G4VPreCompoundFragment:

Public Member Functions

 G4VPreCompoundFragment (const G4ParticleDefinition *, G4VCoulombBarrier *aCoulombBarrier)
 
virtual ~G4VPreCompoundFragment ()
 
void Initialize (const G4Fragment &aFragment)
 
virtual G4double CalcEmissionProbability (const G4Fragment &aFragment)=0
 
virtual G4double SampleKineticEnergy (const G4Fragment &aFragment)=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 &value)
 
void SetOPTxs (G4int)
 
void UseSICB (G4bool)
 

Protected Member Functions

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

Protected Attributes

G4PreCompoundParameterstheParameters
 
G4Powg4pow
 
G4int theA
 
G4int theZ
 
G4int theResA
 
G4int theResZ
 
G4int theFragA
 
G4int theFragZ
 
G4double theResA13
 
G4double theBindingEnergy
 
G4double theMaxKinEnergy
 
G4double theResMass
 
G4double theReducedMass
 
G4double theMass
 
G4double theEmissionProbability
 
G4double theCoulombBarrier
 
G4int OPTxs
 
G4bool useSICB
 

Private Member Functions

 G4VPreCompoundFragment ()
 
 G4VPreCompoundFragment (const G4VPreCompoundFragment &right)
 
const G4VPreCompoundFragmentoperator= (const G4VPreCompoundFragment &right)
 
G4int operator== (const G4VPreCompoundFragment &right) const
 
G4int operator!= (const G4VPreCompoundFragment &right) const
 

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

Constructor & Destructor Documentation

◆ G4VPreCompoundFragment() [1/3]

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

Definition at line 39 of file G4VPreCompoundFragment.cc.

41  : particle(part), theCoulombBarrierPtr(aCoulombBarrier),
42  theMomentum(0.,0.,0.,0.),
46  theMaxKinEnergy(0.0),theResMass(0.0),
47  theReducedMass(0.0),
49  OPTxs(3),useSICB(false)
50 {
54  theResA13 = 0;
55 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
const G4ParticleDefinition * particle
G4VCoulombBarrier * theCoulombBarrierPtr
int G4lrint(double ad)
Definition: templates.hh:163
G4PreCompoundParameters * theParameters
G4double GetPDGCharge() const
Here is the call graph for this function:

◆ ~G4VPreCompoundFragment()

G4VPreCompoundFragment::~G4VPreCompoundFragment ( )
virtual

Definition at line 57 of file G4VPreCompoundFragment.cc.

58 {
59  delete theParameters;
60 }
G4PreCompoundParameters * theParameters

◆ G4VPreCompoundFragment() [2/3]

G4VPreCompoundFragment::G4VPreCompoundFragment ( )
private

◆ G4VPreCompoundFragment() [3/3]

G4VPreCompoundFragment::G4VPreCompoundFragment ( const G4VPreCompoundFragment right)
private

Member Function Documentation

◆ CalcEmissionProbability()

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

Implemented in G4PreCompoundFragment, and G4HETCFragment.

◆ GetA()

G4int G4VPreCompoundFragment::GetA ( ) const
inline
Here is the caller graph for this function:

◆ GetAlpha()

◆ GetBeta()

virtual G4double G4VPreCompoundFragment::GetBeta ( ) const
protectedpure virtual

◆ GetBindingEnergy()

G4double G4VPreCompoundFragment::GetBindingEnergy ( ) const
inline
Here is the caller graph for this function:

◆ GetEmissionProbability()

G4double G4VPreCompoundFragment::GetEmissionProbability ( ) const
inline

◆ GetEnergyThreshold()

G4double G4VPreCompoundFragment::GetEnergyThreshold ( ) const
inline
Here is the caller graph for this function:

◆ GetMomentum()

const G4LorentzVector& G4VPreCompoundFragment::GetMomentum ( ) const
inline

◆ GetNuclearMass()

G4double G4VPreCompoundFragment::GetNuclearMass ( ) const
inline
Here is the caller graph for this function:

◆ GetReactionProduct()

G4ReactionProduct* G4VPreCompoundFragment::GetReactionProduct ( ) const
inline
Here is the caller graph for this function:

◆ GetRestA()

G4int G4VPreCompoundFragment::GetRestA ( ) const
inline
Here is the caller graph for this function:

◆ GetRestNuclearMass()

G4double G4VPreCompoundFragment::GetRestNuclearMass ( ) const
inline

◆ GetRestZ()

G4int G4VPreCompoundFragment::GetRestZ ( ) const
inline
Here is the caller graph for this function:

◆ GetZ()

G4int G4VPreCompoundFragment::GetZ ( ) const
inline
Here is the caller graph for this function:

◆ Initialize()

void G4VPreCompoundFragment::Initialize ( const G4Fragment aFragment)

Definition at line 80 of file G4VPreCompoundFragment.cc.

81 {
82  theFragA = aFragment.GetA_asInt();
83  theFragZ = aFragment.GetZ_asInt();
84  theResA = theFragA - theA;
85  theResZ = theFragZ - theZ;
86 
87  if ((theResA < theResZ) || (theResA < theA) || (theResZ < theZ))
88  {
89  // In order to be sure that emission probability will be 0.
90  theMaxKinEnergy = 0.0;
91  return;
92  }
93 
95 
96  // Calculate Coulomb barrier
98  GetCoulombBarrier(theResA,theResZ,aFragment.GetExcitationEnergy());
99 
100  // Calculate masses
102  theReducedMass = theResMass*theMass/(theResMass + theMass);
103 
104  // Compute Binding Energies for fragments
105  // needed to separate a fragment from the nucleus
106  theBindingEnergy = theResMass + theMass - aFragment.GetGroundStateMass();
107 
108  // Compute Maximal Kinetic Energy which can be carried by fragments
109  // after separation - the true assimptotic value
110  G4double Ecm = aFragment.GetMomentum().m();
111  theMaxKinEnergy = ((Ecm-theResMass)*(Ecm+theResMass) + theMass*theMass)
112  /(2.0*Ecm) - theMass;
113 }
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:273
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4int GetA_asInt() const
Definition: G4Fragment.hh:256
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:278
G4int GetZ_asInt() const
Definition: G4Fragment.hh:261
G4VCoulombBarrier * theCoulombBarrierPtr
G4double Z13(G4int Z) const
Definition: G4Pow.hh:127
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:289
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ IsItPossible()

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

◆ operator!=()

G4int G4VPreCompoundFragment::operator!= ( const G4VPreCompoundFragment right) const
private

◆ operator=()

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

◆ operator==()

G4int G4VPreCompoundFragment::operator== ( const G4VPreCompoundFragment right) const
private

◆ SampleKineticEnergy()

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

Implemented in G4PreCompoundFragment, G4HETCNeutron, and G4HETCChargedFragment.

Here is the caller graph for this function:

◆ SetMomentum()

void G4VPreCompoundFragment::SetMomentum ( const G4LorentzVector value)
inline
Here is the caller graph for this function:

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

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

◆ operator<< [2/2]

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

Definition at line 63 of file G4VPreCompoundFragment.cc.

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

Member Data Documentation

◆ g4pow

G4Pow* G4VPreCompoundFragment::g4pow
protected

Definition at line 151 of file G4VPreCompoundFragment.hh.

◆ OPTxs

G4int G4VPreCompoundFragment::OPTxs
protected

Definition at line 171 of file G4VPreCompoundFragment.hh.

◆ particle

const G4ParticleDefinition* G4VPreCompoundFragment::particle
private

Definition at line 143 of file G4VPreCompoundFragment.hh.

◆ theA

G4int G4VPreCompoundFragment::theA
protected

Definition at line 153 of file G4VPreCompoundFragment.hh.

◆ theBindingEnergy

G4double G4VPreCompoundFragment::theBindingEnergy
protected

Definition at line 161 of file G4VPreCompoundFragment.hh.

◆ theCoulombBarrier

G4double G4VPreCompoundFragment::theCoulombBarrier
protected

Definition at line 168 of file G4VPreCompoundFragment.hh.

◆ theCoulombBarrierPtr

G4VCoulombBarrier* G4VPreCompoundFragment::theCoulombBarrierPtr
private

Definition at line 144 of file G4VPreCompoundFragment.hh.

◆ theEmissionProbability

G4double G4VPreCompoundFragment::theEmissionProbability
protected

Definition at line 167 of file G4VPreCompoundFragment.hh.

◆ theFragA

G4int G4VPreCompoundFragment::theFragA
protected

Definition at line 157 of file G4VPreCompoundFragment.hh.

◆ theFragZ

G4int G4VPreCompoundFragment::theFragZ
protected

Definition at line 158 of file G4VPreCompoundFragment.hh.

◆ theMass

G4double G4VPreCompoundFragment::theMass
protected

Definition at line 165 of file G4VPreCompoundFragment.hh.

◆ theMaxKinEnergy

G4double G4VPreCompoundFragment::theMaxKinEnergy
protected

Definition at line 162 of file G4VPreCompoundFragment.hh.

◆ theMomentum

G4LorentzVector G4VPreCompoundFragment::theMomentum
private

Definition at line 146 of file G4VPreCompoundFragment.hh.

◆ theParameters

G4PreCompoundParameters* G4VPreCompoundFragment::theParameters
protected

Definition at line 150 of file G4VPreCompoundFragment.hh.

◆ theReducedMass

G4double G4VPreCompoundFragment::theReducedMass
protected

Definition at line 164 of file G4VPreCompoundFragment.hh.

◆ theResA

G4int G4VPreCompoundFragment::theResA
protected

Definition at line 155 of file G4VPreCompoundFragment.hh.

◆ theResA13

G4double G4VPreCompoundFragment::theResA13
protected

Definition at line 160 of file G4VPreCompoundFragment.hh.

◆ theResMass

G4double G4VPreCompoundFragment::theResMass
protected

Definition at line 163 of file G4VPreCompoundFragment.hh.

◆ theResZ

G4int G4VPreCompoundFragment::theResZ
protected

Definition at line 156 of file G4VPreCompoundFragment.hh.

◆ theZ

G4int G4VPreCompoundFragment::theZ
protected

Definition at line 154 of file G4VPreCompoundFragment.hh.

◆ useSICB

G4bool G4VPreCompoundFragment::useSICB
protected

Definition at line 173 of file G4VPreCompoundFragment.hh.


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