Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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
68 G4UIparameter* part = new G4UIparameter("partName",'s',false);
69 paiCmd->SetParameter(part);
70
71 G4UIparameter* pregName = new G4UIparameter("regName",'s',false);
72 paiCmd->SetParameter(pregName);
73
74 G4UIparameter* ptype = new G4UIparameter("type",'s',false);
75 paiCmd->SetParameter(ptype);
76 ptype->SetParameterCandidates("pai PAI PAIphoton");
77
78 mscoCmd = new G4UIcommand("/process/em/AddEmRegion",this);
79 mscoCmd->SetGuidance("Add optional EM configuration for a G4Region.");
80 mscoCmd->SetGuidance(" regName : G4Region name");
81 mscoCmd->SetGuidance(" emType : G4EmStandard, G4EmStandard_opt1, ...");
83
84 G4UIparameter* mregName = new G4UIparameter("regName",'s',false);
85 mscoCmd->SetParameter(mregName);
86
87 G4UIparameter* mtype = new G4UIparameter("mscType",'s',false);
88 mscoCmd->SetParameter(mtype);
89 mtype->SetParameterCandidates("G4EmStandard G4EmStandard_opt1 G4EmStandard_opt2 G4EmStandard_opt3 G4EmStandard_opt4 G4EmStandardGS G4EmStandardSS G4EmLivermore G4EmPenelope G4RadioactiveDecay");
90
91 SubSecCmd = new G4UIcommand("/process/eLoss/subsec",this);
92 SubSecCmd->SetGuidance("Switch true/false the subcutoff generation per region.");
93 SubSecCmd->SetGuidance(" subSec : true/false");
94 SubSecCmd->SetGuidance(" Region : region name");
96
97 StepFuncCmd = new G4UIcommand("/process/eLoss/StepFunction",this);
98 StepFuncCmd->SetGuidance("Set the energy loss step limitation parameters for e+-.");
99 StepFuncCmd->SetGuidance(" dRoverR : max Range variation per step");
100 StepFuncCmd->SetGuidance(" finalRange: range for final step");
101 StepFuncCmd->SetGuidance(" unit : unit of finalRange");
103
104 G4UIparameter* dRoverRPrm = new G4UIparameter("dRoverR",'d',false);
105 dRoverRPrm->SetParameterRange("dRoverR>0. && dRoverR<=1.");
106 StepFuncCmd->SetParameter(dRoverRPrm);
107
108 G4UIparameter* finalRangePrm = new G4UIparameter("finalRange",'d',false);
109 finalRangePrm->SetParameterRange("finalRange>0.");
110 StepFuncCmd->SetParameter(finalRangePrm);
111
112 G4UIparameter* unitPrm = new G4UIparameter("unit",'s',true);
113 unitPrm->SetDefaultUnit("mm");
114 StepFuncCmd->SetParameter(unitPrm);
115
116 StepFuncCmd1 = new G4UIcommand("/process/eLoss/StepFunctionMuHad",this);
117 StepFuncCmd1->SetGuidance("Set the energy loss step limitation parameters for muon/hadron.");
118 StepFuncCmd1->SetGuidance(" dRoverR : max Range variation per step");
119 StepFuncCmd1->SetGuidance(" finalRange: range for final step");
121
122 G4UIparameter* dRoverRPrm1 = new G4UIparameter("dRoverRMuHad",'d',false);
123 dRoverRPrm1->SetParameterRange("dRoverRMuHad>0. && dRoverRMuHad<=1.");
124 StepFuncCmd1->SetParameter(dRoverRPrm1);
125
126 G4UIparameter* finalRangePrm1 = new G4UIparameter("finalRangeMuHad",'d',false);
127 finalRangePrm1->SetParameterRange("finalRangeMuHad>0.");
128 StepFuncCmd1->SetParameter(finalRangePrm1);
129
130 G4UIparameter* unitPrm1 = new G4UIparameter("unit",'s',true);
131 unitPrm1->SetDefaultValue("mm");
132 StepFuncCmd1->SetParameter(unitPrm1);
133
134 StepFuncCmd2 = new G4UIcommand("/process/eLoss/StepFunctionLightIons",this);
135 StepFuncCmd2->SetGuidance("Set the energy loss step limitation parameters for light ions.");
136 StepFuncCmd2->SetGuidance(" dRoverR : max Range variation per step");
137 StepFuncCmd2->SetGuidance(" finalRange: range for final step");
139
140 G4UIparameter* dRoverRPrm2 = new G4UIparameter("dRoverRLIons",'d',false);
141 dRoverRPrm2->SetParameterRange("dRoverRLIons>0. && dRoverRLIons<=1.");
142 StepFuncCmd2->SetParameter(dRoverRPrm2);
143
144 G4UIparameter* finalRangePrm2 = new G4UIparameter("finalRangeLIons",'d',false);
145 finalRangePrm2->SetParameterRange("finalRangeLIons>0.");
146 StepFuncCmd2->SetParameter(finalRangePrm2);
147
148 G4UIparameter* unitPrm2 = new G4UIparameter("unit",'s',true);
149 unitPrm2->SetDefaultValue("mm");
150 StepFuncCmd2->SetParameter(unitPrm2);
151
152 StepFuncCmd3 = new G4UIcommand("/process/eLoss/StepFunctionIons",this);
153 StepFuncCmd3->SetGuidance("Set the energy loss step limitation parameters for ions.");
154 StepFuncCmd3->SetGuidance(" dRoverR : max Range variation per step");
155 StepFuncCmd3->SetGuidance(" finalRange: range for final step");
157
158 G4UIparameter* dRoverRPrm3 = new G4UIparameter("dRoverRMuHad",'d',false);
159 dRoverRPrm3->SetParameterRange("dRoverRIons>0. && dRoverRIons<=1.");
160 StepFuncCmd3->SetParameter(dRoverRPrm3);
161
162 G4UIparameter* finalRangePrm3 = new G4UIparameter("finalRangeIons",'d',false);
163 finalRangePrm3->SetParameterRange("finalRangeIons>0.");
164 StepFuncCmd3->SetParameter(finalRangePrm3);
165
166 G4UIparameter* unitPrm3 = new G4UIparameter("unit",'s',true);
167 unitPrm3->SetDefaultValue("mm");
168 StepFuncCmd3->SetParameter(unitPrm3);
169
170 G4UIparameter* subSec = new G4UIparameter("subSec",'s',false);
171 SubSecCmd->SetParameter(subSec);
172
173 G4UIparameter* subSecReg = new G4UIparameter("Region",'s',false);
174 SubSecCmd->SetParameter(subSecReg);
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");
182
183 G4UIparameter* procName = new G4UIparameter("procName",'s',false);
184 bfCmd->SetParameter(procName);
185
186 G4UIparameter* procFact = new G4UIparameter("procFact",'d',false);
187 bfCmd->SetParameter(procFact);
188
189 G4UIparameter* 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
201 G4UIparameter* procNam = new G4UIparameter("procNam",'s',false);
202 fiCmd->SetParameter(procNam);
203
204 G4UIparameter* regNam = new G4UIparameter("regNam",'s',false);
205 fiCmd->SetParameter(regNam);
206
207 G4UIparameter* tlength = new G4UIparameter("tlength",'d',false);
208 tlength->SetParameterRange("tlength>0");
209 fiCmd->SetParameter(tlength);
210
211 G4UIparameter* unitT = new G4UIparameter("unitT",'s',true);
212 unitT->SetDefaultUnit("mm");
213 fiCmd->SetParameter(unitT);
214
215 G4UIparameter* flagT = new G4UIparameter("tflag",'b',true);
216 flagT->SetDefaultValue(true);
217 fiCmd->SetParameter(flagT);
218
219 bsCmd = new G4UIcommand("/process/em/setSecBiasing",this);
220 bsCmd->SetGuidance("Set bremsstrahlung or delta-e- splitting/Russian roulette per region.");
221 bsCmd->SetGuidance(" bProcNam : process name");
222 bsCmd->SetGuidance(" bRegNam : region name");
223 bsCmd->SetGuidance(" bFactor : number of split gamma or probability of Russian roulette");
224 bsCmd->SetGuidance(" bEnergy : max energy of a secondary for this biasing method");
225 bsCmd->SetGuidance(" bUnit : energy unit");
227
228 G4UIparameter* bProcNam = new G4UIparameter("bProcNam",'s',false);
229 bsCmd->SetParameter(bProcNam);
230
231 G4UIparameter* bRegNam = new G4UIparameter("bRegNam",'s',false);
232 bsCmd->SetParameter(bRegNam);
233
234 G4UIparameter* bFactor = new G4UIparameter("bFactor",'d',false);
235 bsCmd->SetParameter(bFactor);
236
237 G4UIparameter* bEnergy = new G4UIparameter("bEnergy",'d',false);
238 bsCmd->SetParameter(bEnergy);
239
240 G4UIparameter* bUnit = new G4UIparameter("bUnit",'s',true);
241 bUnit->SetDefaultUnit("MeV");
242 bsCmd->SetParameter(bUnit);
243
244 dirSplitCmd = new G4UIcmdWithABool("/process/em/setDirectionalSplitting",this);
245 dirSplitCmd->SetGuidance("Enable directional brem splitting");
247
248 qeCmd = new G4UIcmdWithABool("/process/em/QuantumEntanglement",this);
249 qeCmd->SetGuidance("Enable quantum entanglement");
251
252 dirSplitTargetCmd = new G4UIcmdWith3VectorAndUnit("/process/em/setDirectionalSplittingTarget",this);
253 dirSplitTargetCmd->SetGuidance("Position of arget for directional splitting");
255
256 dirSplitRadiusCmd = new G4UIcmdWithADoubleAndUnit("/process/em/setDirectionalSplittingRadius",this);
257 dirSplitRadiusCmd->SetGuidance("Radius of target for directional splitting");
259}
260
261//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
262
264{
265 delete paiCmd;
266 delete mscoCmd;
267 delete SubSecCmd;
268 delete bfCmd;
269 delete fiCmd;
270 delete bsCmd;
271 delete qeCmd;
272 delete StepFuncCmd;
273 delete StepFuncCmd1;
274 delete StepFuncCmd2;
275 delete StepFuncCmd3;
276 delete dirSplitCmd;
277 delete dirSplitTargetCmd;
278 delete dirSplitRadiusCmd;
279}
280
281//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
282
284 G4String newValue)
285{
286 G4bool physicsModified = false;
287
288 if (command == paiCmd) {
289 G4String s1(""),s2(""),s3("");
290 std::istringstream is(newValue);
291 is >> s1 >> s2 >> s3;
292 theParameters->AddPAIModel(s1, s2, s3);
293 } else if (command == mscoCmd) {
294 G4String s1(""),s2("");
295 std::istringstream is(newValue);
296 is >> s1 >> s2;
297 theParameters->AddPhysics(s1, s2);
298 } else if (command == StepFuncCmd || command == StepFuncCmd1 || command == StepFuncCmd2 || command == StepFuncCmd3) {
299 G4double v1,v2;
300 G4String unt;
301 std::istringstream is(newValue);
302 is >> v1 >> v2 >> unt;
303 v2 *= G4UIcommand::ValueOf(unt);
304 if(command == StepFuncCmd) {
305 theParameters->SetStepFunction(v1,v2);
306 } else if(command == StepFuncCmd1) {
307 theParameters->SetStepFunctionMuHad(v1,v2);
308 } else if(command == StepFuncCmd2) {
309 theParameters->SetStepFunctionLightIons(v1,v2);
310 } else {
311 theParameters->SetStepFunctionIons(v1,v2);
312 }
313 physicsModified = true;
314 } else if (command == SubSecCmd) {
315 G4String s1, s2;
316 std::istringstream is(newValue);
317 is >> s1 >> s2;
318 G4bool yes = false;
319 if(s1 == "true") { yes = true; }
320 theParameters->SetSubCutoff(yes,s2);
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 *)
virtual void SetNewValue(G4UIcommand *, G4String) override
void SetQuantumEntanglement(G4bool v)
void SetStepFunctionMuHad(G4double v1, G4double v2)
void SetDirectionalSplittingRadius(G4double r)
void SetSubCutoff(G4bool val, const G4String &region="")
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 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)
static G4double ValueOf(const char *unitName)
Definition: G4UIcommand.cc:348
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:146
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:273
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:485
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
void SetDefaultValue(const char *theDefaultValue)
void SetParameterRange(const char *theRange)
void SetParameterCandidates(const char *theString)
void SetDefaultUnit(const char *theDefaultUnit)