Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4INCLParticleTable.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// INCL++ intra-nuclear cascade model
27// Alain Boudard, CEA-Saclay, France
28// Joseph Cugnon, University of Liege, Belgium
29// Jean-Christophe David, CEA-Saclay, France
30// Pekka Kaitaniemi, CEA-Saclay, France, and Helsinki Institute of Physics, Finland
31// Sylvie Leray, CEA-Saclay, France
32// Davide Mancusi, CEA-Saclay, France
33//
34#define INCLXX_IN_GEANT4_MODE 1
35
36#include "globals.hh"
37
40#include <algorithm>
41// #include <cassert>
42#include <cmath>
43#include <cctype>
44#include <sstream>
45#ifdef INCLXX_IN_GEANT4_MODE
46#include "G4SystemOfUnits.hh"
47#endif
48
49#ifdef INCLXX_IN_GEANT4_MODE
51#include "G4SystemOfUnits.hh"
52#endif
53
54namespace G4INCL {
55
56 namespace ParticleTable {
57
58 namespace {
59
60 /// \brief Static instance of the NaturalIsotopicAbundances class
61 const NaturalIsotopicDistributions *theNaturalIsotopicDistributions = NULL;
62
63 const G4double theINCLNucleonMass = 938.2796;
64 const G4double theINCLPionMass = 138.0;
65 const G4double theINCLLambdaMass = 1115.683;
66 //const G4double theINCLSigmaMass = 1197.45;
67 //const G4double theINCLKaonMass = 497.614;
68 const G4double theINCLEtaMass = 547.862;
69 const G4double theINCLOmegaMass = 782.65;
70 const G4double theINCLEtaPrimeMass = 957.78;
71 const G4double theINCLPhotonMass = 0.0;
72 G4ThreadLocal G4double protonMass = 0.0;
73 G4ThreadLocal G4double neutronMass = 0.0;
74 G4ThreadLocal G4double piPlusMass = 0.0;
75 G4ThreadLocal G4double piMinusMass = 0.0;
76 G4ThreadLocal G4double piZeroMass = 0.0;
77 G4ThreadLocal G4double SigmaPlusMass = 0.0;
78 G4ThreadLocal G4double SigmaZeroMass = 0.0;
79 G4ThreadLocal G4double SigmaMinusMass = 0.0;
80 G4ThreadLocal G4double LambdaMass = 0.0;
81 G4ThreadLocal G4double XiMinusMass = 0.0;
82 G4ThreadLocal G4double XiZeroMass = 0.0;
83 G4ThreadLocal G4double antiProtonMass = 0.0;
84 G4ThreadLocal G4double antiNeutronMass = 0.0;
85 G4ThreadLocal G4double antiSigmaPlusMass = 0.0;
86 G4ThreadLocal G4double antiSigmaZeroMass = 0.0;
87 G4ThreadLocal G4double antiSigmaMinusMass = 0.0;
88 G4ThreadLocal G4double antiLambdaMass = 0.0;
89 G4ThreadLocal G4double antiXiMinusMass = 0.0;
90 G4ThreadLocal G4double antiXiZeroMass = 0.0;
91 G4ThreadLocal G4double KPlusMass = 0.0;
92 G4ThreadLocal G4double KZeroMass = 0.0;
93 G4ThreadLocal G4double KZeroBarMass = 0.0;
94 G4ThreadLocal G4double KShortMass = 0.0;
95 G4ThreadLocal G4double KLongMass = 0.0;
96 G4ThreadLocal G4double KMinusMass = 0.0;
97 G4ThreadLocal G4double etaMass = 0.0;
98 G4ThreadLocal G4double omegaMass = 0.0;
99 G4ThreadLocal G4double etaPrimeMass = 0.0;
100 G4ThreadLocal G4double photonMass = 0.0;
101
102 // Hard-coded values of the real particle masses (MeV/c^2)
103 G4ThreadLocal G4double theRealProtonMass = 938.27203;
104 G4ThreadLocal G4double theRealNeutronMass = 939.56536;
105 G4ThreadLocal G4double theRealChargedPiMass = 139.57018;
106 G4ThreadLocal G4double theRealPiZeroMass = 134.9766;
107 G4ThreadLocal G4double theRealLambdaMass = 1115.683;
108 G4ThreadLocal G4double theRealSigmaPlusMass = 1189.37;
109 G4ThreadLocal G4double theRealSigmaZeroMass = 1192.64;
110 G4ThreadLocal G4double theRealSigmaMinusMass = 1197.45;
111 G4ThreadLocal G4double theRealAntiProtonMass = 938.27203;
112 G4ThreadLocal G4double theRealXiMinusMass = 1321.71;
113 G4ThreadLocal G4double theRealXiZeroMass = 1314.86;
114 G4ThreadLocal G4double theRealAntiNeutronMass = 939.56536;
115 G4ThreadLocal G4double theRealAntiLambdaMass = 1115.683;
116 G4ThreadLocal G4double theRealAntiSigmaPlusMass = 1189.37;
117 G4ThreadLocal G4double theRealAntiSigmaZeroMass = 1192.64;
118 G4ThreadLocal G4double theRealAntiSigmaMinusMass = 1197.45;
119 G4ThreadLocal G4double theRealAntiXiMinusMass = 1321.71;
120 G4ThreadLocal G4double theRealAntiXiZeroMass = 1314.86;
121 G4ThreadLocal G4double theRealChargedKaonMass = 493.677;
122 G4ThreadLocal G4double theRealNeutralKaonMass = 497.614;
123 G4ThreadLocal G4double theRealEtaMass = 547.862;
124 G4ThreadLocal G4double theRealOmegaMass = 782.65;
125 G4ThreadLocal G4double theRealEtaPrimeMass = 957.78;
126 G4ThreadLocal G4double theRealPhotonMass = 0.0;
127
128 // Width (second)
129 const G4double theChargedPiWidth = 2.6033e-08;
130 const G4double thePiZeroWidth = 8.52e-17;
131 const G4double theEtaWidth = 5.025e-19; // 1.31 keV
132 const G4double theOmegaWidth = 7.7528e-23; // 8.49 MeV
133 const G4double theEtaPrimeWidth = 3.3243e-21; // 0.198 MeV
134 const G4double theChargedKaonWidth = 1.238e-08;
135 const G4double theKShortWidth = 8.954e-11;
136 const G4double theKLongWidth = 5.116e-08;
137 const G4double theLambdaWidth = 2.632e-10;
138 const G4double theSigmaPlusWidth = 8.018e-11;
139 const G4double theSigmaZeroWidth = 7.4e-20;
140 const G4double theSigmaMinusWidth = 1.479e-10;
141 //const G4double theXiMinusWidth = 1.639e-10;
142 //const G4double theXiZeroWidth = 2.90e-10;
143 //const G4double theAntiLambdaWidth = 2.632e-10;
144 //const G4double theAntiSigmaPlusWidth = 8.018e-11;
145 //const G4double theAntiSigmaZeroWidth = 7.4e-20;
146 //const G4double theAntiSigmaMinusWidth = 1.479e-10;
147 //const G4double theAntiXiMinusWidth = 1.639e-10;
148 //const G4double theAntiXiZeroWidth = 2.90e-10;
149 G4ThreadLocal G4double piPlusWidth = 0.0;
150 G4ThreadLocal G4double piMinusWidth = 0.0;
151 G4ThreadLocal G4double piZeroWidth = 0.0;
152 G4ThreadLocal G4double etaWidth = 0.0;
153 G4ThreadLocal G4double omegaWidth = 0.0;
154 G4ThreadLocal G4double etaPrimeWidth = 0.0;
155 G4ThreadLocal G4double LambdaWidth = 0.0;
156 G4ThreadLocal G4double SigmaPlusWidth = 0.0;
157 G4ThreadLocal G4double SigmaZeroWidth = 0.0;
158 G4ThreadLocal G4double SigmaMinusWidth = 0.0;
159 G4ThreadLocal G4double KPlusWidth = 0.0;
160 G4ThreadLocal G4double KMinusWidth = 0.0;
161 G4ThreadLocal G4double KShortWidth = 0.0;
162 G4ThreadLocal G4double KLongWidth = 0.0;
163 G4ThreadLocal G4double XiMinusWidth = 0.0;
164 G4ThreadLocal G4double XiZeroWidth = 0.0;
165 G4ThreadLocal G4double antiLambdaWidth = 0.0;
166 G4ThreadLocal G4double antiSigmaZeroWidth = 0.0;
167 G4ThreadLocal G4double antiSigmaMinusWidth = 0.0;
168 G4ThreadLocal G4double antiSigmaPlusWidth = 0.0;
169 G4ThreadLocal G4double antiXiZeroWidth = 0.0;
170 G4ThreadLocal G4double antiXiMinusWidth = 0.0;
171
172 const G4int mediumNucleiTableSize = 30;
173
174 const G4double mediumDiffuseness[mediumNucleiTableSize] =
175 {0.0,0.0,0.0,0.0,0.0,1.78,1.77,1.77,1.69,1.71,
176 1.69,1.72,1.635,1.730,1.81,1.833,1.798,
177 1.93,0.567,0.571, 0.560,0.549,0.550,0.551,
178 0.580,0.575,0.569,0.537,0.0,0.0};
179 const G4double mediumRadius[mediumNucleiTableSize] =
180 {0.0,0.0,0.0,0.0,0.0,0.334,0.327,0.479,0.631,0.838,
181 0.811,0.84,1.403,1.335,1.25,1.544,1.498,1.57,
182 2.58,2.77, 2.775,2.78,2.88,2.98,3.22,3.03,2.84,
183 3.14,0.0,0.0};
184
185 const G4double positionRMS[clusterTableZSize][clusterTableASize] = {
186 /* A= 0 1 2 3 4 5 6 7 8 9 10 11 12 */
187 /* Z=0 */ {-1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0},
188 /* Z=1 */ {-1.0, -1.0, 2.10, 1.80, 1.70, 1.83, 2.60, 2.50, -1.0, -1.0, -1.0, -1.0, -1.0},
189 /* Z=2 */ {-1.0, -1.0, -1.0, 1.80, 1.68, 1.70, 2.60, 2.50, 2.50, 2.50, 2.50, -1.0, -1.0},
190 /* Z=3 */ {-1.0, -1.0, -1.0, -1.0, 1.70, 1.83, 2.56, 2.40, 2.50, 2.50, 2.50, 2.50, 2.50},
191 /* Z=4 */ {-1.0, -1.0, -1.0, -1.0, -1.0, -1.0, 2.60, 2.50, 2.50, 2.51, 2.50, 2.50, 2.50},
192 /* Z=5 */ {-1.0, -1.0, -1.0, -1.0, -1.0, -1.0, 2.50, 2.50, 2.50, 2.50, 2.45, 2.40, 2.50},
193 /* Z=6 */ {-1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, 2.50, 2.50, 2.50, 2.50, 2.47},
194 /* Z=7 */ {-1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, 2.50, 2.50, 2.50},
195 /* Z=8 */ {-1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, 2.50}
196 };
197
198 const G4double momentumRMS[clusterTableZSize][clusterTableASize] = {
199 /* A= 0 1 2 3 4 5 6 7 8 9 10 11 12 */
200 /* Z=0 */ {-1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0},
201 /* Z=1 */ {-1.0, -1.0, 77.0, 110., 153., 100., 100., 100., -1.0, -1.0, -1.0, -1.0, -1.0},
202 /* Z=2 */ {-1.0, -1.0, -1.0, 110., 153., 100., 100., 100., 100., 100., 100., -1.0, -1.0},
203 /* Z=3 */ {-1.0, -1.0, -1.0, -1.0, 153., 100., 100., 100., 100., 100., 100., 100., 100.},
204 /* Z=4 */ {-1.0, -1.0, -1.0, -1.0, -1.0, -1.0, 100., 100., 100., 100., 100., 100., 100.},
205 /* Z=5 */ {-1.0, -1.0, -1.0, -1.0, -1.0, -1.0, 100., 100., 100., 100., 100., 100., 100.},
206 /* Z=6 */ {-1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, 100., 100., 100., 100., 100.},
207 /* Z=7 */ {-1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, 100., 100., 100.},
208 /* Z=8 */ {-1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, -1.0, 100.}
209 };
210
211 const G4int elementTableSize = 113; // up to Cn
212
213 /// \brief Table of chemical element names
214 const std::string elementTable[elementTableSize] = {
215 "",
216 "H",
217 "He",
218 "Li",
219 "Be",
220 "B",
221 "C",
222 "N",
223 "O",
224 "F",
225 "Ne",
226 "Na",
227 "Mg",
228 "Al",
229 "Si",
230 "P",
231 "S",
232 "Cl",
233 "Ar",
234 "K",
235 "Ca",
236 "Sc",
237 "Ti",
238 "V",
239 "Cr",
240 "Mn",
241 "Fe",
242 "Co",
243 "Ni",
244 "Cu",
245 "Zn",
246 "Ga",
247 "Ge",
248 "As",
249 "Se",
250 "Br",
251 "Kr",
252 "Rb",
253 "Sr",
254 "Y",
255 "Zr",
256 "Nb",
257 "Mo",
258 "Tc",
259 "Ru",
260 "Rh",
261 "Pd",
262 "Ag",
263 "Cd",
264 "In",
265 "Sn",
266 "Sb",
267 "Te",
268 "I",
269 "Xe",
270 "Cs",
271 "Ba",
272 "La",
273 "Ce",
274 "Pr",
275 "Nd",
276 "Pm",
277 "Sm",
278 "Eu",
279 "Gd",
280 "Tb",
281 "Dy",
282 "Ho",
283 "Er",
284 "Tm",
285 "Yb",
286 "Lu",
287 "Hf",
288 "Ta",
289 "W",
290 "Re",
291 "Os",
292 "Ir",
293 "Pt",
294 "Au",
295 "Hg",
296 "Tl",
297 "Pb",
298 "Bi",
299 "Po",
300 "At",
301 "Rn",
302 "Fr",
303 "Ra",
304 "Ac",
305 "Th",
306 "Pa",
307 "U",
308 "Np",
309 "Pu",
310 "Am",
311 "Cm",
312 "Bk",
313 "Cf",
314 "Es",
315 "Fm",
316 "Md",
317 "No",
318 "Lr",
319 "Rf",
320 "Db",
321 "Sg",
322 "Bh",
323 "Hs",
324 "Mt",
325 "Ds",
326 "Rg",
327 "Cn"
328 };
329
330 /// \brief Digit names to compose IUPAC element names
331 const std::string elementIUPACDigits = "nubtqphsoe";
332
333#define INCL_DEFAULT_SEPARATION_ENERGY 6.83
334 const G4double theINCLProtonSeparationEnergy = INCL_DEFAULT_SEPARATION_ENERGY;
335 const G4double theINCLNeutronSeparationEnergy = INCL_DEFAULT_SEPARATION_ENERGY;
336 const G4double theINCLLambdaSeparationEnergy = INCL_DEFAULT_SEPARATION_ENERGY;
337 //const G4double theINCLantiProtonSeparationEnergy = INCL_DEFAULT_SEPARATION_ENERGY;
338 const G4double theINCLantiProtonSeparationEnergy = 0.;
340 G4ThreadLocal G4double neutronSeparationEnergy = INCL_DEFAULT_SEPARATION_ENERGY;
342 //G4ThreadLocal G4double antiprotonSeparationEnergy = INCL_DEFAULT_SEPARATION_ENERGY;
343 //G4ThreadLocal G4double antiprotonSeparationEnergy = 0.;
344#undef INCL_DEFAULT_SEPARATION_ENERGY
345
346 G4ThreadLocal G4double rpCorrelationCoefficient[UnknownParticle];
347
348 G4ThreadLocal G4double neutronSkin = 0.0;
349 G4ThreadLocal G4double neutronHalo = 0.0;
350
351#ifdef INCLXX_IN_GEANT4_MODE
352 G4ThreadLocal G4IonTable *theG4IonTable;
353#endif
354
355 /// \brief Default value for constant Fermi momentum
356 G4ThreadLocal G4double constantFermiMomentum = 0.0;
357
358 /// \brief Transform a IUPAC char to an char representing an integer digit
359 char iupacToInt(char c) {
360 return (char)(((G4int)'0')+elementIUPACDigits.find(c));
361 }
362
363 /// \brief Transform an integer digit (represented by a char) to a IUPAC char
364 char intToIUPAC(char n) { return elementIUPACDigits.at(n); }
365
366 /// \brief Get the singleton instance of the natural isotopic distributions
367 const NaturalIsotopicDistributions *getNaturalIsotopicDistributions() {
368 if(!theNaturalIsotopicDistributions)
369 theNaturalIsotopicDistributions = new NaturalIsotopicDistributions;
370 return theNaturalIsotopicDistributions;
371 }
372
373 } // namespace
374
375 void initialize(Config const * const theConfig /*=0*/) {
376 protonMass = theINCLNucleonMass;
377 neutronMass = theINCLNucleonMass;
378 piPlusMass = theINCLPionMass;
379 piMinusMass = theINCLPionMass;
380 piZeroMass = theINCLPionMass;
381
382 etaMass = theINCLEtaMass;
383 omegaMass = theINCLOmegaMass;
384 etaPrimeMass = theINCLEtaPrimeMass;
385 photonMass = theINCLPhotonMass;
386
387 SigmaPlusMass = theRealSigmaPlusMass;
388 SigmaMinusMass = theRealSigmaMinusMass;
389 SigmaZeroMass = theRealSigmaZeroMass;
390 LambdaMass = theINCLLambdaMass;
391 KPlusMass = theRealChargedKaonMass;
392 KZeroMass = theRealNeutralKaonMass;
393 KZeroBarMass = theRealNeutralKaonMass;
394 KShortMass = theRealNeutralKaonMass;
395 KLongMass = theRealNeutralKaonMass;
396 KMinusMass = theRealChargedKaonMass;
397
398 antiProtonMass = theRealAntiProtonMass;
399 XiZeroMass = theRealXiZeroMass;
400 XiMinusMass = theRealXiMinusMass;
401 antiNeutronMass = theRealAntiNeutronMass;
402 antiSigmaPlusMass = theRealAntiSigmaPlusMass;
403 antiSigmaMinusMass = theRealAntiSigmaMinusMass;
404 antiSigmaZeroMass = theRealAntiSigmaZeroMass;
405 antiLambdaMass = theRealAntiLambdaMass; //!
406 antiXiZeroMass = theRealAntiXiZeroMass;
407 antiXiMinusMass = theRealAntiXiMinusMass;
408
409 if(theConfig && theConfig->getUseRealMasses()) {
412 } else {
415 }
416
417#ifndef INCLXX_IN_GEANT4_MODE
418 std::string dataFilePath;
419 if(theConfig)
420 dataFilePath = theConfig->getINCLXXDataFilePath();
421 NuclearMassTable::initialize(dataFilePath, getRealMass(Proton), getRealMass(Neutron));
422#endif
423
424#ifdef INCLXX_IN_GEANT4_MODE
426 theG4IonTable = theG4ParticleTable->GetIonTable();
427 theRealProtonMass = theG4ParticleTable->FindParticle("proton")->GetPDGMass() / MeV;
428 theRealNeutronMass = theG4ParticleTable->FindParticle("neutron")->GetPDGMass() / MeV;
429 theRealChargedPiMass = theG4ParticleTable->FindParticle("pi+")->GetPDGMass() / MeV;
430 theRealPiZeroMass = theG4ParticleTable->FindParticle("pi0")->GetPDGMass() / MeV;
431
432 theRealEtaMass = theG4ParticleTable->FindParticle("eta")->GetPDGMass() / MeV;
433 theRealOmegaMass = theG4ParticleTable->FindParticle("omega")->GetPDGMass() / MeV;
434 theRealEtaPrimeMass = theG4ParticleTable->FindParticle("eta_prime")->GetPDGMass() / MeV;
435 theRealPhotonMass = theG4ParticleTable->FindParticle("gamma")->GetPDGMass() / MeV;
436
437 theRealSigmaPlusMass = theG4ParticleTable->FindParticle("sigma+")->GetPDGMass() / MeV;
438 theRealSigmaZeroMass = theG4ParticleTable->FindParticle("sigma0")->GetPDGMass() / MeV;
439 theRealSigmaMinusMass = theG4ParticleTable->FindParticle("sigma-")->GetPDGMass() / MeV;
440 theRealLambdaMass = theG4ParticleTable->FindParticle("lambda")->GetPDGMass() / MeV;
441 theRealChargedKaonMass = theG4ParticleTable->FindParticle("kaon+")->GetPDGMass() / MeV;
442 theRealNeutralKaonMass = theG4ParticleTable->FindParticle("kaon0")->GetPDGMass() / MeV;
443
444 theRealAntiProtonMass = theG4ParticleTable->FindParticle("anti_proton")->GetPDGMass() / MeV;
445 theRealAntiNeutronMass = theG4ParticleTable->FindParticle("anti_neutron")->GetPDGMass() / MeV;
446 theRealXiZeroMass = theG4ParticleTable->FindParticle("xi0")->GetPDGMass() / MeV;
447 theRealXiMinusMass = theG4ParticleTable->FindParticle("xi-")->GetPDGMass() / MeV;
448 theRealAntiSigmaPlusMass = theG4ParticleTable->FindParticle("anti_sigma+")->GetPDGMass() / MeV;
449 theRealAntiSigmaZeroMass = theG4ParticleTable->FindParticle("anti_sigma0")->GetPDGMass() / MeV;
450 theRealAntiSigmaMinusMass = theG4ParticleTable->FindParticle("anti_sigma-")->GetPDGMass() / MeV;
451 theRealAntiLambdaMass = theG4ParticleTable->FindParticle("anti_lambda")->GetPDGMass() / MeV;
452 theRealAntiXiZeroMass = theG4ParticleTable->FindParticle("anti_xi0")->GetPDGMass() / MeV;
453 theRealAntiXiMinusMass = theG4ParticleTable->FindParticle("anti_xi-")->GetPDGMass() / MeV;
454#endif
455
456 minDeltaMass = theRealNeutronMass + theRealChargedPiMass + 0.5;
459
460 piPlusWidth = theChargedPiWidth;
461 piMinusWidth = theChargedPiWidth;
462 piZeroWidth = thePiZeroWidth;
463 etaWidth = theEtaWidth;
464 omegaWidth = theOmegaWidth;
465 etaPrimeWidth = theEtaPrimeWidth;
466
467 SigmaMinusWidth = theSigmaMinusWidth;
468 SigmaPlusWidth = theSigmaPlusWidth;
469 SigmaZeroWidth = theSigmaZeroWidth;
470 LambdaWidth = theLambdaWidth;
471 KPlusWidth = theChargedKaonWidth;
472 KMinusWidth = theChargedKaonWidth;
473 KShortWidth = theKShortWidth;
474 KLongWidth = theKLongWidth;
475
476 // Initialise HFB tables
477#ifdef INCLXX_IN_GEANT4_MODE
479#else
480 HFB::initialize(dataFilePath);
481#endif
482
483 // Initialise the separation-energy function
484 if(!theConfig || theConfig->getSeparationEnergyType()==INCLSeparationEnergy)
486 else if(theConfig->getSeparationEnergyType()==RealSeparationEnergy)
490 else {
491 INCL_FATAL("Unrecognized separation-energy type in ParticleTable initialization: " << theConfig->getSeparationEnergyType() << '\n');
492 return;
493 }
494
495 // Initialise the Fermi-momentum function
496 if(!theConfig || theConfig->getFermiMomentumType()==ConstantFermiMomentum) {
498 if(theConfig) {
499 const G4double aFermiMomentum = theConfig->getFermiMomentum();
500 if(aFermiMomentum>0.)
501 constantFermiMomentum = aFermiMomentum;
502 else
503 constantFermiMomentum = PhysicalConstants::Pf;
504 } else {
505 constantFermiMomentum = PhysicalConstants::Pf;
506 }
507 } else if(theConfig->getFermiMomentumType()==ConstantLightFermiMomentum)
511 else {
512 INCL_FATAL("Unrecognized Fermi-momentum type in ParticleTable initialization: " << theConfig->getFermiMomentumType() << '\n');
513 return;
514 }
515
516 // Initialise the r-p correlation coefficients
517 std::fill(rpCorrelationCoefficient, rpCorrelationCoefficient + UnknownParticle, 1.);
518 if(theConfig) {
519 rpCorrelationCoefficient[Proton] = theConfig->getRPCorrelationCoefficient(Proton);
520 rpCorrelationCoefficient[Neutron] = theConfig->getRPCorrelationCoefficient(Neutron);
521 }
522
523 // Initialise the neutron-skin parameters
524 if(theConfig) {
525 neutronSkin = theConfig->getNeutronSkin();
526 neutronHalo = theConfig->getNeutronHalo();
527 }
528
529 }
530
532 // Actually this is the 3rd component of isospin (I_z) multiplied by 2!
533 if(t == Proton) {
534 return 1;
535 } else if(t == Neutron) {
536 return -1;
537 } else if(t == PiPlus) {
538 return 2;
539 } else if(t == PiMinus) {
540 return -2;
541 } else if(t == PiZero) {
542 return 0;
543 } else if(t == DeltaPlusPlus) {
544 return 3;
545 } else if(t == DeltaPlus) {
546 return 1;
547 } else if(t == DeltaZero) {
548 return -1;
549 } else if(t == DeltaMinus) {
550 return -3;
551 } else if(t == Lambda) {
552 return 0;
553 } else if(t == SigmaPlus) {
554 return 2;
555 } else if(t == SigmaZero) {
556 return 0;
557 } else if(t == SigmaMinus) {
558 return -2;
559 } else if(t == KPlus) {
560 return 1;
561 } else if(t == KZero) {
562 return -1;
563 } else if(t == KZeroBar) {
564 return 1;
565 } else if(t == KShort) {
566 return 0;
567 } else if(t == KLong) {
568 return 0;
569 } else if(t == KMinus) {
570 return -1;
571 } else if(t == Eta) {
572 return 0;
573 } else if(t == Omega) {
574 return 0;
575 } else if(t == EtaPrime) {
576 return 0;
577 } else if(t == Photon) {
578 return 0;
579 } else if(t == antiProton) {
580 return -1;
581 } else if(t == XiMinus) {
582 return -1;
583 } else if(t == XiZero) {
584 return 1;
585 } else if(t == antiNeutron) {
586 return 1;
587 } else if(t == antiLambda) {
588 return 0;
589 } else if(t == antiSigmaPlus) {
590 return -2;
591 } else if(t == antiSigmaZero) {
592 return 0;
593 } else if(t == antiSigmaMinus) {
594 return 2;
595 } else if(t == antiXiMinus) {
596 return 1;
597 } else if(t == antiXiZero) {
598 return -1;
599 }
600 INCL_ERROR("Requested isospin of an unknown particle!");
601 return -10; // Unknown
602 }
603
604 std::string getShortName(const ParticleSpecies &sp) {
605 if(sp.theType==Composite && sp.theS == 0)
606 return getShortName(sp.theA,sp.theZ);
607 else if(sp.theType==Composite)
608 return getName(sp.theA,sp.theZ,sp.theS);
609 else
610 return getShortName(sp.theType);
611 }
612
613 std::string getName(const ParticleSpecies &sp) {
614 if(sp.theType==Composite && sp.theS == 0)
615 return getName(sp.theA,sp.theZ);
616 else if(sp.theType==Composite)
617 return getName(sp.theA,sp.theZ,sp.theS);
618 else
619 return getName(sp.theType);
620 }
621
622 std::string getName(const G4int A, const G4int Z) {
623 std::stringstream stream;
624 stream << getElementName(Z) << "-" << A;
625 return stream.str();
626 }
627
628 std::string getName(const G4int A, const G4int Z, const G4int S) {
629 std::stringstream stream;
630 if(S >= 0) // S < 0 for hypernuclei
631 return getName(A, Z);
632 else if(S == -1)
633 stream << getElementName(Z) << "-" << A << "_" << "Lambda";
634 else
635 stream << getElementName(Z) << "-" << A << "_" << S << "-Lambda";
636 return stream.str();
637 }
638
639 std::string getShortName(const G4int A, const G4int Z) {
640 std::stringstream stream;
641 stream << getElementName(Z);
642 if(A>0)
643 stream << A;
644 return stream.str();
645 }
646
647 std::string getName(const ParticleType p) {
648 if(p == G4INCL::Proton) {
649 return std::string("proton");
650 } else if(p == G4INCL::Neutron) {
651 return std::string("neutron");
652 } else if(p == G4INCL::DeltaPlusPlus) {
653 return std::string("delta++");
654 } else if(p == G4INCL::DeltaPlus) {
655 return std::string("delta+");
656 } else if(p == G4INCL::DeltaZero) {
657 return std::string("delta0");
658 } else if(p == G4INCL::DeltaMinus) {
659 return std::string("delta-");
660 } else if(p == G4INCL::PiPlus) {
661 return std::string("pi+");
662 } else if(p == G4INCL::PiZero) {
663 return std::string("pi0");
664 } else if(p == G4INCL::PiMinus) {
665 return std::string("pi-");
666 } else if(p == G4INCL::Lambda) {
667 return std::string("lambda");
668 } else if(p == G4INCL::SigmaPlus) {
669 return std::string("sigma+");
670 } else if(p == G4INCL::SigmaZero) {
671 return std::string("sigma0");
672 } else if(p == G4INCL::SigmaMinus) {
673 return std::string("sigma-");
674 } else if(p == G4INCL::antiProton) {
675 return std::string("antiproton");
676 } else if(p == G4INCL::XiMinus) {
677 return std::string("xi-");
678 } else if(p == G4INCL::XiZero) {
679 return std::string("xi0");
680 } else if(p == G4INCL::antiNeutron) {
681 return std::string("antineutron");
682 } else if(p == G4INCL::antiSigmaPlus) {
683 return std::string("antisigma+");
684 } else if(p == G4INCL::antiSigmaZero) {
685 return std::string("antisigma0");
686 } else if(p == G4INCL::antiSigmaMinus) {
687 return std::string("antisigma-");
688 } else if(p == G4INCL::antiLambda) {
689 return std::string("antilambda");
690 } else if(p == G4INCL::antiXiMinus) {
691 return std::string("antixi-");
692 } else if(p == G4INCL::antiXiZero) {
693 return std::string("antixi0");
694 } else if(p == G4INCL::KPlus) {
695 return std::string("kaon+");
696 } else if(p == G4INCL::KZero) {
697 return std::string("kaon0");
698 } else if(p == G4INCL::KZeroBar) {
699 return std::string("kaon0bar");
700 } else if(p == G4INCL::KMinus) {
701 return std::string("kaon-");
702 } else if(p == G4INCL::KShort) {
703 return std::string("kaonshort");
704 } else if(p == G4INCL::KLong) {
705 return std::string("kaonlong");
706 } else if(p == G4INCL::Composite) {
707 return std::string("composite");
708 } else if(p == G4INCL::Eta) {
709 return std::string("eta");
710 } else if(p == G4INCL::Omega) {
711 return std::string("omega");
712 } else if(p == G4INCL::EtaPrime) {
713 return std::string("etaprime");
714 } else if(p == G4INCL::Photon) {
715 return std::string("photon");
716 }
717 return std::string("unknown");
718 }
719
720 std::string getShortName(const ParticleType p) {
721 if(p == G4INCL::Proton) {
722 return std::string("p");
723 } else if(p == G4INCL::Neutron) {
724 return std::string("n");
725 } else if(p == G4INCL::DeltaPlusPlus) {
726 return std::string("d++");
727 } else if(p == G4INCL::DeltaPlus) {
728 return std::string("d+");
729 } else if(p == G4INCL::DeltaZero) {
730 return std::string("d0");
731 } else if(p == G4INCL::DeltaMinus) {
732 return std::string("d-");
733 } else if(p == G4INCL::PiPlus) {
734 return std::string("pi+");
735 } else if(p == G4INCL::PiZero) {
736 return std::string("pi0");
737 } else if(p == G4INCL::PiMinus) {
738 return std::string("pi-");
739 } else if(p == G4INCL::Lambda) {
740 return std::string("l");
741 } else if(p == G4INCL::SigmaPlus) {
742 return std::string("s+");
743 } else if(p == G4INCL::SigmaZero) {
744 return std::string("s0");
745 } else if(p == G4INCL::SigmaMinus) {
746 return std::string("s-");
747 } else if(p == G4INCL::antiProton) {
748 return std::string("pb");
749 } else if(p == G4INCL::XiMinus) {
750 return std::string("x-");
751 } else if(p == G4INCL::XiZero) {
752 return std::string("x0");
753 } else if(p == G4INCL::antiNeutron) {
754 return std::string("nb");
755 } else if(p == G4INCL::antiSigmaPlus) {
756 return std::string("s+b");
757 } else if(p == G4INCL::antiSigmaZero) {
758 return std::string("s0b");
759 } else if(p == G4INCL::antiSigmaMinus) {
760 return std::string("s-b");
761 } else if(p == G4INCL::antiLambda) {
762 return std::string("lb");
763 } else if(p == G4INCL::antiXiMinus) {
764 return std::string("x-b");
765 } else if(p == G4INCL::antiXiZero) {
766 return std::string("x0b");
767 } else if(p == G4INCL::KPlus) {
768 return std::string("k+");
769 } else if(p == G4INCL::KZero) {
770 return std::string("k0");
771 } else if(p == G4INCL::KZeroBar) {
772 return std::string("k0b");
773 } else if(p == G4INCL::KMinus) {
774 return std::string("k-");
775 } else if(p == G4INCL::KShort) {
776 return std::string("ks");
777 } else if(p == G4INCL::KLong) {
778 return std::string("kl");
779 } else if(p == G4INCL::Composite) {
780 return std::string("comp");
781 } else if(p == G4INCL::Eta) {
782 return std::string("eta");
783 } else if(p == G4INCL::Omega) {
784 return std::string("omega");
785 } else if(p == G4INCL::EtaPrime) {
786 return std::string("etap");
787 } else if(p == G4INCL::Photon) {
788 return std::string("photon");
789 }
790 return std::string("unknown");
791 }
792
794 if(pt == Proton) {
795 return protonMass;
796 } else if(pt == Neutron) {
797 return neutronMass;
798 } else if(pt == PiPlus) {
799 return piPlusMass;
800 } else if(pt == PiMinus) {
801 return piMinusMass;
802 } else if(pt == PiZero) {
803 return piZeroMass;
804 } else if(pt == SigmaPlus) {
805 return SigmaPlusMass;
806 } else if(pt == SigmaMinus) {
807 return SigmaMinusMass;
808 } else if(pt == SigmaZero) {
809 return SigmaZeroMass;
810 } else if(pt == Lambda) {
811 return LambdaMass;
812 } else if(pt == antiProton) {
813 return antiProtonMass;
814 } else if(pt == XiMinus) {
815 return XiMinusMass;
816 } else if(pt == XiZero) {
817 return XiZeroMass;
818 } else if(pt == antiNeutron) {
819 return antiNeutronMass;
820 } else if(pt == antiSigmaPlus) {
821 return antiSigmaPlusMass;
822 } else if(pt == antiSigmaMinus) {
823 return antiSigmaMinusMass;
824 } else if(pt == antiSigmaZero) {
825 return antiSigmaZeroMass;
826 } else if(pt == antiLambda) {
827 return antiLambdaMass;
828 } else if(pt == antiXiMinus) {
829 return antiXiMinusMass;
830 } else if(pt == antiXiZero) {
831 return antiXiZeroMass;
832 } else if(pt == KPlus) {
833 return KPlusMass;
834 } else if(pt == KZero) {
835 return KZeroMass;
836 } else if(pt == KZeroBar) {
837 return KZeroBarMass;
838 } else if(pt == KMinus) {
839 return KMinusMass;
840 } else if(pt == KShort) {
841 return KShortMass;
842 } else if(pt == KLong) {
843 return KLongMass;
844 } else if(pt == Eta) {
845 return etaMass;
846 } else if(pt == Omega) {
847 return omegaMass;
848 } else if(pt == EtaPrime) {
849 return etaPrimeMass;
850 } else if(pt == Photon) {
851 return photonMass;
852 } else {
853 INCL_ERROR("getMass : Unknown particle type." << '\n');
854 return 0.0;
855 }
856 }
857
859 switch(t) {
860 case Proton:
861 return theRealProtonMass;
862 break;
863 case Neutron:
864 return theRealNeutronMass;
865 break;
866 case PiPlus:
867 case PiMinus:
868 return theRealChargedPiMass;
869 break;
870 case PiZero:
871 return theRealPiZeroMass;
872 break;
873 case Eta:
874 return theRealEtaMass;
875 break;
876 case Omega:
877 return theRealOmegaMass;
878 break;
879 case EtaPrime:
880 return theRealEtaPrimeMass;
881 break;
882 case Photon:
883 return theRealPhotonMass;
884 break;
885 case Lambda:
886 return theRealLambdaMass;
887 break;
888 case KPlus:
889 case KMinus:
890 return theRealChargedKaonMass;
891 break;
892 case KZero:
893 case KZeroBar:
894 case KShort:
895 case KLong:
896 return theRealNeutralKaonMass;
897 break;
898 case SigmaPlus:
899 return theRealSigmaPlusMass;
900 break;
901 case SigmaZero:
902 return theRealSigmaZeroMass;
903 break;
904 case SigmaMinus:
905 return theRealSigmaMinusMass;
906 break;
907 case antiProton:
908 return theRealAntiProtonMass;
909 break;
910 case XiMinus:
911 return theRealXiMinusMass;
912 break;
913 case XiZero:
914 return theRealXiZeroMass;
915 break;
916 case antiNeutron:
917 return theRealAntiNeutronMass;
918 break;
919 case antiSigmaPlus:
920 return theRealAntiSigmaPlusMass;
921 break;
922 case antiSigmaZero:
923 return theRealAntiSigmaZeroMass;
924 break;
925 case antiSigmaMinus:
926 return theRealAntiSigmaMinusMass;
927 break;
928 case antiXiMinus:
929 return theRealAntiXiMinusMass;
930 break;
931 case antiXiZero:
932 return theRealAntiXiZeroMass;
933 break;
934 case antiLambda:
935 return theRealAntiLambdaMass;
936 break;
937 default:
938 INCL_ERROR("Particle::getRealMass : Unknown particle type." << '\n');
939 return 0.0;
940 break;
941 }
942 }
943
944 G4double getRealMass(const G4int A, const G4int Z, const G4int S) {
945// assert(A>=0);
946 // For nuclei with Z<0 or Z>A, assume that the exotic charge state is due to pions
947 if(Z<0 && S<0)
948 return (A+S)*theRealNeutronMass - S*LambdaMass - Z*getRealMass(PiMinus);
949 else if(Z>A && S<0)
950 return (A+S)*theRealProtonMass - S*LambdaMass + (A+S-Z)*getRealMass(PiPlus);
951 if(Z<0)
952 return (A)*theRealNeutronMass - Z*getRealMass(PiMinus);
953 else if(Z>A)
954 return (A)*theRealProtonMass + (A-Z)*getRealMass(PiPlus);
955 else if(Z==0 && S==0)
956 return A*theRealNeutronMass;
957 else if(A==Z)
958 return A*theRealProtonMass;
959 else if(Z==0 && S<0)
960 return (A+S)*theRealNeutronMass-S*LambdaMass;
961 else if(A>1) {
962#ifndef INCLXX_IN_GEANT4_MODE
963 return ::G4INCL::NuclearMassTable::getMass(A,Z,S);
964#else
965 if(S<0) return theG4IonTable->GetNucleusMass(Z,A,std::abs(S)) / MeV;
966 else return theG4IonTable->GetNucleusMass(Z,A) / MeV;
967#endif
968 } else
969 return 0.;
970 }
971
972 G4double getINCLMass(const G4int A, const G4int Z, const G4int S) {
973// assert(A>=0);
974 // For nuclei with Z<0 or Z>A, assume that the exotic charge state is due to pions
975 // Note that S<0 for lambda
976 if(Z<0 && S<0)
977 return (A+S)*neutronMass - S*LambdaMass - Z*getINCLMass(PiMinus);
978 else if(Z>A && S<0)
979 return (A+S)*protonMass - S*LambdaMass + (A+S-Z)*getINCLMass(PiPlus);
980 else if(Z<0)
981 return (A)*neutronMass - Z*getINCLMass(PiMinus);
982 else if(Z>A)
983 return (A)*protonMass + (A-Z)*getINCLMass(PiPlus);
984 else if(A>1 && S<0)
985 return Z*(protonMass - protonSeparationEnergy) + (A+S-Z)*(neutronMass - neutronSeparationEnergy) + std::abs(S)*(LambdaMass - lambdaSeparationEnergy);
986 else if(A>1)
987 return Z*(protonMass - protonSeparationEnergy) + (A-Z)*(neutronMass - neutronSeparationEnergy);
988 else if(A==1 && Z==0 && S==0)
989 return getINCLMass(Neutron);
990 else if(A==1 && Z==1 && S==0)
991 return getINCLMass(Proton);
992 else if(A==1 && Z==0 && S==-1)
993 return getINCLMass(Lambda);
994 else
995 return 0.;
996 }
997
998 G4double getTableQValue(const G4int A1, const G4int Z1, const G4int S1, const G4int A2, const G4int Z2, const G4int S2) {
999 return getTableMass(A1,Z1,S1) + getTableMass(A2,Z2,S2) - getTableMass(A1+A2,Z1+Z2,S1+S2);
1000 }
1001
1002 G4double getTableQValue(const G4int A1, const G4int Z1, const G4int S1, const G4int A2, const G4int Z2, const G4int S2, const G4int A3, const G4int Z3, const G4int S3) {
1003 return getTableMass(A1,Z1,S1) + getTableMass(A2,Z2,S2) - getTableMass(A3,Z3,S3) - getTableMass(A1+A2-A3,Z1+Z2-Z3,S1+S2-S3);
1004 }
1005
1007 if(p.theType == Composite)
1008 return (*getTableMass)(p.theA, p.theZ, p.theS);
1009 else
1010 return (*getTableParticleMass)(p.theType);
1011 }
1012
1014
1015 switch(t) {
1016 case Proton:
1017 case Neutron:
1018 case DeltaPlusPlus:
1019 case DeltaPlus:
1020 case DeltaZero:
1021 case DeltaMinus:
1022 case SigmaPlus:
1023 case SigmaZero:
1024 case SigmaMinus:
1025 case Lambda:
1026 case XiZero:
1027 case XiMinus:
1028 return 1;
1029 break;
1030 case antiProton:
1031 case antiNeutron:
1032 case antiSigmaPlus:
1033 case antiSigmaZero:
1034 case antiSigmaMinus:
1035 case antiLambda:
1036 case antiXiZero:
1037 case antiXiMinus:
1038 return -1;
1039 break;
1040 case PiPlus:
1041 case PiMinus:
1042 case PiZero:
1043 case KPlus:
1044 case KZero:
1045 case KZeroBar:
1046 case KShort:
1047 case KLong:
1048 case KMinus:
1049 case Eta:
1050 case Omega:
1051 case EtaPrime:
1052 case Photon:
1053 return 0;
1054 break;
1055 default:
1056 return 0;
1057 break;
1058 }
1059 }
1060
1062 switch(t) {
1063 case DeltaPlusPlus:
1064 return 2;
1065 break;
1066 case Proton:
1067 case DeltaPlus:
1068 case PiPlus:
1069 case SigmaPlus:
1070 case KPlus:
1071 case antiSigmaMinus:
1072 case antiXiMinus:
1073 return 1;
1074 break;
1075 case Neutron:
1076 case DeltaZero:
1077 case PiZero:
1078 case SigmaZero:
1079 case Lambda:
1080 case KZero:
1081 case KZeroBar:
1082 case KShort:
1083 case KLong:
1084 case Eta:
1085 case Omega:
1086 case EtaPrime:
1087 case Photon:
1088 case XiZero:
1089 case antiNeutron:
1090 case antiLambda:
1091 case antiSigmaZero:
1092 case antiXiZero:
1093 return 0;
1094 break;
1095 case DeltaMinus:
1096 case PiMinus:
1097 case SigmaMinus:
1098 case KMinus:
1099 case antiProton:
1100 case XiMinus:
1101 case antiSigmaPlus:
1102 return -1;
1103 break;
1104 default:
1105 return 0;
1106 break;
1107 }
1108 }
1109
1111 switch(t) {
1112 case DeltaPlusPlus:
1113 case DeltaPlus:
1114 case DeltaZero:
1115 case DeltaMinus:
1116 case Proton:
1117 case Neutron:
1118 case PiPlus:
1119 case PiZero:
1120 case PiMinus:
1121 case Eta:
1122 case Omega:
1123 case EtaPrime:
1124 case Photon:
1125 case antiProton:
1126 case antiNeutron:
1127 return 0;
1128 break;
1129 case XiMinus:
1130 case XiZero:
1131 case antiXiMinus:
1132 case antiXiZero:
1133 return 2;
1134 break;
1135 case antiLambda:
1136 case antiSigmaPlus:
1137 case antiSigmaZero:
1138 case antiSigmaMinus:
1139 return 1;
1140 break;
1141 case Lambda:
1142 case SigmaPlus:
1143 case SigmaZero:
1144 case SigmaMinus:
1145 case KZeroBar:
1146 case KMinus:
1147 return -1;
1148 break;
1149 case KPlus:
1150 case KZero:
1151 return 1;
1152 break;
1153 case KShort:
1154 return 0;
1155 break;
1156 case KLong:
1157 return 0;
1158 break;
1159 default:
1160 return 0;
1161 break;
1162 }
1163 }
1164
1166// assert(A>=0);
1167 if(A > 19 || (A < 6 && A >= 2)) {
1168 // For large (Woods-Saxon or Modified Harmonic Oscillator) or small
1169 // (Gaussian) nuclei, the radius parameter is just the nuclear radius
1170 return getRadiusParameter(t,A,Z);
1171 } else if(A < clusterTableASize && Z>=0 && Z < clusterTableZSize && A >= 6) {
1172 const G4double thisRMS = positionRMS[Z][A];
1173 if(thisRMS>0.0)
1174 return thisRMS;
1175 else {
1176 INCL_DEBUG("getNuclearRadius: Radius for nucleus A = " << A << " Z = " << Z << " is not available" << '\n'
1177 << "returning radius for C12");
1178 return positionRMS[6][12];
1179 }
1180 } else if(A <= 19) {
1181 const G4double theRadiusParameter = getRadiusParameter(t, A, Z);
1182 const G4double theDiffusenessParameter = getSurfaceDiffuseness(t, A, Z);
1183 // The formula yields the nuclear RMS radius based on the parameters of
1184 // the nuclear-density function
1185 return 1.225*theDiffusenessParameter*
1186 std::sqrt((2.+5.*theRadiusParameter)/(2.+3.*theRadiusParameter));
1187 } else {
1188 INCL_ERROR("getNuclearRadius: No radius for nucleus A = " << A << " Z = " << Z << '\n');
1189 return 0.0;
1190 }
1191 }
1192
1196
1198// assert(A>0);
1199 if(A > 19) {
1200 // radius fit for lambdas
1201 if(t==Lambda){
1202 G4double r0 = (1.128+0.439*std::pow(A,-2./3.)) * std::pow(A, 1.0/3.0);
1203 return r0;
1204 }
1205 // phenomenological radius fit
1206 G4double r0 = (2.745e-4 * A + 1.063) * std::pow(A, 1.0/3.0);
1207 // HFB calculations
1210 if(r0hfb>0.)r0 = r0hfb;
1211 }
1212 //
1213 if(t==Neutron)
1214 r0 += neutronSkin;
1215 return r0;
1216 } else if(A < 6 && A >= 2) {
1217 if(Z<clusterTableZSize && Z>=0) {
1218 const G4double thisRMS = positionRMS[Z][A];
1219 if(thisRMS>0.0)
1220 return thisRMS;
1221 else {
1222 INCL_DEBUG("getRadiusParameter: Radius for nucleus A = " << A << " Z = " << Z << " is not available" << '\n'
1223 << "returning radius for C12");
1224 return positionRMS[6][12];
1225 }
1226 } else {
1227 INCL_DEBUG("getRadiusParameter: Radius for nucleus A = " << A << " Z = " << Z << " is not available" << '\n'
1228 << "returning radius for C12");
1229 return positionRMS[6][12];
1230 }
1231 } else if(A <= 19 && A >= 6) {
1232 if(t==Lambda){
1233 G4double r0 = (1.128+0.439*std::pow(A,-2./3.)) * std::pow(A, 1.0/3.0);
1234 return r0;
1235 }
1236 // HFB calculations
1239 if(r0hfb>0.)return r0hfb;
1240 }
1241 return mediumRadius[A-1];
1242 // return 1.581*mediumDiffuseness[A-1]*(2.+5.*mediumRadius[A-1])/(2.+3.*mediumRadius[A-1]);
1243 } else {
1244 INCL_ERROR("getRadiusParameter: No radius for nucleus A = " << A << " Z = " << Z << '\n');
1245 return 0.0;
1246 }
1247 }
1248
1250 const G4double XFOISA = 8.0;
1251 if(A > 19) {
1252 return getNuclearRadius(t,A,Z) + XFOISA * getSurfaceDiffuseness(t,A,Z);
1253 } else if(A <= 19 && A >= 6) {
1254 return 5.5 + 0.3 * (G4double(A) - 6.0)/12.0;
1255 } else if(A >= 2) {
1256 return getNuclearRadius(t, A, Z) + 4.5;
1257 } else {
1258 INCL_ERROR("getMaximumNuclearRadius : No maximum radius for nucleus A = " << A << " Z = " << Z << '\n');
1259 return 0.0;
1260 }
1261 }
1262
1264 if(A > 19) {
1265 // phenomenological fit
1266 G4double a = 1.63e-4 * A + 0.510;
1267 // HFB calculations
1270 if(ahfb>0.)a=ahfb;
1271 }
1272 //
1273 if(t==Lambda){
1274 // Like for neutrons
1276 if(ahfb>0.)a=ahfb;
1277 }
1278 if(t==Neutron)
1279 a += neutronHalo;
1280 return a;
1281 } else if(A <= 19 && A >= 6) {
1282 // HFB calculations
1285 if(ahfb>0.)return ahfb;
1286 }
1287 return mediumDiffuseness[A-1];
1288 } else if(A < 6 && A >= 2) {
1289 INCL_ERROR("getSurfaceDiffuseness: was called for A = " << A << " Z = " << Z << '\n');
1290 return 0.0;
1291 } else {
1292 INCL_ERROR("getSurfaceDiffuseness: No diffuseness for nucleus A = " << A << " Z = " << Z << '\n');
1293 return 0.0;
1294 }
1295 }
1296
1298// assert(Z>=0 && A>=0 && Z<=A);
1300 }
1301
1302 G4double getSeparationEnergyINCL(const ParticleType t, const G4int /*A*/, const G4int /*Z*/) {
1303 if(t==Proton)
1304 return theINCLProtonSeparationEnergy;
1305 else if(t==Neutron)
1306 return theINCLNeutronSeparationEnergy;
1307 else if(t==Lambda)
1308 return theINCLLambdaSeparationEnergy;
1309 else if(t==antiProton)
1310 return theINCLantiProtonSeparationEnergy;
1311 else {
1312 INCL_ERROR("ParticleTable::getSeparationEnergyINCL : Unknown particle type." << '\n');
1313 return 0.0;
1314 }
1315 }
1316
1318 // Real separation energies for all nuclei
1319 if(t==Proton)
1320 return (*getTableParticleMass)(Proton) + (*getTableMass)(A-1,Z-1,0) - (*getTableMass)(A,Z,0);
1321 else if(t==Neutron)
1322 return (*getTableParticleMass)(Neutron) + (*getTableMass)(A-1,Z,0) - (*getTableMass)(A,Z,0);
1323 else if(t==Lambda)
1324 return (*getTableParticleMass)(Lambda) + (*getTableMass)(A-1,Z,0) - (*getTableMass)(A,Z,-1);
1325 else {
1326 INCL_ERROR("ParticleTable::getSeparationEnergyReal : Unknown particle type." << '\n');
1327 return 0.0;
1328 }
1329 }
1330
1332 // Real separation energies for light nuclei, fixed values for heavy nuclei
1334 return getSeparationEnergyReal(t, A, Z);
1335 else
1336 return getSeparationEnergyINCL(t, A, Z);
1337 }
1338
1339 G4double getProtonSeparationEnergy() { return protonSeparationEnergy; }
1340
1341 G4double getNeutronSeparationEnergy() { return neutronSeparationEnergy; }
1342
1343 G4double getLambdaSeparationEnergy() { return lambdaSeparationEnergy; }
1344
1345 void setProtonSeparationEnergy(const G4double sen) { protonSeparationEnergy = sen; }
1346
1347 void setNeutronSeparationEnergy(const G4double sen) { neutronSeparationEnergy = sen; }
1348
1349 void setLambdaSeparationEnergy(const G4double sen) { lambdaSeparationEnergy = sen; }
1350
1351 std::string getElementName(const G4int Z) {
1352 if(Z<1) {
1353 INCL_WARN("getElementName called with Z<1" << '\n');
1354 return elementTable[0];
1355 } else if(Z<elementTableSize)
1356 return elementTable[Z];
1357 else
1358 return getIUPACElementName(Z);
1359 }
1360
1361 std::string getIUPACElementName(const G4int Z) {
1362 std::stringstream elementStream;
1363 elementStream << Z;
1364 std::string elementName = elementStream.str();
1365 std::transform(elementName.begin(), elementName.end(), elementName.begin(), intToIUPAC);
1366 elementName[0] = (char)std::toupper(elementName.at(0));
1367 return elementName;
1368 }
1369
1370 G4int parseElement(std::string pS) {
1371 // Normalize the element name
1372 std::transform(pS.begin(), pS.end(), pS.begin(), ::tolower);
1373 pS[0] = (char)std::toupper(pS[0]);
1374
1375 const std::string *iter = std::find(elementTable, elementTable+elementTableSize, pS);
1376 if(iter != elementTable+elementTableSize)
1377 return G4int(iter - elementTable);
1378 else
1380 }
1381
1382 G4int parseIUPACElement(std::string const &sel) {
1383 // Normalise to lower case
1384 std::string elementName(sel);
1385 std::transform(elementName.begin(), elementName.end(), elementName.begin(), ::tolower);
1386 // Return 0 if the element name contains anything but IUPAC digits
1387 if(elementName.find_first_not_of(elementIUPACDigits)!=std::string::npos)
1388 return 0;
1389 std::transform(elementName.begin(), elementName.end(), elementName.begin(), iupacToInt);
1390 std::stringstream elementStream(elementName);
1391 G4int Z;
1392 elementStream >> Z;
1393 return Z;
1394 }
1395
1397 return getNaturalIsotopicDistributions()->getIsotopicDistribution(Z);
1398 }
1399
1401 return getNaturalIsotopicDistributions()->drawRandomIsotope(Z);
1402 }
1403
1404 G4double getFermiMomentumConstant(const G4int /*A*/, const G4int /*Z*/) {
1405 return constantFermiMomentum;
1406 }
1407
1409// assert(Z>0 && A>0 && Z<=A);
1411 const G4double rms = momentumRMS[Z][A];
1412 return ((rms>0.) ? rms : momentumRMS[6][12]) * Math::sqrtFiveThirds;
1413 } else
1414 return getFermiMomentumConstant(A,Z);
1415 }
1416
1418// assert(A>0);
1419 static const G4double alphaParam = 259.416; // MeV/c
1420 static const G4double betaParam = 152.824; // MeV/c
1421 static const G4double gammaParam = 9.5157E-2;
1422 return alphaParam - betaParam*std::exp(-gammaParam*((G4double)A));
1423 }
1424
1426// assert(t==Proton || t==Neutron || t==Lambda);
1427 return rpCorrelationCoefficient[t];
1428 }
1429
1430 G4double getNeutronSkin() { return neutronSkin; }
1431
1432 G4double getNeutronHalo() { return neutronHalo; }
1433
1441
1443// assert(isosp == -2 || isosp == 0 || isosp == 2);
1444 if (isosp == -2) {
1445 return PiMinus;
1446 }
1447 else if (isosp == 0) {
1448 return PiZero;
1449 }
1450 else {
1451 return PiPlus;
1452 }
1453 }
1454
1456// assert(isosp == -1 || isosp == 1);
1457 if (isosp == -1) {
1458 return Neutron;
1459 }
1460 else {
1461 return Proton;
1462 }
1463 }
1464
1466// assert(isosp == -3 || isosp == -1 || isosp == 1 || isosp == 3);
1467 if (isosp == -3) {
1468 return DeltaMinus;
1469 }
1470 else if (isosp == -1) {
1471 return DeltaZero;
1472 }
1473 else if (isosp == 1) {
1474 return DeltaPlus;
1475 }
1476 else {
1477 return DeltaPlusPlus;
1478 }
1479 }
1480
1482// assert(isosp == -2 || isosp == 0 || isosp == 2);
1483 if (isosp == -2) {
1484 return SigmaMinus;
1485 }
1486 else if (isosp == 0) {
1487 return SigmaZero;
1488 }
1489 else {
1490 return SigmaPlus;
1491 }
1492 }
1493
1495// assert(isosp == -1 || isosp == 1);
1496 if (isosp == -1) {
1497 return XiMinus;
1498 }
1499 else {
1500 return XiZero;
1501 }
1502 }
1503
1504 /*ParticleType getAntiNucleonType(const G4int isosp) {
1505// assert(isosp == -1); //|| isosp == 1
1506 if (isosp == -1) {
1507 return antiProton;
1508 }
1509 else {
1510 return antiNeutron;
1511 }
1512 }*/
1513
1515// assert(isosp == -2 || isosp == 0 || isosp == 2);
1516 if (isosp == -2) {
1517 return antiSigmaPlus;
1518 }
1519 else if (isosp == 0) {
1520 return antiSigmaZero;
1521 }
1522 else {
1523 return antiSigmaMinus;
1524 }
1525 }
1526
1528// assert(isosp == -1 || isosp == 1);
1529 if (isosp == -1) {
1530 return antiXiZero;
1531 }
1532 else {
1533 return antiXiMinus;
1534 }
1535 }
1536
1538// assert(isosp == -1 || isosp == 1);
1539 if (isosp == -1) {
1540 return KZero;
1541 }
1542 else {
1543 return KPlus;
1544 }
1545 }
1546
1548// assert(isosp == -1 || isosp == 1);
1549 if (isosp == -1) {
1550 return KMinus;
1551 }
1552 else {
1553 return KZeroBar;
1554 }
1555 }
1556
1558// assert(pt == PiPlus || pt == PiMinus || pt == PiZero || pt == Eta || pt == Omega || pt == EtaPrime || pt == KShort || pt == KLong || pt== KPlus || pt == KMinus || pt == Lambda || pt == SigmaPlus || pt == SigmaZero || pt == SigmaMinus || pt == antiLambda || pt == antiSigmaPlus || pt == antiSigmaZero || pt == antiSigmaMinus || pt == XiMinus || pt == XiZero || pt == antiXiZero || pt == antiXiMinus || );
1559 if(pt == PiPlus) {
1560 return piPlusWidth;
1561 } else if(pt == PiMinus) {
1562 return piMinusWidth;
1563 } else if(pt == PiZero) {
1564 return piZeroWidth;
1565 } else if(pt == Eta) {
1566 return etaWidth;
1567 } else if(pt == Omega) {
1568 return omegaWidth;
1569 } else if(pt == EtaPrime) {
1570 return etaPrimeWidth;
1571 } else if(pt == SigmaPlus) {
1572 return SigmaPlusWidth;
1573 } else if(pt == SigmaZero) {
1574 return SigmaZeroWidth;
1575 } else if(pt == SigmaMinus) {
1576 return SigmaMinusWidth;
1577 } else if(pt == KPlus) {
1578 return KPlusWidth;
1579 } else if(pt == KMinus) {
1580 return KMinusWidth;
1581 } else if(pt == KShort) {
1582 return KShortWidth;
1583 } else if(pt == KLong) {
1584 return KLongWidth;
1585 } else if(pt == Lambda) {
1586 return LambdaWidth;
1587 } else if(pt == XiMinus) {
1588 return XiMinusWidth;
1589 } else if(pt == XiZero) {
1590 return XiZeroWidth;
1591 } else if(pt == antiSigmaPlus) {
1592 return antiSigmaPlusWidth;
1593 } else if(pt == antiSigmaZero) {
1594 return antiSigmaZeroWidth;
1595 } else if(pt == antiSigmaMinus) {
1596 return antiSigmaMinusWidth;
1597 } else if(pt == antiLambda) {
1598 return antiLambdaWidth;
1599 } else if(pt == antiXiMinus) {
1600 return antiXiMinusWidth;
1601 } else if(pt == antiXiZero) {
1602 return antiXiZeroWidth;
1603 } else {
1604 INCL_ERROR("getWidth : Unknown particle type." << '\n');
1605 return 0.0;
1606 }
1607 }
1608
1609 } // namespace ParticleTable
1610} // namespace G4INCL
1611
G4double S(G4double temp)
#define INCL_ERROR(x)
#define INCL_WARN(x)
#define INCL_FATAL(x)
#define INCL_DEBUG(x)
Functions that encapsulate a mass table.
#define INCL_DEFAULT_SEPARATION_ENERGY
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4double getNeutronHalo() const
Get the neutron-halo size.
FermiMomentumType getFermiMomentumType() const
Get the Fermi-momentum type.
SeparationEnergyType getSeparationEnergyType() const
Get the separation-energy type.
G4double getRPCorrelationCoefficient(const ParticleType t) const
Get the r-p correlation coefficient.
std::string const & getINCLXXDataFilePath() const
Set the ABLAXX datafile path.
G4double getNeutronSkin() const
Get the neutron-skin thickness.
G4double getFermiMomentum() const
Get the Fermi momentum.
G4bool getUseRealMasses() const
Whether to use real masses.
Class that stores isotopic abundances for a given element.
G4int drawRandomIsotope() const
Draw a random isotope based on the abundance vector.
G4IonTable * GetIonTable() const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
void initialize()
Definition G4INCLHFB.cc:83
G4double getSurfaceDiffusenessHFB(const ParticleType t, const G4int A, const G4int Z)
Definition G4INCLHFB.cc:142
G4double getRadiusParameterHFB(const ParticleType t, const G4int A, const G4int Z)
Get the radius and diffuseness parameters from HFB calculations.
Definition G4INCLHFB.cc:132
T max(const T t1, const T t2)
brief Return the largest of the two arguments
const G4double sqrtFiveThirds
const G4double sqrtThreeFifths
G4int getMassNumber(const ParticleType t)
Get mass number from particle type.
G4ThreadLocal FermiMomentumFn getFermiMomentum
const G4double effectiveDeltaWidth
G4int parseElement(std::string pS)
Get the name of the element from the atomic number.
G4ThreadLocal G4double minDeltaMass2
G4double(* FermiMomentumFn)(const G4int, const G4int)
G4ThreadLocal NuclearMassFn getTableMass
Static pointer to the mass function for nuclei.
G4ThreadLocal SeparationEnergyFn getSeparationEnergy
Static pointer to the separation-energy function.
G4double getTableQValue(const G4int A1, const G4int Z1, const G4int S1, const G4int A2, const G4int Z2, const G4int S2)
Get Q-value (in MeV/c^2)
G4ThreadLocal ParticleMassFn getTableParticleMass
Static pointer to the mass function for particles.
void initialize(Config const *const theConfig=0)
Initialize the particle table.
const G4double effectiveDeltaMass
G4double getFermiMomentumMassDependent(const G4int A, const G4int)
Return the value Fermi momentum from a fit.
G4double getTableSpeciesMass(const ParticleSpecies &p)
G4int drawRandomNaturalIsotope(const G4int Z)
G4double getSeparationEnergyReal(const ParticleType t, const G4int A, const G4int Z)
Return the real separation energy.
G4double getNeutronSeparationEnergy()
Getter for neutronSeparationEnergy.
G4ThreadLocal G4double minDeltaMass
G4double getRadiusParameter(const ParticleType t, const G4int A, const G4int Z)
G4double getLargestNuclearRadius(const G4int A, const G4int Z)
ParticleType getKaonType(const G4int isosp)
Get the type of kaon.
G4double getNeutronHalo()
Get the size of the neutron halo.
G4double getRealMass(const G4INCL::ParticleType t)
Get particle mass (in MeV/c^2)
ParticleType getSigmaType(const G4int isosp)
Get the type of sigma.
G4double getINCLMass(const G4int A, const G4int Z, const G4int S)
Get INCL nuclear mass (in MeV/c^2)
G4double(* ParticleMassFn)(const ParticleType)
G4int getStrangenessNumber(const ParticleType t)
Get strangeness number from particle type.
G4double getMaximumNuclearRadius(const ParticleType t, const G4int A, const G4int Z)
G4double getRPCorrelationCoefficient(const ParticleType t)
Get the value of the r-p correlation coefficient.
G4int parseIUPACElement(std::string const &pS)
Parse a IUPAC element name.
G4double getSeparationEnergyINCL(const ParticleType t, const G4int, const G4int)
Return INCL's default separation energy.
void setNeutronSeparationEnergy(const G4double s)
Setter for protonSeparationEnergy.
G4double getFermiMomentumConstant(const G4int, const G4int)
Return the constant value of the Fermi momentum.
std::string getName(const ParticleType t)
Get the native INCL name of the particle.
G4ThreadLocal G4double minDeltaMassRndm
G4double(* SeparationEnergyFn)(const ParticleType, const G4int, const G4int)
G4double getNeutronSkin()
Get the thickness of the neutron skin.
std::string getIUPACElementName(const G4int Z)
Get the name of an unnamed element from the IUPAC convention.
G4int getIsospin(const ParticleType t)
Get the isospin of a particle.
ParticleType getNucleonType(const G4int isosp)
Get the type of nucleon.
ParticleType getAntiXiType(const G4int isosp)
Get the type of antidelta.
G4double getSurfaceDiffuseness(const ParticleType t, const G4int A, const G4int Z)
G4double getFermiMomentumConstantLight(const G4int A, const G4int Z)
Return the constant value of the Fermi momentum - special for light.
void setProtonSeparationEnergy(const G4double s)
Setter for protonSeparationEnergy.
ParticleType getPionType(const G4int isosp)
Get the type of pion.
ParticleType getDeltaType(const G4int isosp)
Get the type of delta.
void setLambdaSeparationEnergy(const G4double sen)
G4double(* NuclearMassFn)(const G4int, const G4int, const G4int)
G4int getChargeNumber(const ParticleType t)
Get charge number from particle type.
G4double getProtonSeparationEnergy()
Getter for protonSeparationEnergy.
IsotopicDistribution const & getNaturalIsotopicDistribution(const G4int Z)
ParticleType getAntiSigmaType(const G4int isosp)
Get the type of antisigma.
G4double getMomentumRMS(const G4int A, const G4int Z)
Return the RMS of the momentum distribution (light clusters)
ParticleType getAntiKaonType(const G4int isosp)
Get the type of antikaon.
G4double getSeparationEnergyRealForLight(const ParticleType t, const G4int A, const G4int Z)
Return the real separation energy only for light nuclei.
G4double getNuclearRadius(const ParticleType t, const G4int A, const G4int Z)
G4double getWidth(const ParticleType t)
Get particle width (in s)
std::string getShortName(const ParticleType t)
Get the short INCL name of the particle.
std::string getElementName(const G4int Z)
Get the name of the element from the atomic number.
ParticleType getXiType(const G4int isosp)
Get the type of xi.
const G4double Pf
Fermi momentum [MeV/c].
@ MassDependentFermiMomentum
@ ConstantLightFermiMomentum
@ RealForLightSeparationEnergy
#define G4ThreadLocal
Definition tls.hh:77