Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EvaporationChannel Class Reference

#include <G4EvaporationChannel.hh>

Inheritance diagram for G4EvaporationChannel:
Collaboration diagram for G4EvaporationChannel:

Public Member Functions

 G4EvaporationChannel (G4int A, G4int Z, const G4String &aName, G4EvaporationProbability *, G4VCoulombBarrier *)
 
virtual ~G4EvaporationChannel ()
 
void Initialise ()
 
virtual G4double GetEmissionProbability (G4Fragment *fragment)
 
virtual G4FragmentEmittedFragment (G4Fragment *theNucleus)
 
- Public Member Functions inherited from G4VEvaporationChannel
 G4VEvaporationChannel (const G4String &aName="")
 
virtual ~G4VEvaporationChannel ()
 
virtual G4double GetLifeTime (G4Fragment *theNucleus)
 
virtual G4bool BreakUpChain (G4FragmentVector *theResult, G4Fragment *theNucleus)
 
G4FragmentVectorBreakUpFragment (G4Fragment *theNucleus)
 
virtual void Dump () const
 
virtual void SetICM (G4bool)
 
virtual void RDMForced (G4bool)
 
virtual G4double GetFinalLevelEnergy (G4int Z, G4int A, G4double energy)
 
virtual G4double GetUpperLevelEnergy (G4int Z, G4int A)
 
G4double GetMaxLevelEnergy (G4int Z, G4int A)
 
G4double GetNearestLevelEnergy (G4int Z, G4int A, G4double energy)
 
void SetPhotonEvaporation (G4VEvaporationChannel *p)
 
void SetOPTxs (G4int opt)
 
void UseSICB (G4bool use)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEvaporationChannel
G4int OPTxs
 
G4bool useSICB
 

Detailed Description

Definition at line 47 of file G4EvaporationChannel.hh.

Constructor & Destructor Documentation

G4EvaporationChannel::G4EvaporationChannel ( G4int  A,
G4int  Z,
const G4String aName,
G4EvaporationProbability aprob,
G4VCoulombBarrier barrier 
)
explicit

Definition at line 53 of file G4EvaporationChannel.cc.

56  :
57  G4VEvaporationChannel(aName),
58  theA(anA),
59  theZ(aZ),
60  theProbability(aprob),
61  theCoulombBarrier(barrier)
62 {
63  ResA = ResZ = 0;
64  Mass = CoulombBarrier = MinKinEnergy = MaxKinEnergy = EmissionProbability = 0.0;
65  EvapMass = G4NucleiProperties::GetNuclearMass(theA, theZ);
66  pairingCorrection = G4PairingCorrection::GetInstance();
67 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4VEvaporationChannel(const G4String &aName="")
static G4PairingCorrection * GetInstance()

Here is the call graph for this function:

G4EvaporationChannel::~G4EvaporationChannel ( )
virtual

Definition at line 69 of file G4EvaporationChannel.cc.

70 {}

Member Function Documentation

G4Fragment * G4EvaporationChannel::EmittedFragment ( G4Fragment theNucleus)
virtual

Reimplemented from G4VEvaporationChannel.

Definition at line 130 of file G4EvaporationChannel.cc.

131 {
132  G4Fragment* evFragment = nullptr;
133  G4double ekin = 0.0;
134  if(ResA <= 4 &&
135  ((ResA == 4 && ResZ == 2) || (ResA == 3 && ResZ == 2) ||
136  (ResA == 3 && ResZ == 1) || (ResA == 2 && ResZ == 1) ||
137  (ResA == 1 && ResZ == 1) || (ResA == 1 && ResZ == 0) )) {
138  G4double mres = G4NucleiProperties::GetNuclearMass(ResA, ResZ);
139  ekin = 0.5*(Mass*Mass - mres*mres + EvapMass*EvapMass)/Mass - EvapMass;
140  } else {
141  ekin = theProbability->SampleKineticEnergy(MinKinEnergy, MaxKinEnergy,
142  CoulombBarrier);
143  }
144  G4LorentzVector lv0 = theNucleus->GetMomentum();
145  G4LorentzVector lv(std::sqrt(ekin*(ekin + 2.0*EvapMass))*G4RandomDirection(),
146  ekin + EvapMass);
147  lv.boost(lv0.boostVector());
148 
149  evFragment = new G4Fragment(theA, theZ, lv);
150  lv0 -= lv;
151  theNucleus->SetZandA_asInt(ResZ, ResA);
152  theNucleus->SetMomentum(lv0);
153 
154  return evFragment;
155 }
Hep3Vector boostVector() const
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ThreeVector G4RandomDirection()
G4double SampleKineticEnergy(G4double minKineticEnergy, G4double maxKineticEnergy, G4double CoulombBarrier=0.0)
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:307
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:312
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:276
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

G4double G4EvaporationChannel::GetEmissionProbability ( G4Fragment fragment)
virtual

Implements G4VEvaporationChannel.

Definition at line 78 of file G4EvaporationChannel.cc.

79 {
80  G4int FragA = fragment->GetA_asInt();
81  G4int FragZ = fragment->GetZ_asInt();
82  ResA = FragA - theA;
83  ResZ = FragZ - theZ;
84 
85  G4double FragmentMass = fragment->GetGroundStateMass();
86  G4double ExEnergy = fragment->GetExcitationEnergy();
87  Mass = FragmentMass + ExEnergy;
88  //G4cout << "G4EvaporationChannel::Initialize Z= " << theZ << " A= " << theA
89  // << " FragZ= " << FragZ << " FragA= " << FragA << G4endl;
90  EmissionProbability = 0.0;
91 
92  // Only channels which are physically allowed are taken into account
93  if (ResA >= ResZ && ResZ > 0 && ResA >= theA) {
94 
95  //Effective excitation energy
96  G4double ResMass = G4NucleiProperties::GetNuclearMass(ResA, ResZ);
97 
98  CoulombBarrier = (0 == theZ) ? 0.0 :
99  theCoulombBarrier->GetCoulombBarrier(ResA,ResZ,ExEnergy);
100 
101  G4double delta0 =
102  std::max(0.0,pairingCorrection->GetPairingCorrection(FragA,FragZ));
103  G4double delta1 =
104  std::max(0.0,pairingCorrection->GetPairingCorrection(ResA,ResZ));
105  ResMass += delta1;
106  /*
107  G4cout << "ExEnergy= " << ExEnergy << " Ec= " << CoulombBarrier
108  << " delta0= " << delta0 << " delta1= " << delta1
109  << " Free= " << Mass - ResMass - EvapMass
110  << G4endl;
111  */
112  // for OPTxs >0 penetration under the barrier is taken into account
113  G4double elim = (0 == OPTxs) ? CoulombBarrier : CoulombBarrier*0.7;
114  if(ExEnergy >= delta0 && Mass >= ResMass + EvapMass + elim) {
115  G4double xm2 = (Mass - EvapMass)*(Mass - EvapMass);
116  G4double xm = Mass - EvapMass - elim;
117  MinKinEnergy = (0.0 == elim) ? 0.0 : std::max(0.5*(xm2 - xm*xm)/Mass, 0.0);
118  MaxKinEnergy = std::max(0.5*(xm2 - ResMass*ResMass)/Mass, 0.0);
119  //G4cout << "Emin= " << MinKinEnergy << " Emax= " << MaxKinEnergy
120  // << " xm= " << xm << G4endl;
121  EmissionProbability = theProbability->
122  TotalProbability(*fragment, MinKinEnergy, MaxKinEnergy, CoulombBarrier);
123  }
124  }
125  //G4cout << "G4EvaporationChannel:: probability= "
126  // << EmissionProbability << G4endl;
127  return EmissionProbability;
128 }
static G4double GetNuclearMass(const G4double A, const G4double Z)
int G4int
Definition: G4Types.hh:78
G4int GetA_asInt() const
Definition: G4Fragment.hh:266
G4double GetPairingCorrection(G4int A, G4int Z) const
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:288
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4int GetZ_asInt() const
Definition: G4Fragment.hh:271
virtual G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U) const =0
double G4double
Definition: G4Types.hh:76
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:283

Here is the call graph for this function:

void G4EvaporationChannel::Initialise ( )
virtual

Reimplemented from G4VEvaporationChannel.

Definition at line 72 of file G4EvaporationChannel.cc.

73 {
74  theProbability->Initialise();
76 }

Here is the call graph for this function:


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