Geant4-11
G4PreCompoundFragment.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// J. M. Quesada (August 2008).
28// Based on previous work by V. Lara
29//
30// Modified:
31// 06.09.2008 JMQ Also external choice has been added for:
32// - superimposed Coulomb barrier (if useSICB=true)
33// 20.08.2010 V.Ivanchenko cleanup
34//
35
40#include "Randomize.hh"
41
43 G4VCoulombBarrier* aCoulBarrier)
44 : G4VPreCompoundFragment(p, aCoulBarrier)
45{
46 muu = probmax = 0.0;
47 if(0 == theZ) { index = 0; }
48 else if(1 == theZ) { index = theA; }
49 else { index = theA + 1; }
50}
51
53{}
54
56CalcEmissionProbability(const G4Fragment & aFragment)
57{
58 //G4cout << theCoulombBarrier << " " << GetMaximalKineticEnergy() << G4endl;
59 // If theCoulombBarrier effect is included in the emission probabilities
60 // Coulomb barrier is the lower limit of integration over kinetic energy
61
63
64 if (theMaxKinEnergy <= theMinKinEnergy) { return 0.0; }
65
66 // compute power once
67 if(0 < index) {
69 }
70
73 /*
74 G4cout << "## G4PreCompoundFragment::CalcEmisProb "
75 << "Z= " << aFragment.GetZ_asInt()
76 << " A= " << aFragment.GetA_asInt()
77 << " Elow= " << LowerLimit/MeV
78 << " Eup= " << UpperLimit/MeV
79 << " prob= " << theEmissionProbability
80 << G4endl;
81 */
83}
84
87 const G4Fragment & aFragment)
88{
89 static const G4double den = 1.0/CLHEP::MeV;
90 G4double del = (up - low);
91 G4int nbins = del*den;
92 nbins = std::max(nbins, 4);
93 del /= static_cast<G4double>(nbins);
94 G4double e = low + 0.5*del;
96 //G4cout << " 0. e= " << e << " y= " << probmax << G4endl;
97
98 G4double sum = probmax;
99 for (G4int i=1; i<nbins; ++i) {
100 e += del;
101
102 G4double y = ProbabilityDistributionFunction(e, aFragment);
104 sum += y;
105 if(y < sum*0.01) { break; }
106 //G4cout <<" "<<i<<". e= "<<e<<" y= "<<y<<" sum= "<<sum<< G4endl;
107 }
108 sum *= del;
109 //G4cout << "Evap prob: " << sum << " probmax= " << probmax << G4endl;
110 return sum;
111}
112
114{
115 G4double res;
116 if(OPTxs == 0 || (OPTxs == 4 && theMaxKinEnergy < 10.)) {
117 res = GetOpt0(ekin);
118
119 } else if(OPTxs <= 2) {
122 theResA13, muu,
123 index, theZ, theResA);
124
125 } else {
128 theZ, theA, theResA);
129 }
130 return res;
131}
132
134// OPT=0 : Dostrovski's cross section
135{
136 G4double r0 = theParameters->GetR0()*theResA13;
137 // cross section is now given in mb (r0 is in mm) for the sake of consistency
138 //with the rest of the options
139 return 1.e+25*CLHEP::pi*r0*r0*theResA13*GetAlpha()*(1.0 + GetBeta()/ekin);
140}
141
143{
145 static const G4double toler = 1.25;
146 probmax *= toler;
147 G4double prob, T(0.0);
148 CLHEP::HepRandomEngine* rndm = G4Random::getTheEngine();
149 G4int i;
150 for(i=0; i<100; ++i) {
151 T = theMinKinEnergy + delta*rndm->flat();
152 prob = ProbabilityDistributionFunction(T, fragment);
153 /*
154 if(prob > probmax) {
155 G4cout << "G4PreCompoundFragment WARNING: prob= " << prob
156 << " probmax= " << probmax << G4endl;
157 G4cout << "i= " << i << " Z= " << theZ << " A= " << theA
158 << " resZ= " << theResZ << " resA= " << theResA << "\n"
159 << " T= " << T << " Tmax= " << theMaxKinEnergy
160 << " Tmin= " << limit
161 << G4endl;
162 for(G4int i=0; i<N; ++i) { G4cout << " " << probability[i]; }
163 G4cout << G4endl;
164 }
165 */
166 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
167 if(probmax*rndm->flat() <= prob) { break; }
168 }
169 /*
170 G4cout << "G4PreCompoundFragment: i= " << i << " Z= " << theZ
171 << " A= " << theA <<" T(MeV)= " << T << " Emin(MeV)= "
172 << theMinKinEnergy << " Emax= " << theMaxKinEnergy << G4endl;
173 */
174 return T;
175}
176
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
virtual double flat()=0
static G4double ComputeCrossSection(G4double K, G4double cb, G4double resA13, G4double amu1, G4int idx, G4int Z, G4int resA)
static G4double ComputePowerParameter(G4int resA, G4int idx)
static G4double ComputeCrossSection(G4double K, G4double cb, G4double resA13, G4double amu1, G4int idx, G4int Z, G4int A, G4int resA)
G4double CalcEmissionProbability(const G4Fragment &aFragment) override
G4double GetOpt0(G4double ekin) const
G4double IntegrateEmissionProbability(G4double low, G4double up, const G4Fragment &aFragment)
G4PreCompoundFragment(const G4ParticleDefinition *, G4VCoulombBarrier *aCoulombBarrier)
G4double CrossSection(G4double ekin) const
G4double SampleKineticEnergy(const G4Fragment &aFragment) override
virtual G4double ProbabilityDistributionFunction(G4double ekin, const G4Fragment &aFragment)=0
static constexpr double MeV
static constexpr double pi
Definition: SystemOfUnits.h:55
T max(const T t1, const T t2)
brief Return the largest of the two arguments