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

#include <G4EmExtraParametersMessenger.hh>

Inheritance diagram for G4EmExtraParametersMessenger:

Public Member Functions

 G4EmExtraParametersMessenger (G4EmExtraParameters *)
 ~G4EmExtraParametersMessenger () override
void SetNewValue (G4UIcommand *, G4String) override
G4EmExtraParametersMessengeroperator= (const G4EmExtraParametersMessenger &right)=delete
 G4EmExtraParametersMessenger (const G4EmExtraParametersMessenger &)=delete
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

Additional Inherited Members

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 (const 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)
Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir = nullptr
G4String baseDirName = ""
G4bool commandsShouldBeInMaster = false

Detailed Description

Definition at line 64 of file G4EmExtraParametersMessenger.hh.

Constructor & Destructor Documentation

◆ G4EmExtraParametersMessenger() [1/2]

G4EmExtraParametersMessenger::G4EmExtraParametersMessenger ( G4EmExtraParameters * ptr)
explicit

Definition at line 58 of file G4EmExtraParametersMessenger.cc.

59 : theParameters(ptr)
60{
61 paiCmd = new G4UIcommand("/process/em/AddPAIRegion",this);
62 paiCmd->SetGuidance("Activate PAI in the G4Region.");
63 paiCmd->SetGuidance(" partName : particle name (default - all)");
64 paiCmd->SetGuidance(" regName : G4Region name");
65 paiCmd->SetGuidance(" paiType : PAI, PAIphoton");
66 paiCmd->AvailableForStates(G4State_PreInit);
67 paiCmd->SetToBeBroadcasted(false);
68
69 auto part = new G4UIparameter("partName",'s',false);
70 paiCmd->SetParameter(part);
71
72 auto pregName = new G4UIparameter("regName",'s',false);
73 paiCmd->SetParameter(pregName);
74
75 auto ptype = new G4UIparameter("type",'s',false);
76 paiCmd->SetParameter(ptype);
77 ptype->SetParameterCandidates("pai PAI PAIphoton");
78
79 mscoCmd = new G4UIcommand("/process/em/AddEmRegion",this);
80 mscoCmd->SetGuidance("Add optional EM configuration for a G4Region.");
81 mscoCmd->SetGuidance(" regName : G4Region name");
82 mscoCmd->SetGuidance(" emType : G4EmStandard, G4EmStandard_opt1, ...");
83 mscoCmd->AvailableForStates(G4State_PreInit);
84 mscoCmd->SetToBeBroadcasted(false);
85
86 auto mregName = new G4UIparameter("regName",'s',false);
87 mscoCmd->SetParameter(mregName);
88
89 auto mtype = new G4UIparameter("mscType",'s',false);
90 mscoCmd->SetParameter(mtype);
91 mtype->SetParameterCandidates("G4EmStandard G4EmStandard_opt1 G4EmStandard_opt2 G4EmStandard_opt3 G4EmStandard_opt4 G4EmStandardGS G4EmStandardSS G4EmLivermore G4EmPenelope G4RadioactiveDecay");
92
93 SubSecCmd = new G4UIcmdWithAString("/process/eLoss/subsecRegion",this);
94 SubSecCmd->SetGuidance("Enable subcut generation per region.");
95 SubSecCmd->SetGuidance(" Region : region name");
96 SubSecCmd->AvailableForStates(G4State_PreInit);
97 SubSecCmd->SetToBeBroadcasted(false);
98
99 StepFuncCmd = new G4UIcommand("/process/eLoss/StepFunction",this);
100 StepFuncCmd->SetGuidance("Set the energy loss step limitation parameters for e+-.");
101 StepFuncCmd->SetGuidance(" dRoverR : max Range variation per step");
102 StepFuncCmd->SetGuidance(" finalRange: range for final step");
103 StepFuncCmd->SetGuidance(" unit : unit of finalRange");
104 StepFuncCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
105 StepFuncCmd->SetToBeBroadcasted(false);
106
107 auto dRoverRPrm = new G4UIparameter("dRoverR",'d',false);
108 dRoverRPrm->SetParameterRange("dRoverR>0. && dRoverR<=1.");
109 StepFuncCmd->SetParameter(dRoverRPrm);
110
111 auto finalRangePrm = new G4UIparameter("finalRange",'d',false);
112 finalRangePrm->SetParameterRange("finalRange>0.");
113 StepFuncCmd->SetParameter(finalRangePrm);
114
115 auto unitPrm = new G4UIparameter("unit",'s',true);
116 unitPrm->SetDefaultUnit("mm");
117 StepFuncCmd->SetParameter(unitPrm);
118
119 StepFuncCmd1 = new G4UIcommand("/process/eLoss/StepFunctionMuHad",this);
120 StepFuncCmd1->SetGuidance("Set the energy loss step limitation parameters for muon/hadron.");
121 StepFuncCmd1->SetGuidance(" dRoverR : max Range variation per step");
122 StepFuncCmd1->SetGuidance(" finalRange: range for final step");
123 StepFuncCmd1->AvailableForStates(G4State_PreInit,G4State_Idle);
124 StepFuncCmd1->SetToBeBroadcasted(false);
125
126 auto dRoverRPrm1 = new G4UIparameter("dRoverRMuHad",'d',false);
127 dRoverRPrm1->SetParameterRange("dRoverRMuHad>0. && dRoverRMuHad<=1.");
128 StepFuncCmd1->SetParameter(dRoverRPrm1);
129
130 auto finalRangePrm1 = new G4UIparameter("finalRangeMuHad",'d',false);
131 finalRangePrm1->SetParameterRange("finalRangeMuHad>0.");
132 StepFuncCmd1->SetParameter(finalRangePrm1);
133
134 auto unitPrm1 = new G4UIparameter("unit",'s',true);
135 unitPrm1->SetDefaultValue("mm");
136 StepFuncCmd1->SetParameter(unitPrm1);
137
138 StepFuncCmd2 = new G4UIcommand("/process/eLoss/StepFunctionLightIons",this);
139 StepFuncCmd2->SetGuidance("Set the energy loss step limitation parameters for light ions.");
140 StepFuncCmd2->SetGuidance(" dRoverR : max Range variation per step");
141 StepFuncCmd2->SetGuidance(" finalRange: range for final step");
142 StepFuncCmd2->AvailableForStates(G4State_PreInit,G4State_Idle);
143 StepFuncCmd2->SetToBeBroadcasted(false);
144
145 auto dRoverRPrm2 = new G4UIparameter("dRoverRLIons",'d',false);
146 dRoverRPrm2->SetParameterRange("dRoverRLIons>0. && dRoverRLIons<=1.");
147 StepFuncCmd2->SetParameter(dRoverRPrm2);
148
149 auto finalRangePrm2 = new G4UIparameter("finalRangeLIons",'d',false);
150 finalRangePrm2->SetParameterRange("finalRangeLIons>0.");
151 StepFuncCmd2->SetParameter(finalRangePrm2);
152
153 auto unitPrm2 = new G4UIparameter("unit",'s',true);
154 unitPrm2->SetDefaultValue("mm");
155 StepFuncCmd2->SetParameter(unitPrm2);
156
157 StepFuncCmd3 = new G4UIcommand("/process/eLoss/StepFunctionIons",this);
158 StepFuncCmd3->SetGuidance("Set the energy loss step limitation parameters for ions.");
159 StepFuncCmd3->SetGuidance(" dRoverR : max Range variation per step");
160 StepFuncCmd3->SetGuidance(" finalRange: range for final step");
161 StepFuncCmd3->AvailableForStates(G4State_PreInit,G4State_Idle);
162 StepFuncCmd3->SetToBeBroadcasted(false);
163
164 auto dRoverRPrm3 = new G4UIparameter("dRoverRIons",'d',false);
165 dRoverRPrm3->SetParameterRange("dRoverRIons>0. && dRoverRIons<=1.");
166 StepFuncCmd3->SetParameter(dRoverRPrm3);
167
168 auto finalRangePrm3 = new G4UIparameter("finalRangeIons",'d',false);
169 finalRangePrm3->SetParameterRange("finalRangeIons>0.");
170 StepFuncCmd3->SetParameter(finalRangePrm3);
171
172 auto unitPrm3 = new G4UIparameter("unit",'s',true);
173 unitPrm3->SetDefaultValue("mm");
174 StepFuncCmd3->SetParameter(unitPrm3);
175
176 bfCmd = new G4UIcommand("/process/em/setBiasingFactor",this);
177 bfCmd->SetGuidance("Set factor for the process cross section.");
178 bfCmd->SetGuidance(" procName : process name");
179 bfCmd->SetGuidance(" procFact : factor");
180 bfCmd->SetGuidance(" flagFact : flag to change weight");
181 bfCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
182 bfCmd->SetToBeBroadcasted(false);
183
184 auto procName = new G4UIparameter("procName",'s',false);
185 bfCmd->SetParameter(procName);
186
187 auto procFact = new G4UIparameter("procFact",'d',false);
188 bfCmd->SetParameter(procFact);
189
190 auto flagFact = new G4UIparameter("flagFact",'s',false);
191 bfCmd->SetParameter(flagFact);
192
193 fiCmd = new G4UIcommand("/process/em/setForcedInteraction",this);
194 fiCmd->SetGuidance("Set factor for the process cross section.");
195 fiCmd->SetGuidance(" procNam : process name");
196 fiCmd->SetGuidance(" regNam : region name");
197 fiCmd->SetGuidance(" tlength : fixed target length");
198 fiCmd->SetGuidance(" unitT : length unit");
199 fiCmd->SetGuidance(" tflag : flag to change weight");
200 fiCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
201 fiCmd->SetToBeBroadcasted(false);
202
203 auto procNam = new G4UIparameter("procNam",'s',false);
204 fiCmd->SetParameter(procNam);
205
206 auto regNam = new G4UIparameter("regNam",'s',false);
207 fiCmd->SetParameter(regNam);
208
209 auto tlength = new G4UIparameter("tlength",'d',false);
210 tlength->SetParameterRange("tlength>0");
211 fiCmd->SetParameter(tlength);
212
213 auto unitT = new G4UIparameter("unitT",'s',true);
214 unitT->SetDefaultUnit("mm");
215 fiCmd->SetParameter(unitT);
216
217 auto flagT = new G4UIparameter("tflag",'b',true);
218 flagT->SetDefaultValue(true);
219 fiCmd->SetParameter(flagT);
220
221 bsCmd = new G4UIcommand("/process/em/setSecBiasing",this);
222 bsCmd->SetGuidance("Set bremsstrahlung or delta-e- splitting/Russian roulette per region.");
223 bsCmd->SetGuidance(" bProcNam : process name");
224 bsCmd->SetGuidance(" bRegNam : region name");
225 bsCmd->SetGuidance(" bFactor : number of split gamma or probability of Russian roulette");
226 bsCmd->SetGuidance(" bEnergy : max energy of a secondary for this biasing method");
227 bsCmd->SetGuidance(" bUnit : energy unit");
228 bsCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
229 bsCmd->SetToBeBroadcasted(false);
230
231 auto bProcNam = new G4UIparameter("bProcNam",'s',false);
232 bsCmd->SetParameter(bProcNam);
233
234 auto bRegNam = new G4UIparameter("bRegNam",'s',false);
235 bsCmd->SetParameter(bRegNam);
236
237 auto bFactor = new G4UIparameter("bFactor",'d',false);
238 bsCmd->SetParameter(bFactor);
239
240 auto bEnergy = new G4UIparameter("bEnergy",'d',false);
241 bsCmd->SetParameter(bEnergy);
242
243 auto bUnit = new G4UIparameter("bUnit",'s',true);
244 bUnit->SetDefaultUnit("MeV");
245 bsCmd->SetParameter(bUnit);
246
247 dirSplitCmd = new G4UIcmdWithABool("/process/em/setDirectionalSplitting",this);
248 dirSplitCmd->SetGuidance("Enable directional brem splitting");
249 dirSplitCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
250 dirSplitCmd->SetToBeBroadcasted(false);
251
252 qeCmd = new G4UIcmdWithABool("/process/em/QuantumEntanglement",this);
253 qeCmd->SetGuidance("Enable quantum entanglement");
254 qeCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
255 qeCmd->SetToBeBroadcasted(false);
256
257 dirSplitTargetCmd = new G4UIcmdWith3VectorAndUnit("/process/em/setDirectionalSplittingTarget",this);
258 dirSplitTargetCmd->SetGuidance("Position of arget for directional splitting");
259 dirSplitTargetCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
260
261 dirSplitRadiusCmd = new G4UIcmdWithADoubleAndUnit("/process/em/setDirectionalSplittingRadius",this);
262 dirSplitRadiusCmd->SetGuidance("Radius of target for directional splitting");
263 dirSplitRadiusCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
264 dirSplitRadiusCmd->SetToBeBroadcasted(false);
265}
@ G4State_Idle
@ G4State_PreInit

Referenced by G4EmExtraParametersMessenger(), and operator=().

◆ ~G4EmExtraParametersMessenger()

G4EmExtraParametersMessenger::~G4EmExtraParametersMessenger ( )
override

Definition at line 269 of file G4EmExtraParametersMessenger.cc.

270{
271 delete paiCmd;
272 delete mscoCmd;
273 delete SubSecCmd;
274 delete bfCmd;
275 delete fiCmd;
276 delete bsCmd;
277 delete qeCmd;
278 delete StepFuncCmd;
279 delete StepFuncCmd1;
280 delete StepFuncCmd2;
281 delete StepFuncCmd3;
282 delete dirSplitCmd;
283 delete dirSplitTargetCmd;
284 delete dirSplitRadiusCmd;
285}

◆ G4EmExtraParametersMessenger() [2/2]

G4EmExtraParametersMessenger::G4EmExtraParametersMessenger ( const G4EmExtraParametersMessenger & )
delete

Member Function Documentation

◆ operator=()

G4EmExtraParametersMessenger & G4EmExtraParametersMessenger::operator= ( const G4EmExtraParametersMessenger & right)
delete

◆ SetNewValue()

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

Reimplemented from G4UImessenger.

Definition at line 289 of file G4EmExtraParametersMessenger.cc.

291{
292 G4bool physicsModified = false;
293
294 if (command == paiCmd) {
295 G4String s1(""),s2(""),s3("");
296 std::istringstream is(newValue);
297 is >> s1 >> s2 >> s3;
298 theParameters->AddPAIModel(s1, s2, s3);
299 } else if (command == mscoCmd) {
300 G4String s1(""),s2("");
301 std::istringstream is(newValue);
302 is >> s1 >> s2;
303 theParameters->AddPhysics(s1, s2);
304 } else if (command == StepFuncCmd || command == StepFuncCmd1 || command == StepFuncCmd2 || command == StepFuncCmd3) {
305 G4double v1,v2;
306 G4String unt;
307 std::istringstream is(newValue);
308 is >> v1 >> v2 >> unt;
309 v2 *= G4UIcommand::ValueOf(unt);
310 if(command == StepFuncCmd) {
311 theParameters->SetStepFunction(v1,v2);
312 } else if(command == StepFuncCmd1) {
313 theParameters->SetStepFunctionMuHad(v1,v2);
314 } else if(command == StepFuncCmd2) {
315 theParameters->SetStepFunctionLightIons(v1,v2);
316 } else {
317 theParameters->SetStepFunctionIons(v1,v2);
318 }
319 physicsModified = true;
320 } else if (command == SubSecCmd) {
321 theParameters->SetSubCutRegion(newValue);
322 } else if (command == bfCmd) {
323 G4double v1(1.0);
324 G4String s0(""),s1("");
325 std::istringstream is(newValue);
326 is >> s0 >> v1 >> s1;
327 G4bool yes = false;
328 if(s1 == "true") { yes = true; }
329 theParameters->SetProcessBiasingFactor(s0,v1,yes);
330 physicsModified = true;
331 } else if (command == fiCmd) {
332 G4double v1(0.0);
333 G4String s1(""),s2(""),s3(""),unt("mm");
334 std::istringstream is(newValue);
335 is >> s1 >> s2 >> v1 >> unt >> s3;
336 G4bool yes = false;
337 if(s3 == "true") { yes = true; }
338 v1 *= G4UIcommand::ValueOf(unt);
339 theParameters->ActivateForcedInteraction(s1,s2,v1,yes);
340 physicsModified = true;
341 } else if (command == bsCmd) {
342 G4double fb(1.0),en(1.e+30);
343 G4String s1(""),s2(""),unt("MeV");
344 std::istringstream is(newValue);
345 is >> s1 >> s2 >> fb >> en >> unt;
346 en *= G4UIcommand::ValueOf(unt);
347 theParameters->ActivateSecondaryBiasing(s1,s2,fb,en);
348 physicsModified = true;
349 } else if (command == qeCmd) {
350 theParameters->SetQuantumEntanglement(qeCmd->GetNewBoolValue(newValue));
351 } else if (command == dirSplitCmd) {
352 theParameters->SetDirectionalSplitting(
353 dirSplitCmd->GetNewBoolValue(newValue));
354 physicsModified = true;
355 } else if (command == dirSplitTargetCmd) {
356 G4ThreeVector t = dirSplitTargetCmd->GetNew3VectorValue(newValue);
357 theParameters->SetDirectionalSplittingTarget(t);
358 physicsModified = true;
359 } else if (command == dirSplitRadiusCmd) {
360 G4double r = dirSplitRadiusCmd->GetNewDoubleValue(newValue);
361 theParameters->SetDirectionalSplittingRadius(r);
362 physicsModified = true;
363 }
364
365 if(physicsModified) {
366 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
367 }
368}
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
static G4double ValueOf(const char *unitName)
G4int ApplyCommand(const char *aCommand)
static G4UImanager * GetUIpointer()

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