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