Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmParametersMessenger Class Reference

#include <G4EmParametersMessenger.hh>

+ Inheritance diagram for G4EmParametersMessenger:

Public Member Functions

 G4EmParametersMessenger (G4EmParameters *)
 
virtual ~G4EmParametersMessenger ()
 
virtual void SetNewValue (G4UIcommand *, G4String) override
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
virtual void SetNewValue (G4UIcommand *command, G4String newValue)
 
G4bool operator== (const G4UImessenger &messenger) const
 
G4bool operator!= (const G4UImessenger &messenger) const
 
G4bool CommandsShouldBeInMaster () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (G4String s)
 
G4long StoL (G4String s)
 
G4double StoD (G4String s)
 
G4bool StoB (G4String s)
 
void AddUIcommand (G4UIcommand *newCommand)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 
- Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Detailed Description

Definition at line 65 of file G4EmParametersMessenger.hh.

Constructor & Destructor Documentation

◆ G4EmParametersMessenger()

G4EmParametersMessenger::G4EmParametersMessenger ( G4EmParameters ptr)
explicit

Definition at line 61 of file G4EmParametersMessenger.cc.

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);
78 flucCmd->SetDefaultValue(true);
80
81 rangeCmd = new G4UIcmdWithABool("/process/eLoss/CSDARange",this);
82 rangeCmd->SetGuidance("Enable/disable CSDA range calculation");
83 rangeCmd->SetParameterName("range",true);
84 rangeCmd->SetDefaultValue(false);
86
87 lpmCmd = new G4UIcmdWithABool("/process/eLoss/LPM",this);
88 lpmCmd->SetGuidance("Enable/disable LPM effect calculation");
89 lpmCmd->SetParameterName("lpm",true);
90 lpmCmd->SetDefaultValue(true);
92
93 splCmd = new G4UIcmdWithABool("/process/em/spline",this);
94 splCmd->SetGuidance("Enable/disable usage spline for Physics Vectors");
95 splCmd->SetParameterName("spl",true);
96 splCmd->SetDefaultValue(false);
98
99 rsCmd = new G4UIcmdWithABool("/process/eLoss/useCutAsFinalRange",this);
100 rsCmd->SetGuidance("Enable/disable use of cut in range as a final range");
101 rsCmd->SetParameterName("choice",true);
102 rsCmd->SetDefaultValue(false);
104
105 aplCmd = new G4UIcmdWithABool("/process/em/applyCuts",this);
106 aplCmd->SetGuidance("Enable/disable applying cuts for gamma processes");
107 aplCmd->SetParameterName("apl",true);
108 aplCmd->SetDefaultValue(false);
110
111 latCmd = new G4UIcmdWithABool("/process/msc/LateralDisplacement",this);
112 latCmd->SetGuidance("Enable/disable sampling of lateral displacement");
113 latCmd->SetParameterName("lat",true);
114 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);
120 lat96Cmd->SetDefaultValue(false);
122
123 mulatCmd = new G4UIcmdWithABool("/process/msc/MuHadLateralDisplacement",this);
124 mulatCmd->SetGuidance("Enable/disable sampling of lateral displacement for muons and hadrons");
125 mulatCmd->SetParameterName("mulat",true);
126 mulatCmd->SetDefaultValue(true);
128
129 delCmd = new G4UIcmdWithABool("/process/eLoss/UseAngularGenerator",this);
130 delCmd->SetGuidance("Enable usage of angular generator for ionisation");
131 delCmd->SetParameterName("del",true);
132 delCmd->SetDefaultValue(false);
134
135 IntegCmd = new G4UIcmdWithABool("/process/eLoss/integral",this);
136 IntegCmd->SetGuidance("Switch true/false the integral option");
137 IntegCmd->SetParameterName("integ",true);
138 IntegCmd->SetDefaultValue(true);
140
141 mottCmd = new G4UIcmdWithABool("/process/msc/UseMottCorrection",this);
142 mottCmd->SetGuidance("Enable usage of Mott corrections for e- elastic scattering");
143 mottCmd->SetParameterName("mott",true);
144 mottCmd->SetDefaultValue(false);
146
147 birksCmd = new G4UIcmdWithABool("/process/msc/UseG4EmSaturation",this);
148 birksCmd->SetGuidance("Enable usage of built-in Birks saturation");
149 birksCmd->SetParameterName("birks",true);
150 birksCmd->SetDefaultValue(false);
152
153 sharkCmd = new G4UIcmdWithABool("/process/em/UseGeneralProcess",this);
154 sharkCmd->SetGuidance("Enable gamma, e+- general process");
155 sharkCmd->SetParameterName("gen",true);
156 sharkCmd->SetDefaultValue(false);
158
159 sampleTCmd = new G4UIcmdWithABool("/process/em/enableSamplingTable",this);
160 sampleTCmd->SetGuidance("Enable usage of sampling table for secondary generation");
161 sampleTCmd->SetParameterName("sampleT",true);
162 sampleTCmd->SetDefaultValue(false);
164
165 icru90Cmd = new G4UIcmdWithABool("/process/eLoss/UseICRU90",this);
166 icru90Cmd->SetGuidance("Enable usage of ICRU90 stopping powers");
167 icru90Cmd->SetParameterName("icru90",true);
168 icru90Cmd->SetDefaultValue(false);
170
171 mudatCmd = new G4UIcmdWithABool("/process/em/MuDataFromFile",this);
172 mudatCmd->SetGuidance("Enable usage of muon data from file");
173 mudatCmd->SetParameterName("mudat",true);
174 mudatCmd->SetDefaultValue(false);
176
177 minSubSecCmd = new G4UIcmdWithADouble("/process/eLoss/minsubsec",this);
178 minSubSecCmd->SetGuidance("Set the ratio subcut/cut ");
179 minSubSecCmd->SetParameterName("rcmin",true);
180 minSubSecCmd->AvailableForStates(G4State_PreInit);
181
182 minEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/minKinEnergy",this);
183 minEnCmd->SetGuidance("Set the min kinetic energy for EM tables");
184 minEnCmd->SetParameterName("emin",true);
185 minEnCmd->SetUnitCategory("Energy");
187
188 maxEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergy",this);
189 maxEnCmd->SetGuidance("Set the max kinetic energy for EM tables");
190 maxEnCmd->SetParameterName("emax",true);
191 maxEnCmd->SetUnitCategory("Energy");
193
194 cenCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergyCSDA",this);
195 cenCmd->SetGuidance("Set the max kinetic energy for CSDA table");
196 cenCmd->SetParameterName("emaxCSDA",true);
197 cenCmd->SetUnitCategory("Energy");
199
200 max5DCmd = new G4UIcmdWithADoubleAndUnit("/process/em/max5DMuPairEnergy",this);
201 max5DCmd->SetGuidance("Set the max kinetic energy for 5D muon pair production");
202 max5DCmd->SetParameterName("emax5D",true);
203 max5DCmd->SetUnitCategory("Energy");
205
206 lowEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestElectronEnergy",this);
207 lowEnCmd->SetGuidance("Set the lowest kinetic energy for e+-");
208 lowEnCmd->SetParameterName("elow",true);
209 lowEnCmd->SetUnitCategory("Energy");
211
212 lowhEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestMuHadEnergy",this);
213 lowhEnCmd->SetGuidance("Set the lowest kinetic energy for muons and hadrons");
214 lowhEnCmd->SetParameterName("elowh",true);
215 lowhEnCmd->SetUnitCategory("Energy");
217
218 lowEn3Cmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestTripletEnergy",this);
219 lowEn3Cmd->SetGuidance("Set the lowest kinetic energy for triplet production");
220 lowEn3Cmd->SetParameterName("elow3",true);
221 lowEn3Cmd->SetUnitCategory("Energy");
223
224 lllCmd = new G4UIcmdWithADouble("/process/eLoss/linLossLimit",this);
225 lllCmd->SetGuidance("Set linearLossLimit parameter");
226 lllCmd->SetParameterName("linlim",true);
228
229 brCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/bremThreshold",this);
230 brCmd->SetGuidance("Set e+- bremsstrahlung energy threshold");
231 brCmd->SetParameterName("emaxBrem",true);
232 brCmd->SetUnitCategory("Energy");
234
235 br1Cmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/bremMuHadThreshold",this);
236 br1Cmd->SetGuidance("Set muon/hadron bremsstrahlung energy threshold");
237 br1Cmd->SetParameterName("emaxMuHadBrem",true);
238 br1Cmd->SetUnitCategory("Energy");
240
241 labCmd = new G4UIcmdWithADouble("/process/eLoss/LambdaFactor",this);
242 labCmd->SetGuidance("Set lambdaFactor parameter for integral option");
243 labCmd->SetParameterName("Fl",true);
245
246 mscfCmd = new G4UIcmdWithADouble("/process/msc/FactorForAngleLimit",this);
247 mscfCmd->SetGuidance("Set factor for computation of a limit for -t (invariant transfer)");
248 mscfCmd->SetParameterName("Fact",true);
249 mscfCmd->SetRange("Fact>0");
250 mscfCmd->SetDefaultValue(1.);
252
253 angCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/ThetaLimit",this);
254 angCmd->SetGuidance("Set the limit on the polar angle for msc and single scattering");
255 angCmd->SetParameterName("theta",true);
256 angCmd->SetUnitCategory("Angle");
258
259 msceCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/EnergyLimit",this);
260 msceCmd->SetGuidance("Set the upper energy limit for msc");
261 msceCmd->SetParameterName("mscE",true);
262 msceCmd->SetUnitCategory("Energy");
264
265 nielCmd = new G4UIcmdWithADoubleAndUnit("/process/em/MaxEnergyNIEL",this);
266 nielCmd->SetGuidance("Set the upper energy limit for NIEL");
267 nielCmd->SetParameterName("niel",true);
268 nielCmd->SetUnitCategory("Energy");
270
271 frCmd = new G4UIcmdWithADouble("/process/msc/RangeFactor",this);
272 frCmd->SetGuidance("Set RangeFactor for msc processes of e+-");
273 frCmd->SetParameterName("Fr",true);
274 frCmd->SetRange("Fr>0");
275 frCmd->SetDefaultValue(0.04);
277
278 fr1Cmd = new G4UIcmdWithADouble("/process/msc/RangeFactorMuHad",this);
279 fr1Cmd->SetGuidance("Set RangeFactor for msc processes of muons/hadrons");
280 fr1Cmd->SetParameterName("Fr1",true);
281 fr1Cmd->SetRange("Fr1>0");
282 fr1Cmd->SetDefaultValue(0.2);
284
285 fgCmd = new G4UIcmdWithADouble("/process/msc/GeomFactor",this);
286 fgCmd->SetGuidance("Set GeomFactor parameter for msc processes");
287 fgCmd->SetParameterName("Fg",true);
288 fgCmd->SetRange("Fg>0");
289 fgCmd->SetDefaultValue(3.5);
291
292 skinCmd = new G4UIcmdWithADouble("/process/msc/Skin",this);
293 skinCmd->SetGuidance("Set skin parameter for msc processes");
294 skinCmd->SetParameterName("skin",true);
296
297 screCmd = new G4UIcmdWithADouble("/process/msc/ScreeningFactor",this);
298 screCmd->SetGuidance("Set screening factor");
299 screCmd->SetParameterName("screen",true);
301
302 safCmd = new G4UIcmdWithADouble("/process/msc/SafetyFactor",this);
303 safCmd->SetGuidance("Set safety factor");
304 safCmd->SetParameterName("fsafe",true);
306
307 llimCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/LambdaLimit",this);
308 llimCmd->SetGuidance("Set the upper energy limit for NIEL");
309 llimCmd->SetParameterName("ll",true);
310 llimCmd->SetUnitCategory("Length");
312
313 dedxCmd = new G4UIcmdWithAnInteger("/process/eLoss/binsDEDX",this);
314 dedxCmd->SetGuidance("Set number of bins for EM tables");
315 dedxCmd->SetParameterName("binsDEDX",true);
316 dedxCmd->SetDefaultValue(84);
318
319 lamCmd = new G4UIcmdWithAnInteger("/process/eLoss/binsLambda",this);
320 lamCmd->SetGuidance("Set number of bins for EM tables");
321 lamCmd->SetParameterName("binsL",true);
322 lamCmd->SetDefaultValue(84);
324
325 amCmd = new G4UIcmdWithAnInteger("/process/eLoss/binsPerDecade",this);
326 amCmd->SetGuidance("Set number of bins per decade for EM tables");
327 amCmd->SetParameterName("bins",true);
328 amCmd->SetDefaultValue(7);
330
331 verCmd = new G4UIcmdWithAnInteger("/process/eLoss/verbose",this);
332 verCmd->SetGuidance("Set verbose level for EM physics");
333 verCmd->SetParameterName("verb",true);
334 verCmd->SetDefaultValue(1);
336
337 ver1Cmd = new G4UIcmdWithAnInteger("/process/em/verbose",this);
338 ver1Cmd->SetGuidance("Set verbose level for EM physics");
339 ver1Cmd->SetParameterName("verb1",true);
340 ver1Cmd->SetDefaultValue(1);
342
343 ver2Cmd = new G4UIcmdWithAnInteger("/process/em/workerVerbose",this);
344 ver2Cmd->SetGuidance("Set worker verbose level for EM physics");
345 ver2Cmd->SetParameterName("verb2",true);
346 ver2Cmd->SetDefaultValue(1);
348
349 mscCmd = new G4UIcmdWithAString("/process/msc/StepLimit",this);
350 mscCmd->SetGuidance("Set msc step limitation type");
351 mscCmd->SetParameterName("StepLim",true);
352 mscCmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
354
355 msc1Cmd = new G4UIcmdWithAString("/process/msc/StepLimitMuHad",this);
356 msc1Cmd->SetGuidance("Set msc step limitation type for muons/hadrons");
357 msc1Cmd->SetParameterName("StepLim1",true);
358 msc1Cmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
360
361 dumpCmd = new G4UIcommand("/process/em/printParameters",this);
362 dumpCmd->SetGuidance("Print all EM parameters.");
363
364 nffCmd = new G4UIcmdWithAString("/process/em/setNuclearFormFactor",this);
365 nffCmd->SetGuidance("Define type of nuclear form-factor");
366 nffCmd->SetParameterName("NucFF",true);
367 nffCmd->SetCandidates("None Exponential Gaussian Flat");
369
370 ssCmd = new G4UIcmdWithAString("/process/em/setSingleScattering",this);
371 ssCmd->SetGuidance("Define type of e+- single scattering model");
372 ssCmd->SetParameterName("SS",true);
373 ssCmd->SetCandidates("WVI Mott DPWA");
375
376 tripletCmd = new G4UIcmdWithAnInteger("/process/gconv/conversionType",this);
377 tripletCmd->SetGuidance("gamma conversion triplet/nuclear generation type:");
378 tripletCmd->SetGuidance("0 - (default) both triplet and nuclear");
379 tripletCmd->SetGuidance("1 - force nuclear");
380 tripletCmd->SetGuidance("2 - force triplet");
381 tripletCmd->SetParameterName("type",false);
382 tripletCmd->SetRange("type >= 0 && type <= 2");
383 tripletCmd->SetDefaultValue(0);
385
386 onIsolatedCmd = new G4UIcmdWithABool("/process/gconv/onIsolated",this);
387 onIsolatedCmd->SetGuidance("Conversion on isolated charged particles");
388 onIsolatedCmd->SetGuidance("false (default) : atomic electron screening");
389 onIsolatedCmd->SetGuidance("true : conversion on isolated particles.");
390 onIsolatedCmd->SetParameterName("flag",false);
391 onIsolatedCmd->SetDefaultValue(false);
393}
@ G4State_Init
@ G4State_Idle
@ G4State_PreInit
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(G4bool defVal)
void SetUnitCategory(const char *unitCategory)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
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)
void SetDefaultValue(G4int defVal)
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:273

◆ ~G4EmParametersMessenger()

G4EmParametersMessenger::~G4EmParametersMessenger ( )
virtual

Definition at line 397 of file G4EmParametersMessenger.cc.

398{
399 delete gconvDirectory;
400 delete eLossDirectory;
401 delete mscDirectory;
402 delete emDirectory;
403 delete dnaDirectory;
404
405 delete flucCmd;
406 delete rangeCmd;
407 delete lpmCmd;
408 delete splCmd;
409 delete rsCmd;
410 delete aplCmd;
411 delete latCmd;
412 delete lat96Cmd;
413 delete mulatCmd;
414 delete delCmd;
415 delete IntegCmd;
416 delete mottCmd;
417 delete birksCmd;
418 delete sharkCmd;
419 delete onIsolatedCmd;
420 delete sampleTCmd;
421 delete icru90Cmd;
422 delete mudatCmd;
423
424 delete minSubSecCmd;
425 delete minEnCmd;
426 delete maxEnCmd;
427 delete max5DCmd;
428 delete cenCmd;
429 delete lowEnCmd;
430 delete lowhEnCmd;
431 delete lowEn3Cmd;
432 delete lllCmd;
433 delete brCmd;
434 delete br1Cmd;
435 delete labCmd;
436 delete mscfCmd;
437 delete angCmd;
438 delete msceCmd;
439 delete nielCmd;
440 delete frCmd;
441 delete fr1Cmd;
442 delete fgCmd;
443 delete skinCmd;
444 delete safCmd;
445 delete llimCmd;
446 delete screCmd;
447
448 delete dedxCmd;
449 delete lamCmd;
450 delete amCmd;
451 delete verCmd;
452 delete ver1Cmd;
453 delete ver2Cmd;
454 delete tripletCmd;
455
456 delete mscCmd;
457 delete msc1Cmd;
458 delete nffCmd;
459 delete ssCmd;
460
461 delete dumpCmd;
462}

Member Function Documentation

◆ SetNewValue()

void G4EmParametersMessenger::SetNewValue ( G4UIcommand command,
G4String  newValue 
)
overridevirtual

Reimplemented from G4UImessenger.

Definition at line 466 of file G4EmParametersMessenger.cc.

468{
469 G4bool physicsModified = false;
470 if (command == flucCmd) {
471 theParameters->SetLossFluctuations(flucCmd->GetNewBoolValue(newValue));
472 physicsModified = true;
473 } else if (command == rangeCmd) {
474 theParameters->SetBuildCSDARange(rangeCmd->GetNewBoolValue(newValue));
475 } else if (command == lpmCmd) {
476 theParameters->SetLPM(lpmCmd->GetNewBoolValue(newValue));
477 physicsModified = true;
478 } else if (command == splCmd) {
479 theParameters->SetSpline(splCmd->GetNewBoolValue(newValue));
480 } else if (command == rsCmd) {
481 theParameters->SetUseCutAsFinalRange(rsCmd->GetNewBoolValue(newValue));
482 physicsModified = true;
483 } else if (command == aplCmd) {
484 theParameters->SetApplyCuts(aplCmd->GetNewBoolValue(newValue));
485 physicsModified = true;
486 } else if (command == latCmd) {
487 theParameters->SetLateralDisplacement(latCmd->GetNewBoolValue(newValue));
488 physicsModified = true;
489 } else if (command == lat96Cmd) {
490 theParameters->SetLateralDisplacementAlg96(lat96Cmd->GetNewBoolValue(newValue));
491 physicsModified = true;
492 } else if (command == mulatCmd) {
493 theParameters->SetMuHadLateralDisplacement(mulatCmd->GetNewBoolValue(newValue));
494 physicsModified = true;
495 } else if (command == delCmd) {
496 theParameters->ActivateAngularGeneratorForIonisation(delCmd->GetNewBoolValue(newValue));
497 } else if (command == IntegCmd) {
498 theParameters->SetIntegral(IntegCmd->GetNewBoolValue(newValue));
499 physicsModified = true;
500 } else if (command == mottCmd) {
501 theParameters->SetUseMottCorrection(mottCmd->GetNewBoolValue(newValue));
502 } else if (command == birksCmd) {
503 theParameters->SetBirksActive(birksCmd->GetNewBoolValue(newValue));
504 } else if (command == icru90Cmd) {
505 theParameters->SetUseICRU90Data(icru90Cmd->GetNewBoolValue(newValue));
506 } else if (command == sharkCmd) {
507 theParameters->SetGeneralProcessActive(sharkCmd->GetNewBoolValue(newValue));
508 } else if (command == sampleTCmd) {
509 theParameters->SetEnableSamplingTable(sampleTCmd->GetNewBoolValue(newValue));
510 } else if (command == mudatCmd) {
511 theParameters->SetRetrieveMuDataFromFile(mudatCmd->GetNewBoolValue(newValue));
512
513 } else if (command == minSubSecCmd) {
514 theParameters->SetMinSubRange(minSubSecCmd->GetNewDoubleValue(newValue));
515 } else if (command == minEnCmd) {
516 theParameters->SetMinEnergy(minEnCmd->GetNewDoubleValue(newValue));
517 } else if (command == maxEnCmd) {
518 theParameters->SetMaxEnergy(maxEnCmd->GetNewDoubleValue(newValue));
519 } else if (command == max5DCmd) {
520 theParameters->SetMaxEnergyFor5DMuPair(max5DCmd->GetNewDoubleValue(newValue));
521 } else if (command == cenCmd) {
522 theParameters->SetMaxEnergyForCSDARange(cenCmd->GetNewDoubleValue(newValue));
523 physicsModified = true;
524 } else if (command == lowEnCmd) {
525 theParameters->SetLowestElectronEnergy(lowEnCmd->GetNewDoubleValue(newValue));
526 physicsModified = true;
527 } else if (command == lowEn3Cmd) {
528 theParameters->SetLowestTripletEnergy(lowEn3Cmd->GetNewDoubleValue(newValue));
529 physicsModified = true;
530 } else if (command == lowhEnCmd) {
531 theParameters->SetLowestMuHadEnergy(lowhEnCmd->GetNewDoubleValue(newValue));
532 physicsModified = true;
533 } else if (command == lllCmd) {
534 theParameters->SetLinearLossLimit(lllCmd->GetNewDoubleValue(newValue));
535 physicsModified = true;
536 } else if (command == brCmd) {
537 theParameters->SetBremsstrahlungTh(brCmd->GetNewDoubleValue(newValue));
538 physicsModified = true;
539 } else if (command == br1Cmd) {
540 theParameters->SetMuHadBremsstrahlungTh(br1Cmd->GetNewDoubleValue(newValue));
541 physicsModified = true;
542 } else if (command == labCmd) {
543 theParameters->SetLambdaFactor(labCmd->GetNewDoubleValue(newValue));
544 physicsModified = true;
545 } else if (command == mscfCmd) {
546 theParameters->SetFactorForAngleLimit(mscfCmd->GetNewDoubleValue(newValue));
547 } else if (command == angCmd) {
548 theParameters->SetMscThetaLimit(angCmd->GetNewDoubleValue(newValue));
549 physicsModified = true;
550 } else if (command == msceCmd) {
551 theParameters->SetMscEnergyLimit(msceCmd->GetNewDoubleValue(newValue));
552 } else if (command == nielCmd) {
553 theParameters->SetMaxNIELEnergy(nielCmd->GetNewDoubleValue(newValue));
554 } else if (command == frCmd) {
555 theParameters->SetMscRangeFactor(frCmd->GetNewDoubleValue(newValue));
556 physicsModified = true;
557 } else if (command == fr1Cmd) {
558 theParameters->SetMscMuHadRangeFactor(fr1Cmd->GetNewDoubleValue(newValue));
559 physicsModified = true;
560 } else if (command == fgCmd) {
561 theParameters->SetMscGeomFactor(fgCmd->GetNewDoubleValue(newValue));
562 physicsModified = true;
563 } else if (command == skinCmd) {
564 theParameters->SetMscSkin(skinCmd->GetNewDoubleValue(newValue));
565 physicsModified = true;
566 } else if (command == safCmd) {
567 theParameters->SetMscSafetyFactor(safCmd->GetNewDoubleValue(newValue));
568 physicsModified = true;
569 } else if (command == llimCmd) {
570 theParameters->SetMscLambdaLimit(llimCmd->GetNewDoubleValue(newValue));
571 physicsModified = true;
572 } else if (command == screCmd) {
573 theParameters->SetScreeningFactor(screCmd->GetNewDoubleValue(newValue));
574 } else if (command == dedxCmd) {
575 theParameters->SetNumberOfBins(dedxCmd->GetNewIntValue(newValue));
576 } else if (command == lamCmd) {
577 theParameters->SetNumberOfBins(lamCmd->GetNewIntValue(newValue));
578 } else if (command == amCmd) {
579 theParameters->SetNumberOfBinsPerDecade(amCmd->GetNewIntValue(newValue));
580 } else if (command == verCmd) {
581 theParameters->SetVerbose(verCmd->GetNewIntValue(newValue));
582 } else if (command == ver1Cmd) {
583 theParameters->SetVerbose(ver1Cmd->GetNewIntValue(newValue));
584 physicsModified = true;
585 } else if (command == ver2Cmd) {
586 theParameters->SetWorkerVerbose(ver2Cmd->GetNewIntValue(newValue));
587 physicsModified = true;
588 } else if (command == dumpCmd) {
589 theParameters->Dump();
590 } else if (command == mscCmd || command == msc1Cmd) {
592 if(newValue == "Minimal") {
593 msctype = fMinimal;
594 } else if(newValue == "UseDistanceToBoundary") {
595 msctype = fUseDistanceToBoundary;
596 } else if(newValue == "UseSafety") {
597 msctype = fUseSafety;
598 } else if(newValue == "UseSafetyPlus") {
599 msctype = fUseSafetyPlus;
600 } else {
602 ed << " StepLimit type <" << newValue << "> unknown!";
603 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
604 return;
605 }
606 if (command == mscCmd) {
607 theParameters->SetMscStepLimitType(msctype);
608 } else {
609 theParameters->SetMscMuHadStepLimitType(msctype);
610 }
611 physicsModified = true;
612 } else if (command == nffCmd) {
614 if(newValue == "Exponential") { x = fExponentialNF; }
615 else if(newValue == "Gaussian") { x = fGaussianNF; }
616 else if(newValue == "Flat") { x = fFlatNF; }
617 else if(newValue != "None") {
619 ed << " NuclearFormFactor type <" << newValue << "> unknown!";
620 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
621 return;
622 }
623 theParameters->SetNuclearFormfactorType(x);
624 } else if (command == ssCmd) {
626 if(newValue == "DPWA") { x = fDPWA; }
627 else if(newValue == "Mott") { x = fMott; }
628 else if(newValue != "WVI") {
630 ed << " G4eSingleScatteringType type <" << newValue << "> unknown!";
631 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
632 return;
633 }
634 theParameters->SetSingleScatteringType(x);
635 } else if ( command==tripletCmd ) {
636 theParameters->SetConversionType(tripletCmd->GetNewIntValue(newValue));
637 physicsModified = true;
638 } else if ( command==onIsolatedCmd ) {
639 theParameters->SetOnIsolated(onIsolatedCmd->GetNewBoolValue(newValue));
640 physicsModified = true;
641 }
642
643 if(physicsModified) {
644 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
645 }
646}
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
void SetLambdaFactor(G4double val)
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetBuildCSDARange(G4bool val)
void Dump() const
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 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 SetConversionType(G4int val)
void SetUseICRU90Data(G4bool val)
void SetOnIsolated(G4bool val)
void SetIntegral(G4bool val)
void SetUseMottCorrection(G4bool val)
void SetLowestMuHadEnergy(G4double val)
void SetMaxEnergy(G4double val)
void SetNumberOfBins(G4int val)
void SetSpline(G4bool val)
void SetMinSubRange(G4double val)
void SetSingleScatteringType(G4eSingleScatteringType val)
void SetMscRangeFactor(G4double val)
static G4bool GetNewBoolValue(const char *paramString)
static G4double GetNewDoubleValue(const char *paramString)
static G4double GetNewDoubleValue(const char *paramString)
static G4int GetNewIntValue(const char *paramString)
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:485
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77

The documentation for this class was generated from the following files: