Geant4  10.02
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: G4PreCompoundEmission.cc 90337 2015-05-26 08:34:27Z gcosmo $
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 
45 #include "G4PreCompoundEmission.hh"
46 #include "G4PhysicalConstants.hh"
47 #include "G4SystemOfUnits.hh"
48 #include "G4Pow.hh"
49 #include "G4Exp.hh"
50 #include "G4Log.hh"
51 #include "Randomize.hh"
54 #include "G4HETCEmissionFactory.hh"
55 #include "G4HadronicException.hh"
56 
58 {
65  fFermiEnergy = param.GetFermiEnergy();
66 }
67 
69 {
70  delete theFragmentsFactory;
71  delete theFragmentsVector;
72 }
73 
75 {
78  if (theFragmentsVector) {
80  } else {
83  }
84 }
85 
87 {
90  if (theFragmentsVector) {
92  } else {
95  }
96 }
97 
100 {
101  // Choose a Fragment for emission
102  G4VPreCompoundFragment * thePreFragment =
104  if (thePreFragment == 0)
105  {
106  G4cout << "G4PreCompoundEmission::PerformEmission : "
107  << "I couldn't choose a fragment\n"
108  << "while trying to de-excite\n"
109  << aFragment << G4endl;
110  throw G4HadronicException(__FILE__, __LINE__, "");
111  }
112 
113  //G4cout << "Chosen fragment: " << G4endl;
114  //G4cout << *thePreFragment << G4endl;
115 
116  // Kinetic Energy of emitted fragment
117  G4double kinEnergy = thePreFragment->SampleKineticEnergy(aFragment);
118  // if(kinEnergy < MeV) {
119  // G4cout << "Chosen fragment: " << G4endl;
120  // G4cout << *thePreFragment << G4endl;
121  // G4cout << "Ekin= " << kinEnergy << G4endl;
122  // }
123  if(kinEnergy < 0.0) { kinEnergy = 0.0; }
124 
125  // Calculate the fragment momentum (three vector)
126  AngularDistribution(thePreFragment,aFragment,kinEnergy);
127 
128  // Mass of emittef fragment
129  G4double EmittedMass = thePreFragment->GetNuclearMass();
130  // Now we can calculate the four momentum
131  // both options are valid and give the same result but 2nd one is faster
132  G4LorentzVector Emitted4Momentum(theFinalMomentum,EmittedMass + kinEnergy);
133 
134  // Perform Lorentz boost
135  G4LorentzVector Rest4Momentum = aFragment.GetMomentum();
136  Emitted4Momentum.boost(Rest4Momentum.boostVector());
137 
138  // Set emitted fragment momentum
139  thePreFragment->SetMomentum(Emitted4Momentum);
140 
141  // NOW THE RESIDUAL NUCLEUS
142  // ------------------------
143 
144  Rest4Momentum -= Emitted4Momentum;
145 
146  // Update nucleus parameters:
147  // --------------------------
148 
149  // Z and A
150  aFragment.SetZandA_asInt(thePreFragment->GetRestZ(),
151  thePreFragment->GetRestA());
152 
153  // Number of excitons
154  aFragment.SetNumberOfParticles(aFragment.GetNumberOfParticles()-
155  thePreFragment->GetA());
156  // Number of charges
157  aFragment.SetNumberOfCharged(aFragment.GetNumberOfCharged()-
158  thePreFragment->GetZ());
159 
160  // Update nucleus momentum
161  // A check on consistence of Z, A, and mass will be performed
162  aFragment.SetMomentum(Rest4Momentum);
163 
164  // Create a G4ReactionProduct
165  G4ReactionProduct * MyRP = thePreFragment->GetReactionProduct();
166 
167  // if(kinEnergy < MeV) {
168  // G4cout << "G4PreCompoundEmission::Fragment emitted" << G4endl;
169  // G4cout << *thePreFragment << G4endl;
170  // }
171  return MyRP;
172 }
173 
175  G4VPreCompoundFragment* thePreFragment,
176  const G4Fragment& aFragment,
177  G4double ekin)
178 {
179  G4int p = aFragment.GetNumberOfParticles();
180  G4int h = aFragment.GetNumberOfHoles();
181  G4double U = aFragment.GetExcitationEnergy();
182 
183  // Emission particle separation energy
184  G4double Bemission = thePreFragment->GetBindingEnergy();
185 
186  G4double gg = (6.0/pi2)*aFragment.GetA_asInt()*fLevelDensity;
187 
188  // Average exciton energy relative to bottom of nuclear well
189  G4double Eav = 2*p*(p+1)/((p+h)*gg);
190 
191  // Excitation energy relative to the Fermi Level
192  G4double Uf = std::max(U - (p - h)*fFermiEnergy , 0.0);
193  // G4double Uf = U - KineticEnergyOfEmittedFragment - Bemission;
194 
195  G4double w_num = rho(p+1, h, gg, Uf, fFermiEnergy);
196  G4double w_den = rho(p, h, gg, Uf, fFermiEnergy);
197  if (w_num > 0.0 && w_den > 0.0)
198  {
199  Eav *= (w_num/w_den);
200  Eav += - Uf/(p+h) + fFermiEnergy;
201  }
202  else
203  {
204  Eav = fFermiEnergy;
205  }
206 
207  // VI + JMQ 19/01/2010 update computation of the parameter an
208  //
209  G4double an = 0.0;
210  G4double Eeff = ekin + Bemission + fFermiEnergy;
211  if(ekin > DBL_MIN && Eeff > DBL_MIN) {
212 
213  G4double zeta = std::max(1.0,9.3/std::sqrt(ekin/MeV));
214 
215  // This should be the projectile energy. If I would know which is
216  // the projectile (proton, neutron) I could remove the binding energy.
217  // But, what happens if INC precedes precompound? This approximation
218  // seems to work well enough
219  G4double ProjEnergy = aFragment.GetExcitationEnergy();
220 
221  an = 3*std::sqrt((ProjEnergy+fFermiEnergy)*Eeff)/(zeta*Eav);
222 
223  G4int ne = aFragment.GetNumberOfExcitons() - 1;
224  if ( ne > 1 ) { an /= (G4double)ne; }
225 
226  // protection of exponent
227  if ( an > 10. ) { an = 10.; }
228  }
229 
230  // sample cosine of theta and not theta as in old versions
231  G4double random = G4UniformRand();
232  G4double cost;
233 
234  if(an < 0.1) { cost = 1. - 2*random; }
235  else {
236  G4double exp2an = G4Exp(-2*an);
237  cost = 1. + G4Log(1-random*(1-exp2an))/an;
238  if(cost > 1.) { cost = 1.; }
239  else if(cost < -1.) {cost = -1.; }
240  }
241 
243 
244  // Calculate the momentum magnitude of emitted fragment
245  G4double pmag =
246  std::sqrt(ekin*(ekin + 2.0*thePreFragment->GetNuclearMass()));
247 
248  G4double sint = std::sqrt((1.0-cost)*(1.0+cost));
249 
250  theFinalMomentum.set(pmag*std::cos(phi)*sint,pmag*std::sin(phi)*sint,
251  pmag*cost);
252 
253  // theta is the angle wrt the incident direction
254  G4ThreeVector theIncidentDirection = aFragment.GetMomentum().vect().unit();
255  theFinalMomentum.rotateUz(theIncidentDirection);
256 }
257 
259  G4double E, G4double Ef) const
260 {
261  // 25.02.2010 V.Ivanchenko added more protections
262  G4double Aph = (p*p + h*h + p - 3.0*h)/(4.0*gg);
263  // G4double alpha = (p*p + h*h)/(2.0*gg);
264 
265  if ( E - Aph < 0.0) { return 0.0; }
266 
267  G4double logConst = (p+h)*G4Log(gg)
268  - g4pow->logfactorial(p+h-1) - g4pow->logfactorial(p)
269  - g4pow->logfactorial(h);
270 
271  // initialise values using j=0
272 
273  G4double t1=1;
274  G4double t2=1;
275  G4double logt3 = (p+h-1) * G4Log(E-Aph) + logConst;
276  const G4double logmax = 200.;
277  if(logt3 > logmax) { logt3 = logmax; }
278  G4double tot = G4Exp( logt3 );
279 
280  // and now sum rest of terms
281  // 25.02.2010 V.Ivanchenko change while to for loop and cleanup
282  G4double Eeff = E - Aph;
283  for(G4int j=1; j<=h; ++j)
284  {
285  Eeff -= Ef;
286  if(Eeff < 0.0) { break; }
287  t1 *= -1.;
288  t2 *= (G4double)(h+1-j)/(G4double)j;
289  logt3 = (p+h-1) * G4Log( Eeff) + logConst;
290  if(logt3 > logmax) { logt3 = logmax; }
291  tot += t1*t2*G4Exp(logt3);
292  }
293 
294  return tot;
295 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4ReactionProduct * GetReactionProduct() const
static const double MeV
Definition: G4SIunits.hh:211
CLHEP::Hep3Vector G4ThreeVector
std::vector< G4VPreCompoundFragment * > * GetFragmentVector()
void AngularDistribution(G4VPreCompoundFragment *theFragment, const G4Fragment &aFragment, G4double KineticEnergy)
G4int GetA() const
G4int GetNumberOfParticles() const
Definition: G4Fragment.hh:322
G4double rho(G4int p, G4int h, G4double gg, G4double E, G4double Ef) const
static const double pi2
Definition: G4SIunits.hh:77
int G4int
Definition: G4Types.hh:78
G4double GetBindingEnergy() const
G4int GetNumberOfHoles() const
Definition: G4Fragment.hh:342
void SetMomentum(const G4LorentzVector &value)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4int GetA_asInt() const
Definition: G4Fragment.hh:251
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:284
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:289
static const double twopi
Definition: G4SIunits.hh:75
G4int GetRestZ() const
G4PreCompoundFragmentVector * theFragmentsVector
void SetNumberOfParticles(G4int value)
Definition: G4Fragment.hh:361
G4VPreCompoundEmissionFactory * theFragmentsFactory
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4VPreCompoundFragment * ChooseFragment()
G4double GetNuclearMass() const
G4ReactionProduct * PerformEmission(G4Fragment &aFragment)
G4int GetNumberOfExcitons() const
Definition: G4Fragment.hh:317
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double logfactorial(G4int Z) const
Definition: G4Pow.hh:269
void SetZandA_asInt(G4int Znew, G4int Anew)
Definition: G4Fragment.hh:261
virtual G4double SampleKineticEnergy(const G4Fragment &aFragment)=0
#define DBL_MIN
Definition: templates.hh:75
G4int GetRestA() const
#define G4endl
Definition: G4ios.hh:61
G4int GetZ() const
void SetNumberOfCharged(G4int value)
Definition: G4Fragment.hh:366
double G4double
Definition: G4Types.hh:76
G4int GetNumberOfCharged() const
Definition: G4Fragment.hh:327
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:268
CLHEP::HepLorentzVector G4LorentzVector