Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GNASHTransitions Class Reference

#include <G4GNASHTransitions.hh>

+ Inheritance diagram for G4GNASHTransitions:

Public Member Functions

 G4GNASHTransitions ()
 
virtual ~G4GNASHTransitions ()
 
virtual G4double CalculateProbability (const G4Fragment &aFragment)
 
virtual void PerformTransition (G4Fragment &aFragment)
 
- Public Member Functions inherited from G4VPreCompoundTransitions
 G4VPreCompoundTransitions ()
 
virtual ~G4VPreCompoundTransitions ()
 
virtual G4double CalculateProbability (const G4Fragment &aFragment)=0
 
virtual void PerformTransition (G4Fragment &aFragment)=0
 
G4double GetTransitionProb1 () const
 
G4double GetTransitionProb2 () const
 
G4double GetTransitionProb3 () const
 
void UseNGB (G4bool use)
 
void UseCEMtr (G4bool use)
 
 G4VPreCompoundTransitions (const G4VPreCompoundTransitions &)=delete
 
const G4VPreCompoundTransitionsoperator= (const G4VPreCompoundTransitions &right)=delete
 
G4bool operator== (const G4VPreCompoundTransitions &right) const =delete
 
G4bool operator!= (const G4VPreCompoundTransitions &right) const =delete
 

Additional Inherited Members

- Protected Attributes inherited from G4VPreCompoundTransitions
G4bool useNGB
 
G4bool useCEMtr
 
G4double TransitionProb1
 
G4double TransitionProb2
 
G4double TransitionProb3
 

Detailed Description

Definition at line 35 of file G4GNASHTransitions.hh.

Constructor & Destructor Documentation

◆ G4GNASHTransitions()

G4GNASHTransitions::G4GNASHTransitions ( )
explicit

Definition at line 37 of file G4GNASHTransitions.cc.

38{
40}
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()

◆ ~G4GNASHTransitions()

G4GNASHTransitions::~G4GNASHTransitions ( )
virtual

Definition at line 42 of file G4GNASHTransitions.cc.

43{}

Member Function Documentation

◆ CalculateProbability()

G4double G4GNASHTransitions::CalculateProbability ( const G4Fragment aFragment)
virtual

Implements G4VPreCompoundTransitions.

Definition at line 45 of file G4GNASHTransitions.cc.

47{
48 static const G4double k = 135.0 *CLHEP::MeV*CLHEP::MeV*CLHEP::MeV;
49 G4double E = aFragment.GetExcitationEnergy();
50 G4double P = aFragment.GetNumberOfParticles();
51 G4double H = aFragment.GetNumberOfHoles();
52 G4double N = P + H;
53 G4int Z = aFragment.GetZ_asInt();
54 G4int A = aFragment.GetA_asInt();
55
56 G4double theMatrixElement(k*N/((A*A*A)*E));
57 G4double x = E/(N*CLHEP::MeV);
58 static const G4double xf = std::sqrt(2.0/7.0);
59 if ( x < 2.0) { x *= xf; }
60 else if ( x < 7.0) { x *= std::sqrt(x/7.0); }
61 else if ( x > 15.0){ x *= std::sqrt(15.0/x); }
62 theMatrixElement *= x;
63
65
66 G4double Epauli = ((P+1.0)*(P+1.0) + (H+1.0)*(H+1.0) + (P+1.0) - 3.0*(H-1.0))*0.25;
67
68 G4double Probability = gg*gg*gg *(E-Epauli)*(E-Epauli);
69 Probability *= theMatrixElement/(2.0*(N+1.0)*CLHEP::h_Planck);
70
71 return Probability;
72}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
G4int GetNumberOfParticles() const
Definition: G4Fragment.hh:366
G4int GetNumberOfHoles() const
Definition: G4Fragment.hh:386
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:312
G4int GetZ_asInt() const
Definition: G4Fragment.hh:289
G4int GetA_asInt() const
Definition: G4Fragment.hh:284
G4double GetLevelDensity(G4int Z, G4int A, G4double U)
#define N
Definition: crc32.c:56

◆ PerformTransition()

void G4GNASHTransitions::PerformTransition ( G4Fragment aFragment)
virtual

Implements G4VPreCompoundTransitions.

Definition at line 74 of file G4GNASHTransitions.cc.

75{
76 result.SetNumberOfParticles(result.GetNumberOfParticles()+1);
77 result.SetNumberOfHoles(result.GetNumberOfHoles()+1);
78 if (G4UniformRand()*result.GetA_asInt() <= G4double(result.GetZ_asInt()))
79 {
80 result.SetNumberOfCharged(result.GetNumberOfCharged()+1);
81 }
82
83 if (result.GetNumberOfParticles() < result.GetNumberOfCharged())
84 {
85 result.SetNumberOfCharged(result.GetNumberOfParticles());
86 }
87}
#define G4UniformRand()
Definition: Randomize.hh:52

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