Geant4-11
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 gconvDirectory = new G4UIdirectory("/process/gconv/");
65 gconvDirectory->SetGuidance("Commands for EM gamma conversion BH5D model.");
66 eLossDirectory = new G4UIdirectory("/process/eLoss/");
67 eLossDirectory->SetGuidance("Commands for EM processes.");
68 mscDirectory = new G4UIdirectory("/process/msc/");
69 mscDirectory->SetGuidance("Commands for EM scattering processes.");
70 emDirectory = new G4UIdirectory("/process/em/");
71 emDirectory->SetGuidance("General commands for EM processes.");
72 dnaDirectory = new G4UIdirectory("/process/dna/");
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);
81
82 rangeCmd = new G4UIcmdWithABool("/process/eLoss/CSDARange",this);
83 rangeCmd->SetGuidance("Enable/disable CSDA range calculation");
84 rangeCmd->SetParameterName("range",true);
88
89 lpmCmd = new G4UIcmdWithABool("/process/eLoss/LPM",this);
90 lpmCmd->SetGuidance("Enable/disable LPM effect calculation");
91 lpmCmd->SetParameterName("lpm",true);
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);
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);
109
110 latCmd = new G4UIcmdWithABool("/process/msc/LateralDisplacement",this);
111 latCmd->SetGuidance("Enable/disable sampling of lateral displacement");
112 latCmd->SetParameterName("lat",true);
113 latCmd->SetDefaultValue(true);
116
117 lat96Cmd = new G4UIcmdWithABool("/process/msc/LateralDisplacementAlg96",this);
118 lat96Cmd->SetGuidance("Enable/disable sampling of lateral displacement");
119 lat96Cmd->SetParameterName("lat96",true);
123
124 mulatCmd = new G4UIcmdWithABool("/process/msc/MuHadLateralDisplacement",this);
125 mulatCmd->SetGuidance("Enable/disable sampling of lateral displacement for muons and hadrons");
126 mulatCmd->SetParameterName("mulat",true);
130
131 delCmd = new G4UIcmdWithABool("/process/eLoss/UseAngularGenerator",this);
132 delCmd->SetGuidance("Enable usage of angular generator for ionisation");
133 delCmd->SetParameterName("del",true);
134 delCmd->SetDefaultValue(false);
137
138 mottCmd = new G4UIcmdWithABool("/process/msc/UseMottCorrection",this);
139 mottCmd->SetGuidance("Enable usage of Mott corrections for e- elastic scattering");
140 mottCmd->SetParameterName("mott",true);
141 mottCmd->SetDefaultValue(false);
144
145 birksCmd = new G4UIcmdWithABool("/process/msc/UseG4EmSaturation",this);
146 birksCmd->SetGuidance("Enable usage of built-in Birks saturation");
147 birksCmd->SetParameterName("birks",true);
151
152 sharkCmd = new G4UIcmdWithABool("/process/em/UseGeneralProcess",this);
153 sharkCmd->SetGuidance("Enable gamma, e+- general process");
154 sharkCmd->SetParameterName("gen",true);
158
159 poCmd = new G4UIcmdWithABool("/process/em/Polarisation",this);
160 poCmd->SetGuidance("Enable polarisation");
163
164 sampleTCmd = new G4UIcmdWithABool("/process/em/enableSamplingTable",this);
165 sampleTCmd->SetGuidance("Enable usage of sampling table for secondary generation");
166 sampleTCmd->SetParameterName("sampleT",true);
170
171 icru90Cmd = new G4UIcmdWithABool("/process/eLoss/UseICRU90",this);
172 icru90Cmd->SetGuidance("Enable usage of ICRU90 stopping powers");
173 icru90Cmd->SetParameterName("icru90",true);
177
178 mudatCmd = new G4UIcmdWithABool("/process/em/MuDataFromFile",this);
179 mudatCmd->SetGuidance("Enable usage of muon data from file");
180 mudatCmd->SetParameterName("mudat",true);
184
185 minEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/minKinEnergy",this);
186 minEnCmd->SetGuidance("Set the min kinetic energy for EM tables");
187 minEnCmd->SetParameterName("emin",true);
188 minEnCmd->SetUnitCategory("Energy");
191
192 maxEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergy",this);
193 maxEnCmd->SetGuidance("Set the max kinetic energy for EM tables");
194 maxEnCmd->SetParameterName("emax",true);
195 maxEnCmd->SetUnitCategory("Energy");
198
199 cenCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergyCSDA",this);
200 cenCmd->SetGuidance("Set the max kinetic energy for CSDA table");
201 cenCmd->SetParameterName("emaxCSDA",true);
202 cenCmd->SetUnitCategory("Energy");
205
206 max5DCmd = new G4UIcmdWithADoubleAndUnit("/process/em/max5DMuPairEnergy",this);
207 max5DCmd->SetGuidance("Set the max kinetic energy for 5D muon pair production");
208 max5DCmd->SetParameterName("emax5D",true);
209 max5DCmd->SetUnitCategory("Energy");
212
213 lowEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestElectronEnergy",this);
214 lowEnCmd->SetGuidance("Set the lowest kinetic energy for e+-");
215 lowEnCmd->SetParameterName("elow",true);
216 lowEnCmd->SetUnitCategory("Energy");
219
220 lowhEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestMuHadEnergy",this);
221 lowhEnCmd->SetGuidance("Set the lowest kinetic energy for muons and hadrons");
222 lowhEnCmd->SetParameterName("elowh",true);
223 lowhEnCmd->SetUnitCategory("Energy");
226
227 lowEn3Cmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestTripletEnergy",this);
228 lowEn3Cmd->SetGuidance("Set the lowest kinetic energy for triplet production");
229 lowEn3Cmd->SetParameterName("elow3",true);
230 lowEn3Cmd->SetUnitCategory("Energy");
233
234 lllCmd = new G4UIcmdWithADouble("/process/eLoss/linLossLimit",this);
235 lllCmd->SetGuidance("Set linearLossLimit parameter");
236 lllCmd->SetParameterName("linlim",true);
239
240 brCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/bremThreshold",this);
241 brCmd->SetGuidance("Set e+- bremsstrahlung energy threshold");
242 brCmd->SetParameterName("emaxBrem",true);
243 brCmd->SetUnitCategory("Energy");
246
247 br1Cmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/bremMuHadThreshold",this);
248 br1Cmd->SetGuidance("Set muon/hadron bremsstrahlung energy threshold");
249 br1Cmd->SetParameterName("emaxMuHadBrem",true);
250 br1Cmd->SetUnitCategory("Energy");
253
254 labCmd = new G4UIcmdWithADouble("/process/eLoss/LambdaFactor",this);
255 labCmd->SetGuidance("Set lambdaFactor parameter for integral option");
256 labCmd->SetParameterName("Fl",true);
259
260 mscfCmd = new G4UIcmdWithADouble("/process/msc/FactorForAngleLimit",this);
261 mscfCmd->SetGuidance("Set factor for computation of a limit for -t (invariant transfer)");
262 mscfCmd->SetParameterName("Fact",true);
263 mscfCmd->SetRange("Fact>0");
267
268 angCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/ThetaLimit",this);
269 angCmd->SetGuidance("Set the limit on the polar angle for msc and single scattering");
270 angCmd->SetParameterName("theta",true);
271 angCmd->SetUnitCategory("Angle");
274
275 msceCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/EnergyLimit",this);
276 msceCmd->SetGuidance("Set the upper energy limit for msc");
277 msceCmd->SetParameterName("mscE",true);
278 msceCmd->SetUnitCategory("Energy");
281
282 nielCmd = new G4UIcmdWithADoubleAndUnit("/process/em/MaxEnergyNIEL",this);
283 nielCmd->SetGuidance("Set the upper energy limit for NIEL");
284 nielCmd->SetParameterName("niel",true);
285 nielCmd->SetUnitCategory("Energy");
288
289 frCmd = new G4UIcmdWithADouble("/process/msc/RangeFactor",this);
290 frCmd->SetGuidance("Set RangeFactor for msc processes of e+-");
291 frCmd->SetParameterName("Fr",true);
292 frCmd->SetRange("Fr>0");
293 frCmd->SetDefaultValue(0.04);
296
297 fr1Cmd = new G4UIcmdWithADouble("/process/msc/RangeFactorMuHad",this);
298 fr1Cmd->SetGuidance("Set RangeFactor for msc processes of muons/hadrons");
299 fr1Cmd->SetParameterName("Fr1",true);
300 fr1Cmd->SetRange("Fr1>0");
304
305 fgCmd = new G4UIcmdWithADouble("/process/msc/GeomFactor",this);
306 fgCmd->SetGuidance("Set GeomFactor parameter for msc processes");
307 fgCmd->SetParameterName("Fg",true);
308 fgCmd->SetRange("Fg>0");
312
313 skinCmd = new G4UIcmdWithADouble("/process/msc/Skin",this);
314 skinCmd->SetGuidance("Set skin parameter for msc processes");
315 skinCmd->SetParameterName("skin",true);
318
319 screCmd = new G4UIcmdWithADouble("/process/msc/ScreeningFactor",this);
320 screCmd->SetGuidance("Set screening factor");
321 screCmd->SetParameterName("screen",true);
324
325 safCmd = new G4UIcmdWithADouble("/process/msc/SafetyFactor",this);
326 safCmd->SetGuidance("Set safety factor");
327 safCmd->SetParameterName("fsafe",true);
330
331 llimCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/LambdaLimit",this);
332 llimCmd->SetGuidance("Set the upper energy limit for NIEL");
333 llimCmd->SetParameterName("ll",true);
334 llimCmd->SetUnitCategory("Length");
337
338 amCmd = new G4UIcmdWithAnInteger("/process/eLoss/binsPerDecade",this);
339 amCmd->SetGuidance("Set number of bins per decade for EM tables");
340 amCmd->SetParameterName("bins",true);
344
345 verCmd = new G4UIcmdWithAnInteger("/process/eLoss/verbose",this);
346 verCmd->SetGuidance("Set verbose level for EM physics");
347 verCmd->SetParameterName("verb",true);
351
352 ver1Cmd = new G4UIcmdWithAnInteger("/process/em/verbose",this);
353 ver1Cmd->SetGuidance("Set verbose level for EM physics");
354 ver1Cmd->SetParameterName("verb1",true);
358
359 ver2Cmd = new G4UIcmdWithAnInteger("/process/em/workerVerbose",this);
360 ver2Cmd->SetGuidance("Set worker verbose level for EM physics");
361 ver2Cmd->SetParameterName("verb2",true);
365
366 mscCmd = new G4UIcmdWithAString("/process/msc/StepLimit",this);
367 mscCmd->SetGuidance("Set msc step limitation type");
368 mscCmd->SetParameterName("StepLim",true);
369 mscCmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
372
373 msc1Cmd = new G4UIcmdWithAString("/process/msc/StepLimitMuHad",this);
374 msc1Cmd->SetGuidance("Set msc step limitation type for muons/hadrons");
375 msc1Cmd->SetParameterName("StepLim1",true);
376 msc1Cmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
379
380 dumpCmd = new G4UIcommand("/process/em/printParameters",this);
381 dumpCmd->SetGuidance("Print all EM parameters.");
384
385 nffCmd = new G4UIcmdWithAString("/process/em/setNuclearFormFactor",this);
386 nffCmd->SetGuidance("Define type of nuclear form-factor");
387 nffCmd->SetParameterName("NucFF",true);
388 nffCmd->SetCandidates("None Exponential Gaussian Flat");
391
392 ssCmd = new G4UIcmdWithAString("/process/em/setSingleScattering",this);
393 ssCmd->SetGuidance("Define type of e+- single scattering model");
394 ssCmd->SetParameterName("SS",true);
395 ssCmd->SetCandidates("WVI Mott DPWA");
398
399 tripletCmd = new G4UIcmdWithAnInteger("/process/gconv/conversionType",this);
400 tripletCmd->SetGuidance("gamma conversion triplet/nuclear generation type:");
401 tripletCmd->SetGuidance("0 - (default) both triplet and nuclear");
402 tripletCmd->SetGuidance("1 - force nuclear");
403 tripletCmd->SetGuidance("2 - force triplet");
404 tripletCmd->SetParameterName("type",false);
405 tripletCmd->SetRange("type >= 0 && type <= 2");
409
410 onIsolatedCmd = new G4UIcmdWithABool("/process/gconv/onIsolated",this);
411 onIsolatedCmd->SetGuidance("Conversion on isolated charged particles");
412 onIsolatedCmd->SetGuidance("false (default) : atomic electron screening");
413 onIsolatedCmd->SetGuidance("true : conversion on isolated particles.");
414 onIsolatedCmd->SetParameterName("flag",false);
418}
419
420//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
421
423{
424 delete gconvDirectory;
425 delete eLossDirectory;
426 delete mscDirectory;
427 delete emDirectory;
428 delete dnaDirectory;
429
430 delete flucCmd;
431 delete rangeCmd;
432 delete lpmCmd;
433 delete rsCmd;
434 delete aplCmd;
435 delete latCmd;
436 delete lat96Cmd;
437 delete mulatCmd;
438 delete delCmd;
439 delete mottCmd;
440 delete birksCmd;
441 delete sharkCmd;
442 delete onIsolatedCmd;
443 delete sampleTCmd;
444 delete poCmd;
445 delete icru90Cmd;
446 delete mudatCmd;
447
448 delete minEnCmd;
449 delete maxEnCmd;
450 delete max5DCmd;
451 delete cenCmd;
452 delete lowEnCmd;
453 delete lowhEnCmd;
454 delete lowEn3Cmd;
455 delete lllCmd;
456 delete brCmd;
457 delete br1Cmd;
458 delete labCmd;
459 delete mscfCmd;
460 delete angCmd;
461 delete msceCmd;
462 delete nielCmd;
463 delete frCmd;
464 delete fr1Cmd;
465 delete fgCmd;
466 delete skinCmd;
467 delete safCmd;
468 delete llimCmd;
469 delete screCmd;
470
471 delete amCmd;
472 delete verCmd;
473 delete ver1Cmd;
474 delete ver2Cmd;
475 delete tripletCmd;
476
477 delete mscCmd;
478 delete msc1Cmd;
479 delete nffCmd;
480 delete ssCmd;
481
482 delete dumpCmd;
483}
484
485//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
486
488 G4String newValue)
489{
490 G4bool physicsModified = false;
491 if (command == flucCmd) {
493 physicsModified = true;
494 } else if (command == rangeCmd) {
496 } else if (command == lpmCmd) {
498 physicsModified = true;
499 } else if (command == rsCmd) {
501 physicsModified = true;
502 } else if (command == aplCmd) {
504 physicsModified = true;
505 } else if (command == latCmd) {
507 physicsModified = true;
508 } else if (command == lat96Cmd) {
510 physicsModified = true;
511 } else if (command == mulatCmd) {
513 physicsModified = true;
514 } else if (command == delCmd) {
516 } else if (command == mottCmd) {
518 } else if (command == birksCmd) {
520 } else if (command == icru90Cmd) {
522 } else if (command == sharkCmd) {
524 } else if (command == poCmd) {
526 } else if (command == sampleTCmd) {
528 } else if (command == mudatCmd) {
530
531 } else if (command == minEnCmd) {
533 } else if (command == maxEnCmd) {
535 } else if (command == max5DCmd) {
537 } else if (command == cenCmd) {
539 physicsModified = true;
540 } else if (command == lowEnCmd) {
542 physicsModified = true;
543 } else if (command == lowEn3Cmd) {
545 physicsModified = true;
546 } else if (command == lowhEnCmd) {
548 physicsModified = true;
549 } else if (command == lllCmd) {
551 physicsModified = true;
552 } else if (command == brCmd) {
554 physicsModified = true;
555 } else if (command == br1Cmd) {
557 physicsModified = true;
558 } else if (command == labCmd) {
560 physicsModified = true;
561 } else if (command == mscfCmd) {
563 } else if (command == angCmd) {
565 } else if (command == msceCmd) {
567 } else if (command == nielCmd) {
569 } else if (command == frCmd) {
571 physicsModified = true;
572 } else if (command == fr1Cmd) {
574 physicsModified = true;
575 } else if (command == fgCmd) {
577 physicsModified = true;
578 } else if (command == skinCmd) {
580 physicsModified = true;
581 } else if (command == safCmd) {
583 } else if (command == llimCmd) {
585 } else if (command == screCmd) {
587 } else if (command == verCmd) {
589 } else if (command == ver1Cmd) {
591 } else if (command == ver2Cmd) {
593 } else if (command == dumpCmd) {
596 } else if (command == mscCmd || command == msc1Cmd) {
598 if(newValue == "Minimal") {
599 msctype = fMinimal;
600 } else if(newValue == "UseDistanceToBoundary") {
601 msctype = fUseDistanceToBoundary;
602 } else if(newValue == "UseSafety") {
603 msctype = fUseSafety;
604 } else if(newValue == "UseSafetyPlus") {
605 msctype = fUseSafetyPlus;
606 } else {
608 ed << " StepLimit type <" << newValue << "> unknown!";
609 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
610 return;
611 }
612 if (command == mscCmd) {
614 } else {
616 }
617 physicsModified = true;
618 } else if (command == nffCmd) {
620 if(newValue == "Exponential") { x = fExponentialNF; }
621 else if(newValue == "Gaussian") { x = fGaussianNF; }
622 else if(newValue == "Flat") { x = fFlatNF; }
623 else if(newValue != "None") {
625 ed << " NuclearFormFactor type <" << newValue << "> unknown!";
626 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
627 return;
628 }
630 } else if (command == ssCmd) {
632 if(newValue == "DPWA") { x = fDPWA; }
633 else if(newValue == "Mott") { x = fMott; }
634 else if(newValue != "WVI") {
636 ed << " G4eSingleScatteringType type <" << newValue << "> unknown!";
637 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
638 return;
639 }
641 } else if ( command==tripletCmd ) {
643 } else if ( command==onIsolatedCmd ) {
645 physicsModified = true;
646 }
647
648 if(physicsModified) {
649 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
650 }
651}
652
653//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
@ G4State_Init
@ G4State_Idle
@ G4State_PreInit
G4eSingleScatteringType
@ fWVI
@ fMott
@ fDPWA
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
G4MscStepLimitType
@ fMinimal
@ fUseSafety
@ fUseSafetyPlus
@ fUseDistanceToBoundary
G4NuclearFormfactorType
bool G4bool
Definition: G4Types.hh:86
G4UIcmdWithADoubleAndUnit * lowEn3Cmd
G4UIcmdWithADoubleAndUnit * angCmd
G4UIcmdWithADoubleAndUnit * brCmd
G4UIcmdWithADoubleAndUnit * max5DCmd
G4UIcmdWithADoubleAndUnit * llimCmd
G4UIcmdWithAnInteger * ver1Cmd
G4UIcmdWithAnInteger * verCmd
G4UIcmdWithADoubleAndUnit * cenCmd
void SetNewValue(G4UIcommand *, G4String) override
G4UIcmdWithAnInteger * amCmd
G4EmParametersMessenger(G4EmParameters *)
G4UIcmdWithADoubleAndUnit * msceCmd
G4UIcmdWithADoubleAndUnit * maxEnCmd
G4UIcmdWithADoubleAndUnit * nielCmd
G4UIcmdWithADoubleAndUnit * lowhEnCmd
G4UIcmdWithADoubleAndUnit * br1Cmd
G4UIcmdWithADoubleAndUnit * minEnCmd
G4UIcmdWithADoubleAndUnit * lowEnCmd
G4UIcmdWithAnInteger * tripletCmd
G4UIcmdWithAnInteger * ver2Cmd
void SetLambdaFactor(G4double val)
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetBuildCSDARange(G4bool val)
void SetEnablePolarisation(G4bool 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 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 SetUseMottCorrection(G4bool val)
void SetLowestMuHadEnergy(G4double val)
void SetMaxEnergy(G4double val)
void SetSingleScatteringType(G4eSingleScatteringType val)
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:156
void SetRange(const char *rs)
Definition: G4UIcommand.hh:120
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:288
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:485
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77