Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4PreCompoundEmission.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//
27// -------------------------------------------------------------------
28//
29// GEANT4 Class file
30//
31//
32// File name: G4PreCompoundEmission
33//
34// Author: V.Lara
35//
36// Modified:
37// 15.01.2010 J.M.Quesada added protection against unphysical values of parameter an
38// 19.01.2010 V.Ivanchenko simplified computation of parameter an, sample cosTheta
39// instead of theta; protect all calls to sqrt
40// 20.08.2010 V.Ivanchenko added G4Pow and G4PreCompoundParameters pointers
41// use int Z and A and cleanup
42//
43
46#include "G4SystemOfUnits.hh"
47#include "G4Pow.hh"
48#include "G4Exp.hh"
49#include "G4Log.hh"
50#include "Randomize.hh"
51#include "G4RandomDirection.hh"
55#include "G4NuclearLevelData.hh"
58
59
61{
62 theFragmentsFactory = new G4PreCompoundEmissionFactory();
63 theFragmentsVector =
64 new G4PreCompoundFragmentVector(theFragmentsFactory->GetFragmentVector());
65 g4calc = G4Pow::GetInstance();
67 G4DeexPrecoParameters* param = fNuclData->GetParameters();
68 fFermiEnergy = param->GetFermiEnergy();
69 fUseAngularGenerator = param->UseAngularGen();
70 fModelID = G4PhysicsModelCatalog::GetModelID("model_PRECO");
71}
72
74{
75 delete theFragmentsFactory;
76 delete theFragmentsVector;
77}
78
80{
81 if (theFragmentsFactory) { delete theFragmentsFactory; }
82 theFragmentsFactory = new G4PreCompoundEmissionFactory();
83 if (theFragmentsVector) {
84 theFragmentsVector->SetVector(theFragmentsFactory->GetFragmentVector());
85 } else {
86 theFragmentsVector =
87 new G4PreCompoundFragmentVector(theFragmentsFactory->GetFragmentVector());
88 }
89}
90
92{
93 if (theFragmentsFactory) delete theFragmentsFactory;
94 theFragmentsFactory = new G4HETCEmissionFactory();
95 if (theFragmentsVector) {
96 theFragmentsVector->SetVector(theFragmentsFactory->GetFragmentVector());
97 } else {
98 theFragmentsVector =
99 new G4PreCompoundFragmentVector(theFragmentsFactory->GetFragmentVector());
100 }
101}
102
105{
106 G4ReactionProduct* res = nullptr;
107 // Choose a Fragment for emission
108 G4VPreCompoundFragment * thePreFragment =
109 theFragmentsVector->ChooseFragment();
110 if (thePreFragment == nullptr)
111 {
112 G4cout << "G4PreCompoundEmission::PerformEmission : "
113 << "I couldn't choose a fragment\n"
114 << "while trying to de-excite\n"
115 << aFragment << G4endl;
116 throw G4HadronicException(__FILE__, __LINE__, "");
117 return res;
118 }
119
120 // Kinetic Energy of emitted fragment
121 G4double kinEnergy = thePreFragment->SampleKineticEnergy(aFragment);
122 kinEnergy = std::max(kinEnergy, 0.0);
123
124 // Calculate the fragment momentum (three vector)
125 if(fUseAngularGenerator) {
126 AngularDistribution(thePreFragment,aFragment,kinEnergy);
127 } else {
128 G4double pmag =
129 std::sqrt(kinEnergy*(kinEnergy + 2.0*thePreFragment->GetNuclearMass()));
130 theFinalMomentum = pmag*G4RandomDirection();
131 }
132
133 // Mass of emittef fragment
134 G4double EmittedMass = thePreFragment->GetNuclearMass();
135 // Now we can calculate the four momentum
136 // both options are valid and give the same result but 2nd one is faster
137 G4LorentzVector Emitted4Momentum(theFinalMomentum,EmittedMass + kinEnergy);
138
139 if (2 < fVerbose) {
140 G4cout << " Emitted Z="
141 << thePreFragment->GetZ() << " A=" << thePreFragment->GetA()
142 << " Ekin(MeV)=" << kinEnergy << " 4-mom C.M.S.: "
143 << Emitted4Momentum << G4endl;
144 }
145
146 // Perform Lorentz boost
147 G4LorentzVector Rest4Momentum = aFragment.GetMomentum();
148 Emitted4Momentum.boost(Rest4Momentum.boostVector());
149
150 // Set emitted fragment momentum
151 thePreFragment->SetMomentum(Emitted4Momentum);
152
153 // NOW THE RESIDUAL NUCLEUS
154 // ------------------------
155
156 Rest4Momentum -= Emitted4Momentum;
157
158 // Update nucleus parameters:
159 // --------------------------
160
161 // Z and A
162 aFragment.SetZandA_asInt(thePreFragment->GetRestZ(),
163 thePreFragment->GetRestA());
164
165 // Number of excitons
166 aFragment.SetNumberOfParticles(aFragment.GetNumberOfParticles()-
167 thePreFragment->GetA());
168 // Number of charges
169 aFragment.SetNumberOfCharged(aFragment.GetNumberOfCharged()-
170 thePreFragment->GetZ());
171
172 // Update nucleus momentum
173 // A check on consistence of Z, A, and mass will be performed
174 aFragment.SetMomentum(Rest4Momentum);
175
176 // Create a G4ReactionProduct
177 res = thePreFragment->GetReactionProduct();
178
179 // Set the creator model ID
180 aFragment.SetCreatorModelID(fModelID);
181 if (res != nullptr) { res->SetCreatorModelID(fModelID); }
182
183 return res;
184}
185
186void G4PreCompoundEmission::AngularDistribution(
187 G4VPreCompoundFragment* thePreFragment,
188 const G4Fragment& aFragment,
189 G4double ekin)
190{
191 G4int p = aFragment.GetNumberOfParticles();
192 G4int h = aFragment.GetNumberOfHoles();
193 G4double U = aFragment.GetExcitationEnergy();
194
195 // Emission particle separation energy
196 G4double Bemission = thePreFragment->GetBindingEnergy();
197
198 G4double gg = (6.0/pi2)*fNuclData->GetLevelDensity(aFragment.GetZ_asInt(),
199 aFragment.GetA_asInt(),U);
200
201 // Average exciton energy relative to bottom of nuclear well
202 G4double Eav = 2*p*(p+1)/((p+h)*gg);
203
204 // Excitation energy relative to the Fermi Level
205 G4double Uf = std::max(U - (p - h)*fFermiEnergy , 0.0);
206 // G4double Uf = U - KineticEnergyOfEmittedFragment - Bemission;
207
208 G4double w_num = rho(p+1, h, gg, Uf, fFermiEnergy);
209 G4double w_den = rho(p, h, gg, Uf, fFermiEnergy);
210 if (w_num > 0.0 && w_den > 0.0)
211 {
212 Eav *= (w_num/w_den);
213 Eav += - Uf/(p+h) + fFermiEnergy;
214 }
215 else
216 {
217 Eav = fFermiEnergy;
218 }
219
220 // VI + JMQ 19/01/2010 update computation of the parameter an
221 //
222 G4double an = 0.0;
223 G4double Eeff = ekin + Bemission + fFermiEnergy;
224 if(ekin > DBL_MIN && Eeff > DBL_MIN) {
225
226 G4double zeta = std::max(1.0,9.3/std::sqrt(ekin/CLHEP::MeV));
227
228 // This should be the projectile energy. If I would know which is
229 // the projectile (proton, neutron) I could remove the binding energy.
230 // But, what happens if INC precedes precompound? This approximation
231 // seems to work well enough
232 G4double ProjEnergy = aFragment.GetExcitationEnergy();
233
234 an = 3*std::sqrt((ProjEnergy+fFermiEnergy)*Eeff)/(zeta*Eav);
235
236 G4int ne = aFragment.GetNumberOfExcitons() - 1;
237 if ( ne > 1 ) { an /= static_cast<G4double>(ne); }
238
239 // protection of exponent
240 if ( an > 10. ) { an = 10.; }
241 }
242
243 // sample cosine of theta and not theta as in old versions
244 G4double random = G4UniformRand();
245 G4double cost;
246
247 if(an < 0.1) { cost = 1. - 2*random; }
248 else {
249 G4double exp2an = G4Exp(-2*an);
250 cost = 1. + G4Log(1-random*(1-exp2an))/an;
251 if(cost > 1.) { cost = 1.; }
252 else if(cost < -1.) {cost = -1.; }
253 }
254
255 G4double phi = CLHEP::twopi*G4UniformRand();
256
257 // Calculate the momentum magnitude of emitted fragment
258 G4double pmag =
259 std::sqrt(ekin*(ekin + 2.0*thePreFragment->GetNuclearMass()));
260
261 G4double sint = std::sqrt((1.0-cost)*(1.0+cost));
262
263 theFinalMomentum.set(pmag*std::cos(phi)*sint,pmag*std::sin(phi)*sint,
264 pmag*cost);
265
266 // theta is the angle wrt the incident direction
267 G4ThreeVector theIncidentDirection = aFragment.GetMomentum().vect().unit();
268 theFinalMomentum.rotateUz(theIncidentDirection);
269}
270
271G4double G4PreCompoundEmission::rho(G4int p, G4int h, G4double gg,
272 G4double E, G4double Ef) const
273{
274 // 25.02.2010 V.Ivanchenko added more protections
275 G4double Aph = (p*p + h*h + p - 3.0*h)/(4.0*gg);
276
277 if ( E - Aph < 0.0) { return 0.0; }
278
279 G4double logConst = (p+h)*G4Log(gg)
280 - g4calc->logfactorial(p+h-1) - g4calc->logfactorial(p)
281 - g4calc->logfactorial(h);
282
283 // initialise values using j=0
284
285 G4double t1=1;
286 G4double t2=1;
287 G4double logt3 = (p+h-1) * G4Log(E-Aph) + logConst;
288 const G4double logmax = 200.;
289 if(logt3 > logmax) { logt3 = logmax; }
290 G4double tot = G4Exp( logt3 );
291
292 // and now sum rest of terms
293 // 25.02.2010 V.Ivanchenko change while to for loop and cleanup
294 G4double Eeff = E - Aph;
295 for(G4int j=1; j<=h; ++j)
296 {
297 Eeff -= Ef;
298 if(Eeff < 0.0) { break; }
299 t1 *= -1.;
300 t2 *= static_cast<G4double>(h+1-j)/static_cast<G4double>(j);
301 logt3 = (p+h-1) * G4Log( Eeff) + logConst;
302 if(logt3 > logmax) { logt3 = logmax; }
303 tot += t1*t2*G4Exp(logt3);
304 }
305
306 return tot;
307}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:132
G4double G4Log(G4double x)
Definition G4Log.hh:169
CLHEP::HepLorentzVector G4LorentzVector
G4ThreeVector G4RandomDirection()
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
G4int GetNumberOfParticles() const
G4int GetNumberOfHoles() const
void SetZandA_asInt(G4int Znew, G4int Anew, G4int Lnew=0)
void SetNumberOfCharged(G4int value)
G4double GetExcitationEnergy() const
const G4LorentzVector & GetMomentum() const
void SetCreatorModelID(G4int value)
G4int GetZ_asInt() const
G4int GetNumberOfExcitons() const
void SetMomentum(const G4LorentzVector &value)
void SetNumberOfParticles(G4int value)
G4int GetNumberOfCharged() const
G4int GetA_asInt() const
G4double GetLevelDensity(G4int Z, G4int A, G4double U)
static G4NuclearLevelData * GetInstance()
static G4int GetModelID(const G4int modelIndex)
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4ReactionProduct * PerformEmission(G4Fragment &aFragment)
void SetCreatorModelID(const G4int mod)
virtual G4double SampleKineticEnergy(const G4Fragment &)
void SetMomentum(const G4LorentzVector &lv)
G4ReactionProduct * GetReactionProduct() const
#define DBL_MIN
Definition templates.hh:54