Geant4 11.2.2
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 *)
 
 ~G4EmParametersMessenger () override
 
void SetNewValue (G4UIcommand *, G4String) override
 
G4EmParametersMessengeroperator= (const G4EmParametersMessenger &right)=delete
 
 G4EmParametersMessenger (const G4EmParametersMessenger &)=delete
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()=default
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
G4bool CommandsShouldBeInMaster () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String LtoS (G4long l)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (const G4String &s)
 
G4long StoL (const G4String &s)
 
G4double StoD (const 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() [1/2]

G4EmParametersMessenger::G4EmParametersMessenger ( G4EmParameters * ptr)
explicit

Definition at line 61 of file G4EmParametersMessenger.cc.

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 nFreeCmd = new G4UIcmdWithAnInteger("/process/em/nForFreeVector",this);
388 nFreeCmd->SetGuidance("Set number for logarithmic bin search algorithm");
389 nFreeCmd->SetParameterName("nFree",true);
390 nFreeCmd->SetDefaultValue(2);
392 nFreeCmd->SetToBeBroadcasted(false);
393
394 transWithMscCmd = new G4UIcmdWithAString("/process/em/transportationWithMsc",this);
395 transWithMscCmd->SetGuidance("Enable/disable the G4TransportationWithMsc process");
396 transWithMscCmd->SetParameterName("trans",true);
397 transWithMscCmd->SetCandidates("Disabled Enabled MultipleSteps");
398 transWithMscCmd->AvailableForStates(G4State_PreInit);
399 transWithMscCmd->SetToBeBroadcasted(false);
400
401 mscCmd = new G4UIcmdWithAString("/process/msc/StepLimit",this);
402 mscCmd->SetGuidance("Set msc step limitation type");
403 mscCmd->SetParameterName("StepLim",true);
404 mscCmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
406 mscCmd->SetToBeBroadcasted(false);
407
408 msc1Cmd = new G4UIcmdWithAString("/process/msc/StepLimitMuHad",this);
409 msc1Cmd->SetGuidance("Set msc step limitation type for muons/hadrons");
410 msc1Cmd->SetParameterName("StepLim1",true);
411 msc1Cmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
413 msc1Cmd->SetToBeBroadcasted(false);
414
415 dumpCmd = new G4UIcommand("/process/em/printParameters",this);
416 dumpCmd->SetGuidance("Print all EM parameters.");
418 dumpCmd->SetToBeBroadcasted(false);
419
420 nffCmd = new G4UIcmdWithAString("/process/em/setNuclearFormFactor",this);
421 nffCmd->SetGuidance("Define type of nuclear form-factor");
422 nffCmd->SetParameterName("NucFF",true);
423 nffCmd->SetCandidates("None Exponential Gaussian Flat");
425 nffCmd->SetToBeBroadcasted(false);
426
427 ssCmd = new G4UIcmdWithAString("/process/em/setSingleScattering",this);
428 ssCmd->SetGuidance("Define type of e+- single scattering model");
429 ssCmd->SetParameterName("SS",true);
430 ssCmd->SetCandidates("WVI Mott DPWA");
432 ssCmd->SetToBeBroadcasted(false);
433
434 fluc1Cmd = new G4UIcmdWithAString("/process/eLoss/setFluctModel",this);
435 fluc1Cmd->SetGuidance("Define type of energy loss fluctuation model");
436 fluc1Cmd->SetParameterName("Fluc1",true);
437 fluc1Cmd->SetCandidates("Dummy Universal Urban");
439 fluc1Cmd->SetToBeBroadcasted(false);
440
441 tripletCmd = new G4UIcmdWithAnInteger("/process/gconv/conversionType",this);
442 tripletCmd->SetGuidance("gamma conversion triplet/nuclear generation type:");
443 tripletCmd->SetGuidance("0 - (default) both triplet and nuclear");
444 tripletCmd->SetGuidance("1 - force nuclear");
445 tripletCmd->SetGuidance("2 - force triplet");
446 tripletCmd->SetParameterName("type",false);
447 tripletCmd->SetRange("type >= 0 && type <= 2");
448 tripletCmd->SetDefaultValue(0);
450 tripletCmd->SetToBeBroadcasted(false);
451
452 onIsolatedCmd = new G4UIcmdWithABool("/process/gconv/onIsolated",this);
453 onIsolatedCmd->SetGuidance("Conversion on isolated charged particles");
454 onIsolatedCmd->SetGuidance("false (default) : atomic electron screening");
455 onIsolatedCmd->SetGuidance("true : conversion on isolated particles.");
456 onIsolatedCmd->SetParameterName("flag",false);
457 onIsolatedCmd->SetDefaultValue(false);
459 onIsolatedCmd->SetToBeBroadcasted(false);
460}
@ 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 SetToBeBroadcasted(G4bool val)
void SetGuidance(const char *aGuidance)
void SetRange(const char *rs)
void AvailableForStates(G4ApplicationState s1)

◆ ~G4EmParametersMessenger()

G4EmParametersMessenger::~G4EmParametersMessenger ( )
override

Definition at line 464 of file G4EmParametersMessenger.cc.

465{
466 delete gconvDirectory;
467 delete eLossDirectory;
468 delete mscDirectory;
469 delete emDirectory;
470 delete dnaDirectory;
471
472 delete flucCmd;
473 delete rangeCmd;
474 delete lpmCmd;
475 delete rsCmd;
476 delete aplCmd;
477 delete intCmd;
478 delete latCmd;
479 delete lat96Cmd;
480 delete mulatCmd;
481 delete delCmd;
482 delete mottCmd;
483 delete birksCmd;
484 delete sharkCmd;
485 delete onIsolatedCmd;
486 delete sampleTCmd;
487 delete poCmd;
488 delete icru90Cmd;
489 delete mudatCmd;
490 delete peKCmd;
491 delete mscPCmd;
492
493 delete minEnCmd;
494 delete maxEnCmd;
495 delete max5DCmd;
496 delete cenCmd;
497 delete lowEnCmd;
498 delete lowhEnCmd;
499 delete lowEn3Cmd;
500 delete lllCmd;
501 delete brCmd;
502 delete br1Cmd;
503 delete labCmd;
504 delete mscfCmd;
505 delete angCmd;
506 delete msceCmd;
507 delete nielCmd;
508 delete frCmd;
509 delete fr1Cmd;
510 delete fgCmd;
511 delete skinCmd;
512 delete safCmd;
513 delete llimCmd;
514 delete screCmd;
515
516 delete amCmd;
517 delete verCmd;
518 delete ver1Cmd;
519 delete ver2Cmd;
520 delete transWithMscCmd;
521 delete nFreeCmd;
522 delete tripletCmd;
523
524 delete mscCmd;
525 delete msc1Cmd;
526 delete nffCmd;
527 delete ssCmd;
528 delete fluc1Cmd;
529
530 delete dumpCmd;
531}

◆ G4EmParametersMessenger() [2/2]

G4EmParametersMessenger::G4EmParametersMessenger ( const G4EmParametersMessenger & )
delete

Member Function Documentation

◆ operator=()

G4EmParametersMessenger & G4EmParametersMessenger::operator= ( const G4EmParametersMessenger & right)
delete

◆ SetNewValue()

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

Reimplemented from G4UImessenger.

Definition at line 535 of file G4EmParametersMessenger.cc.

537{
538 G4bool physicsModified = false;
539 if (command == flucCmd) {
540 theParameters->SetLossFluctuations(flucCmd->GetNewBoolValue(newValue));
541 physicsModified = true;
542 } else if (command == rangeCmd) {
543 theParameters->SetBuildCSDARange(rangeCmd->GetNewBoolValue(newValue));
544 } else if (command == lpmCmd) {
545 theParameters->SetLPM(lpmCmd->GetNewBoolValue(newValue));
546 physicsModified = true;
547 } else if (command == rsCmd) {
548 theParameters->SetUseCutAsFinalRange(rsCmd->GetNewBoolValue(newValue));
549 physicsModified = true;
550 } else if (command == aplCmd) {
551 theParameters->SetApplyCuts(aplCmd->GetNewBoolValue(newValue));
552 physicsModified = true;
553 } else if (command == intCmd) {
554 theParameters->SetIntegral(intCmd->GetNewBoolValue(newValue));
555 } else if (command == latCmd) {
556 theParameters->SetLateralDisplacement(latCmd->GetNewBoolValue(newValue));
557 physicsModified = true;
558 } else if (command == lat96Cmd) {
559 theParameters->SetLateralDisplacementAlg96(lat96Cmd->GetNewBoolValue(newValue));
560 physicsModified = true;
561 } else if (command == mulatCmd) {
562 theParameters->SetMuHadLateralDisplacement(mulatCmd->GetNewBoolValue(newValue));
563 physicsModified = true;
564 } else if (command == delCmd) {
565 theParameters->ActivateAngularGeneratorForIonisation(delCmd->GetNewBoolValue(newValue));
566 } else if (command == mottCmd) {
567 theParameters->SetUseMottCorrection(mottCmd->GetNewBoolValue(newValue));
568 } else if (command == birksCmd) {
569 theParameters->SetBirksActive(birksCmd->GetNewBoolValue(newValue));
570 } else if (command == icru90Cmd) {
571 theParameters->SetUseICRU90Data(icru90Cmd->GetNewBoolValue(newValue));
572 } else if (command == sharkCmd) {
573 theParameters->SetGeneralProcessActive(sharkCmd->GetNewBoolValue(newValue));
574 } else if (command == poCmd) {
575 theParameters->SetEnablePolarisation(poCmd->GetNewBoolValue(newValue));
576 } else if (command == sampleTCmd) {
577 theParameters->SetEnableSamplingTable(sampleTCmd->GetNewBoolValue(newValue));
578 } else if (command == mudatCmd) {
579 theParameters->SetRetrieveMuDataFromFile(mudatCmd->GetNewBoolValue(newValue));
580 } else if (command == peKCmd) {
581 theParameters->SetPhotoeffectBelowKShell(peKCmd->GetNewBoolValue(newValue));
582 } else if (command == mscPCmd) {
583 theParameters->SetMscPositronCorrection(mscPCmd->GetNewBoolValue(newValue));
584
585 } else if (command == minEnCmd) {
586 theParameters->SetMinEnergy(minEnCmd->GetNewDoubleValue(newValue));
587 } else if (command == maxEnCmd) {
588 theParameters->SetMaxEnergy(maxEnCmd->GetNewDoubleValue(newValue));
589 } else if (command == max5DCmd) {
590 theParameters->SetMaxEnergyFor5DMuPair(max5DCmd->GetNewDoubleValue(newValue));
591 } else if (command == cenCmd) {
592 theParameters->SetMaxEnergyForCSDARange(cenCmd->GetNewDoubleValue(newValue));
593 physicsModified = true;
594 } else if (command == lowEnCmd) {
595 theParameters->SetLowestElectronEnergy(lowEnCmd->GetNewDoubleValue(newValue));
596 physicsModified = true;
597 } else if (command == lowEn3Cmd) {
598 theParameters->SetLowestTripletEnergy(lowEn3Cmd->GetNewDoubleValue(newValue));
599 physicsModified = true;
600 } else if (command == lowhEnCmd) {
601 theParameters->SetLowestMuHadEnergy(lowhEnCmd->GetNewDoubleValue(newValue));
602 physicsModified = true;
603 } else if (command == lllCmd) {
604 theParameters->SetLinearLossLimit(lllCmd->GetNewDoubleValue(newValue));
605 physicsModified = true;
606 } else if (command == brCmd) {
607 theParameters->SetBremsstrahlungTh(brCmd->GetNewDoubleValue(newValue));
608 physicsModified = true;
609 } else if (command == br1Cmd) {
610 theParameters->SetMuHadBremsstrahlungTh(br1Cmd->GetNewDoubleValue(newValue));
611 physicsModified = true;
612 } else if (command == labCmd) {
613 theParameters->SetLambdaFactor(labCmd->GetNewDoubleValue(newValue));
614 physicsModified = true;
615 } else if (command == mscfCmd) {
616 theParameters->SetFactorForAngleLimit(mscfCmd->GetNewDoubleValue(newValue));
617 } else if (command == angCmd) {
618 theParameters->SetMscThetaLimit(angCmd->GetNewDoubleValue(newValue));
619 } else if (command == msceCmd) {
620 theParameters->SetMscEnergyLimit(msceCmd->GetNewDoubleValue(newValue));
621 } else if (command == nielCmd) {
622 theParameters->SetMaxNIELEnergy(nielCmd->GetNewDoubleValue(newValue));
623 } else if (command == frCmd) {
624 theParameters->SetMscRangeFactor(frCmd->GetNewDoubleValue(newValue));
625 physicsModified = true;
626 } else if (command == fr1Cmd) {
627 theParameters->SetMscMuHadRangeFactor(fr1Cmd->GetNewDoubleValue(newValue));
628 physicsModified = true;
629 } else if (command == fgCmd) {
630 theParameters->SetMscGeomFactor(fgCmd->GetNewDoubleValue(newValue));
631 physicsModified = true;
632 } else if (command == skinCmd) {
633 theParameters->SetMscSkin(skinCmd->GetNewDoubleValue(newValue));
634 physicsModified = true;
635 } else if (command == safCmd) {
636 theParameters->SetMscSafetyFactor(safCmd->GetNewDoubleValue(newValue));
637 } else if (command == llimCmd) {
638 theParameters->SetMscLambdaLimit(llimCmd->GetNewDoubleValue(newValue));
639 } else if (command == screCmd) {
640 theParameters->SetScreeningFactor(screCmd->GetNewDoubleValue(newValue));
641 } else if (command == amCmd) {
642 theParameters->SetNumberOfBinsPerDecade(amCmd->GetNewIntValue(newValue));
643 } else if (command == verCmd) {
644 theParameters->SetVerbose(verCmd->GetNewIntValue(newValue));
645 } else if (command == ver1Cmd) {
646 theParameters->SetVerbose(ver1Cmd->GetNewIntValue(newValue));
647 } else if (command == ver2Cmd) {
648 theParameters->SetWorkerVerbose(ver2Cmd->GetNewIntValue(newValue));
649 } else if (command == nFreeCmd) {
650 theParameters->SetNumberForFreeVector(nFreeCmd->GetNewIntValue(newValue));
651 } else if (command == dumpCmd) {
652 theParameters->SetIsPrintedFlag(false);
653 theParameters->Dump();
654 } else if (command == transWithMscCmd) {
656 if(newValue == "Disabled") {
658 } else if(newValue == "Enabled") {
660 } else if(newValue == "MultipleSteps") {
662 } else {
664 ed << " TransportationWithMsc type <" << newValue << "> unknown!";
665 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
666 }
667 theParameters->SetTransportationWithMsc(type);
668 } else if (command == mscCmd || command == msc1Cmd) {
670 if(newValue == "Minimal") {
671 msctype = fMinimal;
672 } else if(newValue == "UseDistanceToBoundary") {
673 msctype = fUseDistanceToBoundary;
674 } else if(newValue == "UseSafety") {
675 msctype = fUseSafety;
676 } else if(newValue == "UseSafetyPlus") {
677 msctype = fUseSafetyPlus;
678 } else {
680 ed << " StepLimit type <" << newValue << "> unknown!";
681 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
682 return;
683 }
684 if (command == mscCmd) {
685 theParameters->SetMscStepLimitType(msctype);
686 } else {
687 theParameters->SetMscMuHadStepLimitType(msctype);
688 }
689 physicsModified = true;
690 } else if (command == nffCmd) {
692 if(newValue == "Exponential") { x = fExponentialNF; }
693 else if(newValue == "Gaussian") { x = fGaussianNF; }
694 else if(newValue == "Flat") { x = fFlatNF; }
695 else if(newValue != "None") {
697 ed << " NuclearFormFactor type <" << newValue << "> unknown!";
698 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
699 return;
700 }
701 theParameters->SetNuclearFormfactorType(x);
702 } else if (command == ssCmd) {
704 if(newValue == "DPWA") { x = fDPWA; }
705 else if(newValue == "Mott") { x = fMott; }
706 else if(newValue != "WVI") {
708 ed << " G4eSingleScatteringType type <" << newValue << "> unknown!";
709 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed);
710 return;
711 }
712 theParameters->SetSingleScatteringType(x);
713 } else if (command == fluc1Cmd) {
715 if(newValue == "Dummy") { x = fDummyFluctuation; }
716 else if(newValue == "Urban") { x = fUrbanFluctuation; }
717 theParameters->SetFluctuationType(x);
718 } else if ( command==tripletCmd ) {
719 theParameters->SetConversionType(tripletCmd->GetNewIntValue(newValue));
720 } else if ( command==onIsolatedCmd ) {
721 theParameters->SetOnIsolated(onIsolatedCmd->GetNewBoolValue(newValue));
722 physicsModified = true;
723 }
724
725 if(physicsModified) {
726 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
727 }
728}
G4EmFluctuationType
@ fUniversalFluctuation
@ fUrbanFluctuation
@ fDummyFluctuation
G4eSingleScatteringType
@ fWVI
@ fMott
@ fDPWA
G4TransportationWithMscType
@ 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 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 SetNumberForFreeVector(G4int 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)
static G4double GetNewDoubleValue(const char *paramString)
static G4double GetNewDoubleValue(const char *paramString)
static G4int GetNewIntValue(const char *paramString)
G4int ApplyCommand(const char *aCommand)
static G4UImanager * GetUIpointer()

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