Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4GEMChannel.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// Hadronic Process: Nuclear De-excitations
29// by V. Lara (Oct 1998)
30//
31// J. M. Quesada (September 2009) bugs fixed in probability distribution for kinetic
32// energy sampling:
33// -hbarc instead of hbar_Planck (BIG BUG)
34// -quantities for initial nucleus and residual are calculated separately
35// V.Ivanchenko (September 2009) Added proper protection for unphysical final state
36// J. M. Quesada (October 2009) fixed bug in CoulombBarrier calculation
37
38#include "G4GEMChannel.hh"
41#include "G4SystemOfUnits.hh"
42#include "G4Pow.hh"
43
45 G4GEMProbability * aEmissionStrategy,
46 G4VCoulombBarrier * aCoulombBarrier) :
48 A(theA),
49 Z(theZ),
50 theEvaporationProbabilityPtr(aEmissionStrategy),
51 theCoulombBarrierPtr(aCoulombBarrier),
52 EmissionProbability(0.0),
53 MaximalKineticEnergy(-CLHEP::GeV)
54{
55 theLevelDensityPtr = new G4EvaporationLevelDensityParameter;
56 MyOwnLevelDensity = true;
57 EvaporatedMass = G4NucleiProperties::GetNuclearMass(A, Z);
58 ResidualMass = CoulombBarrier = 0.0;
59 fG4pow = G4Pow::GetInstance();
60 ResidualZ = ResidualA = 0;
61}
62
64{
65 if (MyOwnLevelDensity) { delete theLevelDensityPtr; }
66}
67
69{
70 G4int anA = fragment->GetA_asInt();
71 G4int aZ = fragment->GetZ_asInt();
72 ResidualA = anA - A;
73 ResidualZ = aZ - Z;
74 //G4cout << "G4GEMChannel::Initialize: Z= " << aZ << " A= " << anA
75 // << " Zres= " << ResidualZ << " Ares= " << ResidualA << G4endl;
76
77 // We only take into account channels which are physically allowed
78 if (ResidualA <= 0 || ResidualZ <= 0 || ResidualA < ResidualZ ||
79 (ResidualA == ResidualZ && ResidualA > 1))
80 {
81 CoulombBarrier = 0.0;
82 MaximalKineticEnergy = -CLHEP::GeV;
83 EmissionProbability = 0.0;
84 }
85 else
86 {
87 // Effective excitation energy
88 // JMQ 071009: pairing in ExEnergy should be the one of parent compound nucleus
89 // FIXED the bug causing reported crash by VI (negative Probabilities
90 // due to inconsistency in Coulomb barrier calculation (CoulombBarrier and -Beta
91 // param for protons must be the same)
92 // G4double ExEnergy = fragment.GetExcitationEnergy() -
93 // G4PairingCorrection::GetInstance()->GetPairingCorrection(ResidualA,ResidualZ);
94 G4double ExEnergy = fragment->GetExcitationEnergy() -
96
97 //G4cout << "Eexc(MeV)= " << ExEnergy/MeV << G4endl;
98
99 if( ExEnergy <= 0.0) {
100 CoulombBarrier = 0.0;
101 MaximalKineticEnergy = -1000.0*MeV;
102 EmissionProbability = 0.0;
103
104 } else {
105
106 ResidualMass = G4NucleiProperties::GetNuclearMass(ResidualA, ResidualZ);
107
108 // Coulomb Barrier calculation
109 CoulombBarrier = theCoulombBarrierPtr->GetCoulombBarrier(ResidualA,ResidualZ,ExEnergy);
110 //G4cout << "CBarrier(MeV)= " << CoulombBarrier/MeV << G4endl;
111
112 //Maximal kinetic energy (JMQ : at the Coulomb barrier)
113 MaximalKineticEnergy =
114 CalcMaximalKineticEnergy(fragment->GetGroundStateMass()+ExEnergy);
115 //G4cout << "MaxE(MeV)= " << MaximalKineticEnergy/MeV << G4endl;
116
117 // Emission probability
118 if (MaximalKineticEnergy <= 0.0)
119 {
120 EmissionProbability = 0.0;
121 }
122 else
123 {
124 // Total emission probability for this channel
125 EmissionProbability =
126 theEvaporationProbabilityPtr->EmissionProbability(*fragment,
127 MaximalKineticEnergy);
128 }
129 }
130 }
131 //G4cout << "Prob= " << EmissionProbability << G4endl;
132 return EmissionProbability;
133}
134
136{
137 G4double EvaporatedKineticEnergy = CalcKineticEnergy(theNucleus);
138 G4double EvaporatedEnergy = EvaporatedKineticEnergy + EvaporatedMass;
139
140 G4ThreeVector momentum(IsotropicVector(std::sqrt(EvaporatedKineticEnergy*
141 (EvaporatedKineticEnergy+2.0*EvaporatedMass))));
142
143 momentum.rotateUz(theNucleus.GetMomentum().vect().unit());
144
145 G4LorentzVector EvaporatedMomentum(momentum,EvaporatedEnergy);
146 EvaporatedMomentum.boost(theNucleus.GetMomentum().boostVector());
147 G4Fragment * EvaporatedFragment = new G4Fragment(A,Z,EvaporatedMomentum);
148 // ** And now the residual nucleus **
149 G4double theExEnergy = theNucleus.GetExcitationEnergy();
150 G4double theMass = theNucleus.GetGroundStateMass();
151 G4double ResidualEnergy =
152 theMass + (theExEnergy - EvaporatedKineticEnergy) - EvaporatedMass;
153
154 G4LorentzVector ResidualMomentum(-momentum,ResidualEnergy);
155 ResidualMomentum.boost(theNucleus.GetMomentum().boostVector());
156
157 G4Fragment * ResidualFragment = new G4Fragment( ResidualA, ResidualZ, ResidualMomentum );
158
159 G4FragmentVector * theResult = new G4FragmentVector;
160
161 theResult->push_back(EvaporatedFragment);
162 theResult->push_back(ResidualFragment);
163 return theResult;
164}
165
166G4double G4GEMChannel::CalcMaximalKineticEnergy(const G4double NucleusTotalE)
167// Calculate maximal kinetic energy that can be carried by fragment.
168//JMQ this is not the assimptotic kinetic energy but the one at the Coulomb barrier
169{
170 G4double T = (NucleusTotalE*NucleusTotalE + EvaporatedMass*EvaporatedMass - ResidualMass*ResidualMass)/
171 (2.0*NucleusTotalE) - EvaporatedMass - CoulombBarrier;
172
173 return T;
174}
175
176G4double G4GEMChannel::CalcKineticEnergy(const G4Fragment & fragment)
177// Samples fragment kinetic energy.
178{
179 G4double U = fragment.GetExcitationEnergy();
180
181 G4double Alpha = theEvaporationProbabilityPtr->CalcAlphaParam(fragment);
182 G4double Beta = theEvaporationProbabilityPtr->CalcBetaParam(fragment);
183
184 G4double Normalization = theEvaporationProbabilityPtr->GetNormalization();
185
186 // ***RESIDUAL***
187 //JMQ (September 2009) the following quantities refer to the RESIDUAL:
188 G4double delta0 = G4PairingCorrection::GetInstance()->GetPairingCorrection(ResidualA,ResidualZ);
189 G4double Ux = (2.5 + 150.0/ResidualA)*MeV;
190 G4double Ex = Ux + delta0;
191 G4double InitialLevelDensity;
192 // ***end RESIDUAL ***
193
194 // ***PARENT***
195 //JMQ (September 2009) the following quantities refer to the PARENT:
196
197 G4double deltaCN =
199 fragment.GetZ_asInt());
200 G4double aCN = theLevelDensityPtr->LevelDensityParameter(fragment.GetA_asInt(),
201 fragment.GetZ_asInt(),U-deltaCN);
202 G4double UxCN = (2.5 + 150.0/fragment.GetA())*MeV;
203 G4double ExCN = UxCN + deltaCN;
204 G4double TCN = 1.0/(std::sqrt(aCN/UxCN) - 1.5/UxCN);
205 // ***end PARENT***
206
207 //JMQ quantities calculated for CN in InitialLevelDensity
208 if ( U < ExCN )
209 {
210 G4double E0CN = ExCN - TCN*(std::log(TCN/MeV) - std::log(aCN*MeV)/4.0
211 - 1.25*std::log(UxCN/MeV) + 2.0*std::sqrt(aCN*UxCN));
212 InitialLevelDensity = (pi/12.0)*std::exp((U-E0CN)/TCN)/TCN;
213 }
214 else
215 {
216 G4double x = U-deltaCN;
217 G4double x1 = std::sqrt(aCN*x);
218 InitialLevelDensity = (pi/12.0)*std::exp(2*x1)/(x*std::sqrt(x1));
219 //InitialLevelDensity =
220 //(pi/12.0)*std::exp(2*std::sqrt(aCN*(U-deltaCN)))/std::pow(aCN*std::pow(U-deltaCN,5.0),1.0/4.0);
221 }
222
223 const G4double Spin = theEvaporationProbabilityPtr->GetSpin();
224 //JMQ BIG BUG fixed: hbarc instead of hbar_Planck !!!!
225 // it was fixed in total probability (for this channel) but remained still here!!
226 // G4double g = (2.0*Spin+1.0)*NuclearMass/(pi2* hbar_Planck*hbar_Planck);
227 G4double gg = (2.0*Spin+1.0)*EvaporatedMass/(pi2* hbarc*hbarc);
228 //
229 //JMQ fix on Rb and geometrical cross sections according to Furihata's paper
230 // (JAERI-Data/Code 2001-105, p6)
231 G4double Rb = 0.0;
232 if (A > 4)
233 {
234 G4double Ad = fG4pow->Z13(ResidualA);
235 G4double Aj = fG4pow->Z13(A);
236 // RN = 1.12*(R1 + R2) - 0.86*((R1+R2)/(R1*R2));
237 Rb = 1.12*(Aj + Ad) - 0.86*((Aj+Ad)/(Aj*Ad))+2.85;
238 Rb *= fermi;
239 }
240 else if (A>1)
241 {
242 G4double Ad = fG4pow->Z13(ResidualA);
243 G4double Aj = fG4pow->Z13(A);
244 Rb=1.5*(Aj+Ad)*fermi;
245 }
246 else
247 {
248 G4double Ad = fG4pow->Z13(ResidualA);
249 Rb = 1.5*Ad*fermi;
250 }
251 // G4double GeometricalXS = pi*RN*RN*std::pow(G4double(fragment.GetA()-A),2./3.);
252 G4double GeometricalXS = pi*Rb*Rb;
253
254 G4double ConstantFactor = gg*GeometricalXS*Alpha/InitialLevelDensity;
255 ConstantFactor *= pi/12.0;
256 //JMQ : this is the assimptotic maximal kinetic energy of the ejectile
257 // (obtained by energy-momentum conservation).
258 // In general smaller than U-theSeparationEnergy
259 G4double theEnergy = MaximalKineticEnergy + CoulombBarrier;
260 G4double KineticEnergy;
261 G4double Probability;
262
263 do
264 {
265 KineticEnergy = CoulombBarrier + G4UniformRand()*MaximalKineticEnergy;
266 Probability = ConstantFactor*(KineticEnergy + Beta);
267 G4double a =
268 theLevelDensityPtr->LevelDensityParameter(ResidualA,ResidualZ,theEnergy-KineticEnergy-delta0);
269 G4double T = 1.0/(std::sqrt(a/Ux) - 1.5/Ux);
270 //JMQ fix in units
271
272 if ( theEnergy-KineticEnergy < Ex)
273 {
274 G4double E0 = Ex - T*(std::log(T/MeV) - std::log(a*MeV)/4.0
275 - 1.25*std::log(Ux/MeV) + 2.0*std::sqrt(a*Ux));
276 Probability *= std::exp((theEnergy-KineticEnergy-E0)/T)/T;
277 }
278 else
279 {
280 Probability *= std::exp(2*std::sqrt(a*(theEnergy-KineticEnergy-delta0)))/
281 std::pow(a*fG4pow->powN(theEnergy-KineticEnergy-delta0,5), 0.25);
282 }
283 }
284 while (Normalization*G4UniformRand() > Probability);
285
286 return KineticEnergy;
287}
288
289
290G4ThreeVector G4GEMChannel::IsotropicVector(const G4double Magnitude)
291 // Samples a isotropic random vectorwith a magnitud given by Magnitude.
292 // By default Magnitude = 1.0
293{
294 G4double CosTheta = 1.0 - 2.0*G4UniformRand();
295 G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta);
296 G4double Phi = twopi*G4UniformRand();
297 G4ThreeVector Vector(Magnitude*std::cos(Phi)*SinTheta,
298 Magnitude*std::sin(Phi)*SinTheta,
299 Magnitude*CosTheta);
300 return Vector;
301}
302
303
304
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
Hep3Vector boostVector() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:240
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:235
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:251
G4int GetZ_asInt() const
Definition: G4Fragment.hh:223
G4double GetA() const
Definition: G4Fragment.hh:283
G4int GetA_asInt() const
Definition: G4Fragment.hh:218
G4GEMChannel(const G4int theA, const G4int theZ, const G4String &aName, G4GEMProbability *aEmissionStrategy, G4VCoulombBarrier *aCoulombBarrier)
Definition: G4GEMChannel.cc:44
virtual ~G4GEMChannel()
Definition: G4GEMChannel.cc:63
virtual G4double GetEmissionProbability(G4Fragment *theNucleus)
Definition: G4GEMChannel.cc:68
virtual G4FragmentVector * BreakUp(const G4Fragment &theNucleus)
G4double GetNormalization(void) const
G4double CalcAlphaParam(const G4Fragment &) const
G4double CalcBetaParam(const G4Fragment &) const
G4double GetSpin(void) const
G4double EmissionProbability(const G4Fragment &fragment, G4double anEnergy)
static G4double GetNuclearMass(const G4double A, const G4double Z)
static G4PairingCorrection * GetInstance()
G4double GetPairingCorrection(G4int A, G4int Z) const
static G4Pow * GetInstance()
Definition: G4Pow.cc:50
G4double powN(G4double x, G4int n)
Definition: G4Pow.cc:98
G4double Z13(G4int Z)
Definition: G4Pow.hh:110
virtual G4double GetCoulombBarrier(G4int ARes, G4int ZRes, G4double U) const =0
virtual G4double LevelDensityParameter(G4int A, G4int Z, G4double U) const =0
Definition: DoubConv.h:17
const G4double pi