Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NucleiModel.hh
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// GEANT4 tag: $Name: $
28//
29// 20100319 M. Kelsey -- Remove "using" directory and unnecessary #includes,
30// move ctor to .cc file
31// 20100407 M. Kelsey -- Create "partners thePartners" data member to act
32// as buffer between ::generateInteractionPartners() and
33// ::generateParticleFate(), and make "outgoing_cparticles" a
34// data member returned from the latter by const-ref. Replace
35// return-by-value of initializeCascad() with an input buffer.
36// 20100409 M. Kelsey -- Add function to sort list of partnerts by pathlen,
37// move non-inlinable code to .cc.
38// 20100421 M. Kelsey -- Move getFermiKinetic() to .cc, no hardwired masses.
39// 20100517 M. Kelsey -- Change cross-section tables to static arrays. Move
40// absorptionCrossSection() from SpecialFunc.
41// 20100520 M. Kelsey -- Add function to separate momentum from nucleon
42// 20100617 M. Kelsey -- Add setVerboseLevel() function, add generateModel()
43// with particle input, and ctor with A/Z input.
44// 20100715 M. Kelsey -- Add G4InuclNuclei object for use with balance checks
45// 20100723 M. Kelsey -- Move G4CollisionOutput buffer, along with all
46// std::vector<> buffers, here for reuse
47// 20100914 M. Kelsey -- Migrate to integer A and Z
48// 20101004 M. Kelsey -- Rename and create functions used to generate model
49// 20101019 M. Kelsey -- CoVerity report: dtor leak; move dtor to .cc file
50// 20110223 M. Kelsey -- Add static parameters for radius and cross-section
51// scaling factors.
52// 20110303 M. Kelsey -- Add accessors for model parameters and units
53// 20110304 M. Kelsey -- Extend reset() to fill neutron and proton counts
54// 20110324 D. Wright -- Add list of nucleon interaction points, and nucleon
55// effective radius, for trailing effect.
56// 20110324 M. Kelsey -- Extend reset() to pass list of points; move
57// implementation to .cc file.
58// 20110405 M. Kelsey -- Add "passTrailing()" to encapsulate trailing effect
59// 20110808 M. Kelsey -- Pass buffer into generateParticleFate instead of
60// returning a vector to be copied.
61// 20110823 M. Kelsey -- Remove local cross-section tables entirely
62// 20120125 M. Kelsey -- Add special case for photons to have zero potential.
63// 20120307 M. Kelsey -- Add zone volume array for use with quasideuteron
64// density, functions to identify projectile, compute interaction
65// distance.
66
67#ifndef G4NUCLEI_MODEL_HH
68#define G4NUCLEI_MODEL_HH
69
70#include <algorithm>
71#include <vector>
73
75#include "G4CascadParticle.hh"
76#include "G4CollisionOutput.hh"
77#include "G4LorentzConvertor.hh"
78
79class G4InuclNuclei;
81
83public:
86 explicit G4NucleiModel(G4InuclNuclei* nuclei);
87
89
90 void setVerboseLevel(G4int verbose) { verboseLevel = verbose; }
91
92 void generateModel(G4InuclNuclei* nuclei);
93 void generateModel(G4int a, G4int z);
94
95 // Arguments here are used for rescattering (::Propagate)
96 void reset(G4int nHitNeutrons=0, G4int nHitProtons=0,
97 const std::vector<G4ThreeVector>* hitPoints=0);
98
99 void printModel() const;
100
101 G4double getDensity(G4int ip, G4int izone) const {
102 return nucleon_densities[ip - 1][izone];
103 }
104
106 return fermi_momenta[ip - 1][izone];
107 }
108
109 G4double getFermiKinetic(G4int ip, G4int izone) const;
110
111 G4double getPotential(G4int ip, G4int izone) const {
112 if (ip == 9) return 0.0; // Special case for photons
113 G4int ip0 = ip < 3 ? ip - 1 : 2;
114 if (ip > 10 && ip < 18) ip0 = 3;
115 if (ip > 20) ip0 = 4;
116 return izone < number_of_zones ? zone_potentials[ip0][izone] : 0.0;
117 }
118
119 // Factor to convert GEANT4 lengths to internal units
120 G4double getRadiusUnits() const { return radiusUnits*CLHEP::fermi; }
121
122 G4double getRadius() const { return nuclei_radius; }
123 G4double getRadius(G4int izone) const {
124 return ( (izone<0) ? 0.
125 : (izone<number_of_zones) ? zone_radii[izone] : nuclei_radius);
126 }
127 G4double getVolume(G4int izone) const {
128 return ( (izone<0) ? 0.
129 : (izone<number_of_zones) ? zone_volumes[izone] : nuclei_volume);
130 }
131
132 G4int getNumberOfZones() const { return number_of_zones; }
134 for (G4int iz=0; iz<number_of_zones; iz++) if (r<zone_radii[iz]) return iz;
135 return number_of_zones;
136 }
137
138 G4int getNumberOfNeutrons() const { return neutronNumberCurrent; }
139 G4int getNumberOfProtons() const { return protonNumberCurrent; }
140
141 G4bool empty() const {
142 return neutronNumberCurrent < 1 && protonNumberCurrent < 1;
143 }
144
146 return cparticle.getCurrentZone() < number_of_zones;
147 }
148
149
151
152 typedef std::pair<std::vector<G4CascadParticle>, std::vector<G4InuclElementaryParticle> > modelLists;
153
154 void initializeCascad(G4InuclNuclei* bullet, G4InuclNuclei* target,
155 modelLists& output);
156
157
158 std::pair<G4int, G4int> getTypesOfNucleonsInvolved() const {
159 return std::pair<G4int, G4int>(current_nucl1, current_nucl2);
160 }
161
163 G4ElementaryParticleCollider* theEPCollider,
164 std::vector<G4CascadParticle>& cascade);
165
166 G4bool isProjectile(const G4CascadParticle& cparticle) const;
167 G4bool worthToPropagate(const G4CascadParticle& cparticle) const;
168
170
172
174 G4double totalCrossSection(G4double ke, G4int rtype) const;
175
176private:
177 G4int verboseLevel;
178
179 G4bool passFermi(const std::vector<G4InuclElementaryParticle>& particles,
180 G4int zone);
181
182 G4bool passTrailing(const G4ThreeVector& hit_position);
183
184 void boundaryTransition(G4CascadParticle& cparticle);
185
186 G4InuclElementaryParticle generateQuasiDeutron(G4int type1,
187 G4int type2,
188 G4int zone) const;
189
190 typedef std::pair<G4InuclElementaryParticle, G4double> partner;
191
192 std::vector<partner> thePartners; // Buffer for output below
193 void generateInteractionPartners(G4CascadParticle& cparticle);
194
195 // Function for std::sort() to use in organizing partners by path length
196 static G4bool sortPartners(const partner& p1, const partner& p2) {
197 return (p2.second > p1.second);
198 }
199
200 // Functions used to generate model nuclear structure
201 void fillBindingEnergies();
202
203 void fillZoneRadii(G4double nuclearRadius);
204
205 G4double fillZoneVolumes(G4double nuclearRadius);
206
207 void fillPotentials(G4int type, G4double tot_vol);
208
209 G4double zoneIntegralWoodsSaxon(G4double ur1, G4double ur2,
210 G4double nuclearRadius) const;
211
212 G4double zoneIntegralGaussian(G4double ur1, G4double ur2,
213 G4double nuclearRadius) const;
214
215 G4double getRatio(G4int ip) const; // Fraction of nucleons still available
216
217 // Scale nuclear density with interactions so far
218 G4double getCurrentDensity(G4int ip, G4int izone) const;
219
220 // Average path length for any interaction of projectile and target
221 // = 1. / (density * cross-section)
222 G4double inverseMeanFreePath(const G4CascadParticle& cparticle,
223 const G4InuclElementaryParticle& target);
224 // NOTE: Function is non-const in order to use dummy_converter
225
226 // Use path-length and MFP (above) to throw random distance to collision
227 G4double generateInteractionLength(const G4CascadParticle& cparticle,
228 G4double path, G4double invmfp) const;
229
230 // Buffers for processing interactions on each cycle
231 G4LorentzConvertor dummy_convertor; // For getting collision frame
232 G4CollisionOutput EPCoutput; // For N-body inelastic collisions
233
234 std::vector<G4InuclElementaryParticle> qdeutrons; // For h+(NN) trials
235 std::vector<G4double> acsecs;
236
237 std::vector<G4ThreeVector> coordinates; // for initializeCascad()
238 std::vector<G4LorentzVector> momentums;
239 std::vector<G4InuclElementaryParticle> raw_particles;
240
241 std::vector<G4ThreeVector> collisionPts;
242
243 // Temporary buffers for computing nuclear model
244 G4double ur[7]; // Number of skin depths for each zone
245 G4double v[6]; // Density integrals by zone
246 G4double v1[6]; // Pseudo-volume (delta r^3) by zone
247 std::vector<G4double> rod; // Nucleon density
248 std::vector<G4double> pf; // Fermi momentum
249 std::vector<G4double> vz; // Nucleo potential
250
251 // Nuclear model configuration
252 std::vector<std::vector<G4double> > nucleon_densities;
253 std::vector<std::vector<G4double> > zone_potentials;
254 std::vector<std::vector<G4double> > fermi_momenta;
255 std::vector<G4double> zone_radii;
256 std::vector<G4double> zone_volumes;
257 std::vector<G4double> binding_energies;
258 G4double nuclei_radius;
259 G4double nuclei_volume;
260 G4int number_of_zones;
261
262 G4int A;
263 G4int Z;
264 G4InuclNuclei* theNucleus;
265
266 G4int neutronNumber;
267 G4int protonNumber;
268
269 G4int neutronNumberCurrent;
270 G4int protonNumberCurrent;
271
272 G4int current_nucl1;
273 G4int current_nucl2;
274
275 // Symbolic names for nuclear potentials
276 enum PotentialType { WoodsSaxon=0, Gaussian=1 };
277
278 // Parameters for nuclear structure
279 static const G4double skinDepth;
280 static const G4double radiusScale; // Coefficients for two-parameter fit
281 static const G4double radiusScale2; // R = 1.16*cbrt(A) - 1.3456/cbrt(A)
282 static const G4double radiusForSmall; // Average radius of light A<5 nuclei
283 static const G4double radScaleAlpha; // Scaling factor R_alpha/R_small
284 static const G4double fermiMomentum;
285 static const G4double R_nucleon;
286 static const G4double alfa3[3], alfa6[6];
287 static const G4double pion_vp;
288 static const G4double pion_vp_small;
289 static const G4double kaon_vp;
290 static const G4double hyperon_vp;
291
292 // FIXME: We should not be using this!
293 static const G4double piTimes4thirds;
294 static const G4double crossSectionUnits;
295 static const G4double radiusUnits;
296};
297
298#endif // G4NUCLEI_MODEL_HH
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4int getCurrentZone() const
G4int getZone(G4double r) const
G4int getNumberOfNeutrons() const
G4double getDensity(G4int ip, G4int izone) const
G4CascadParticle initializeCascad(G4InuclElementaryParticle *particle)
G4bool worthToPropagate(const G4CascadParticle &cparticle) const
G4int getNumberOfProtons() const
void generateParticleFate(G4CascadParticle &cparticle, G4ElementaryParticleCollider *theEPCollider, std::vector< G4CascadParticle > &cascade)
void printModel() const
G4double getRadius(G4int izone) const
G4bool empty() const
G4InuclElementaryParticle generateNucleon(G4int type, G4int zone) const
G4double getVolume(G4int izone) const
G4double totalCrossSection(G4double ke, G4int rtype) const
G4double getPotential(G4int ip, G4int izone) const
G4double getFermiMomentum(G4int ip, G4int izone) const
void reset(G4int nHitNeutrons=0, G4int nHitProtons=0, const std::vector< G4ThreeVector > *hitPoints=0)
G4double absorptionCrossSection(G4double e, G4int type) const
G4bool stillInside(const G4CascadParticle &cparticle)
void generateModel(G4InuclNuclei *nuclei)
std::pair< std::vector< G4CascadParticle >, std::vector< G4InuclElementaryParticle > > modelLists
G4double getFermiKinetic(G4int ip, G4int izone) const
std::pair< G4int, G4int > getTypesOfNucleonsInvolved() const
G4LorentzVector generateNucleonMomentum(G4int type, G4int zone) const
G4double getRadiusUnits() const
G4bool isProjectile(const G4CascadParticle &cparticle) const
G4double getRadius() const
G4int getNumberOfZones() const
void setVerboseLevel(G4int verbose)