Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CascadeParameters.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// Encapsulate all user-configurable parameters with associated envvars
27//
28// 20120912 M. Kelsey -- Add interface to support UI commands
29// 20130304 M. Kelsey -- Add flag to collect and display cascade structure
30// 20130308 M. Kelsey -- Add flag to use separate 3-body momentum generators
31// 20130421 M. Kelsey -- Add flag for CHECK_ECONS, replacing #ifdef's
32// 20130702 M. Kelsey -- Add flag to use N-body phase-space generator
33// 20140311 G. Cosmo -- Implement standard (non-const) singleton pattern
34// 20140929 M. Kelsey -- Enable some parameters as default true (must be set
35// '0' for false): PreCompound, phase-space, clustering,
36// trailing effect
37// 20141111 M. Kelsey -- Revert defaults for PreCompound, phase-space,
38// and trailing effect.
39// 20141121 Use G4AutoDelete to avoid end-of-thread memory leaks
40// 20141211 M. Kelsey -- Change PIN_ABSORPTION flag to double, for energy cut
41
44#include "G4AutoDelete.hh"
45#include <stdlib.h>
46#include <iostream>
47using std::endl;
48
49#define OLD_RADIUS_UNITS (3.3836/1.2) // Used with NucModel params
50
52namespace {
54 class BERTParameters {
55 public:
56 BERTParameters(){
57 //HDP.SetDefault("NAME",VALUE,LOWER_LIMIT(default=-DBL_MAX),UPPER_LIMIT(default=DBL_MAX));
58 HDP.SetDefault( "BERT_FERMI_SCALE" , 1.932/OLD_RADIUS_UNITS , 1.932/OLD_RADIUS_UNITS/2 , 1.932/OLD_RADIUS_UNITS*2 );
59 HDP.SetDefault( "BERT_RADIUS_SCALE" , OLD_RADIUS_UNITS , OLD_RADIUS_UNITS/2 , OLD_RADIUS_UNITS*2 );
60 HDP.SetDefault( "BERT_RAD_TRAILING" , 0. , 0. , 2. );
61 HDP.SetDefault( "BERT_XSEC_SCALE" , 1.0 , 0.5 , 2. );
62 }
63 };
64 BERTParameters BP;
65}
66
67
68// Singleton accessor
69
71 static auto _instance = []() {
72 auto _ptr = new G4CascadeParameters{};
74 return _ptr;
75 }();
76 return _instance;
77}
78
79
80// Constructor initializes everything once
81
82//#define OLD_RADIUS_UNITS (3.3836/1.2) // Used with NucModel params
83
84G4CascadeParameters::G4CascadeParameters()
85 : G4CASCADE_VERBOSE(std::getenv("G4CASCADE_VERBOSE")),
86 G4CASCADE_CHECK_ECONS(std::getenv("G4CASCADE_CHECK_ECONS")),
87 G4CASCADE_USE_PRECOMPOUND(std::getenv("G4CASCADE_USE_PRECOMPOUND")),
88 G4CASCADE_USE_ABLA(std::getenv("G4CASCADE_USE_ABLA")),
89 G4CASCADE_DO_COALESCENCE(std::getenv("G4CASCADE_DO_COALESCENCE")),
90 G4CASCADE_SHOW_HISTORY(std::getenv("G4CASCADE_SHOW_HISTORY")),
91 G4CASCADE_USE_3BODYMOM(std::getenv("G4CASCADE_USE_3BODYMOM")),
92 G4CASCADE_USE_PHASESPACE(std::getenv("G4CASCADE_USE_PHASESPACE")),
93 G4CASCADE_PIN_ABSORPTION(std::getenv("G4CASCADE_PIN_ABSORPTION")),
94 G4CASCADE_RANDOM_FILE(std::getenv("G4CASCADE_RANDOM_FILE")),
95 G4NUCMODEL_USE_BEST(std::getenv("G4NUCMODEL_USE_BEST")),
96 G4NUCMODEL_RAD_2PAR(std::getenv("G4NUCMODEL_RAD_2PAR")),
97 G4NUCMODEL_RAD_SCALE(std::getenv("G4NUCMODEL_RAD_SCALE")),
98 G4NUCMODEL_RAD_SMALL(std::getenv("G4NUCMODEL_RAD_SMALL")),
99 G4NUCMODEL_RAD_ALPHA(std::getenv("G4NUCMODEL_RAD_ALPHA")),
100 G4NUCMODEL_RAD_TRAILING(std::getenv("G4NUCMODEL_RAD_TRAILING")),
101 G4NUCMODEL_FERMI_SCALE(std::getenv("G4NUCMODEL_FERMI_SCALE")),
102 G4NUCMODEL_XSEC_SCALE(std::getenv("G4NUCMODEL_XSEC_SCALE")),
103 G4NUCMODEL_GAMMAQD(std::getenv("G4NUCMODEL_GAMMAQD")),
104 DPMAX_2CLUSTER(std::getenv("DPMAX_2CLUSTER")),
105 DPMAX_3CLUSTER(std::getenv("DPMAX_3CLUSTER")),
106 DPMAX_4CLUSTER(std::getenv("DPMAX_4CLUSTER")),
107 messenger(0) {
108 messenger = new G4CascadeParamMessenger(this);
109 Initialize();
110}
111
112void G4CascadeParameters::Initialize() {
113 VERBOSE_LEVEL = (G4CASCADE_VERBOSE ? atoi(G4CASCADE_VERBOSE) : 0);
114 CHECK_ECONS = (0!=G4CASCADE_CHECK_ECONS);
115 USE_PRECOMPOUND = (0!=G4CASCADE_USE_PRECOMPOUND &&
116 G4CASCADE_USE_PRECOMPOUND[0]!='0');
117 USE_ABLA = (0!=G4CASCADE_USE_ABLA && G4CASCADE_USE_ABLA[0]!='0');
118 DO_COALESCENCE = (0==G4CASCADE_DO_COALESCENCE ||
119 G4CASCADE_DO_COALESCENCE[0]!='0');
120 SHOW_HISTORY = (0!=G4CASCADE_SHOW_HISTORY);
121 USE_3BODYMOM = (0!=G4CASCADE_USE_3BODYMOM);
122 USE_PHASESPACE = (0!=G4CASCADE_USE_PHASESPACE &&
123 G4CASCADE_USE_PHASESPACE[0]!='0');
124 PIN_ABSORPTION = (G4CASCADE_PIN_ABSORPTION ? strtod(G4CASCADE_PIN_ABSORPTION,0)
125 : 0.);
126 RANDOM_FILE = (G4CASCADE_RANDOM_FILE ? G4CASCADE_RANDOM_FILE : "");
127 BEST_PAR = (0!=G4NUCMODEL_USE_BEST);
128 TWOPARAM_RADIUS = (0!=G4NUCMODEL_RAD_2PAR); // && G4NUCMODEL_RAD_2PAR[0]!='0');
129 RADIUS_SCALE = (G4NUCMODEL_RAD_SCALE ? strtod(G4NUCMODEL_RAD_SCALE,0)
130 : (BEST_PAR?1.0:OLD_RADIUS_UNITS));
131 if ( G4NUCMODEL_RAD_SCALE == 0 && BEST_PAR == 0 ) HDP.DeveloperGet("BERT_RADIUS_SCALE",RADIUS_SCALE);
132 RADIUS_SMALL = ((G4NUCMODEL_RAD_SMALL ? strtod(G4NUCMODEL_RAD_SMALL,0)
133 : (BEST_PAR?1.992:(8.0/OLD_RADIUS_UNITS))) * RADIUS_SCALE);
134 RADIUS_ALPHA = (G4NUCMODEL_RAD_ALPHA ? strtod(G4NUCMODEL_RAD_ALPHA,0)
135 : (BEST_PAR?0.84:0.70));
136 RADIUS_TRAILING = ((G4NUCMODEL_RAD_TRAILING ? strtod(G4NUCMODEL_RAD_TRAILING,0)
137 : 0.) * RADIUS_SCALE);
138 if ( G4NUCMODEL_RAD_TRAILING == 0 ) HDP.DeveloperGet("BERT_RAD_TRAILING",RADIUS_TRAILING),RADIUS_TRAILING*=RADIUS_SCALE;
139 FERMI_SCALE = ((G4NUCMODEL_FERMI_SCALE ? strtod(G4NUCMODEL_FERMI_SCALE,0)
140 : (BEST_PAR?0.685:(1.932/OLD_RADIUS_UNITS))) * RADIUS_SCALE);
141 if ( G4NUCMODEL_FERMI_SCALE == 0 && BEST_PAR == 0 ) HDP.DeveloperGet("BERT_FERMI_SCALE",FERMI_SCALE),FERMI_SCALE*=RADIUS_SCALE;
142 XSEC_SCALE = (G4NUCMODEL_XSEC_SCALE ? strtod(G4NUCMODEL_XSEC_SCALE,0)
143 : (BEST_PAR?0.1:1.0) );
144 if ( G4NUCMODEL_XSEC_SCALE == 0 && BEST_PAR == 0 ) HDP.DeveloperGet("BERT_XSEC_SCALE",XSEC_SCALE);
145 GAMMAQD_SCALE = (G4NUCMODEL_GAMMAQD?strtod(G4NUCMODEL_GAMMAQD,0):1.);
146 DPMAX_DOUBLET = (DPMAX_2CLUSTER ? strtod(DPMAX_2CLUSTER,0) : 0.090);
147 DPMAX_TRIPLET = (DPMAX_3CLUSTER ? strtod(DPMAX_3CLUSTER,0) : 0.108);
148 DPMAX_ALPHA = (DPMAX_4CLUSTER ? strtod(DPMAX_4CLUSTER,0) : 0.115);
149}
150
152 delete messenger;
153}
154
155
156// Report any non-default parameters (used by G4CascadeInterface)
157
158void G4CascadeParameters::DumpConfig(std::ostream& os) const {
159 if (G4CASCADE_VERBOSE)
160 os << "G4CASCADE_VERBOSE = " << G4CASCADE_VERBOSE << endl;
161 if (G4CASCADE_CHECK_ECONS)
162 os << "G4CASCADE_CHECK_ECONS = " << G4CASCADE_CHECK_ECONS << endl;
163 if (G4CASCADE_USE_PRECOMPOUND)
164 os << "G4CASCADE_USE_PRECOMPOUND = " << G4CASCADE_USE_PRECOMPOUND << endl;
165 if (G4CASCADE_USE_ABLA)
166 os << "G4CASCADE_USE_ABLA = " << G4CASCADE_USE_ABLA << endl;
167 if (G4CASCADE_DO_COALESCENCE)
168 os << "G4CASCADE_DO_COALESCENCE = " << G4CASCADE_DO_COALESCENCE << endl;
169 if (G4CASCADE_PIN_ABSORPTION)
170 os << "G4CASCADE_PIN_ABSORPTION = " << G4CASCADE_PIN_ABSORPTION << endl;
171 if (G4CASCADE_SHOW_HISTORY)
172 os << "G4CASCADE_SHOW_HISTORY = " << G4CASCADE_SHOW_HISTORY << endl;
173 if (G4CASCADE_USE_3BODYMOM)
174 os << "G4CASCADE_USE_3BODYMOM = " << G4CASCADE_USE_3BODYMOM << endl;
175 if (G4CASCADE_USE_PHASESPACE)
176 os << "G4CASCADE_USE_PHASESPACE = " << G4CASCADE_USE_PHASESPACE << endl;
177 if (G4CASCADE_RANDOM_FILE)
178 os << "G4CASCADE_RANDOM_FILE = " << G4CASCADE_RANDOM_FILE << endl;
179 if (G4NUCMODEL_USE_BEST)
180 os << "G4NUCMODEL_USE_BEST = " << G4NUCMODEL_USE_BEST << endl;
181 if (G4NUCMODEL_RAD_2PAR)
182 os << "G4NUCMODEL_RAD_2PAR = " << G4NUCMODEL_RAD_2PAR << endl;
183 if (G4NUCMODEL_RAD_SCALE)
184 os << "G4NUCMODEL_RAD_SCALE = " << G4NUCMODEL_RAD_SCALE << endl;
185 if (G4NUCMODEL_RAD_SMALL)
186 os << "G4NUCMODEL_RAD_SMALL = " << G4NUCMODEL_RAD_SMALL << endl;
187 if (G4NUCMODEL_RAD_ALPHA)
188 os << "G4NUCMODEL_RAD_ALPHA = " << G4NUCMODEL_RAD_ALPHA << endl;
189 if (G4NUCMODEL_RAD_TRAILING)
190 os << "G4NUCMODEL_RAD_TRAILING = " << G4NUCMODEL_RAD_TRAILING << endl;
191 if (G4NUCMODEL_FERMI_SCALE)
192 os << "G4NUCMODEL_FERMI_SCALE = " << G4NUCMODEL_FERMI_SCALE << endl;
193 if (G4NUCMODEL_XSEC_SCALE)
194 os << "G4NUCMODEL_XSEC_SCALE = " << G4NUCMODEL_XSEC_SCALE << endl;
195 if (G4NUCMODEL_GAMMAQD)
196 os << "G4NUCMODEL_GAMMAQD = " << G4NUCMODEL_GAMMAQD << endl;
197 if (DPMAX_2CLUSTER)
198 os << "DPMAX_2CLUSTER = " << DPMAX_2CLUSTER << endl;
199 if (DPMAX_3CLUSTER)
200 os << "DPMAX_3CLUSTER = " << DPMAX_3CLUSTER << endl;
201 if (DPMAX_4CLUSTER)
202 os << "DPMAX_4CLUSTER = " << DPMAX_4CLUSTER << endl;
203}
#define OLD_RADIUS_UNITS
G4HadronicDeveloperParameters & HDP
static const G4CascadeParameters * Instance()
friend class G4CascadeParamMessenger
G4bool DeveloperGet(const std::string name, G4bool &value)
static G4HadronicDeveloperParameters & GetInstance()
G4bool SetDefault(const std::string name, const G4bool value)
void Register(T *inst)