Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4INCLConfig.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// 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
38#ifndef G4INCLConfig_hh
39#define G4INCLConfig_hh 1
40
42#include "G4INCLConfigEnums.hh"
44#include <iostream>
45#include <string>
46#include <sstream>
47// #include <cassert>
48
49class ConfigParser;
50
51namespace G4INCL {
52
53 /**
54 * The INCL configuration object
55 *
56 * The Config object keeps track of various INCL physics options
57 * (e.g. which Pauli blocking scheme to use, whether to use local
58 * energy option or not, etc.
59 */
60 class Config {
61 public:
62 /// \brief Default constructor
63 Config();
64
65 /// \brief Default destructor
66 ~Config();
67
68 /// \brief Initialise the members
69 void init();
70
71 /// \brief Return a summary of the run configuration.
72 std::string summary();
73
74 /// \brief Get the verbosity.
75 G4int getVerbosity() const { return verbosity; }
76
77 /// \brief Get the run title.
78 std::string const &getCalculationTitle() const { return title; }
79
80 /// \brief Get the output file root.
81 std::string const &getOutputFileRoot() const { return outputFileRoot; }
82
83 /// \brief Get the number of shots.
84 G4int getNumberOfShots() const { return nShots; }
85
86 /// \brief Natural targets.
87 G4bool isNaturalTarget() const { return naturalTarget; }
88
89 /** \brief Get the target mass number.
90 *
91 * Note that A==0 means natural target. You should first check the
92 * isNaturalTarget() method.
93 */
94 G4int getTargetA() const { return targetSpecies.theA; }
95
96 /// \brief Get the target charge number.
97 G4int getTargetZ() const { return targetSpecies.theZ; }
98
99 /// \brief Get the target strangess number.
100 G4int getTargetS() const { return targetSpecies.theS; }
101
102 /// \brief Set target mass number
103 void setTargetA(G4int A) { targetSpecies.theA = A; }
104
105 /// \brief Set target charge number
106 void setTargetZ(G4int Z) { targetSpecies.theZ = Z; }
107
108 /// \brief Set target strangess number
109 void setTargetS(G4int S) { targetSpecies.theS = S; }
110
111 /// \brief Get the projectile type
112 ParticleType getProjectileType() const { return projectileSpecies.theType; }
113
114 /// \brief Get the projectile species
115 ParticleSpecies getProjectileSpecies() const { return projectileSpecies; }
116
117 /// \brief Set the projectile species
118 void setProjectileSpecies(ParticleSpecies const &pars) { projectileSpecies=pars; }
119
120 /// \brief Get the projectile kinetic energy.
121 G4double getProjectileKineticEnergy() const { return projectileKineticEnergy; }
122
123 /// \brief Set the projectile kinetic energy.
124 void setProjectileKineticEnergy(G4double const kinE) { projectileKineticEnergy=kinE; }
125
126 /// \brief Get the number of the verbose event.
127 G4int getVerboseEvent() const { return verboseEvent; }
128
129 /// \brief Get the INCL version ID.
130 static std::string const getVersionID();
131
132 /// \brief Get the INCL version hash.
133 static std::string const getVersionHash();
134
135 /// \brief Get the INCL version string.
136 static std::string const getVersionString() {
137 std::stringstream ss;
138 ss << getVersionID() << "-" << getVersionHash();
139 return ss.str();
140 }
141
142 /// \brief Get the seeds for the random-number generator.
144 return randomSeedVector;
145 }
146
147 /// \brief Get the Pauli-blocking algorithm.
148 PauliType getPauliType() const { return pauliType; }
149
150 /// \brief Do we want CDPP?
151 G4bool getCDPP() const { return CDPP; }
152
153 /// \brief Get the Coulomb-distortion algorithm.
154 CoulombType getCoulombType() const { return coulombType; }
155
156 /// \brief Set the Coulomb-distortion algorithm.
157 void setCoulombType(CoulombType const c) { coulombType = c; }
158
159 /// \brief Get the type of the potential for nucleons.
160 PotentialType getPotentialType() const { return potentialType; }
161
162 /// \brief Set the type of the potential for nucleons.
163 void setPotentialType(PotentialType type) { potentialType = type; }
164
165 /// \brief Do we want the pion potential?
166 G4bool getPionPotential() const { return pionPotential; }
167
168 /// \brief Set the type of the potential for nucleons.
169 void setPionPotential(const G4bool pionPot) { pionPotential = pionPot; }
170
171 /// \brief Get the type of local energy for N-N avatars.
172 LocalEnergyType getLocalEnergyBBType() const { return localEnergyBBType; }
173
174 /// \brief Set the type of local energy for N-N avatars.
175 void setLocalEnergyBBType(const LocalEnergyType t) { localEnergyBBType=t; }
176
177 /// \brief Get the type of local energy for pi-N and decay avatars.
178 LocalEnergyType getLocalEnergyPiType() const { return localEnergyPiType; }
179
180 /// \brief Set the type of local energy for N-N avatars.
181 void setLocalEnergyPiType(const LocalEnergyType t) { localEnergyPiType=t; }
182
183 /// \brief Get the log file name.
184 std::string const &getLogFileName() const { return logFileName; }
185
186 /// \brief Get the de-excitation model.
187 DeExcitationType getDeExcitationType() const { return deExcitationType; }
188
189 /// \brief Get the de-excitation string.
190 std::string getDeExcitationString() const { return deExcitationString; }
191
192 /// \brief Get the clustering algorithm.
193 ClusterAlgorithmType getClusterAlgorithm() const { return clusterAlgorithmType; }
194
195 /// \brief Set the clustering algorithm.
196 void setClusterAlgorithm(ClusterAlgorithmType const c) { clusterAlgorithmType = c; }
197
198 /// \brief Get the maximum mass for production of clusters.
199 G4int getClusterMaxMass() const { return clusterMaxMass; }
200
201 /// \brief Set the maximum mass for production of clusters.
202 void setClusterMaxMass(const G4int clm){ clusterMaxMass=clm; }
203
204 /// \brief Get back-to-spectator
205 G4bool getBackToSpectator() const { return backToSpectator; }
206
207 /// \brief Set back-to-spectator
208 void setBackToSpectator(const G4bool b) { backToSpectator = b; }
209
210 /// \brief Whether to use real masses
211 G4bool getUseRealMasses() const { return useRealMasses; }
212
213 /// \brief Set whether to use real masses
214 void setUseRealMasses(G4bool use) { useRealMasses = use; }
215
216 /// \brief Set the INCLXX datafile path
217 void setINCLXXDataFilePath(std::string const &path) { INCLXXDataFilePath=path; }
218
219 /// \brief Set the ABLAXX datafile path
220#ifdef INCL_DEEXCITATION_ABLAXX
221 void setABLAXXDataFilePath(std::string const &path) { ablaxxDataFilePath=path; }
222#endif
223
224 std::string const &getINCLXXDataFilePath() const {
225 return INCLXXDataFilePath;
226 }
227
228#ifdef INCL_DEEXCITATION_ABLAXX
229 std::string const &getABLAXXDataFilePath() const {
230 return ablaxxDataFilePath;
231 }
232#endif
233
234#ifdef INCL_DEEXCITATION_ABLA07
235 std::string const &getABLA07DataFilePath() const {
236 return abla07DataFilePath;
237 }
238#endif
239#ifdef INCL_DEEXCITATION_GEMINIXX
240 std::string const &getGEMINIXXDataFilePath() const {
241 return geminixxDataFilePath;
242 }
243#endif
244
245 G4double getImpactParameter() const { return impactParameter; }
246
247 /// \brief Get the separation-energy type
248 SeparationEnergyType getSeparationEnergyType() const { return separationEnergyType; }
249
250 /// \brief Get the Fermi-momentum type
251 FermiMomentumType getFermiMomentumType() const { return fermiMomentumType; }
252
253 /// \brief Set the Fermi-momentum type
254 void setFermiMomentumType(FermiMomentumType const f) { fermiMomentumType=f; }
255
256 /// \brief Get the Fermi momentum
257 G4double getFermiMomentum() const { return fermiMomentum; }
258
259 /// \brief Set the Fermi momentum
260 void setFermiMomentum(const G4double p) { fermiMomentum = p; }
261
262 G4double getCutNN() const { return cutNN; }
263
264#ifdef INCL_ROOT_USE
265 std::string const &getROOTSelectionString() const {
266 return rootSelectionString;
267 }
268#endif
269
270#ifdef INCL_DEEXCITATION_FERMI_BREAKUP
271 G4int getMaxMassFermiBreakUp() const {
272 return maxMassFermiBreakUp;
273 }
274
275 G4int getMaxChargeFermiBreakUp() const {
276 return maxChargeFermiBreakUp;
277 }
278#endif
279
280 /// \brief Get the r-p correlation coefficient
282// assert(t==Proton || t==Neutron);
283 return ((t==Proton) ? rpCorrelationCoefficientProton : rpCorrelationCoefficientNeutron);
284 }
285
286 /// \brief Set the r-p correlation coefficient
287 void setRPCorrelationCoefficient(const ParticleType t, const G4double corrCoeff) {
288// assert(t==Proton || t==Neutron);
289 if(t==Proton)
290 rpCorrelationCoefficientProton=corrCoeff;
291 else
292 rpCorrelationCoefficientNeutron=corrCoeff;
293 }
294
295 /// \brief Set the r-p correlation coefficient
300
301 /// \brief Get the neutron-skin thickness
302 G4double getNeutronSkin() const { return neutronSkin; }
303
304 /// \brief Set the neutron-skin thickness
305 void setNeutronSkin(const G4double d) { neutronSkin=d; }
306
307 /// \brief Get the neutron-halo size
308 G4double getNeutronHalo() const { return neutronHalo; }
309
310 /// \brief Set the neutron-halo size
311 void setNeutronHalo(const G4double d) { neutronHalo=d; }
312
313 /// \brief True if we should use refraction
314 G4bool getRefraction() const { return refraction; }
315
316 /// \brief Set the refraction variable
317 void setRefraction(const G4bool r) { refraction = r; }
318
319 /// \brief Get the RNG type
320 RNGType getRNGType() const { return rngType; }
321
322 /// \brief Set the RNG type
323 void setRNGType(RNGType const r) { rngType=r; }
324
325 /// \brief Get the phase-space-generator type
326 PhaseSpaceGeneratorType getPhaseSpaceGeneratorType() const { return phaseSpaceGeneratorType; }
327
328 /// \brief Set the phase-space-generator type
329 void setPhaseSpaceGeneratorType(PhaseSpaceGeneratorType const p) { phaseSpaceGeneratorType=p; }
330
331 /// \brief Get the cascade-action type
332 CascadeActionType getCascadeActionType() const { return cascadeActionType; }
333
334 /// \brief Set the cascade-action type
335 void setCascadeActionType(CascadeActionType const c) { cascadeActionType=c; }
336
337 /// \brief Get the autosave frequency
338 unsigned int getAutosaveFrequency() const { return autosaveFrequency; }
339
340 /// \brief Set the autosave frequency
341 void setAutosaveFrequency(const unsigned int f) { autosaveFrequency=f; }
342
343 /// \brief Get the Cross Section type
344 CrossSectionsType getCrossSectionsType() const { return crossSectionsType; }
345
346 /// \brief Get the maximum number of pions for multipion collisions
347 G4int getMaxNumberMultipions() const { return maxNumberMultipions; }
348
349 /// \brief Set the maximum number of pions for multipion collisions
350 void setMaxNumberMultipions(const G4int n) { maxNumberMultipions=n; }
351
352 /// \brief Set the Cross Section type
353 void setCrossSectionsType(CrossSectionsType const c) { crossSectionsType=c; }
354
355 /// \brief Get the hadronization time
356 G4double getHadronizationTime() const { return hadronizationTime; }
357
358 /// \brief Set the hadronization time
359 void setHadronizationTime(const G4double t) { hadronizationTime=t; }
360
361#ifdef INCL_ROOT_USE
362 G4bool getConciseROOTTree() const { return conciseROOTTree; }
363#endif
364
365 G4bool getInverseKinematics() const { return inverseKinematics; }
366
367 G4bool getsrcPairConfig() const { return srcPairCorrelations; }
368
369 G4float getsrcPairDist() const { return srcPairDistance; }
370
371 /// \brief Get the decay time threshold time
372 G4double getDecayTimeThreshold() const { return decayTimeThreshold; }
373
374 /// \brief Set decay time threshold time
375 void setDecayTimeThreshold(const G4double t) { decayTimeThreshold=t; }
376
377 /// \brief Get the bias
378 G4double getBias() const { return bias; }
379
380 /// \brief Get the pbar at rest annihilation threshold
381 G4double getAtrestThreshold() const { return atrestThreshold; }
382
383 /// \brief Set the pbar at rest annihilation threshold
384 void setAtrestThreshold(const G4double t) { atrestThreshold=t; }
385
386 private:
387
388 G4int verbosity;
389 std::string inputFileName;
390 std::string title;
391 std::string outputFileRoot;
392 std::string fileSuffix;
393 std::string logFileName;
394
395 G4int nShots;
396
397 std::string targetString;
398 ParticleSpecies targetSpecies;
399 G4bool naturalTarget;
400
401 std::string projectileString;
402 ParticleSpecies projectileSpecies;
403 G4double projectileKineticEnergy;
404
405 G4int verboseEvent;
406
407 std::string randomSeeds;
408 Random::SeedVector randomSeedVector;
409
410 std::string pauliString;
411 PauliType pauliType;
412 G4bool CDPP;
413
414 std::string coulombString;
415 CoulombType coulombType;
416
417 std::string potentialString;
418 PotentialType potentialType;
419 G4bool pionPotential;
420
421 std::string localEnergyBBString;
422 LocalEnergyType localEnergyBBType;
423
424 std::string localEnergyPiString;
425 LocalEnergyType localEnergyPiType;
426
427 std::string deExcitationModelList;
428 std::string deExcitationOptionDescription;
429 std::string deExcitationString;
430 DeExcitationType deExcitationType;
431#ifdef INCL_DEEXCITATION_ABLAXX
432 std::string ablaxxDataFilePath;
433#endif
434#ifdef INCL_DEEXCITATION_ABLA07
435 std::string abla07DataFilePath;
436#endif
437#ifdef INCL_DEEXCITATION_GEMINIXX
438 std::string geminixxDataFilePath;
439#endif
440 std::string INCLXXDataFilePath;
441
442 std::string clusterAlgorithmString;
443 ClusterAlgorithmType clusterAlgorithmType;
444
445 G4int clusterMaxMass;
446
447 G4bool backToSpectator;
448
449 G4bool useRealMasses;
450
451 G4double impactParameter;
452
453 std::string separationEnergyString;
454 SeparationEnergyType separationEnergyType;
455
456 std::string fermiMomentumString;
457 FermiMomentumType fermiMomentumType;
458
459 G4double fermiMomentum;
460
461 G4double cutNN;
462
463 G4bool ann;
464
465 G4double bias;
466
467 G4double atrestThreshold;
468
469#ifdef INCL_ROOT_USE
470 std::string rootSelectionString;
471#endif
472
473#ifdef INCL_DEEXCITATION_FERMI_BREAKUP
474 G4int maxMassFermiBreakUp;
475 G4int maxChargeFermiBreakUp;
476#endif
477
478 G4double rpCorrelationCoefficient;
479 G4double rpCorrelationCoefficientProton;
480 G4double rpCorrelationCoefficientNeutron;
481
482 G4double neutronSkin;
483 G4double neutronHalo;
484
485 G4bool refraction;
486
487 std::string randomNumberGenerator;
488 RNGType rngType;
489
490 std::string phaseSpaceGenerator;
491 PhaseSpaceGeneratorType phaseSpaceGeneratorType;
492
493 unsigned int autosaveFrequency;
494
495 std::string crossSectionsString;
496 CrossSectionsType crossSectionsType;
497 G4int maxNumberMultipions;
498
499 std::string cascadeAction;
500 CascadeActionType cascadeActionType;
501
502 G4double hadronizationTime;
503
504#ifdef INCL_ROOT_USE
505 G4bool conciseROOTTree;
506#endif
507
508 G4bool inverseKinematics;
509
510 G4bool srcPairCorrelations;
511
512 G4float srcPairDistance;
513
514 G4double decayTimeThreshold;
515
516 friend class ::ConfigParser;
517 };
518
519}
520
521#endif
G4double S(G4double temp)
float G4float
Definition G4Types.hh:84
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4double getNeutronHalo() const
Get the neutron-halo size.
static std::string const getVersionString()
Get the INCL version string.
void setFermiMomentumType(FermiMomentumType const f)
Set the Fermi-momentum type.
PauliType getPauliType() const
Get the Pauli-blocking algorithm.
void setClusterAlgorithm(ClusterAlgorithmType const c)
Set the clustering algorithm.
G4bool getPionPotential() const
Do we want the pion potential?
FermiMomentumType getFermiMomentumType() const
Get the Fermi-momentum type.
void setClusterMaxMass(const G4int clm)
Set the maximum mass for production of clusters.
PotentialType getPotentialType() const
Get the type of the potential for nucleons.
static std::string const getVersionID()
Get the INCL version ID.
void setDecayTimeThreshold(const G4double t)
Set decay time threshold time.
G4double getImpactParameter() const
void setLocalEnergyPiType(const LocalEnergyType t)
Set the type of local energy for N-N avatars.
void setBackToSpectator(const G4bool b)
Set back-to-spectator.
SeparationEnergyType getSeparationEnergyType() const
Get the separation-energy type.
CascadeActionType getCascadeActionType() const
Get the cascade-action type.
void setINCLXXDataFilePath(std::string const &path)
Set the INCLXX datafile path.
G4double getRPCorrelationCoefficient(const ParticleType t) const
Get the r-p correlation coefficient.
void setRNGType(RNGType const r)
Set the RNG type.
G4bool getRefraction() const
True if we should use refraction.
void setAtrestThreshold(const G4double t)
Set the pbar at rest annihilation threshold.
std::string const & getINCLXXDataFilePath() const
Set the ABLAXX datafile path.
G4int getClusterMaxMass() const
Get the maximum mass for production of clusters.
G4int getVerboseEvent() const
Get the number of the verbose event.
G4float getsrcPairDist() const
void setCoulombType(CoulombType const c)
Set the Coulomb-distortion algorithm.
G4int getTargetA() const
Get the target mass number.
G4double getDecayTimeThreshold() const
Get the decay time threshold time.
G4double getAtrestThreshold() const
Get the pbar at rest annihilation threshold.
std::string const & getCalculationTitle() const
Get the run title.
G4double getProjectileKineticEnergy() const
Get the projectile kinetic energy.
G4int getVerbosity() const
Get the verbosity.
void setTargetS(G4int S)
Set target strangess number.
void setRPCorrelationCoefficient(const ParticleType t, const G4double corrCoeff)
Set the r-p correlation coefficient.
G4int getMaxNumberMultipions() const
Get the maximum number of pions for multipion collisions.
std::string const & getLogFileName() const
Get the log file name.
unsigned int getAutosaveFrequency() const
Get the autosave frequency.
PhaseSpaceGeneratorType getPhaseSpaceGeneratorType() const
Get the phase-space-generator type.
void setCascadeActionType(CascadeActionType const c)
Set the cascade-action type.
DeExcitationType getDeExcitationType() const
Get the de-excitation model.
Random::SeedVector getRandomSeeds() const
Get the seeds for the random-number generator.
void setTargetA(G4int A)
Set target mass number.
G4double getHadronizationTime() const
Get the hadronization time.
ParticleSpecies getProjectileSpecies() const
Get the projectile species.
ClusterAlgorithmType getClusterAlgorithm() const
Get the clustering algorithm.
void setTargetZ(G4int Z)
Set target charge number.
void setPionPotential(const G4bool pionPot)
Set the type of the potential for nucleons.
~Config()
Default destructor.
void setNeutronHalo(const G4double d)
Set the neutron-halo size.
void setPhaseSpaceGeneratorType(PhaseSpaceGeneratorType const p)
Set the phase-space-generator type.
void setPotentialType(PotentialType type)
Set the type of the potential for nucleons.
LocalEnergyType getLocalEnergyPiType() const
Get the type of local energy for pi-N and decay avatars.
void setProjectileSpecies(ParticleSpecies const &pars)
Set the projectile species.
G4bool getInverseKinematics() const
void setCrossSectionsType(CrossSectionsType const c)
Set the Cross Section type.
G4bool isNaturalTarget() const
Natural targets.
void setAutosaveFrequency(const unsigned int f)
Set the autosave frequency.
G4int getTargetS() const
Get the target strangess number.
G4bool getBackToSpectator() const
Get back-to-spectator.
void setProjectileKineticEnergy(G4double const kinE)
Set the projectile kinetic energy.
void setMaxNumberMultipions(const G4int n)
Set the maximum number of pions for multipion collisions.
Config()
Default constructor.
G4bool getCDPP() const
Do we want CDPP?
G4double getNeutronSkin() const
Get the neutron-skin thickness.
CrossSectionsType getCrossSectionsType() const
Get the Cross Section type.
void setHadronizationTime(const G4double t)
Set the hadronization time.
G4double getCutNN() const
G4double getFermiMomentum() const
Get the Fermi momentum.
LocalEnergyType getLocalEnergyBBType() const
Get the type of local energy for N-N avatars.
G4double getBias() const
Get the bias.
void setNeutronSkin(const G4double d)
Set the neutron-skin thickness.
G4bool getsrcPairConfig() const
std::string summary()
Return a summary of the run configuration.
void init()
Initialise the members.
void setFermiMomentum(const G4double p)
Set the Fermi momentum.
G4int getTargetZ() const
Get the target charge number.
void setLocalEnergyBBType(const LocalEnergyType t)
Set the type of local energy for N-N avatars.
RNGType getRNGType() const
Get the RNG type.
void setRefraction(const G4bool r)
Set the refraction variable.
void setUseRealMasses(G4bool use)
Set whether to use real masses.
ParticleType getProjectileType() const
Get the projectile type.
G4bool getUseRealMasses() const
Whether to use real masses.
CoulombType getCoulombType() const
Get the Coulomb-distortion algorithm.
G4int getNumberOfShots() const
Get the number of shots.
std::string getDeExcitationString() const
Get the de-excitation string.
std::string const & getOutputFileRoot() const
Get the output file root.
void setRPCorrelationCoefficient(const G4double corrCoeff)
Set the r-p correlation coefficient.
static std::string const getVersionHash()
Get the INCL version hash.