Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
38#include "G4UIcmdWithABool.hh"
39#include "G4UIcmdWithADouble.hh"
40#include "G4UIcmdWithAString.hh"
43#include "G4UIcommand.hh"
44#include "G4UIcommandTree.hh"
45#include "G4UIdirectory.hh"
46#include "G4UImanager.hh"
47
48
49// Constructor and destructor
50
52 : G4UImessenger(), theParams(params), cmdDir(0), localCmdDir(false) {
53 // NOTE: Put under same top-level tree as EM
54 CreateDirectory("/process/had/","Hadronic processes"); cmdDir=0;
55 CreateDirectory("/process/had/cascade/","Bertini-esque cascade parameters");
56
57 verboseCmd = CreateCommand<G4UIcmdWithAnInteger>("verbose",
58 "Enable information messages");
59 balanceCmd = CreateCommand<G4UIcmdWithABool>("checkBalance",
60 "Enable internal conservation checking");
61 reportCmd = CreateCommand<G4UIcmdWithoutParameter>("report",
62 "Dump all non-default parameter settings");
63 usePreCoCmd = CreateCommand<G4UIcmdWithABool>("usePreCompound",
64 "Use PreCompoundModel for nuclear de-excitation");
65 doCoalCmd = CreateCommand<G4UIcmdWithABool>("doCoalescence",
66 "Apply final-state nucleon clustering");
67 piNAbsCmd = CreateCommand<G4UIcmdWithADouble>("piNAbsorption",
68 "Probability for pion absorption on single nucleon");
69 historyCmd = CreateCommand<G4UIcmdWithABool>("showHistory",
70 "Collect and report full structure of cascade");
71 use3BodyCmd = CreateCommand<G4UIcmdWithABool>("use3BodyMom",
72 "Use three-body momentum parametrizations");
73 usePSCmd = CreateCommand<G4UIcmdWithABool>("usePhaseSpace",
74 "Use Kopylov N-body momentum generator");
75 randomFileCmd = CreateCommand<G4UIcmdWithAString>("randomFile",
76 "Save random-engine to file at each interaction");
77 nucUseBestCmd = CreateCommand<G4UIcmdWithABool>("useBestNuclearModel",
78 "Use all physical-units for nuclear structure");
79 nucRad2parCmd = CreateCommand<G4UIcmdWithADouble>("useTwoParamNuclearRadius",
80 "Use R = C1*cbrt(A) + C2/cbrt(A)");
81 nucRadScaleCmd = CreateCommand<G4UIcmdWithADouble>("nuclearRadiusScale",
82 "Set length scale for nuclear model");
83 nucRadSmallCmd = CreateCommand<G4UIcmdWithADouble>("smallNucleusRadius",
84 "Set radius of A<4 nuclei");
85 nucRadAlphaCmd = CreateCommand<G4UIcmdWithADouble>("alphaRadiusScale",
86 "Fraction of small-radius for He-4");
87 nucRadTrailingCmd = CreateCommand<G4UIcmdWithADouble>("shadowningRadius",
88 "Effective nucleon radius for trailing effect");
89 nucFermiScaleCmd = CreateCommand<G4UIcmdWithADouble>("fermiScale",
90 "Scale factor for fermi momentum");
91 nucXsecScaleCmd = CreateCommand<G4UIcmdWithADouble>("crossSectionScale",
92 "Scale fator for total cross-sections");
93 nucGammaQDCmd = CreateCommand<G4UIcmdWithADouble>("gammaQuasiDeutScale",
94 "Scale factor for gamma-quasideutron cross-sections");
95 coalDPmax2Cmd = CreateCommand<G4UIcmdWithADouble>("cluster2DPmax",
96 "Maximum momentum for p-n clusters");
97 coalDPmax3Cmd = CreateCommand<G4UIcmdWithADouble>("cluster3DPmax",
98 "Maximum momentum for ppn/pnn clusters");
99 coalDPmax4Cmd = CreateCommand<G4UIcmdWithADouble>("cluster4DPmax",
100 "Maximum momentum for alpha clusters");
101}
102
104 delete verboseCmd;
105 delete balanceCmd;
106 delete reportCmd;
107 delete usePreCoCmd;
108 delete doCoalCmd;
109 delete piNAbsCmd;
110 delete historyCmd;
111 delete use3BodyCmd;
112 delete usePSCmd;
113 delete randomFileCmd;
114 delete nucUseBestCmd;
115 delete nucRad2parCmd;
116 delete nucRadScaleCmd;
117 delete nucRadSmallCmd;
118 delete nucRadAlphaCmd;
119 delete nucRadTrailingCmd;
120 delete nucFermiScaleCmd;
121 delete nucXsecScaleCmd;
122 delete nucGammaQDCmd;
123 delete coalDPmax2Cmd;
124 delete coalDPmax3Cmd;
125 delete coalDPmax4Cmd;
126 if (localCmdDir) delete cmdDir;
127}
128
129
130// Create or reuse existing UIdirectory path
131
133 const char* desc) {
135 if (!UIman) return;
136
137 // Directory path must be absolute, prepend "/" if ncessary
138 G4String fullPath = path;
139 if (fullPath(0) != '/') fullPath.prepend("/");
140 if (fullPath(fullPath.length()-1) != '/') fullPath.append("/");
141
142 // See if input path has already been registered
143 G4UIcommand* foundPath = UIman->GetTree()->FindPath(fullPath);
144 if (foundPath) cmdDir = dynamic_cast<G4UIdirectory*>(foundPath);
145
146 if (!cmdDir) { // Create local deletable directory
147 localCmdDir = true;
148 cmdDir = new G4UIdirectory(fullPath.c_str());
149 cmdDir->SetGuidance(desc);
150 }
151}
152
153
154// Use command argument (literal string) to set envvar maps in container
155
157 if (cmd == reportCmd) theParams->DumpConfig(G4cout);
158
159 if (cmd == verboseCmd)
160 theParams->G4CASCADE_VERBOSE = strdup(arg.c_str());
161
162 if (cmd == balanceCmd)
163 theParams->G4CASCADE_CHECK_ECONS = StoB(arg) ? strdup(arg.c_str()) : 0;
164
165 if (cmd == usePreCoCmd)
166 theParams->G4CASCADE_USE_PRECOMPOUND = StoB(arg) ? strdup(arg.c_str()) : 0;
167
168 if (cmd == doCoalCmd)
169 theParams->G4CASCADE_DO_COALESCENCE = StoB(arg) ? strdup(arg.c_str()) : 0;
170
171 if (cmd == piNAbsCmd)
172 theParams->G4CASCADE_PIN_ABSORPTION = strdup(arg.c_str());
173
174 if (cmd == historyCmd)
175 theParams->G4CASCADE_SHOW_HISTORY = StoB(arg) ? strdup(arg.c_str()) : 0;
176
177 if (cmd == use3BodyCmd)
178 theParams->G4CASCADE_USE_3BODYMOM = StoB(arg) ? strdup(arg.c_str()) : 0;
179
180 if (cmd == usePSCmd)
181 theParams->G4CASCADE_USE_PHASESPACE = StoB(arg) ? strdup(arg.c_str()) : 0;
182
183 if (cmd == randomFileCmd)
184 theParams->G4CASCADE_RANDOM_FILE = arg.empty() ? 0 : strdup(arg.c_str());
185
186 if (cmd == nucUseBestCmd)
187 theParams->G4NUCMODEL_USE_BEST = StoB(arg) ? strdup(arg.c_str()) : 0;
188
189// if (cmd == nucRad2parCmd)
190// theParams->G4NUCMODEL_RAD_2PAR = strdup(arg.c_str());
191 if (cmd == nucRad2parCmd)
192 theParams->G4NUCMODEL_RAD_2PAR = StoB(arg) ? strdup(arg.c_str()) : 0;
193
194 if (cmd == nucRadScaleCmd)
195 theParams->G4NUCMODEL_RAD_SCALE = strdup(arg.c_str());
196
197 if (cmd == nucRadSmallCmd)
198 theParams->G4NUCMODEL_RAD_SMALL = strdup(arg.c_str());
199
200 if (cmd == nucRadAlphaCmd)
201 theParams->G4NUCMODEL_RAD_ALPHA = strdup(arg.c_str());
202
203 if (cmd == nucRadTrailingCmd)
204 theParams->G4NUCMODEL_RAD_TRAILING = strdup(arg.c_str());
205
206 if (cmd == nucFermiScaleCmd)
207 theParams->G4NUCMODEL_FERMI_SCALE = strdup(arg.c_str());
208
209 if (cmd == nucXsecScaleCmd)
210 theParams->G4NUCMODEL_XSEC_SCALE = strdup(arg.c_str());
211
212 if (cmd == nucGammaQDCmd)
213 theParams->G4NUCMODEL_GAMMAQD = strdup(arg.c_str());
214
215 if (cmd == coalDPmax2Cmd)
216 theParams->DPMAX_2CLUSTER = strdup(arg.c_str());
217
218 if (cmd == coalDPmax3Cmd)
219 theParams->DPMAX_3CLUSTER = strdup(arg.c_str());
220
221 if (cmd == coalDPmax4Cmd)
222 theParams->DPMAX_4CLUSTER = strdup(arg.c_str());
223
224 theParams->Initialize(); // Update numerical values from settings
225}
G4GLOB_DLL std::ostream G4cout
virtual void SetNewValue(G4UIcommand *command, G4String newValue)
G4CascadeParamMessenger(G4CascadeParameters *params)
void CreateDirectory(const char *path, const char *desc)
G4String & append(const G4String &)
G4String & prepend(const char *)
G4UIcommand * FindPath(const char *commandPath) const
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
G4UIcommandTree * GetTree() const
Definition: G4UImanager.hh:179
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
G4bool StoB(G4String s)