Geant4 11.1.1
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);
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);
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);
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);
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);
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);
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);
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);
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);
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);
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);
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);
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);
164 sharkCmd->SetToBeBroadcasted(false);
165
166 poCmd = new G4UIcmdWithABool("/process/em/Polarisation",this);
167 poCmd->SetGuidance("Enable polarisation");
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);
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);
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);
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);
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);
204 mscPCmd->SetToBeBroadcasted(false);
205
206 minEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/minKinEnergy",this);
207 minEnCmd->SetGuidance("Set the min kinetic energy for EM tables");
208 minEnCmd->SetParameterName("emin",true);
209 minEnCmd->SetUnitCategory("Energy");
211 minEnCmd->SetToBeBroadcasted(false);
212
213 maxEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergy",this);
214 maxEnCmd->SetGuidance("Set the max kinetic energy for EM tables");
215 maxEnCmd->SetParameterName("emax",true);
216 maxEnCmd->SetUnitCategory("Energy");
218 maxEnCmd->SetToBeBroadcasted(false);
219
220 cenCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergyCSDA",this);
221 cenCmd->SetGuidance("Set the max kinetic energy for CSDA table");
222 cenCmd->SetParameterName("emaxCSDA",true);
223 cenCmd->SetUnitCategory("Energy");
225 cenCmd->SetToBeBroadcasted(false);
226
227 max5DCmd = new G4UIcmdWithADoubleAndUnit("/process/em/max5DMuPairEnergy",this);
228 max5DCmd->SetGuidance("Set the max kinetic energy for 5D muon pair production");
229 max5DCmd->SetParameterName("emax5D",true);
230 max5DCmd->SetUnitCategory("Energy");
232 max5DCmd->SetToBeBroadcasted(false);
233
234 lowEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestElectronEnergy",this);
235 lowEnCmd->SetGuidance("Set the lowest kinetic energy for e+-");
236 lowEnCmd->SetParameterName("elow",true);
237 lowEnCmd->SetUnitCategory("Energy");
239 lowEnCmd->SetToBeBroadcasted(false);
240
241 lowhEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestMuHadEnergy",this);
242 lowhEnCmd->SetGuidance("Set the lowest kinetic energy for muons and hadrons");
243 lowhEnCmd->SetParameterName("elowh",true);
244 lowhEnCmd->SetUnitCategory("Energy");
246 lowhEnCmd->SetToBeBroadcasted(false);
247
248 lowEn3Cmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestTripletEnergy",this);
249 lowEn3Cmd->SetGuidance("Set the lowest kinetic energy for triplet production");
250 lowEn3Cmd->SetParameterName("elow3",true);
251 lowEn3Cmd->SetUnitCategory("Energy");
253 lowEn3Cmd->SetToBeBroadcasted(false);
254
255 lllCmd = new G4UIcmdWithADouble("/process/eLoss/linLossLimit",this);
256 lllCmd->SetGuidance("Set linearLossLimit parameter");
257 lllCmd->SetParameterName("linlim",true);
259 lllCmd->SetToBeBroadcasted(false);
260
261 brCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/bremThreshold",this);
262 brCmd->SetGuidance("Set e+- bremsstrahlung energy threshold");
263 brCmd->SetParameterName("emaxBrem",true);
264 brCmd->SetUnitCategory("Energy");
266 brCmd->SetToBeBroadcasted(false);
267
268 br1Cmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/bremMuHadThreshold",this);
269 br1Cmd->SetGuidance("Set muon/hadron bremsstrahlung energy threshold");
270 br1Cmd->SetParameterName("emaxMuHadBrem",true);
271 br1Cmd->SetUnitCategory("Energy");
273 br1Cmd->SetToBeBroadcasted(false);
274
275 labCmd = new G4UIcmdWithADouble("/process/eLoss/LambdaFactor",this);
276 labCmd->SetGuidance("Set lambdaFactor parameter for integral option");
277 labCmd->SetParameterName("Fl",true);
279 labCmd->SetToBeBroadcasted(false);
280
281 mscfCmd = new G4UIcmdWithADouble("/process/msc/FactorForAngleLimit",this);
282 mscfCmd->SetGuidance("Set factor for computation of a limit for -t (invariant transfer)");
283 mscfCmd->SetParameterName("Fact",true);
284 mscfCmd->SetRange("Fact>0");
285 mscfCmd->SetDefaultValue(1.);
287 mscfCmd->SetToBeBroadcasted(false);
288
289 angCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/ThetaLimit",this);
290 angCmd->SetGuidance("Set the limit on the polar angle for msc and single scattering");
291 angCmd->SetParameterName("theta",true);
292 angCmd->SetUnitCategory("Angle");
294 angCmd->SetToBeBroadcasted(false);
295
296 msceCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/EnergyLimit",this);
297 msceCmd->SetGuidance("Set the upper energy limit for msc");
298 msceCmd->SetParameterName("mscE",true);
299 msceCmd->SetUnitCategory("Energy");
301 msceCmd->SetToBeBroadcasted(false);
302
303 nielCmd = new G4UIcmdWithADoubleAndUnit("/process/em/MaxEnergyNIEL",this);
304 nielCmd->SetGuidance("Set the upper energy limit for NIEL");
305 nielCmd->SetParameterName("niel",true);
306 nielCmd->SetUnitCategory("Energy");
308 nielCmd->SetToBeBroadcasted(false);
309
310 frCmd = new G4UIcmdWithADouble("/process/msc/RangeFactor",this);
311 frCmd->SetGuidance("Set RangeFactor for msc processes of e+-");
312 frCmd->SetParameterName("Fr",true);
313 frCmd->SetRange("Fr>0");
314 frCmd->SetDefaultValue(0.04);
316 frCmd->SetToBeBroadcasted(false);
317
318 fr1Cmd = new G4UIcmdWithADouble("/process/msc/RangeFactorMuHad",this);
319 fr1Cmd->SetGuidance("Set RangeFactor for msc processes of muons/hadrons");
320 fr1Cmd->SetParameterName("Fr1",true);
321 fr1Cmd->SetRange("Fr1>0");
322 fr1Cmd->SetDefaultValue(0.2);
324 fr1Cmd->SetToBeBroadcasted(false);
325
326 fgCmd = new G4UIcmdWithADouble("/process/msc/GeomFactor",this);
327 fgCmd->SetGuidance("Set GeomFactor parameter for msc processes");
328 fgCmd->SetParameterName("Fg",true);
329 fgCmd->SetRange("Fg>0");
330 fgCmd->SetDefaultValue(2.5);
332 fgCmd->SetToBeBroadcasted(false);
333
334 skinCmd = new G4UIcmdWithADouble("/process/msc/Skin",this);
335 skinCmd->SetGuidance("Set skin parameter for msc processes");
336 skinCmd->SetParameterName("skin",true);
338 skinCmd->SetToBeBroadcasted(false);
339
340 screCmd = new G4UIcmdWithADouble("/process/msc/ScreeningFactor",this);
341 screCmd->SetGuidance("Set screening factor");
342 screCmd->SetParameterName("screen",true);
344 screCmd->SetToBeBroadcasted(false);
345
346 safCmd = new G4UIcmdWithADouble("/process/msc/SafetyFactor",this);
347 safCmd->SetGuidance("Set safety factor");
348 safCmd->SetParameterName("fsafe",true);
350 safCmd->SetToBeBroadcasted(false);
351
352 llimCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/LambdaLimit",this);
353 llimCmd->SetGuidance("Set the upper energy limit for NIEL");
354 llimCmd->SetParameterName("ll",true);
355 llimCmd->SetUnitCategory("Length");
357 llimCmd->SetToBeBroadcasted(false);
358
359 amCmd = new G4UIcmdWithAnInteger("/process/em/binsPerDecade",this);
360 amCmd->SetGuidance("Set number of bins per decade for EM tables");
361 amCmd->SetParameterName("bins",true);
362 amCmd->SetDefaultValue(7);
364 amCmd->SetToBeBroadcasted(false);
365
366 verCmd = new G4UIcmdWithAnInteger("/process/eLoss/verbose",this);
367 verCmd->SetGuidance("Set verbose level for EM physics");
368 verCmd->SetParameterName("verb",true);
369 verCmd->SetDefaultValue(1);
371 verCmd->SetToBeBroadcasted(false);
372
373 ver1Cmd = new G4UIcmdWithAnInteger("/process/em/verbose",this);
374 ver1Cmd->SetGuidance("Set verbose level for EM physics");
375 ver1Cmd->SetParameterName("verb1",true);
376 ver1Cmd->SetDefaultValue(1);
378 ver1Cmd->SetToBeBroadcasted(false);
379
380 ver2Cmd = new G4UIcmdWithAnInteger("/process/em/workerVerbose",this);
381 ver2Cmd->SetGuidance("Set worker verbose level for EM physics");
382 ver2Cmd->SetParameterName("verb2",true);
383 ver2Cmd->SetDefaultValue(0);
385 ver2Cmd->SetToBeBroadcasted(false);
386
387 transWithMscCmd = new G4UIcmdWithAString("/process/em/transportationWithMsc",this);
388 transWithMscCmd->SetGuidance("Enable/disable the G4TransportationWithMsc process");
389 transWithMscCmd->SetParameterName("trans",true);
390 transWithMscCmd->SetCandidates("Disabled Enabled MultipleSteps");
391 transWithMscCmd->AvailableForStates(G4State_PreInit);
392 transWithMscCmd->SetToBeBroadcasted(false);
393
394 mscCmd = new G4UIcmdWithAString("/process/msc/StepLimit",this);
395 mscCmd->SetGuidance("Set msc step limitation type");
396 mscCmd->SetParameterName("StepLim",true);
397 mscCmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
399 mscCmd->SetToBeBroadcasted(false);
400
401 msc1Cmd = new G4UIcmdWithAString("/process/msc/StepLimitMuHad",this);
402 msc1Cmd->SetGuidance("Set msc step limitation type for muons/hadrons");
403 msc1Cmd->SetParameterName("StepLim1",true);
404 msc1Cmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
406 msc1Cmd->SetToBeBroadcasted(false);
407
408 dumpCmd = new G4UIcommand("/process/em/printParameters",this);
409 dumpCmd->SetGuidance("Print all EM parameters.");
411 dumpCmd->SetToBeBroadcasted(false);
412
413 nffCmd = new G4UIcmdWithAString("/process/em/setNuclearFormFactor",this);
414 nffCmd->SetGuidance("Define type of nuclear form-factor");
415 nffCmd->SetParameterName("NucFF",true);
416 nffCmd->SetCandidates("None Exponential Gaussian Flat");
418 nffCmd->SetToBeBroadcasted(false);
419
420 ssCmd = new G4UIcmdWithAString("/process/em/setSingleScattering",this);
421 ssCmd->SetGuidance("Define type of e+- single scattering model");
422 ssCmd->SetParameterName("SS",true);
423 ssCmd->SetCandidates("WVI Mott DPWA");
425 ssCmd->SetToBeBroadcasted(false);
426
427 fluc1Cmd = new G4UIcmdWithAString("/process/eloss/setFluctModel",this);
428 fluc1Cmd->SetGuidance("Define type of energy loss fluctuation model");
429 fluc1Cmd->SetParameterName("Fluc1",true);
430 fluc1Cmd->SetCandidates("Dummy Universal Urban");
432 fluc1Cmd->SetToBeBroadcasted(false);
433
434 tripletCmd = new G4UIcmdWithAnInteger("/process/gconv/conversionType",this);
435 tripletCmd->SetGuidance("gamma conversion triplet/nuclear generation type:");
436 tripletCmd->SetGuidance("0 - (default) both triplet and nuclear");
437 tripletCmd->SetGuidance("1 - force nuclear");
438 tripletCmd->SetGuidance("2 - force triplet");
439 tripletCmd->SetParameterName("type",false);
440 tripletCmd->SetRange("type >= 0 && type <= 2");
441 tripletCmd->SetDefaultValue(0);
443 tripletCmd->SetToBeBroadcasted(false);
444
445 onIsolatedCmd = new G4UIcmdWithABool("/process/gconv/onIsolated",this);
446 onIsolatedCmd->SetGuidance("Conversion on isolated charged particles");
447 onIsolatedCmd->SetGuidance("false (default) : atomic electron screening");
448 onIsolatedCmd->SetGuidance("true : conversion on isolated particles.");
449 onIsolatedCmd->SetParameterName("flag",false);
450 onIsolatedCmd->SetDefaultValue(false);
452 onIsolatedCmd->SetToBeBroadcasted(false);
453}
454
455//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
456
458{
459 delete gconvDirectory;
460 delete eLossDirectory;
461 delete mscDirectory;
462 delete emDirectory;
463 delete dnaDirectory;
464
465 delete flucCmd;
466 delete rangeCmd;
467 delete lpmCmd;
468 delete rsCmd;
469 delete aplCmd;
470 delete intCmd;
471 delete latCmd;
472 delete lat96Cmd;
473 delete mulatCmd;
474 delete delCmd;
475 delete mottCmd;
476 delete birksCmd;
477 delete sharkCmd;
478 delete onIsolatedCmd;
479 delete sampleTCmd;
480 delete poCmd;
481 delete icru90Cmd;
482 delete mudatCmd;
483 delete peKCmd;
484 delete mscPCmd;
485
486 delete minEnCmd;
487 delete maxEnCmd;
488 delete max5DCmd;
489 delete cenCmd;
490 delete lowEnCmd;
491 delete lowhEnCmd;
492 delete lowEn3Cmd;
493 delete lllCmd;
494 delete brCmd;
495 delete br1Cmd;
496 delete labCmd;
497 delete mscfCmd;
498 delete angCmd;
499 delete msceCmd;
500 delete nielCmd;
501 delete frCmd;
502 delete fr1Cmd;
503 delete fgCmd;
504 delete skinCmd;
505 delete safCmd;
506 delete llimCmd;
507 delete screCmd;
508
509 delete amCmd;
510 delete verCmd;
511 delete ver1Cmd;
512 delete ver2Cmd;
513 delete transWithMscCmd;
514 delete tripletCmd;
515
516 delete mscCmd;
517 delete msc1Cmd;
518 delete nffCmd;
519 delete ssCmd;
520 delete fluc1Cmd;
521
522 delete dumpCmd;
523}
524
525//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
526
528 G4String newValue)
529{
530 G4bool physicsModified = false;
531 if (command == flucCmd) {
532 theParameters->SetLossFluctuations(flucCmd->GetNewBoolValue(newValue));
533 physicsModified = true;
534 } else if (command == rangeCmd) {
535 theParameters->SetBuildCSDARange(rangeCmd->GetNewBoolValue(newValue));
536 } else if (command == lpmCmd) {
537 theParameters->SetLPM(lpmCmd->GetNewBoolValue(newValue));
538 physicsModified = true;
539 } else if (command == rsCmd) {
540 theParameters->SetUseCutAsFinalRange(rsCmd->GetNewBoolValue(newValue));
541 physicsModified = true;
542 } else if (command == aplCmd) {
543 theParameters->SetApplyCuts(aplCmd->GetNewBoolValue(newValue));
544 physicsModified = true;
545 } else if (command == intCmd) {
546 theParameters->SetIntegral(intCmd->GetNewBoolValue(newValue));
547 } else if (command == latCmd) {
548 theParameters->SetLateralDisplacement(latCmd->GetNewBoolValue(newValue));
549 physicsModified = true;
550 } else if (command == lat96Cmd) {
551 theParameters->SetLateralDisplacementAlg96(lat96Cmd->GetNewBoolValue(newValue));
552 physicsModified = true;
553 } else if (command == mulatCmd) {
554 theParameters->SetMuHadLateralDisplacement(mulatCmd->GetNewBoolValue(newValue));
555 physicsModified = true;
556 } else if (command == delCmd) {
557 theParameters->ActivateAngularGeneratorForIonisation(delCmd->GetNewBoolValue(newValue));
558 } else if (command == mottCmd) {
559 theParameters->SetUseMottCorrection(mottCmd->GetNewBoolValue(newValue));
560 } else if (command == birksCmd) {
561 theParameters->SetBirksActive(birksCmd->GetNewBoolValue(newValue));
562 } else if (command == icru90Cmd) {
563 theParameters->SetUseICRU90Data(icru90Cmd->GetNewBoolValue(newValue));
564 } else if (command == sharkCmd) {
565 theParameters->SetGeneralProcessActive(sharkCmd->GetNewBoolValue(newValue));
566 } else if (command == poCmd) {
567 theParameters->SetEnablePolarisation(poCmd->GetNewBoolValue(newValue));
568 } else if (command == sampleTCmd) {
569 theParameters->SetEnableSamplingTable(sampleTCmd->GetNewBoolValue(newValue));
570 } else if (command == mudatCmd) {
571 theParameters->SetRetrieveMuDataFromFile(mudatCmd->GetNewBoolValue(newValue));
572 } else if (command == peKCmd) {
573 theParameters->SetPhotoeffectBelowKShell(peKCmd->GetNewBoolValue(newValue));
574 } else if (command == mscPCmd) {
575 theParameters->SetMscPositronCorrection(mscPCmd->GetNewBoolValue(newValue));
576
577 } else if (command == minEnCmd) {
578 theParameters->SetMinEnergy(minEnCmd->GetNewDoubleValue(newValue));
579 } else if (command == maxEnCmd) {
580 theParameters->SetMaxEnergy(maxEnCmd->GetNewDoubleValue(newValue));
581 } else if (command == max5DCmd) {
582 theParameters->SetMaxEnergyFor5DMuPair(max5DCmd->GetNewDoubleValue(newValue));
583 } else if (command == cenCmd) {
584 theParameters->SetMaxEnergyForCSDARange(cenCmd->GetNewDoubleValue(newValue));
585 physicsModified = true;
586 } else if (command == lowEnCmd) {
587 theParameters->SetLowestElectronEnergy(lowEnCmd->GetNewDoubleValue(newValue));
588 physicsModified = true;
589 } else if (command == lowEn3Cmd) {
590 theParameters->SetLowestTripletEnergy(lowEn3Cmd->GetNewDoubleValue(newValue));
591 physicsModified = true;
592 } else if (command == lowhEnCmd) {
593 theParameters->SetLowestMuHadEnergy(lowhEnCmd->GetNewDoubleValue(newValue));
594 physicsModified = true;
595 } else if (command == lllCmd) {
596 theParameters->SetLinearLossLimit(lllCmd->GetNewDoubleValue(newValue));
597 physicsModified = true;
598 } else if (command == brCmd) {
599 theParameters->SetBremsstrahlungTh(brCmd->GetNewDoubleValue(newValue));
600 physicsModified = true;
601 } else if (command == br1Cmd) {
602 theParameters->SetMuHadBremsstrahlungTh(br1Cmd->GetNewDoubleValue(newValue));
603 physicsModified = true;
604 } else if (command == labCmd) {
605 theParameters->SetLambdaFactor(labCmd->GetNewDoubleValue(newValue));
606 physicsModified = true;
607 } else if (command == mscfCmd) {
608 theParameters->SetFactorForAngleLimit(mscfCmd->GetNewDoubleValue(newValue));
609 } else if (command == angCmd) {
610 theParameters->SetMscThetaLimit(angCmd->GetNewDoubleValue(newValue));
611 } else if (command == msceCmd) {
612 theParameters->SetMscEnergyLimit(msceCmd->GetNewDoubleValue(newValue));
613 } else if (command == nielCmd) {
614 theParameters->SetMaxNIELEnergy(nielCmd->GetNewDoubleValue(newValue));
615 } else if (command == frCmd) {
616 theParameters->SetMscRangeFactor(frCmd->GetNewDoubleValue(newValue));
617 physicsModified = true;
618 } else if (command == fr1Cmd) {
619 theParameters->SetMscMuHadRangeFactor(fr1Cmd->GetNewDoubleValue(newValue));
620 physicsModified = true;
621 } else if (command == fgCmd) {
622 theParameters->SetMscGeomFactor(fgCmd->GetNewDoubleValue(newValue));
623 physicsModified = true;
624 } else if (command == skinCmd) {
625 theParameters->SetMscSkin(skinCmd->GetNewDoubleValue(newValue));
626 physicsModified = true;
627 } else if (command == safCmd) {
628 theParameters->SetMscSafetyFactor(safCmd->GetNewDoubleValue(newValue));
629 } else if (command == llimCmd) {
630 theParameters->SetMscLambdaLimit(llimCmd->GetNewDoubleValue(newValue));
631 } else if (command == screCmd) {
632 theParameters->SetScreeningFactor(screCmd->GetNewDoubleValue(newValue));
633 } else if (command == amCmd) {
634 theParameters->SetNumberOfBinsPerDecade(amCmd->GetNewIntValue(newValue));
635 } else if (command == verCmd) {
636 theParameters->SetVerbose(verCmd->GetNewIntValue(newValue));
637 } else if (command == ver1Cmd) {
638 theParameters->SetVerbose(ver1Cmd->GetNewIntValue(newValue));
639 } else if (command == ver2Cmd) {
640 theParameters->SetWorkerVerbose(ver2Cmd->GetNewIntValue(newValue));
641 } else if (command == dumpCmd) {
642 theParameters->SetIsPrintedFlag(false);
643 theParameters->Dump();
644 } else if (command == transWithMscCmd) {
645 G4TransportationWithMscType type = G4TransportationWithMscType::fDisabled;
646 if(newValue == "Disabled") {
647 type = G4TransportationWithMscType::fDisabled;
648 } else if(newValue == "Enabled") {
649 type = G4TransportationWithMscType::fEnabled;
650 } else if(newValue == "MultipleSteps") {
651 type = G4TransportationWithMscType::fMultipleSteps;
652 } else {
654 ed << " TransportationWithMsc type <" << newValue << "> unknown!";
655 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
656 }
657 theParameters->SetTransportationWithMsc(type);
658 } else if (command == mscCmd || command == msc1Cmd) {
660 if(newValue == "Minimal") {
661 msctype = fMinimal;
662 } else if(newValue == "UseDistanceToBoundary") {
663 msctype = fUseDistanceToBoundary;
664 } else if(newValue == "UseSafety") {
665 msctype = fUseSafety;
666 } else if(newValue == "UseSafetyPlus") {
667 msctype = fUseSafetyPlus;
668 } else {
670 ed << " StepLimit type <" << newValue << "> unknown!";
671 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
672 return;
673 }
674 if (command == mscCmd) {
675 theParameters->SetMscStepLimitType(msctype);
676 } else {
677 theParameters->SetMscMuHadStepLimitType(msctype);
678 }
679 physicsModified = true;
680 } else if (command == nffCmd) {
682 if(newValue == "Exponential") { x = fExponentialNF; }
683 else if(newValue == "Gaussian") { x = fGaussianNF; }
684 else if(newValue == "Flat") { x = fFlatNF; }
685 else if(newValue != "None") {
687 ed << " NuclearFormFactor type <" << newValue << "> unknown!";
688 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
689 return;
690 }
691 theParameters->SetNuclearFormfactorType(x);
692 } else if (command == ssCmd) {
694 if(newValue == "DPWA") { x = fDPWA; }
695 else if(newValue == "Mott") { x = fMott; }
696 else if(newValue != "WVI") {
698 ed << " G4eSingleScatteringType type <" << newValue << "> unknown!";
699 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
700 return;
701 }
702 theParameters->SetSingleScatteringType(x);
703 } else if (command == fluc1Cmd) {
705 if(newValue == "Dummy") { x = fDummyFluctuation; }
706 else if(newValue == "Urban") { x = fUrbanFluctuation; }
707 theParameters->SetFluctuationType(x);
708 } else if ( command==tripletCmd ) {
709 theParameters->SetConversionType(tripletCmd->GetNewIntValue(newValue));
710 } else if ( command==onIsolatedCmd ) {
711 theParameters->SetOnIsolated(onIsolatedCmd->GetNewBoolValue(newValue));
712 physicsModified = true;
713 }
714
715 if(physicsModified) {
716 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
717 }
718}
719
720//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ G4State_Init
@ G4State_Idle
@ G4State_PreInit
G4EmFluctuationType
@ fUniversalFluctuation
@ fUrbanFluctuation
@ fDummyFluctuation
G4eSingleScatteringType
@ fWVI
@ fMott
@ fDPWA
G4TransportationWithMscType
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4MscStepLimitType
@ fMinimal
@ fUseSafety
@ fUseSafetyPlus
@ fUseDistanceToBoundary
G4NuclearFormfactorType
bool G4bool
Definition: G4Types.hh:86
void SetNewValue(G4UIcommand *, G4String) override
G4EmParametersMessenger(G4EmParameters *)
void SetLambdaFactor(G4double val)
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetBuildCSDARange(G4bool val)
void SetEnablePolarisation(G4bool val)
void SetNumberOfBinsPerDecade(G4int val)
void SetGeneralProcessActive(G4bool val)
void SetMscSafetyFactor(G4double val)
void SetLateralDisplacementAlg96(G4bool val)
void SetFactorForAngleLimit(G4double val)
void SetRetrieveMuDataFromFile(G4bool v)
void SetMscMuHadRangeFactor(G4double val)
void SetLPM(G4bool val)
void SetMaxEnergyFor5DMuPair(G4double val)
void SetLinearLossLimit(G4double val)
void SetMscThetaLimit(G4double val)
void SetLossFluctuations(G4bool val)
void SetLowestTripletEnergy(G4double val)
void SetMuHadLateralDisplacement(G4bool val)
void ActivateAngularGeneratorForIonisation(G4bool val)
void SetScreeningFactor(G4double val)
void SetNuclearFormfactorType(G4NuclearFormfactorType val)
void SetLateralDisplacement(G4bool val)
void SetWorkerVerbose(G4int val)
void SetUseCutAsFinalRange(G4bool val)
void SetBirksActive(G4bool val)
void SetMuHadBremsstrahlungTh(G4double val)
void SetFluctuationType(G4EmFluctuationType val)
void SetVerbose(G4int val)
void SetMscGeomFactor(G4double val)
void SetMscLambdaLimit(G4double val)
void SetMscSkin(G4double val)
void SetApplyCuts(G4bool val)
void SetEnableSamplingTable(G4bool val)
void SetMaxNIELEnergy(G4double val)
void SetMaxEnergyForCSDARange(G4double val)
void SetMscMuHadStepLimitType(G4MscStepLimitType val)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetMscEnergyLimit(G4double val)
void SetBremsstrahlungTh(G4double val)
void SetIsPrintedFlag(G4bool val)
void SetConversionType(G4int val)
void SetUseICRU90Data(G4bool val)
void SetOnIsolated(G4bool val)
void SetTransportationWithMsc(G4TransportationWithMscType val)
void SetIntegral(G4bool val)
void SetUseMottCorrection(G4bool val)
void SetMscPositronCorrection(G4bool v)
void SetLowestMuHadEnergy(G4double val)
void SetMaxEnergy(G4double val)
void SetSingleScatteringType(G4eSingleScatteringType val)
void SetPhotoeffectBelowKShell(G4bool v)
void SetMscRangeFactor(G4double val)
static G4bool GetNewBoolValue(const char *paramString)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(G4bool defVal)
void SetUnitCategory(const char *unitCategory)
static G4double GetNewDoubleValue(const char *paramString)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
static G4double GetNewDoubleValue(const char *paramString)
void SetDefaultValue(G4double defVal)
void SetCandidates(const char *candidateList)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
static G4int GetNewIntValue(const char *paramString)
void SetDefaultValue(G4int defVal)
void SetToBeBroadcasted(G4bool val)
Definition: G4UIcommand.hh:172
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:157
void SetRange(const char *rs)
Definition: G4UIcommand.hh:121
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:287
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:495
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77