Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EvaporationChannel.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//J.M. Quesada (August2008). Based on:
27//
28// Hadronic Process: Nuclear De-excitations
29// by V. Lara (Oct 1998)
30//
31// Modified:
32// 03-09-2008 J.M. Quesada for external choice of inverse cross section option
33// 06-09-2008 J.M. Quesada Also external choices have been added for superimposed
34// Coulomb barrier (if useSICB is set true, by default is false)
35// 17-11-2010 V.Ivanchenko in constructor replace G4VEmissionProbability by
36// G4EvaporationProbability and do not new and delete probability
37// object at each call; use G4Pow
38
41#include "G4CoulombBarrier.hh"
42#include "G4NuclearLevelData.hh"
43#include "G4NucleiProperties.hh"
44#include "G4Pow.hh"
45#include "G4Log.hh"
46#include "G4Exp.hh"
48#include "G4SystemOfUnits.hh"
49#include "Randomize.hh"
50#include "G4RandomDirection.hh"
52
56 theProbability(aprob),
57 theCoulombBarrier(new G4CoulombBarrier(anA, aZ)),
58 theA(anA), theZ(aZ)
59{
60 secID = G4PhysicsModelCatalog::GetModelID("model_G4EvaporationChannel");
61 evapMass = G4NucleiProperties::GetNuclearMass(theA, theZ);
62 evapMass2 = evapMass*evapMass;
63 theLevelData = G4NuclearLevelData::GetInstance();
64}
65
67{
68 delete theCoulombBarrier;
69}
70
76
78{
79 theProbability->ResetProbability();
80 G4int fragA = fragment->GetA_asInt();
81 G4int fragZ = fragment->GetZ_asInt();
82 resA = fragA - theA;
83 resZ = fragZ - theZ;
84
85 // Only channels which are physically allowed are taken into account
86 if(resA < theA || resA < resZ || resZ < 0 || (resA == theA && resZ < theZ)
87 || ((resA > 1) && (resA == resZ || resZ == 0)))
88 { return 0.0; }
89
90 G4double exEnergy = fragment->GetExcitationEnergy();
91 G4double delta0 = theLevelData->GetPairingCorrection(fragZ,fragA);
92 /*
93 G4cout << "G4EvaporationChannel::Initialize Z= "<<theZ<<" A= "<<theA
94 << " FragZ= " << fragZ << " FragA= " << fragA
95 << " exEnergy= " << exEnergy << " d0= " << delta0 << G4endl;
96 */
97 if(exEnergy < delta0) { return 0.0; }
98
99 G4double fragMass = fragment->GetGroundStateMass();
100 mass = fragMass + exEnergy;
101
102 resMass = G4NucleiProperties::GetNuclearMass(resA, resZ);
103 ekinmax = 0.5*((mass-resMass)*(mass+resMass) + evapMass2)/mass - evapMass;
104
105 G4double elim = 0.0;
106 if(theZ > 0) {
107 bCoulomb = theCoulombBarrier->GetCoulombBarrier(resA, resZ, 0.0);
108
109 // for OPTxs >0 penetration under the barrier is taken into account
110 elim = (0 != OPTxs) ? bCoulomb*0.6 : bCoulomb;
111 }
112 /*
113 G4cout << "exEnergy= " << exEnergy << " Ec= " << bCoulomb
114 << " elim=" << elim << " d0= " << delta0
115 << " Free= " << mass - resMass - evapMass
116 << G4endl;
117 */
118 if(mass <= resMass + evapMass + elim) { return 0.0; }
119
120 G4double ekinmin = 0.0;
121 if(elim > 0.0) {
122 G4double resM = mass - evapMass - elim;
123 ekinmin =
124 std::max(0.5*((mass-resM)*(mass+resM) + evapMass2)/mass - evapMass, 0.0);
125 }
126 /*
127 G4cout << "Emin= " <<ekinmin<<" Emax= "<<ekinmax
128 << " mass= " << mass << " resM= " << resMass
129 << " evapM= " << evapMass << G4endl;
130 */
131 if(ekinmax <= ekinmin) { return 0.0; }
132
133 theProbability->SetDecayKinematics(resZ, resA, resMass, mass);
134 G4double prob = theProbability->TotalProbability(*fragment, ekinmin,
135 ekinmax, bCoulomb,
136 exEnergy - delta0);
137 return prob;
138}
139
141{
142 G4double ekin = ekinmax;
143 // assumed, that TotalProbability(...) was already called
144 // if value iz zero no possiblity to sample final state
145 if(resA > 4 && theProbability->GetProbability() > 0.0) {
146 ekin = theProbability->SampleEnergy();
147 }
148 ekin = std::max(ekin, 0.0);
149 G4LorentzVector lv0 = theNucleus->GetMomentum();
150 G4LorentzVector lv(std::sqrt(ekin*(ekin + 2.0*evapMass))*G4RandomDirection(),
151 ekin + evapMass);
152 lv.boost(lv0.boostVector());
153
154 G4Fragment* evFragment = new G4Fragment(theA, theZ, lv);
155 evFragment->SetCreatorModelID(secID);
156 lv0 -= lv;
157 theNucleus->SetZAandMomentum(lv0, resZ, resA);
158 theNucleus->SetCreatorModelID(secID);
159 return evFragment;
160}
161
163 G4double kinEnergy)
164{
165 ComputeProbability(frag, kinEnergy);
166 return theProbability->CrossSection(kinEnergy, bCoulomb);
167}
168
170 G4double kinEnergy)
171{
172 G4int fragA = frag->GetA_asInt();
173 G4int fragZ = frag->GetZ_asInt();
174 resA = fragA - theA;
175 resZ = fragZ - theZ;
176 bCoulomb = theCoulombBarrier->GetCoulombBarrier(resA, resZ, 0.0);
177 return theProbability->ComputeProbability(kinEnergy, bCoulomb);
178}
G4ThreeVector G4RandomDirection()
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U) const override
G4EvaporationChannel(G4int A, G4int Z, G4EvaporationProbability *)
G4double ComputeInverseXSection(G4Fragment *, G4double kinEnergy) override
G4double GetEmissionProbability(G4Fragment *fragment) override
G4Fragment * EmittedFragment(G4Fragment *theNucleus) override
G4double ComputeProbability(G4Fragment *, G4double kinEnergy) override
G4double CrossSection(G4double K, G4double CB)
G4double ComputeProbability(G4double K, G4double CB) override
virtual G4double TotalProbability(const G4Fragment &fragment, G4double minKinEnergy, G4double maxKinEnergy, G4double CB, G4double exEnergy)
G4double GetGroundStateMass() const
G4double GetExcitationEnergy() const
const G4LorentzVector & GetMomentum() const
void SetCreatorModelID(G4int value)
G4int GetZ_asInt() const
G4int GetA_asInt() const
void SetZAandMomentum(const G4LorentzVector &, G4int Z, G4int A, G4int nLambdas=0)
G4PairingCorrection * GetPairingCorrection()
static G4NuclearLevelData * GetInstance()
static G4double GetNuclearMass(const G4double A, const G4double Z)
static G4int GetModelID(const G4int modelIndex)
void SetDecayKinematics(G4int rZ, G4int rA, G4double rmass, G4double fmass)