Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmParametersMessenger.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: G4EmParametersMessenger
31//
32// Author: Vladimir Ivanchenko created from G4EnergyLossMessenger
33//
34// Creation date: 22-05-2013
35//
36// -------------------------------------------------------------------
37//
38
39//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
40//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
41
43#include "G4UIdirectory.hh"
44#include "G4UIcommand.hh"
45#include "G4UIparameter.hh"
46#include "G4UIcmdWithABool.hh"
48#include "G4UIcmdWithADouble.hh"
50#include "G4UIcmdWithAString.hh"
52#include "G4UImanager.hh"
53#include "G4MscStepLimitType.hh"
55#include "G4EmParameters.hh"
56
57#include <sstream>
58
59//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
60
62 : theParameters(ptr)
63{
64 emDirectory = new G4UIdirectory("/process/em/", false);
65 emDirectory->SetGuidance("General commands for EM processes.");
66 eLossDirectory = new G4UIdirectory("/process/eLoss/", false);
67 eLossDirectory->SetGuidance("Commands for energy loss processes.");
68 mscDirectory = new G4UIdirectory("/process/msc/", false);
69 mscDirectory->SetGuidance("Commands for EM scattering processes.");
70 gconvDirectory = new G4UIdirectory("/process/gconv/", false);
71 gconvDirectory->SetGuidance("Commands for EM gamma conversion BH5D model.");
72 dnaDirectory = new G4UIdirectory("/process/dna/", false);
73 dnaDirectory->SetGuidance("Commands for DNA processes.");
74
75 flucCmd = new G4UIcmdWithABool("/process/eLoss/fluct",this);
76 flucCmd->SetGuidance("Enable/disable energy loss fluctuations.");
77 flucCmd->SetParameterName("choice",true);
78 flucCmd->SetDefaultValue(true);
79 flucCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
80 flucCmd->SetToBeBroadcasted(false);
81
82 rangeCmd = new G4UIcmdWithABool("/process/eLoss/CSDARange",this);
83 rangeCmd->SetGuidance("Enable/disable CSDA range calculation");
84 rangeCmd->SetParameterName("range",true);
85 rangeCmd->SetDefaultValue(false);
86 rangeCmd->AvailableForStates(G4State_PreInit);
87 rangeCmd->SetToBeBroadcasted(false);
88
89 lpmCmd = new G4UIcmdWithABool("/process/eLoss/LPM",this);
90 lpmCmd->SetGuidance("Enable/disable LPM effect calculation");
91 lpmCmd->SetParameterName("lpm",true);
92 lpmCmd->SetDefaultValue(true);
93 lpmCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
94 lpmCmd->SetToBeBroadcasted(false);
95
96 rsCmd = new G4UIcmdWithABool("/process/eLoss/useCutAsFinalRange",this);
97 rsCmd->SetGuidance("Enable/disable use of cut in range as a final range");
98 rsCmd->SetParameterName("choice",true);
99 rsCmd->SetDefaultValue(false);
100 rsCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
101 rsCmd->SetToBeBroadcasted(false);
102
103 aplCmd = new G4UIcmdWithABool("/process/em/applyCuts",this);
104 aplCmd->SetGuidance("Enable/disable applying cuts for gamma processes");
105 aplCmd->SetParameterName("apl",true);
106 aplCmd->SetDefaultValue(false);
107 aplCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
108 aplCmd->SetToBeBroadcasted(false);
109
110 intCmd = new G4UIcmdWithABool("/process/em/integral",this);
111 intCmd->SetGuidance("Enable/disable integral method.");
112 intCmd->SetParameterName("choice",true);
113 intCmd->SetDefaultValue(true);
114 intCmd->AvailableForStates(G4State_PreInit);
115 intCmd->SetToBeBroadcasted(false);
116
117 latCmd = new G4UIcmdWithABool("/process/msc/LateralDisplacement",this);
118 latCmd->SetGuidance("Enable/disable sampling of lateral displacement");
119 latCmd->SetParameterName("lat",true);
120 latCmd->SetDefaultValue(true);
121 latCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
122 latCmd->SetToBeBroadcasted(false);
123
124 lat96Cmd = new G4UIcmdWithABool("/process/msc/LateralDisplacementAlg96",this);
125 lat96Cmd->SetGuidance("Enable/disable sampling of lateral displacement");
126 lat96Cmd->SetParameterName("lat96",true);
127 lat96Cmd->SetDefaultValue(false);
128 lat96Cmd->AvailableForStates(G4State_PreInit,G4State_Idle);
129 lat96Cmd->SetToBeBroadcasted(false);
130
131 mulatCmd = new G4UIcmdWithABool("/process/msc/MuHadLateralDisplacement",this);
132 mulatCmd->SetGuidance("Enable/disable sampling of lateral displacement for muons and hadrons");
133 mulatCmd->SetParameterName("mulat",true);
134 mulatCmd->SetDefaultValue(true);
135 mulatCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
136 mulatCmd->SetToBeBroadcasted(false);
137
138 delCmd = new G4UIcmdWithABool("/process/eLoss/UseAngularGenerator",this);
139 delCmd->SetGuidance("Enable usage of angular generator for ionisation");
140 delCmd->SetParameterName("del",true);
141 delCmd->SetDefaultValue(false);
142 delCmd->AvailableForStates(G4State_PreInit);
143 delCmd->SetToBeBroadcasted(false);
144
145 mottCmd = new G4UIcmdWithABool("/process/msc/UseMottCorrection",this);
146 mottCmd->SetGuidance("Enable usage of Mott corrections for e- elastic scattering");
147 mottCmd->SetParameterName("mott",true);
148 mottCmd->SetDefaultValue(false);
149 mottCmd->AvailableForStates(G4State_PreInit);
150 mottCmd->SetToBeBroadcasted(false);
151
152 birksCmd = new G4UIcmdWithABool("/process/em/UseG4EmSaturation",this);
153 birksCmd->SetGuidance("Enable usage of built-in Birks saturation");
154 birksCmd->SetParameterName("birks",true);
155 birksCmd->SetDefaultValue(false);
156 birksCmd->AvailableForStates(G4State_PreInit,G4State_Init);
157 birksCmd->SetToBeBroadcasted(false);
158
159 sharkCmd = new G4UIcmdWithABool("/process/em/UseGeneralProcess",this);
160 sharkCmd->SetGuidance("Enable gamma, e+- general process");
161 sharkCmd->SetParameterName("gen",true);
162 sharkCmd->SetDefaultValue(false);
163 sharkCmd->AvailableForStates(G4State_PreInit);
164 sharkCmd->SetToBeBroadcasted(false);
165
166 poCmd = new G4UIcmdWithABool("/process/em/Polarisation",this);
167 poCmd->SetGuidance("Enable polarisation");
168 poCmd->AvailableForStates(G4State_PreInit);
169 poCmd->SetToBeBroadcasted(false);
170
171 sampleTCmd = new G4UIcmdWithABool("/process/em/enableSamplingTable",this);
172 sampleTCmd->SetGuidance("Enable usage of sampling table for secondary generation");
173 sampleTCmd->SetParameterName("sampleT",true);
174 sampleTCmd->SetDefaultValue(false);
175 sampleTCmd->AvailableForStates(G4State_PreInit);
176 sampleTCmd->SetToBeBroadcasted(false);
177
178 icru90Cmd = new G4UIcmdWithABool("/process/eLoss/UseICRU90",this);
179 icru90Cmd->SetGuidance("Enable usage of ICRU90 stopping powers");
180 icru90Cmd->SetParameterName("icru90",true);
181 icru90Cmd->SetDefaultValue(false);
182 icru90Cmd->AvailableForStates(G4State_PreInit);
183 icru90Cmd->SetToBeBroadcasted(false);
184
185 mudatCmd = new G4UIcmdWithABool("/process/em/MuDataFromFile",this);
186 mudatCmd->SetGuidance("Enable usage of muon data from file");
187 mudatCmd->SetParameterName("mudat",true);
188 mudatCmd->SetDefaultValue(false);
189 mudatCmd->AvailableForStates(G4State_PreInit);
190 mudatCmd->SetToBeBroadcasted(false);
191
192 peKCmd = new G4UIcmdWithABool("/process/em/PhotoeffectBelowKShell",this);
193 peKCmd->SetGuidance("Enable sampling of photoeffect below K-shell");
194 peKCmd->SetParameterName("peK",true);
195 peKCmd->SetDefaultValue(true);
196 peKCmd->AvailableForStates(G4State_PreInit);
197 peKCmd->SetToBeBroadcasted(false);
198
199 mscPCmd = new G4UIcmdWithABool("/process/msc/PositronCorrection",this);
200 mscPCmd->SetGuidance("Enable msc positron correction");
201 mscPCmd->SetParameterName("mscPC",true);
202 mscPCmd->SetDefaultValue(true);
203 mscPCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
204 mscPCmd->SetToBeBroadcasted(false);
205
206 pepicsCmd = new G4UIcmdWithABool("/process/em/UseEPICS2017XS",this);
207 pepicsCmd->SetGuidance("Use EPICS2017 data for gamma x-ections");
208 pepicsCmd->SetParameterName("pepics",true);
209 pepicsCmd->SetDefaultValue(false);
210 pepicsCmd->AvailableForStates(G4State_PreInit);
211 pepicsCmd->SetToBeBroadcasted(false);
212
213 f3gCmd = new G4UIcmdWithABool("/process/em/3GammaAnnihilationOnFly",this);
214 f3gCmd->SetGuidance("Enable/disable 3 gamma annihilation on fly");
215 f3gCmd->SetParameterName("f3gamma",true);
216 f3gCmd->SetDefaultValue(false);
217 f3gCmd->AvailableForStates(G4State_PreInit);
218 f3gCmd->SetToBeBroadcasted(false);
219
220 fRiGeCmd = new G4UIcmdWithABool("/process/em/PairProd5D",this);
221 fRiGeCmd->SetGuidance("Enable/disable 5D model for e+e- pair production by muons");
222 fRiGeCmd->SetParameterName("ee5D",true);
223 fRiGeCmd->SetDefaultValue(false);
224 fRiGeCmd->AvailableForStates(G4State_PreInit);
225 fRiGeCmd->SetToBeBroadcasted(false);
226
227 minEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/minKinEnergy",this);
228 minEnCmd->SetGuidance("Set the min kinetic energy for EM tables");
229 minEnCmd->SetParameterName("emin",true);
230 minEnCmd->SetUnitCategory("Energy");
231 minEnCmd->AvailableForStates(G4State_PreInit);
232 minEnCmd->SetToBeBroadcasted(false);
233
234 maxEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergy",this);
235 maxEnCmd->SetGuidance("Set the max kinetic energy for EM tables");
236 maxEnCmd->SetParameterName("emax",true);
237 maxEnCmd->SetUnitCategory("Energy");
238 maxEnCmd->AvailableForStates(G4State_PreInit);
239 maxEnCmd->SetToBeBroadcasted(false);
240
241 cenCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergyCSDA",this);
242 cenCmd->SetGuidance("Set the max kinetic energy for CSDA table");
243 cenCmd->SetParameterName("emaxCSDA",true);
244 cenCmd->SetUnitCategory("Energy");
245 cenCmd->AvailableForStates(G4State_PreInit);
246 cenCmd->SetToBeBroadcasted(false);
247
248 max5DCmd = new G4UIcmdWithADoubleAndUnit("/process/em/max5DMuPairEnergy",this);
249 max5DCmd->SetGuidance("Set the max kinetic energy for 5D muon pair production");
250 max5DCmd->SetParameterName("emax5D",true);
251 max5DCmd->SetUnitCategory("Energy");
252 max5DCmd->AvailableForStates(G4State_PreInit);
253 max5DCmd->SetToBeBroadcasted(false);
254
255 lowEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestElectronEnergy",this);
256 lowEnCmd->SetGuidance("Set the lowest kinetic energy for e+-");
257 lowEnCmd->SetParameterName("elow",true);
258 lowEnCmd->SetUnitCategory("Energy");
259 lowEnCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
260 lowEnCmd->SetToBeBroadcasted(false);
261
262 lowhEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestMuHadEnergy",this);
263 lowhEnCmd->SetGuidance("Set the lowest kinetic energy for muons and hadrons");
264 lowhEnCmd->SetParameterName("elowh",true);
265 lowhEnCmd->SetUnitCategory("Energy");
266 lowhEnCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
267 lowhEnCmd->SetToBeBroadcasted(false);
268
269 lowEn3Cmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestTripletEnergy",this);
270 lowEn3Cmd->SetGuidance("Set the lowest kinetic energy for triplet production");
271 lowEn3Cmd->SetParameterName("elow3",true);
272 lowEn3Cmd->SetUnitCategory("Energy");
273 lowEn3Cmd->AvailableForStates(G4State_PreInit,G4State_Idle);
274 lowEn3Cmd->SetToBeBroadcasted(false);
275
276 lllCmd = new G4UIcmdWithADouble("/process/eLoss/linLossLimit",this);
277 lllCmd->SetGuidance("Set linearLossLimit parameter");
278 lllCmd->SetParameterName("linlim",true);
279 lllCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
280 lllCmd->SetToBeBroadcasted(false);
281
282 brCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/bremThreshold",this);
283 brCmd->SetGuidance("Set e+- bremsstrahlung energy threshold");
284 brCmd->SetParameterName("emaxBrem",true);
285 brCmd->SetUnitCategory("Energy");
286 brCmd->AvailableForStates(G4State_PreInit);
287 brCmd->SetToBeBroadcasted(false);
288
289 br1Cmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/bremMuHadThreshold",this);
290 br1Cmd->SetGuidance("Set muon/hadron bremsstrahlung energy threshold");
291 br1Cmd->SetParameterName("emaxMuHadBrem",true);
292 br1Cmd->SetUnitCategory("Energy");
293 br1Cmd->AvailableForStates(G4State_PreInit);
294 br1Cmd->SetToBeBroadcasted(false);
295
296 labCmd = new G4UIcmdWithADouble("/process/eLoss/LambdaFactor",this);
297 labCmd->SetGuidance("Set lambdaFactor parameter for integral option");
298 labCmd->SetParameterName("Fl",true);
299 labCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
300 labCmd->SetToBeBroadcasted(false);
301
302 mscfCmd = new G4UIcmdWithADouble("/process/msc/FactorForAngleLimit",this);
303 mscfCmd->SetGuidance("Set factor for computation of a limit for -t (invariant transfer)");
304 mscfCmd->SetParameterName("Fact",true);
305 mscfCmd->SetRange("Fact>0");
306 mscfCmd->SetDefaultValue(1.);
307 mscfCmd->AvailableForStates(G4State_PreInit);
308 mscfCmd->SetToBeBroadcasted(false);
309
310 angCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/ThetaLimit",this);
311 angCmd->SetGuidance("Set the limit on the polar angle for msc and single scattering");
312 angCmd->SetParameterName("theta",true);
313 angCmd->SetUnitCategory("Angle");
314 angCmd->AvailableForStates(G4State_PreInit);
315 angCmd->SetToBeBroadcasted(false);
316
317 msceCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/EnergyLimit",this);
318 msceCmd->SetGuidance("Set the upper energy limit for msc");
319 msceCmd->SetParameterName("mscE",true);
320 msceCmd->SetUnitCategory("Energy");
321 msceCmd->AvailableForStates(G4State_PreInit);
322 msceCmd->SetToBeBroadcasted(false);
323
324 nielCmd = new G4UIcmdWithADoubleAndUnit("/process/em/MaxEnergyNIEL",this);
325 nielCmd->SetGuidance("Set the upper energy limit for NIEL");
326 nielCmd->SetParameterName("niel",true);
327 nielCmd->SetUnitCategory("Energy");
328 nielCmd->AvailableForStates(G4State_PreInit);
329 nielCmd->SetToBeBroadcasted(false);
330
331 frCmd = new G4UIcmdWithADouble("/process/msc/RangeFactor",this);
332 frCmd->SetGuidance("Set RangeFactor for msc processes of e+-");
333 frCmd->SetParameterName("Fr",true);
334 frCmd->SetRange("Fr>0");
335 frCmd->SetDefaultValue(0.04);
336 frCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
337 frCmd->SetToBeBroadcasted(false);
338
339 fr1Cmd = new G4UIcmdWithADouble("/process/msc/RangeFactorMuHad",this);
340 fr1Cmd->SetGuidance("Set RangeFactor for msc processes of muons/hadrons");
341 fr1Cmd->SetParameterName("Fr1",true);
342 fr1Cmd->SetRange("Fr1>0");
343 fr1Cmd->SetDefaultValue(0.2);
344 fr1Cmd->AvailableForStates(G4State_PreInit,G4State_Idle);
345 fr1Cmd->SetToBeBroadcasted(false);
346
347 fgCmd = new G4UIcmdWithADouble("/process/msc/GeomFactor",this);
348 fgCmd->SetGuidance("Set GeomFactor parameter for msc processes");
349 fgCmd->SetParameterName("Fg",true);
350 fgCmd->SetRange("Fg>0");
351 fgCmd->SetDefaultValue(2.5);
352 fgCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
353 fgCmd->SetToBeBroadcasted(false);
354
355 skinCmd = new G4UIcmdWithADouble("/process/msc/Skin",this);
356 skinCmd->SetGuidance("Set skin parameter for msc processes");
357 skinCmd->SetParameterName("skin",true);
358 skinCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
359 skinCmd->SetToBeBroadcasted(false);
360
361 screCmd = new G4UIcmdWithADouble("/process/msc/ScreeningFactor",this);
362 screCmd->SetGuidance("Set screening factor");
363 screCmd->SetParameterName("screen",true);
364 screCmd->AvailableForStates(G4State_PreInit);
365 screCmd->SetToBeBroadcasted(false);
366
367 safCmd = new G4UIcmdWithADouble("/process/msc/SafetyFactor",this);
368 safCmd->SetGuidance("Set safety factor");
369 safCmd->SetParameterName("fsafe",true);
370 safCmd->AvailableForStates(G4State_PreInit);
371 safCmd->SetToBeBroadcasted(false);
372
373 llimCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/LambdaLimit",this);
374 llimCmd->SetGuidance("Set the upper energy limit for NIEL");
375 llimCmd->SetParameterName("ll",true);
376 llimCmd->SetUnitCategory("Length");
377 llimCmd->AvailableForStates(G4State_PreInit);
378 llimCmd->SetToBeBroadcasted(false);
379
380 amCmd = new G4UIcmdWithAnInteger("/process/em/binsPerDecade",this);
381 amCmd->SetGuidance("Set number of bins per decade for EM tables");
382 amCmd->SetParameterName("bins",true);
383 amCmd->SetDefaultValue(7);
384 amCmd->AvailableForStates(G4State_PreInit);
385 amCmd->SetToBeBroadcasted(false);
386
387 verCmd = new G4UIcmdWithAnInteger("/process/eLoss/verbose",this);
388 verCmd->SetGuidance("Set verbose level for EM physics");
389 verCmd->SetParameterName("verb",true);
390 verCmd->SetDefaultValue(1);
391 verCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
392 verCmd->SetToBeBroadcasted(false);
393
394 ver1Cmd = new G4UIcmdWithAnInteger("/process/em/verbose",this);
395 ver1Cmd->SetGuidance("Set verbose level for EM physics");
396 ver1Cmd->SetParameterName("verb1",true);
397 ver1Cmd->SetDefaultValue(1);
398 ver1Cmd->AvailableForStates(G4State_PreInit,G4State_Idle);
399 ver1Cmd->SetToBeBroadcasted(false);
400
401 ver2Cmd = new G4UIcmdWithAnInteger("/process/em/workerVerbose",this);
402 ver2Cmd->SetGuidance("Set worker verbose level for EM physics");
403 ver2Cmd->SetParameterName("verb2",true);
404 ver2Cmd->SetDefaultValue(0);
405 ver2Cmd->AvailableForStates(G4State_PreInit,G4State_Idle);
406 ver2Cmd->SetToBeBroadcasted(false);
407
408 nFreeCmd = new G4UIcmdWithAnInteger("/process/em/nForFreeVector",this);
409 nFreeCmd->SetGuidance("Set number for logarithmic bin search algorithm");
410 nFreeCmd->SetParameterName("nFree",true);
411 nFreeCmd->SetDefaultValue(2);
412 nFreeCmd->AvailableForStates(G4State_PreInit);
413 nFreeCmd->SetToBeBroadcasted(false);
414
415 transWithMscCmd = new G4UIcmdWithAString("/process/em/transportationWithMsc",this);
416 transWithMscCmd->SetGuidance("Enable/disable the G4TransportationWithMsc process");
417 transWithMscCmd->SetParameterName("trans",true);
418 transWithMscCmd->SetCandidates("Disabled Enabled MultipleSteps");
419 transWithMscCmd->AvailableForStates(G4State_PreInit);
420 transWithMscCmd->SetToBeBroadcasted(false);
421
422 mscCmd = new G4UIcmdWithAString("/process/msc/StepLimit",this);
423 mscCmd->SetGuidance("Set msc step limitation type");
424 mscCmd->SetParameterName("StepLim",true);
425 mscCmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
426 mscCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
427 mscCmd->SetToBeBroadcasted(false);
428
429 msc1Cmd = new G4UIcmdWithAString("/process/msc/StepLimitMuHad",this);
430 msc1Cmd->SetGuidance("Set msc step limitation type for muons/hadrons");
431 msc1Cmd->SetParameterName("StepLim1",true);
432 msc1Cmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
433 msc1Cmd->AvailableForStates(G4State_PreInit,G4State_Idle);
434 msc1Cmd->SetToBeBroadcasted(false);
435
436 dumpCmd = new G4UIcommand("/process/em/printParameters",this);
437 dumpCmd->SetGuidance("Print all EM parameters.");
438 dumpCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
439 dumpCmd->SetToBeBroadcasted(false);
440
441 nffCmd = new G4UIcmdWithAString("/process/em/setNuclearFormFactor",this);
442 nffCmd->SetGuidance("Define type of nuclear form-factor");
443 nffCmd->SetParameterName("NucFF",true);
444 nffCmd->SetCandidates("None Exponential Gaussian Flat");
445 nffCmd->AvailableForStates(G4State_PreInit);
446 nffCmd->SetToBeBroadcasted(false);
447
448 ssCmd = new G4UIcmdWithAString("/process/em/setSingleScattering",this);
449 ssCmd->SetGuidance("Define type of e+- single scattering model");
450 ssCmd->SetParameterName("SS",true);
451 ssCmd->SetCandidates("WVI Mott DPWA");
452 ssCmd->AvailableForStates(G4State_PreInit);
453 ssCmd->SetToBeBroadcasted(false);
454
455 fluc1Cmd = new G4UIcmdWithAString("/process/eLoss/setFluctModel",this);
456 fluc1Cmd->SetGuidance("Define type of energy loss fluctuation model");
457 fluc1Cmd->SetParameterName("Fluc1",true);
458 fluc1Cmd->SetCandidates("Dummy Universal Urban");
459 fluc1Cmd->AvailableForStates(G4State_PreInit);
460 fluc1Cmd->SetToBeBroadcasted(false);
461
462 fluc2Cmd = new G4UIcmdWithAString("/process/eLoss/enableFluctForRegion",this);
463 fluc2Cmd->SetGuidance("Enable dEdx fluctuations for G4Region");
464 fluc2Cmd->SetParameterName("Fluc2",true);
465 fluc2Cmd->AvailableForStates(G4State_PreInit);
466 fluc2Cmd->SetToBeBroadcasted(false);
467
468 fluc3Cmd = new G4UIcmdWithAString("/process/eLoss/disableFluctForRegion",this);
469 fluc3Cmd->SetGuidance("Disable dEdx fluctuations for G4Region");
470 fluc3Cmd->SetParameterName("Fluc3",true);
471 fluc3Cmd->AvailableForStates(G4State_PreInit);
472 fluc3Cmd->SetToBeBroadcasted(false);
473
474 posiCmd = new G4UIcmdWithAString("/process/em/setPositronAtRestModel",this);
475 posiCmd->SetGuidance("Define model of positron annihilation at rest");
476 posiCmd->SetParameterName("Posi",true);
477 posiCmd->SetCandidates("Simple Allison OrePawell OrePowellPolar");
478 posiCmd->AvailableForStates(G4State_PreInit);
479 posiCmd->SetToBeBroadcasted(false);
480
481 tripletCmd = new G4UIcmdWithAnInteger("/process/gconv/conversionType",this);
482 tripletCmd->SetGuidance("gamma conversion triplet/nuclear generation type:");
483 tripletCmd->SetGuidance("0 - (default) both triplet and nuclear");
484 tripletCmd->SetGuidance("1 - force nuclear");
485 tripletCmd->SetGuidance("2 - force triplet");
486 tripletCmd->SetParameterName("type",false);
487 tripletCmd->SetRange("type >= 0 && type <= 2");
488 tripletCmd->SetDefaultValue(0);
489 tripletCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
490 tripletCmd->SetToBeBroadcasted(false);
491
492 onIsolatedCmd = new G4UIcmdWithABool("/process/gconv/onIsolated",this);
493 onIsolatedCmd->SetGuidance("Conversion on isolated charged particles");
494 onIsolatedCmd->SetGuidance("false (default) : atomic electron screening");
495 onIsolatedCmd->SetGuidance("true : conversion on isolated particles.");
496 onIsolatedCmd->SetParameterName("flag",false);
497 onIsolatedCmd->SetDefaultValue(false);
498 onIsolatedCmd->AvailableForStates(G4State_PreInit,G4State_Idle);
499 onIsolatedCmd->SetToBeBroadcasted(false);
500}
501
502//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
503
505{
506 delete gconvDirectory;
507 delete eLossDirectory;
508 delete mscDirectory;
509 delete emDirectory;
510 delete dnaDirectory;
511
512 delete flucCmd;
513 delete rangeCmd;
514 delete lpmCmd;
515 delete rsCmd;
516 delete aplCmd;
517 delete intCmd;
518 delete latCmd;
519 delete lat96Cmd;
520 delete mulatCmd;
521 delete delCmd;
522 delete mottCmd;
523 delete birksCmd;
524 delete sharkCmd;
525 delete onIsolatedCmd;
526 delete sampleTCmd;
527 delete poCmd;
528 delete icru90Cmd;
529 delete mudatCmd;
530 delete peKCmd;
531 delete f3gCmd;
532 delete fRiGeCmd;
533 delete mscPCmd;
534 delete pepicsCmd;
535
536 delete minEnCmd;
537 delete maxEnCmd;
538 delete max5DCmd;
539 delete cenCmd;
540 delete lowEnCmd;
541 delete lowhEnCmd;
542 delete lowEn3Cmd;
543 delete lllCmd;
544 delete brCmd;
545 delete br1Cmd;
546 delete labCmd;
547 delete mscfCmd;
548 delete angCmd;
549 delete msceCmd;
550 delete nielCmd;
551 delete frCmd;
552 delete fr1Cmd;
553 delete fgCmd;
554 delete skinCmd;
555 delete safCmd;
556 delete llimCmd;
557 delete screCmd;
558
559 delete amCmd;
560 delete verCmd;
561 delete ver1Cmd;
562 delete ver2Cmd;
563 delete transWithMscCmd;
564 delete nFreeCmd;
565 delete tripletCmd;
566
567 delete mscCmd;
568 delete msc1Cmd;
569 delete nffCmd;
570 delete ssCmd;
571 delete fluc1Cmd;
572 delete fluc2Cmd;
573 delete fluc3Cmd;
574 delete posiCmd;
575
576 delete dumpCmd;
577}
578
579//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
580
582 G4String newValue)
583{
584 G4bool physicsModified = false;
585 if (command == flucCmd) {
586 theParameters->SetLossFluctuations(flucCmd->GetNewBoolValue(newValue));
587 physicsModified = true;
588 } else if (command == rangeCmd) {
589 theParameters->SetBuildCSDARange(rangeCmd->GetNewBoolValue(newValue));
590 } else if (command == lpmCmd) {
591 theParameters->SetLPM(lpmCmd->GetNewBoolValue(newValue));
592 physicsModified = true;
593 } else if (command == rsCmd) {
594 theParameters->SetUseCutAsFinalRange(rsCmd->GetNewBoolValue(newValue));
595 physicsModified = true;
596 } else if (command == aplCmd) {
597 theParameters->SetApplyCuts(aplCmd->GetNewBoolValue(newValue));
598 physicsModified = true;
599 } else if (command == intCmd) {
600 theParameters->SetIntegral(intCmd->GetNewBoolValue(newValue));
601 } else if (command == latCmd) {
602 theParameters->SetLateralDisplacement(latCmd->GetNewBoolValue(newValue));
603 physicsModified = true;
604 } else if (command == lat96Cmd) {
605 theParameters->SetLateralDisplacementAlg96(lat96Cmd->GetNewBoolValue(newValue));
606 physicsModified = true;
607 } else if (command == mulatCmd) {
608 theParameters->SetMuHadLateralDisplacement(mulatCmd->GetNewBoolValue(newValue));
609 physicsModified = true;
610 } else if (command == delCmd) {
611 theParameters->ActivateAngularGeneratorForIonisation(delCmd->GetNewBoolValue(newValue));
612 } else if (command == mottCmd) {
613 theParameters->SetUseMottCorrection(mottCmd->GetNewBoolValue(newValue));
614 } else if (command == birksCmd) {
615 theParameters->SetBirksActive(birksCmd->GetNewBoolValue(newValue));
616 } else if (command == icru90Cmd) {
617 theParameters->SetUseICRU90Data(icru90Cmd->GetNewBoolValue(newValue));
618 } else if (command == sharkCmd) {
619 theParameters->SetGeneralProcessActive(sharkCmd->GetNewBoolValue(newValue));
620 } else if (command == poCmd) {
621 theParameters->SetEnablePolarisation(poCmd->GetNewBoolValue(newValue));
622 } else if (command == sampleTCmd) {
623 theParameters->SetEnableSamplingTable(sampleTCmd->GetNewBoolValue(newValue));
624 } else if (command == mudatCmd) {
625 theParameters->SetRetrieveMuDataFromFile(mudatCmd->GetNewBoolValue(newValue));
626 } else if (command == peKCmd) {
627 theParameters->SetPhotoeffectBelowKShell(peKCmd->GetNewBoolValue(newValue));
628 } else if (command == f3gCmd) {
629 theParameters->Set3GammaAnnihilationOnFly(f3gCmd->GetNewBoolValue(newValue));
630 } else if (command == fRiGeCmd) {
631 theParameters->SetUseRiGePairProductionModel(fRiGeCmd->GetNewBoolValue(newValue));
632 } else if (command == mscPCmd) {
633 theParameters->SetMscPositronCorrection(mscPCmd->GetNewBoolValue(newValue));
634 } else if (command == pepicsCmd) {
635 theParameters->SetUseEPICS2017XS(pepicsCmd->GetNewBoolValue(newValue));
636
637 } else if (command == minEnCmd) {
638 theParameters->SetMinEnergy(minEnCmd->GetNewDoubleValue(newValue));
639 } else if (command == maxEnCmd) {
640 theParameters->SetMaxEnergy(maxEnCmd->GetNewDoubleValue(newValue));
641 } else if (command == max5DCmd) {
642 theParameters->SetMaxEnergyFor5DMuPair(max5DCmd->GetNewDoubleValue(newValue));
643 } else if (command == cenCmd) {
644 theParameters->SetMaxEnergyForCSDARange(cenCmd->GetNewDoubleValue(newValue));
645 physicsModified = true;
646 } else if (command == lowEnCmd) {
647 theParameters->SetLowestElectronEnergy(lowEnCmd->GetNewDoubleValue(newValue));
648 physicsModified = true;
649 } else if (command == lowEn3Cmd) {
650 theParameters->SetLowestTripletEnergy(lowEn3Cmd->GetNewDoubleValue(newValue));
651 physicsModified = true;
652 } else if (command == lowhEnCmd) {
653 theParameters->SetLowestMuHadEnergy(lowhEnCmd->GetNewDoubleValue(newValue));
654 physicsModified = true;
655 } else if (command == lllCmd) {
656 theParameters->SetLinearLossLimit(lllCmd->GetNewDoubleValue(newValue));
657 physicsModified = true;
658 } else if (command == brCmd) {
659 theParameters->SetBremsstrahlungTh(brCmd->GetNewDoubleValue(newValue));
660 physicsModified = true;
661 } else if (command == br1Cmd) {
662 theParameters->SetMuHadBremsstrahlungTh(br1Cmd->GetNewDoubleValue(newValue));
663 physicsModified = true;
664 } else if (command == labCmd) {
665 theParameters->SetLambdaFactor(labCmd->GetNewDoubleValue(newValue));
666 physicsModified = true;
667 } else if (command == mscfCmd) {
668 theParameters->SetFactorForAngleLimit(mscfCmd->GetNewDoubleValue(newValue));
669 } else if (command == angCmd) {
670 theParameters->SetMscThetaLimit(angCmd->GetNewDoubleValue(newValue));
671 } else if (command == msceCmd) {
672 theParameters->SetMscEnergyLimit(msceCmd->GetNewDoubleValue(newValue));
673 } else if (command == nielCmd) {
674 theParameters->SetMaxNIELEnergy(nielCmd->GetNewDoubleValue(newValue));
675 } else if (command == frCmd) {
676 theParameters->SetMscRangeFactor(frCmd->GetNewDoubleValue(newValue));
677 physicsModified = true;
678 } else if (command == fr1Cmd) {
679 theParameters->SetMscMuHadRangeFactor(fr1Cmd->GetNewDoubleValue(newValue));
680 physicsModified = true;
681 } else if (command == fgCmd) {
682 theParameters->SetMscGeomFactor(fgCmd->GetNewDoubleValue(newValue));
683 physicsModified = true;
684 } else if (command == skinCmd) {
685 theParameters->SetMscSkin(skinCmd->GetNewDoubleValue(newValue));
686 physicsModified = true;
687 } else if (command == safCmd) {
688 theParameters->SetMscSafetyFactor(safCmd->GetNewDoubleValue(newValue));
689 } else if (command == llimCmd) {
690 theParameters->SetMscLambdaLimit(llimCmd->GetNewDoubleValue(newValue));
691 } else if (command == screCmd) {
692 theParameters->SetScreeningFactor(screCmd->GetNewDoubleValue(newValue));
693 } else if (command == amCmd) {
694 theParameters->SetNumberOfBinsPerDecade(amCmd->GetNewIntValue(newValue));
695 } else if (command == verCmd) {
696 theParameters->SetVerbose(verCmd->GetNewIntValue(newValue));
697 } else if (command == ver1Cmd) {
698 theParameters->SetVerbose(ver1Cmd->GetNewIntValue(newValue));
699 } else if (command == ver2Cmd) {
700 theParameters->SetWorkerVerbose(ver2Cmd->GetNewIntValue(newValue));
701 } else if (command == nFreeCmd) {
702 theParameters->SetNumberForFreeVector(nFreeCmd->GetNewIntValue(newValue));
703 } else if (command == dumpCmd) {
704 theParameters->SetIsPrintedFlag(false);
705 theParameters->Dump();
706 } else if (command == transWithMscCmd) {
708 if(newValue == "Disabled") {
710 } else if(newValue == "Enabled") {
712 } else if(newValue == "MultipleSteps") {
714 } else {
716 ed << " TransportationWithMsc type <" << newValue << "> unknown!";
717 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
718 }
719 theParameters->SetTransportationWithMsc(type);
720 } else if (command == mscCmd || command == msc1Cmd) {
722 if(newValue == "Minimal") {
723 msctype = fMinimal;
724 } else if(newValue == "UseDistanceToBoundary") {
725 msctype = fUseDistanceToBoundary;
726 } else if(newValue == "UseSafety") {
727 msctype = fUseSafety;
728 } else if(newValue == "UseSafetyPlus") {
729 msctype = fUseSafetyPlus;
730 } else {
732 ed << " StepLimit type <" << newValue << "> unknown!";
733 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
734 return;
735 }
736 if (command == mscCmd) {
737 theParameters->SetMscStepLimitType(msctype);
738 } else {
739 theParameters->SetMscMuHadStepLimitType(msctype);
740 }
741 physicsModified = true;
742 } else if (command == nffCmd) {
744 if(newValue == "Exponential") { x = fExponentialNF; }
745 else if(newValue == "Gaussian") { x = fGaussianNF; }
746 else if(newValue == "Flat") { x = fFlatNF; }
747 else if(newValue != "None") {
749 ed << " NuclearFormFactor type <" << newValue << "> unknown!";
750 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
751 return;
752 }
753 theParameters->SetNuclearFormfactorType(x);
754 } else if (command == ssCmd) {
756 if(newValue == "DPWA") { x = fDPWA; }
757 else if(newValue == "Mott") { x = fMott; }
758 else if(newValue != "WVI") {
760 ed << " G4eSingleScatteringType type <" << newValue << "> unknown!";
761 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
762 return;
763 }
764 theParameters->SetSingleScatteringType(x);
765 } else if (command == fluc1Cmd) {
767 if(newValue == "Dummy") { x = fDummyFluctuation; }
768 else if(newValue == "Urban") { x = fUrbanFluctuation; }
769 theParameters->SetFluctuationType(x);
770 } else if (command == fluc2Cmd) {
771 theParameters->SetFluctuationsForRegion(newValue, true);
772 } else if (command == fluc3Cmd) {
773 theParameters->SetFluctuationsForRegion(newValue, false);
774 } else if (command == posiCmd) {
776 if (newValue == "Allison") { x = fAllisonPositronium; }
777 else if (newValue == "OrePowell") { x = fOrePowell; }
778 else if (newValue == "OrePowellPolar") { x = fOrePowellPolar; }
779 theParameters->SetPositronAtRestModelType(x);
780 } else if ( command==tripletCmd ) {
781 theParameters->SetConversionType(tripletCmd->GetNewIntValue(newValue));
782 } else if ( command==onIsolatedCmd ) {
783 theParameters->SetOnIsolated(onIsolatedCmd->GetNewBoolValue(newValue));
784 physicsModified = true;
785 }
786
787 if(physicsModified) {
788 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
789 }
790}
791
792//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ G4State_Init
@ G4State_Idle
@ G4State_PreInit
G4EmFluctuationType
@ fUniversalFluctuation
@ fUrbanFluctuation
@ fDummyFluctuation
G4eSingleScatteringType
@ fWVI
@ fMott
@ fDPWA
G4TransportationWithMscType
G4PositronAtRestModelType
@ fSimplePositronium
@ fAllisonPositronium
@ fOrePowell
@ fOrePowellPolar
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4MscStepLimitType
@ fUseSafety
@ fUseSafetyPlus
@ fUseDistanceToBoundary
bool G4bool
Definition G4Types.hh:86
void SetNewValue(G4UIcommand *, G4String) override
G4EmParametersMessenger(G4EmParameters *)
G4int ApplyCommand(const char *aCommand)
static G4UImanager * GetUIpointer()