Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4EmExtraParametersMessenger.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// -------------------------------------------------------------------
27//
28// GEANT4 Class file
29//
30// File name: G4EmExtraParametersMessenger
31//
32// Author: Vladimir Ivanchenko
33//
34// Creation date: 07-05-2019
35//
36// -------------------------------------------------------------------
37//
38
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41
43#include "G4UIcommand.hh"
44#include "G4UIparameter.hh"
45#include "G4UIcmdWithABool.hh"
47#include "G4UIcmdWithADouble.hh"
49#include "G4UIcmdWithAString.hh"
51#include "G4UImanager.hh"
53
54#include <sstream>
55
56//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
57
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");
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, ...");
84
85 auto mregName = new G4UIparameter("regName",'s',false);
86 mscoCmd->SetParameter(mregName);
87
88 auto mtype = new G4UIparameter("mscType",'s',false);
89 mscoCmd->SetParameter(mtype);
90 mtype->SetParameterCandidates("G4EmStandard G4EmStandard_opt1 G4EmStandard_opt2 G4EmStandard_opt3 G4EmStandard_opt4 G4EmStandardGS G4EmStandardSS G4EmLivermore G4EmPenelope G4RadioactiveDecay");
91
92 SubSecCmd = new G4UIcmdWithAString("/process/eLoss/subsecRegion",this);
93 SubSecCmd->SetGuidance("Enable subcut generation per region.");
94 SubSecCmd->SetGuidance(" Region : region name");
96 SubSecCmd->SetToBeBroadcasted(false);
97
98 StepFuncCmd = new G4UIcommand("/process/eLoss/StepFunction",this);
99 StepFuncCmd->SetGuidance("Set the energy loss step limitation parameters for e+-.");
100 StepFuncCmd->SetGuidance(" dRoverR : max Range variation per step");
101 StepFuncCmd->SetGuidance(" finalRange: range for final step");
102 StepFuncCmd->SetGuidance(" unit : unit of finalRange");
104 StepFuncCmd->SetToBeBroadcasted(false);
105
106 auto dRoverRPrm = new G4UIparameter("dRoverR",'d',false);
107 dRoverRPrm->SetParameterRange("dRoverR>0. && dRoverR<=1.");
108 StepFuncCmd->SetParameter(dRoverRPrm);
109
110 auto finalRangePrm = new G4UIparameter("finalRange",'d',false);
111 finalRangePrm->SetParameterRange("finalRange>0.");
112 StepFuncCmd->SetParameter(finalRangePrm);
113
114 auto unitPrm = new G4UIparameter("unit",'s',true);
115 unitPrm->SetDefaultUnit("mm");
116 StepFuncCmd->SetParameter(unitPrm);
117
118 StepFuncCmd1 = new G4UIcommand("/process/eLoss/StepFunctionMuHad",this);
119 StepFuncCmd1->SetGuidance("Set the energy loss step limitation parameters for muon/hadron.");
120 StepFuncCmd1->SetGuidance(" dRoverR : max Range variation per step");
121 StepFuncCmd1->SetGuidance(" finalRange: range for final step");
123 StepFuncCmd1->SetToBeBroadcasted(false);
124
125 auto dRoverRPrm1 = new G4UIparameter("dRoverRMuHad",'d',false);
126 dRoverRPrm1->SetParameterRange("dRoverRMuHad>0. && dRoverRMuHad<=1.");
127 StepFuncCmd1->SetParameter(dRoverRPrm1);
128
129 auto finalRangePrm1 = new G4UIparameter("finalRangeMuHad",'d',false);
130 finalRangePrm1->SetParameterRange("finalRangeMuHad>0.");
131 StepFuncCmd1->SetParameter(finalRangePrm1);
132
133 auto unitPrm1 = new G4UIparameter("unit",'s',true);
134 unitPrm1->SetDefaultValue("mm");
135 StepFuncCmd1->SetParameter(unitPrm1);
136
137 StepFuncCmd2 = new G4UIcommand("/process/eLoss/StepFunctionLightIons",this);
138 StepFuncCmd2->SetGuidance("Set the energy loss step limitation parameters for light ions.");
139 StepFuncCmd2->SetGuidance(" dRoverR : max Range variation per step");
140 StepFuncCmd2->SetGuidance(" finalRange: range for final step");
142 StepFuncCmd2->SetToBeBroadcasted(false);
143
144 auto dRoverRPrm2 = new G4UIparameter("dRoverRLIons",'d',false);
145 dRoverRPrm2->SetParameterRange("dRoverRLIons>0. && dRoverRLIons<=1.");
146 StepFuncCmd2->SetParameter(dRoverRPrm2);
147
148 auto finalRangePrm2 = new G4UIparameter("finalRangeLIons",'d',false);
149 finalRangePrm2->SetParameterRange("finalRangeLIons>0.");
150 StepFuncCmd2->SetParameter(finalRangePrm2);
151
152 auto unitPrm2 = new G4UIparameter("unit",'s',true);
153 unitPrm2->SetDefaultValue("mm");
154 StepFuncCmd2->SetParameter(unitPrm2);
155
156 StepFuncCmd3 = new G4UIcommand("/process/eLoss/StepFunctionIons",this);
157 StepFuncCmd3->SetGuidance("Set the energy loss step limitation parameters for ions.");
158 StepFuncCmd3->SetGuidance(" dRoverR : max Range variation per step");
159 StepFuncCmd3->SetGuidance(" finalRange: range for final step");
161 StepFuncCmd3->SetToBeBroadcasted(false);
162
163 auto dRoverRPrm3 = new G4UIparameter("dRoverRIons",'d',false);
164 dRoverRPrm3->SetParameterRange("dRoverRIons>0. && dRoverRIons<=1.");
165 StepFuncCmd3->SetParameter(dRoverRPrm3);
166
167 auto finalRangePrm3 = new G4UIparameter("finalRangeIons",'d',false);
168 finalRangePrm3->SetParameterRange("finalRangeIons>0.");
169 StepFuncCmd3->SetParameter(finalRangePrm3);
170
171 auto unitPrm3 = new G4UIparameter("unit",'s',true);
172 unitPrm3->SetDefaultValue("mm");
173 StepFuncCmd3->SetParameter(unitPrm3);
174
175 bfCmd = new G4UIcommand("/process/em/setBiasingFactor",this);
176 bfCmd->SetGuidance("Set factor for the process cross section.");
177 bfCmd->SetGuidance(" procName : process name");
178 bfCmd->SetGuidance(" procFact : factor");
179 bfCmd->SetGuidance(" flagFact : flag to change weight");
181 bfCmd->SetToBeBroadcasted(false);
182
183 auto procName = new G4UIparameter("procName",'s',false);
184 bfCmd->SetParameter(procName);
185
186 auto procFact = new G4UIparameter("procFact",'d',false);
187 bfCmd->SetParameter(procFact);
188
189 auto flagFact = new G4UIparameter("flagFact",'s',false);
190 bfCmd->SetParameter(flagFact);
191
192 fiCmd = new G4UIcommand("/process/em/setForcedInteraction",this);
193 fiCmd->SetGuidance("Set factor for the process cross section.");
194 fiCmd->SetGuidance(" procNam : process name");
195 fiCmd->SetGuidance(" regNam : region name");
196 fiCmd->SetGuidance(" tlength : fixed target length");
197 fiCmd->SetGuidance(" unitT : length unit");
198 fiCmd->SetGuidance(" tflag : flag to change weight");
200 fiCmd->SetToBeBroadcasted(false);
201
202 auto procNam = new G4UIparameter("procNam",'s',false);
203 fiCmd->SetParameter(procNam);
204
205 auto regNam = new G4UIparameter("regNam",'s',false);
206 fiCmd->SetParameter(regNam);
207
208 auto tlength = new G4UIparameter("tlength",'d',false);
209 tlength->SetParameterRange("tlength>0");
210 fiCmd->SetParameter(tlength);
211
212 auto unitT = new G4UIparameter("unitT",'s',true);
213 unitT->SetDefaultUnit("mm");
214 fiCmd->SetParameter(unitT);
215
216 auto flagT = new G4UIparameter("tflag",'b',true);
217 flagT->SetDefaultValue(true);
218 fiCmd->SetParameter(flagT);
219
220 bsCmd = new G4UIcommand("/process/em/setSecBiasing",this);
221 bsCmd->SetGuidance("Set bremsstrahlung or delta-e- splitting/Russian roulette per region.");
222 bsCmd->SetGuidance(" bProcNam : process name");
223 bsCmd->SetGuidance(" bRegNam : region name");
224 bsCmd->SetGuidance(" bFactor : number of split gamma or probability of Russian roulette");
225 bsCmd->SetGuidance(" bEnergy : max energy of a secondary for this biasing method");
226 bsCmd->SetGuidance(" bUnit : energy unit");
228 bsCmd->SetToBeBroadcasted(false);
229
230 auto bProcNam = new G4UIparameter("bProcNam",'s',false);
231 bsCmd->SetParameter(bProcNam);
232
233 auto bRegNam = new G4UIparameter("bRegNam",'s',false);
234 bsCmd->SetParameter(bRegNam);
235
236 auto bFactor = new G4UIparameter("bFactor",'d',false);
237 bsCmd->SetParameter(bFactor);
238
239 auto bEnergy = new G4UIparameter("bEnergy",'d',false);
240 bsCmd->SetParameter(bEnergy);
241
242 auto bUnit = new G4UIparameter("bUnit",'s',true);
243 bUnit->SetDefaultUnit("MeV");
244 bsCmd->SetParameter(bUnit);
245
246 dirSplitCmd = new G4UIcmdWithABool("/process/em/setDirectionalSplitting",this);
247 dirSplitCmd->SetGuidance("Enable directional brem splitting");
249 dirSplitCmd->SetToBeBroadcasted(false);
250
251 qeCmd = new G4UIcmdWithABool("/process/em/QuantumEntanglement",this);
252 qeCmd->SetGuidance("Enable quantum entanglement");
254 qeCmd->SetToBeBroadcasted(false);
255
256 dirSplitTargetCmd = new G4UIcmdWith3VectorAndUnit("/process/em/setDirectionalSplittingTarget",this);
257 dirSplitTargetCmd->SetGuidance("Position of arget for directional splitting");
259
260 dirSplitRadiusCmd = new G4UIcmdWithADoubleAndUnit("/process/em/setDirectionalSplittingRadius",this);
261 dirSplitRadiusCmd->SetGuidance("Radius of target for directional splitting");
263 dirSplitRadiusCmd->SetToBeBroadcasted(false);
264}
265
266//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
267
269{
270 delete paiCmd;
271 delete mscoCmd;
272 delete SubSecCmd;
273 delete bfCmd;
274 delete fiCmd;
275 delete bsCmd;
276 delete qeCmd;
277 delete StepFuncCmd;
278 delete StepFuncCmd1;
279 delete StepFuncCmd2;
280 delete StepFuncCmd3;
281 delete dirSplitCmd;
282 delete dirSplitTargetCmd;
283 delete dirSplitRadiusCmd;
284}
285
286//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
287
289 G4String newValue)
290{
291 G4bool physicsModified = false;
292
293 if (command == paiCmd) {
294 G4String s1(""),s2(""),s3("");
295 std::istringstream is(newValue);
296 is >> s1 >> s2 >> s3;
297 theParameters->AddPAIModel(s1, s2, s3);
298 } else if (command == mscoCmd) {
299 G4String s1(""),s2("");
300 std::istringstream is(newValue);
301 is >> s1 >> s2;
302 theParameters->AddPhysics(s1, s2);
303 } else if (command == StepFuncCmd || command == StepFuncCmd1 || command == StepFuncCmd2 || command == StepFuncCmd3) {
304 G4double v1,v2;
305 G4String unt;
306 std::istringstream is(newValue);
307 is >> v1 >> v2 >> unt;
308 v2 *= G4UIcommand::ValueOf(unt);
309 if(command == StepFuncCmd) {
310 theParameters->SetStepFunction(v1,v2);
311 } else if(command == StepFuncCmd1) {
312 theParameters->SetStepFunctionMuHad(v1,v2);
313 } else if(command == StepFuncCmd2) {
314 theParameters->SetStepFunctionLightIons(v1,v2);
315 } else {
316 theParameters->SetStepFunctionIons(v1,v2);
317 }
318 physicsModified = true;
319 } else if (command == SubSecCmd) {
320 theParameters->SetSubCutRegion(newValue);
321 } else if (command == bfCmd) {
322 G4double v1(1.0);
323 G4String s0(""),s1("");
324 std::istringstream is(newValue);
325 is >> s0 >> v1 >> s1;
326 G4bool yes = false;
327 if(s1 == "true") { yes = true; }
328 theParameters->SetProcessBiasingFactor(s0,v1,yes);
329 physicsModified = true;
330 } else if (command == fiCmd) {
331 G4double v1(0.0);
332 G4String s1(""),s2(""),s3(""),unt("mm");
333 std::istringstream is(newValue);
334 is >> s1 >> s2 >> v1 >> unt >> s3;
335 G4bool yes = false;
336 if(s3 == "true") { yes = true; }
337 v1 *= G4UIcommand::ValueOf(unt);
338 theParameters->ActivateForcedInteraction(s1,s2,v1,yes);
339 physicsModified = true;
340 } else if (command == bsCmd) {
341 G4double fb(1.0),en(1.e+30);
342 G4String s1(""),s2(""),unt("MeV");
343 std::istringstream is(newValue);
344 is >> s1 >> s2 >> fb >> en >> unt;
345 en *= G4UIcommand::ValueOf(unt);
346 theParameters->ActivateSecondaryBiasing(s1,s2,fb,en);
347 physicsModified = true;
348 } else if (command == qeCmd) {
349 theParameters->SetQuantumEntanglement(qeCmd->GetNewBoolValue(newValue));
350 } else if (command == dirSplitCmd) {
351 theParameters->SetDirectionalSplitting(
352 dirSplitCmd->GetNewBoolValue(newValue));
353 physicsModified = true;
354 } else if (command == dirSplitTargetCmd) {
355 G4ThreeVector t = dirSplitTargetCmd->GetNew3VectorValue(newValue);
356 theParameters->SetDirectionalSplittingTarget(t);
357 physicsModified = true;
358 } else if (command == dirSplitRadiusCmd) {
359 G4double r = dirSplitRadiusCmd->GetNewDoubleValue(newValue);
360 theParameters->SetDirectionalSplittingRadius(r);
361 physicsModified = true;
362 }
363
364 if(physicsModified) {
365 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
366 }
367}
368
369//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ G4State_Idle
@ G4State_PreInit
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4EmExtraParametersMessenger(G4EmExtraParameters *)
void SetNewValue(G4UIcommand *, G4String) override
void SetQuantumEntanglement(G4bool v)
void SetStepFunctionMuHad(G4double v1, G4double v2)
void SetDirectionalSplittingRadius(G4double r)
void SetStepFunctionIons(G4double v1, G4double v2)
void ActivateForcedInteraction(const G4String &procname, const G4String &region, G4double length, G4bool wflag)
void ActivateSecondaryBiasing(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
void SetDirectionalSplitting(G4bool v)
void SetDirectionalSplittingTarget(const G4ThreeVector &v)
void SetSubCutRegion(const G4String &region)
void SetProcessBiasingFactor(const G4String &procname, G4double val, G4bool wflag)
void SetStepFunction(G4double v1, G4double v2)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void AddPAIModel(const G4String &particle, const G4String &region, const G4String &type)
void AddPhysics(const G4String &region, const G4String &type)
static G4ThreeVector GetNew3VectorValue(const char *paramString)
static G4bool GetNewBoolValue(const char *paramString)
static G4double GetNewDoubleValue(const char *paramString)
void SetToBeBroadcasted(G4bool val)
Definition: G4UIcommand.hh:172
static G4double ValueOf(const char *unitName)
Definition: G4UIcommand.cc:362
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:147
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:157
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:287
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:495
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77