Geant4 11.2.2
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)
 
 ~G4CascadeParamMessenger () override
 
void SetNewValue (G4UIcommand *command, G4String newValue) override
 
- 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)
 
G4bool CommandsShouldBeInMaster () const
 

Protected Member Functions

template<class T >
T * CreateCommand (const G4String &cmd, const G4String &desc)
 
- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String LtoS (G4long l)
 
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 52 of file G4CascadeParamMessenger.cc.

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}
T * CreateCommand(const G4String &cmd, const G4String &desc)
static G4HadronicParameters * Instance()
void SetGuidance(const char *aGuidance)
G4UImessenger()=default

◆ ~G4CascadeParamMessenger()

G4CascadeParamMessenger::~G4CascadeParamMessenger ( )
override

Definition at line 104 of file G4CascadeParamMessenger.cc.

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

Member Function Documentation

◆ CreateCommand()

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

Referenced by G4CascadeParamMessenger().

◆ SetNewValue()

void G4CascadeParamMessenger::SetNewValue ( G4UIcommand * command,
G4String newValue )
overridevirtual

Reimplemented from G4UImessenger.

Definition at line 132 of file G4CascadeParamMessenger.cc.

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

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