Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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 // $Id$
27 //
28 //J.M. Quesada (August2008). Based on:
29 //
30 // Hadronic Process: Nuclear De-excitations
31 // by V. Lara (Oct 1998)
32 //
33 // Modified:
34 // 03-09-2008 J.M. Quesada for external choice of inverse cross section option
35 // 06-09-2008 J.M. Quesada Also external choices have been added for superimposed
36 // Coulomb barrier (if useSICB is set true, by default is false)
37 // 17-11-2010 V.Ivanchenko in constructor replace G4VEmissionProbability by
38 // G4EvaporationProbability and do not new and delete probability
39 // object at each call; use G4Pow
40 
41 #include "G4EvaporationChannel.hh"
42 #include "G4PairingCorrection.hh"
43 #include "G4NucleiProperties.hh"
44 #include "G4Pow.hh"
46 #include "G4PhysicalConstants.hh"
47 #include "G4SystemOfUnits.hh"
48 #include "Randomize.hh"
49 #include "G4Alpha.hh"
50 
52  const G4String & aName,
53  G4EvaporationProbability* aEmissionStrategy,
54  G4VCoulombBarrier* aCoulombBarrier):
55  G4VEvaporationChannel(aName),
56  theA(anA),
57  theZ(aZ),
58  theEvaporationProbabilityPtr(aEmissionStrategy),
59  theCoulombBarrierPtr(aCoulombBarrier),
60  EmissionProbability(0.0),
61  MaximalKineticEnergy(-1000.0)
62 {
63  ResidualA = 0;
64  ResidualZ = 0;
65  ResidualMass = CoulombBarrier=0.0;
66  EvaporatedMass = G4NucleiProperties::GetNuclearMass(theA, theZ);
67  theLevelDensityPtr = new G4EvaporationLevelDensityParameter;
68 }
69 
72  theA(0),
73  theZ(0),
74  theEvaporationProbabilityPtr(0),
75  theCoulombBarrierPtr(0),
76  EmissionProbability(0.0),
77  MaximalKineticEnergy(-1000.0)
78 {
79  ResidualA = 0;
80  ResidualZ = 0;
81  EvaporatedMass = ResidualMass = CoulombBarrier = 0.0;
82  theLevelDensityPtr = new G4EvaporationLevelDensityParameter;
83 }
84 
86 {
87  delete theLevelDensityPtr;
88 }
89 
90 //void G4EvaporationChannel::Initialize(const G4Fragment &)
91 //{}
92 
94 {
95  //for inverse cross section choice
96  theEvaporationProbabilityPtr->SetOPTxs(OPTxs);
97  // for superimposed Coulomb Barrier for inverse cross sections
98  theEvaporationProbabilityPtr->UseSICB(useSICB);
99 
100  G4int FragmentA = fragment->GetA_asInt();
101  G4int FragmentZ = fragment->GetZ_asInt();
102  ResidualA = FragmentA - theA;
103  ResidualZ = FragmentZ - theZ;
104  //G4cout << "G4EvaporationChannel::Initialize Z= " << theZ << " A= " << theA
105  // << " FragZ= " << FragmentZ << " FragA= " << FragmentA << G4endl;
106 
107  //Effective excitation energy
108  G4double ExEnergy = fragment->GetExcitationEnergy() -
110 
111  // Only channels which are physically allowed are taken into account
112  if (ResidualA <= 0 || ResidualZ <= 0 || ResidualA < ResidualZ ||
113  (ResidualA == ResidualZ && ResidualA > 1) || ExEnergy <= 0.0) {
114  CoulombBarrier = ResidualMass = 0.0;
115  MaximalKineticEnergy = -1000.0*MeV;
116  EmissionProbability = 0.0;
117  } else {
118  ResidualMass = G4NucleiProperties::GetNuclearMass(ResidualA, ResidualZ);
119  G4double FragmentMass = fragment->GetGroundStateMass();
120  CoulombBarrier = theCoulombBarrierPtr->GetCoulombBarrier(ResidualA,ResidualZ,ExEnergy);
121  // Maximal Kinetic Energy
122  MaximalKineticEnergy = CalcMaximalKineticEnergy(FragmentMass + ExEnergy);
123  //MaximalKineticEnergy = ExEnergy + fragment.GetGroundStateMass()
124  // - EvaporatedMass - ResidualMass;
125 
126  // Emission probability
127  // Protection for the case Tmax<V. If not set in this way we could end up in an
128  // infinite loop in the method GetKineticEnergy if OPTxs!=0 && useSICB=true.
129  // Of course for OPTxs=0 we have the Coulomb barrier
130 
131  G4double limit;
132  if (OPTxs==0 || (OPTxs!=0 && useSICB))
133  limit= CoulombBarrier;
134  else limit=0.;
135  limit = std::max(limit, FragmentMass - ResidualMass - EvaporatedMass);
136 
137  // The threshold for charged particle emission must be set to 0 if Coulomb
138  //cutoff is included in the cross sections
139  if (MaximalKineticEnergy <= limit) EmissionProbability = 0.0;
140  else {
141  // Total emission probability for this channel
142  EmissionProbability = theEvaporationProbabilityPtr->
143  EmissionProbability(*fragment, MaximalKineticEnergy);
144  }
145  }
146  //G4cout << "G4EvaporationChannel:: probability= " << EmissionProbability << G4endl;
147 
148  return EmissionProbability;
149 }
150 
152 {
153  /*
154  G4double Ecm = GetKineticEnergy(theNucleus) + ResidualMass + EvaporatedMass;
155 
156  G4double EvaporatedEnergy =
157  ((Ecm-ResidualMass)*(Ecm+ResidualMass) + EvaporatedMass*EvaporatedMass)/(2*Ecm);
158  */
159  G4double EvaporatedEnergy = GetKineticEnergy(theNucleus) + EvaporatedMass;
160 
161  G4ThreeVector momentum(IsotropicVector
162  (std::sqrt((EvaporatedEnergy - EvaporatedMass)*
163  (EvaporatedEnergy + EvaporatedMass))));
164 
165  G4LorentzVector EvaporatedMomentum(momentum,EvaporatedEnergy);
166  G4LorentzVector ResidualMomentum = theNucleus.GetMomentum();
167  EvaporatedMomentum.boost(ResidualMomentum.boostVector());
168 
169  G4Fragment * EvaporatedFragment = new G4Fragment(theA,theZ,EvaporatedMomentum);
170  ResidualMomentum -= EvaporatedMomentum;
171 
172  G4Fragment * ResidualFragment = new G4Fragment(ResidualA, ResidualZ, ResidualMomentum);
173 
174  G4FragmentVector * theResult = new G4FragmentVector;
175 
176 #ifdef debug
177  G4double Efinal = ResidualMomentum.e() + EvaporatedMomentum.e();
178  G4ThreeVector Pfinal = ResidualMomentum.vect() + EvaporatedMomentum.vect();
179  if (std::abs(Efinal-theNucleus.GetMomentum().e()) > 1.0*keV) {
180  G4cout << "@@@@@@@@@@@@@@@@@@@@@ G4Evaporation Chanel: ENERGY @@@@@@@@@@@@@@@@" << G4endl;
181  G4cout << "Initial : " << theNucleus.GetMomentum().e()/MeV << " MeV Final :"
182  <<Efinal/MeV << " MeV Delta: " << (Efinal-theNucleus.GetMomentum().e())/MeV
183  << " MeV" << G4endl;
184  }
185  if (std::abs(Pfinal.x()-theNucleus.GetMomentum().x()) > 1.0*keV ||
186  std::abs(Pfinal.y()-theNucleus.GetMomentum().y()) > 1.0*keV ||
187  std::abs(Pfinal.z()-theNucleus.GetMomentum().z()) > 1.0*keV ) {
188  G4cout << "@@@@@@@@@@@@@@@@@@@@@ G4Evaporation Chanel: MOMENTUM @@@@@@@@@@@@@@@@" << G4endl;
189  G4cout << "Initial : " << theNucleus.GetMomentum().vect() << " MeV Final :"
190  <<Pfinal/MeV << " MeV Delta: " << Pfinal-theNucleus.GetMomentum().vect()
191  << " MeV" << G4endl;
192 
193  }
194 #endif
195  theResult->push_back(EvaporatedFragment);
196  theResult->push_back(ResidualFragment);
197  return theResult;
198 }
199 
201 // Calculates the maximal kinetic energy that can be carried by fragment.
202 G4double G4EvaporationChannel::CalcMaximalKineticEnergy(G4double NucleusTotalE)
203 {
204  // This is the "true" assimptotic kinetic energy (from energy conservation)
205  G4double Tmax =
206  ((NucleusTotalE-ResidualMass)*(NucleusTotalE+ResidualMass) + EvaporatedMass*EvaporatedMass)
207  /(2.0*NucleusTotalE) - EvaporatedMass;
208 
209  //JMQ (13-09-08) bug fixed: in the original version the Tmax is calculated
210  //at the Coulomb barrier
211  //IMPORTANT: meaning of Tmax differs in OPTxs=0 and OPTxs!=0
212  //When OPTxs!=0 Tmax is the TRUE (assimptotic) maximal kinetic energy
213 
214  if(OPTxs==0)
215  Tmax=Tmax- CoulombBarrier;
216 
217  return Tmax;
218 }
219 
221 //JMQ: New method for MC sampling of kinetic energy. Substitutes old CalcKineticEnergy
222 G4double G4EvaporationChannel::GetKineticEnergy(const G4Fragment & aFragment)
223 {
224 
225  if (OPTxs==0) {
226  // It uses Dostrovsky's approximation for the inverse reaction cross
227  // in the probability for fragment emission
228  // MaximalKineticEnergy energy in the original version (V.Lara) was calculated at
229  //the Coulomb barrier.
230 
231 
232  if (MaximalKineticEnergy < 0.0) {
233  throw G4HadronicException(__FILE__, __LINE__,
234  "G4EvaporationChannel::CalcKineticEnergy: maximal kinetic at the Coulomb barrier is less than 0");
235  }
236  G4double Rb = 4.0*theLevelDensityPtr->
237  LevelDensityParameter(ResidualA+theA,ResidualZ+theZ,MaximalKineticEnergy)*
238  MaximalKineticEnergy;
239  G4double RbSqrt = std::sqrt(Rb);
240  G4double PEX1 = 0.0;
241  if (RbSqrt < 160.0) PEX1 = std::exp(-RbSqrt);
242  G4double Rk = 0.0;
243  G4double FRk = 0.0;
244  do {
245  G4double RandNumber = G4UniformRand();
246  Rk = 1.0 + (1./RbSqrt)*std::log(RandNumber + (1.0-RandNumber)*PEX1);
247  G4double Q1 = 1.0;
248  G4double Q2 = 1.0;
249  if (theZ == 0) { // for emitted neutron
250  G4double Beta = (2.12/G4Pow::GetInstance()->Z23(ResidualA) - 0.05)*MeV/
251  (0.76 + 2.2/G4Pow::GetInstance()->Z13(ResidualA));
252  Q1 = 1.0 + Beta/(MaximalKineticEnergy);
253  Q2 = Q1*std::sqrt(Q1);
254  }
255 
256  FRk = (3.0*std::sqrt(3.0)/2.0)/Q2 * Rk * (Q1 - Rk*Rk);
257 
258  } while (FRk < G4UniformRand());
259 
260  G4double result = MaximalKineticEnergy * (1.0-Rk*Rk) + CoulombBarrier;
261  return result;
262  } else if (OPTxs==1 || OPTxs==2 || OPTxs==3 || OPTxs==4) {
263 
264  // Coulomb barrier is just included in the cross sections
265  G4double V = 0;
266  if(useSICB) { V= CoulombBarrier; }
267 
268  V = std::max(V, aFragment.GetGroundStateMass()-EvaporatedMass-ResidualMass);
269 
270  G4double Tmax=MaximalKineticEnergy;
271  G4double T(0.0);
272  G4double NormalizedProbability(1.0);
273 
274  // VI: This is very ineffective - create new objects at each call to the method
275  /*
276  // A pointer is created in order to access the distribution function.
277  G4EvaporationProbability * G4EPtemp = 0;
278 
279  if (theA==1 && theZ==0) G4EPtemp=new G4NeutronEvaporationProbability();
280  else if (theA==1 && theZ==1) G4EPtemp=new G4ProtonEvaporationProbability();
281  else if (theA==2 && theZ==1 ) G4EPtemp=new G4DeuteronEvaporationProbability();
282  else if (theA==3 && theZ==1 ) G4EPtemp=new G4TritonEvaporationProbability();
283  else if (theA==3 && theZ==2 ) G4EPtemp=new G4He3EvaporationProbability();
284  else if (theA==4 && theZ==2) G4EPtemp=new G4AlphaEvaporationProbability();
285  else {
286  std::ostringstream errOs;
287  errOs << "ejected particle out of range in G4EvaporationChannel" << G4endl;
288  throw G4HadronicException(__FILE__, __LINE__, errOs.str());
289  }
290 
291  //for cross section selection and superimposed Coulom Barrier for xs
292  G4EPtemp->SetOPTxs(OPTxs);
293  G4EPtemp->UseSICB(useSICB);
294  */
295 
296  // use local pointer and not create a new one
297  do
298  {
299  T=V+G4UniformRand()*(Tmax-V);
300  NormalizedProbability =
301  theEvaporationProbabilityPtr->ProbabilityDistributionFunction(aFragment,T)/
302  GetEmissionProbability(const_cast<G4Fragment*>(&aFragment));
303 
304  }
305  while (G4UniformRand() > NormalizedProbability);
306  // delete G4EPtemp;
307  return T;
308  } else{
309  std::ostringstream errOs;
310  errOs << "Bad option for energy sampling in evaporation" <<G4endl;
311  throw G4HadronicException(__FILE__, __LINE__, errOs.str());
312  }
313 }
314 
315 G4ThreeVector G4EvaporationChannel::IsotropicVector(G4double Magnitude)
316  // Samples a isotropic random vectorwith a magnitud given by Magnitude.
317  // By default Magnitude = 1.0
318 {
319  G4double CosTheta = 1.0 - 2.0*G4UniformRand();
320  G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta);
321  G4double Phi = twopi*G4UniformRand();
322  G4ThreeVector Vector(Magnitude*std::cos(Phi)*SinTheta,
323  Magnitude*std::sin(Phi)*SinTheta,
324  Magnitude*CosTheta);
325  return Vector;
326 }
327 
328 
329 
330 
331 
332 
333 
334 
335 
336