Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CascadeParamMessenger Class Reference

#include <G4CascadeParamMessenger.hh>

+ Inheritance diagram for G4CascadeParamMessenger:

Public Member Functions

 G4CascadeParamMessenger (G4CascadeParameters *params)
 
virtual ~G4CascadeParamMessenger ()
 
virtual void SetNewValue (G4UIcommand *command, G4String newValue)
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()=default
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
virtual void SetNewValue (G4UIcommand *command, G4String newValue)
 
G4bool operator== (const G4UImessenger &messenger) const
 
G4bool operator!= (const G4UImessenger &messenger) const
 
G4bool CommandsShouldBeInMaster () const
 

Protected Member Functions

void CreateDirectory (const char *path, const char *desc)
 
template<class T >
T * CreateCommand (const G4String &cmd, const G4String &desc)
 
- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (const G4String &s)
 
G4long StoL (const G4String &s)
 
G4double StoD (const G4String &s)
 
G4bool StoB (G4String s)
 
void AddUIcommand (G4UIcommand *newCommand)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 

Additional Inherited Members

- Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Detailed Description

Definition at line 51 of file G4CascadeParamMessenger.hh.

Constructor & Destructor Documentation

◆ G4CascadeParamMessenger()

G4CascadeParamMessenger::G4CascadeParamMessenger ( G4CascadeParameters params)

Definition at line 51 of file G4CascadeParamMessenger.cc.

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}
void CreateDirectory(const char *path, const char *desc)
G4UImessenger()=default

◆ ~G4CascadeParamMessenger()

G4CascadeParamMessenger::~G4CascadeParamMessenger ( )
virtual

Definition at line 103 of file G4CascadeParamMessenger.cc.

103 {
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}

Member Function Documentation

◆ CreateCommand()

template<class T >
T * G4CascadeParamMessenger::CreateCommand ( const G4String cmd,
const G4String desc 
)
protected

◆ CreateDirectory()

void G4CascadeParamMessenger::CreateDirectory ( const char *  path,
const char *  desc 
)
protected

Definition at line 132 of file G4CascadeParamMessenger.cc.

133 {
135 if (!UIman) return;
136
137 // Directory path must be absolute, prepend "/" if ncessary
138 G4String fullPath = path;
139 if (fullPath[0] != '/') fullPath.insert(0, '/', 1);
140 if (fullPath.back() != '/') fullPath.append('/', 1);
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}
G4UIcommand * FindPath(const char *commandPath) const
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:157
G4UIcommandTree * GetTree() const
Definition: G4UImanager.hh:186
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77

Referenced by G4CascadeParamMessenger().

◆ SetNewValue()

void G4CascadeParamMessenger::SetNewValue ( G4UIcommand command,
G4String  newValue 
)
virtual

Reimplemented from G4UImessenger.

Definition at line 156 of file G4CascadeParamMessenger.cc.

156 {
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
G4bool StoB(G4String s)

The documentation for this class was generated from the following files: