Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4CascadeParamMessenger.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// Define simple UI commands as alternative to environment variables
27//
28// 20130304 M. Kelsey -- Add flag to collect and display cascade structure
29// 20130621 M. Kelsey -- Add flag for CHECK_ECONS, replacing #ifdef's; add
30// flag to use three-body momentum parametrizations
31// 20130703 M. Kelsey -- Add flag for USE_PHASESPACE
32// 20141030 M. Kelsey -- Add flag to enable direct pi-N absorption
33// 20141211 M. Kelsey -- Change PIN_ABSORPTION flag to double, for energy cut
34// 20200110 M. Kelsey -- Reset cmdDir to 0 before .../cascade/ directory.
35
39#include "G4UIcmdWithABool.hh"
40#include "G4UIcmdWithADouble.hh"
41#include "G4UIcmdWithAString.hh"
44#include "G4UIcommand.hh"
45#include "G4UIcommandTree.hh"
46#include "G4UIdirectory.hh"
47#include "G4UImanager.hh"
48
49
50// Constructor and destructor
51
53 : G4UImessenger(), theParams(params) {
55 cmdDir = new G4UIdirectory( "/process/had/cascade/" );
56 cmdDir->SetGuidance( "Bertini-esque cascade parameters" );
57
58 verboseCmd = CreateCommand<G4UIcmdWithAnInteger>("verbose",
59 "Enable information messages");
60 balanceCmd = CreateCommand<G4UIcmdWithABool>("checkBalance",
61 "Enable internal conservation checking");
62 reportCmd = CreateCommand<G4UIcmdWithoutParameter>("report",
63 "Dump all non-default parameter settings");
64 usePreCoCmd = CreateCommand<G4UIcmdWithABool>("usePreCompound",
65 "Use PreCompoundModel for nuclear de-excitation");
66 doCoalCmd = CreateCommand<G4UIcmdWithABool>("doCoalescence",
67 "Apply final-state nucleon clustering");
68 piNAbsCmd = CreateCommand<G4UIcmdWithADouble>("piNAbsorption",
69 "Probability for pion absorption on single nucleon");
70 historyCmd = CreateCommand<G4UIcmdWithABool>("showHistory",
71 "Collect and report full structure of cascade");
72 use3BodyCmd = CreateCommand<G4UIcmdWithABool>("use3BodyMom",
73 "Use three-body momentum parametrizations");
74 usePSCmd = CreateCommand<G4UIcmdWithABool>("usePhaseSpace",
75 "Use Kopylov N-body momentum generator");
76 randomFileCmd = CreateCommand<G4UIcmdWithAString>("randomFile",
77 "Save random-engine to file at each interaction");
78 nucUseBestCmd = CreateCommand<G4UIcmdWithABool>("useBestNuclearModel",
79 "Use all physical-units for nuclear structure");
80 nucRad2parCmd = CreateCommand<G4UIcmdWithADouble>("useTwoParamNuclearRadius",
81 "Use R = C1*cbrt(A) + C2/cbrt(A)");
82 nucRadScaleCmd = CreateCommand<G4UIcmdWithADouble>("nuclearRadiusScale",
83 "Set length scale for nuclear model");
84 nucRadSmallCmd = CreateCommand<G4UIcmdWithADouble>("smallNucleusRadius",
85 "Set radius of A<4 nuclei");
86 nucRadAlphaCmd = CreateCommand<G4UIcmdWithADouble>("alphaRadiusScale",
87 "Fraction of small-radius for He-4");
88 nucRadTrailingCmd = CreateCommand<G4UIcmdWithADouble>("shadowningRadius",
89 "Effective nucleon radius for trailing effect");
90 nucFermiScaleCmd = CreateCommand<G4UIcmdWithADouble>("fermiScale",
91 "Scale factor for fermi momentum");
92 nucXsecScaleCmd = CreateCommand<G4UIcmdWithADouble>("crossSectionScale",
93 "Scale fator for total cross-sections");
94 nucGammaQDCmd = CreateCommand<G4UIcmdWithADouble>("gammaQuasiDeutScale",
95 "Scale factor for gamma-quasideutron cross-sections");
96 coalDPmax2Cmd = CreateCommand<G4UIcmdWithADouble>("cluster2DPmax",
97 "Maximum momentum for p-n clusters");
98 coalDPmax3Cmd = CreateCommand<G4UIcmdWithADouble>("cluster3DPmax",
99 "Maximum momentum for ppn/pnn clusters");
100 coalDPmax4Cmd = CreateCommand<G4UIcmdWithADouble>("cluster4DPmax",
101 "Maximum momentum for alpha clusters");
102}
103
105 delete verboseCmd;
106 delete balanceCmd;
107 delete reportCmd;
108 delete usePreCoCmd;
109 delete doCoalCmd;
110 delete piNAbsCmd;
111 delete historyCmd;
112 delete use3BodyCmd;
113 delete usePSCmd;
114 delete randomFileCmd;
115 delete nucUseBestCmd;
116 delete nucRad2parCmd;
117 delete nucRadScaleCmd;
118 delete nucRadSmallCmd;
119 delete nucRadAlphaCmd;
120 delete nucRadTrailingCmd;
121 delete nucFermiScaleCmd;
122 delete nucXsecScaleCmd;
123 delete nucGammaQDCmd;
124 delete coalDPmax2Cmd;
125 delete coalDPmax3Cmd;
126 delete coalDPmax4Cmd;
127 delete cmdDir;
128}
129
130// Use command argument (literal string) to set envvar maps in container
131
133 if (cmd == reportCmd) theParams->DumpConfig(G4cout);
134
135 if (cmd == verboseCmd)
136 theParams->G4CASCADE_VERBOSE = strdup(arg.c_str());
137
138 if (cmd == balanceCmd)
139 theParams->G4CASCADE_CHECK_ECONS = StoB(arg) ? strdup(arg.c_str()) : 0;
140
141 if (cmd == usePreCoCmd)
142 theParams->G4CASCADE_USE_PRECOMPOUND = StoB(arg) ? strdup(arg.c_str()) : 0;
143
144 if (cmd == doCoalCmd)
145 theParams->G4CASCADE_DO_COALESCENCE = StoB(arg) ? strdup(arg.c_str()) : 0;
146
147 if (cmd == piNAbsCmd)
148 theParams->G4CASCADE_PIN_ABSORPTION = strdup(arg.c_str());
149
150 if (cmd == historyCmd)
151 theParams->G4CASCADE_SHOW_HISTORY = StoB(arg) ? strdup(arg.c_str()) : 0;
152
153 if (cmd == use3BodyCmd)
154 theParams->G4CASCADE_USE_3BODYMOM = StoB(arg) ? strdup(arg.c_str()) : 0;
155
156 if (cmd == usePSCmd)
157 theParams->G4CASCADE_USE_PHASESPACE = StoB(arg) ? strdup(arg.c_str()) : 0;
158
159 if (cmd == randomFileCmd)
160 theParams->G4CASCADE_RANDOM_FILE = arg.empty() ? 0 : strdup(arg.c_str());
161
162 if (cmd == nucUseBestCmd)
163 theParams->G4NUCMODEL_USE_BEST = StoB(arg) ? strdup(arg.c_str()) : 0;
164
165// if (cmd == nucRad2parCmd)
166// theParams->G4NUCMODEL_RAD_2PAR = strdup(arg.c_str());
167 if (cmd == nucRad2parCmd)
168 theParams->G4NUCMODEL_RAD_2PAR = StoB(arg) ? strdup(arg.c_str()) : 0;
169
170 if (cmd == nucRadScaleCmd)
171 theParams->G4NUCMODEL_RAD_SCALE = strdup(arg.c_str());
172
173 if (cmd == nucRadSmallCmd)
174 theParams->G4NUCMODEL_RAD_SMALL = strdup(arg.c_str());
175
176 if (cmd == nucRadAlphaCmd)
177 theParams->G4NUCMODEL_RAD_ALPHA = strdup(arg.c_str());
178
179 if (cmd == nucRadTrailingCmd)
180 theParams->G4NUCMODEL_RAD_TRAILING = strdup(arg.c_str());
181
182 if (cmd == nucFermiScaleCmd)
183 theParams->G4NUCMODEL_FERMI_SCALE = strdup(arg.c_str());
184
185 if (cmd == nucXsecScaleCmd)
186 theParams->G4NUCMODEL_XSEC_SCALE = strdup(arg.c_str());
187
188 if (cmd == nucGammaQDCmd)
189 theParams->G4NUCMODEL_GAMMAQD = strdup(arg.c_str());
190
191 if (cmd == coalDPmax2Cmd)
192 theParams->DPMAX_2CLUSTER = strdup(arg.c_str());
193
194 if (cmd == coalDPmax3Cmd)
195 theParams->DPMAX_3CLUSTER = strdup(arg.c_str());
196
197 if (cmd == coalDPmax4Cmd)
198 theParams->DPMAX_4CLUSTER = strdup(arg.c_str());
199
200 theParams->Initialize(); // Update numerical values from settings
201}
G4GLOB_DLL std::ostream G4cout
void SetNewValue(G4UIcommand *command, G4String newValue) override
T * CreateCommand(const G4String &cmd, const G4String &desc)
G4CascadeParamMessenger(G4CascadeParameters *params)
static G4HadronicParameters * Instance()
G4UImessenger()=default
G4bool StoB(const G4String &s)