Geant4  10.02.p02
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 // $Id: G4EmParametersMessenger.cc 66241 2012-12-13 18:34:42Z gunter $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 // File name: G4EmParametersMessenger
33 //
34 // Author: Vladimir Ivanchenko created from G4EnergyLossMessenger
35 //
36 // Creation date: 22-05-2013
37 //
38 // Modifications:
39 //
40 // -------------------------------------------------------------------
41 //
42 
43 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45 
47 
48 #include "G4UIdirectory.hh"
49 #include "G4UIcommand.hh"
50 #include "G4UIparameter.hh"
51 #include "G4UIcmdWithABool.hh"
52 #include "G4UIcmdWithAnInteger.hh"
53 #include "G4UIcmdWithADouble.hh"
55 #include "G4UIcmdWithAString.hh"
56 #include "G4UImanager.hh"
57 #include "G4MscStepLimitType.hh"
58 #include "G4EmParameters.hh"
59 
60 #include <sstream>
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
63 
65  : theParameters(ptr)
66 {
67  eLossDirectory = new G4UIdirectory("/process/eLoss/");
68  eLossDirectory->SetGuidance("Commands for EM processes.");
69  mscDirectory = new G4UIdirectory("/process/msc/");
70  mscDirectory->SetGuidance("Commands for EM scattering processes.");
71  emDirectory = new G4UIdirectory("/process/em/");
72  emDirectory->SetGuidance("General commands for EM processes.");
73 
74  flucCmd = new G4UIcmdWithABool("/process/eLoss/fluct",this);
75  flucCmd->SetGuidance("Enable/disable energy loss fluctuations.");
76  flucCmd->SetParameterName("choice",true);
77  flucCmd->SetDefaultValue(true);
79 
80  rangeCmd = new G4UIcmdWithABool("/process/eLoss/CSDARange",this);
81  rangeCmd->SetGuidance("Enable/disable CSDA range calculation");
82  rangeCmd->SetParameterName("range",true);
83  rangeCmd->SetDefaultValue(false);
85 
86  lpmCmd = new G4UIcmdWithABool("/process/eLoss/LPM",this);
87  lpmCmd->SetGuidance("Enable/disable LPM effect calculation");
88  lpmCmd->SetParameterName("lpm",true);
89  lpmCmd->SetDefaultValue(true);
91 
92  splCmd = new G4UIcmdWithABool("/process/em/spline",this);
93  splCmd->SetGuidance("Enable/disable usage spline for Physics Vectors");
94  splCmd->SetParameterName("spl",true);
95  splCmd->SetDefaultValue(false);
97 
98  rsCmd = new G4UIcmdWithABool("/process/eLoss/useCutAsFinalRange",this);
99  rsCmd->SetGuidance("Enable?disable use of cut in range as a final range");
100  rsCmd->SetParameterName("choice",true);
101  rsCmd->SetDefaultValue(false);
103 
104  aplCmd = new G4UIcmdWithABool("/process/em/applyCuts",this);
105  aplCmd->SetGuidance("Enable/disable applying cuts for gamma processes");
106  aplCmd->SetParameterName("apl",true);
107  aplCmd->SetDefaultValue(false);
109 
110  deCmd = new G4UIcmdWithABool("/process/em/fluo",this);
111  deCmd->SetGuidance("Enable/disable atomic deexcitation");
112  deCmd->SetParameterName("fluoFlag",true);
113  deCmd->SetDefaultValue(false);
115 
116  dirFluoCmd = new G4UIcmdWithABool("/process/em/fluoBearden",this);
117  dirFluoCmd->SetGuidance("Enable/disable usage of Bearden fluorescence files");
118  dirFluoCmd->SetParameterName("fluoBeardenFlag",true);
119  dirFluoCmd->SetDefaultValue(false);
121 
122  auCmd = new G4UIcmdWithABool("/process/em/auger",this);
123  auCmd->SetGuidance("Enable/disable Auger electrons production");
124  auCmd->SetParameterName("augerFlag",true);
125  auCmd->SetDefaultValue(false);
127 
128  auCascadeCmd = new G4UIcmdWithABool("/process/em/augerCascade",this);
129  auCascadeCmd->SetGuidance("Enable/disable simulation of cascade of Auger electrons");
130  auCascadeCmd->SetParameterName("augerCascadeFlag",true);
133 
134  pixeCmd = new G4UIcmdWithABool("/process/em/pixe",this);
135  pixeCmd->SetGuidance("Enable/disable PIXE simulation");
136  pixeCmd->SetParameterName("pixeFlag",true);
137  pixeCmd->SetDefaultValue(false);
139 
140  dcutCmd = new G4UIcmdWithABool("/process/em/deexcitationIgnoreCut",this);
141  dcutCmd->SetGuidance("Enable/Disable usage of cuts in de-excitation module");
142  dcutCmd->SetParameterName("deexcut",true);
143  dcutCmd->SetDefaultValue(false);
145 
146  latCmd = new G4UIcmdWithABool("/process/msc/LateralDisplacement",this);
147  latCmd->SetGuidance("Enable/disable sampling of lateral displacement");
148  latCmd->SetParameterName("lat",true);
149  latCmd->SetDefaultValue(true);
151 
152  mulatCmd = new G4UIcmdWithABool("/process/msc/MuHadLateralDisplacement",this);
153  mulatCmd->SetGuidance("Enable/disable sampling of lateral displacement for muons and hadrons");
154  mulatCmd->SetParameterName("mulat",true);
155  mulatCmd->SetDefaultValue(true);
157 
158  catCmd = new G4UIcmdWithABool("/process/msc/DisplacementBeyondSafety",this);
159  catCmd->SetGuidance("Enable/disable displacement at geometry boundary");
160  catCmd->SetParameterName("cat",true);
161  catCmd->SetDefaultValue(false);
163 
164  delCmd = new G4UIcmdWithABool("/process/eLoss/UseAngularGenerator",this);
165  delCmd->SetGuidance("Enable usage of angular generator");
166  delCmd->SetParameterName("del",true);
167  delCmd->SetDefaultValue(false);
169 
170  mottCmd = new G4UIcmdWithABool("/process/msc/UseMottCorrection",this);
171  mottCmd->SetGuidance("Enable usage of Mott corrections for e- elastic scattering");
172  mottCmd->SetParameterName("mott",true);
173  mottCmd->SetDefaultValue(false);
175 
176  minSubSecCmd = new G4UIcmdWithADouble("/process/eLoss/minsubsec",this);
177  minSubSecCmd->SetGuidance("Set the ratio subcut/cut ");
178  minSubSecCmd->SetParameterName("rcmin",true);
180 
181  minEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/minKinEnergy",this);
182  minEnCmd->SetGuidance("Set the min kinetic energy for EM tables");
183  minEnCmd->SetParameterName("emin",true);
184  minEnCmd->SetUnitCategory("Energy");
186 
187  maxEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergy",this);
188  maxEnCmd->SetGuidance("Set the max kinetic energy for EM tables");
189  maxEnCmd->SetParameterName("emax",true);
190  maxEnCmd->SetUnitCategory("Energy");
192 
193  cenCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergyCSDA",this);
194  cenCmd->SetGuidance("Set the max kinetic energy for CSDA table");
195  cenCmd->SetParameterName("emaxCSDA",true);
196  cenCmd->SetUnitCategory("Energy");
198 
199  lowEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestElectronEnergy",this);
200  lowEnCmd->SetGuidance("Set the lowest kinetic energy for e+-");
201  lowEnCmd->SetParameterName("elow",true);
202  lowEnCmd->SetUnitCategory("Energy");
204 
205  lowhEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestMuHadEnergy",this);
206  lowhEnCmd->SetGuidance("Set the lowest kinetic energy for muons and hadrons");
207  lowhEnCmd->SetParameterName("elowh",true);
208  lowhEnCmd->SetUnitCategory("Energy");
210 
211  lllCmd = new G4UIcmdWithADouble("/process/eLoss/linLossLimit",this);
212  lllCmd->SetGuidance("Set linearLossLimit parameter");
213  lllCmd->SetParameterName("linlim",true);
215 
216  brCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/bremThreshold",this);
217  brCmd->SetGuidance("Set bremsstrahlung energy threshold");
218  brCmd->SetParameterName("emaxBrem",true);
219  brCmd->SetUnitCategory("Energy");
221 
222  labCmd = new G4UIcmdWithADouble("/process/eLoss/LambdaFactor",this);
223  labCmd->SetGuidance("Set lambdaFactor parameter for integral option");
224  labCmd->SetParameterName("Fl",true);
226 
227  mscfCmd = new G4UIcmdWithADouble("/process/msc/FactorForAngleLimit",this);
228  mscfCmd->SetGuidance("Set factor for computation of a limit for -t (invariant trasfer)");
229  mscfCmd->SetParameterName("Fact",true);
230  mscfCmd->SetRange("Fact>0");
233 
234  angCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/ThetaLimit",this);
235  angCmd->SetGuidance("Set the limit on the polar angle for msc and single scattering");
236  angCmd->SetParameterName("theta",true);
237  angCmd->SetUnitCategory("Angle");
239 
240  frCmd = new G4UIcmdWithADouble("/process/msc/RangeFactor",this);
241  frCmd->SetGuidance("Set RangeFactor for msc processes of e+-");
242  frCmd->SetParameterName("Fr",true);
243  frCmd->SetRange("Fr>0");
244  frCmd->SetDefaultValue(0.04);
246 
247  fr1Cmd = new G4UIcmdWithADouble("/process/msc/RangeFactorMuHad",this);
248  fr1Cmd->SetGuidance("Set RangeFactor for msc processes of muons/hadrons");
249  fr1Cmd->SetParameterName("Fr1",true);
250  fr1Cmd->SetRange("Fr>0");
251  fr1Cmd->SetDefaultValue(0.2);
253 
254  fgCmd = new G4UIcmdWithADouble("/process/msc/GeomFactor",this);
255  fgCmd->SetGuidance("Set GeomFactor parameter for msc processes");
256  fgCmd->SetParameterName("Fg",true);
257  fgCmd->SetRange("Fg>0");
258  fgCmd->SetDefaultValue(3.5);
260 
261  skinCmd = new G4UIcmdWithADouble("/process/msc/Skin",this);
262  skinCmd->SetGuidance("Set skin parameter for msc processes");
263  skinCmd->SetParameterName("skin",true);
265 
266  dedxCmd = new G4UIcmdWithAnInteger("/process/eLoss/binsDEDX",this);
267  dedxCmd->SetGuidance("Set number of bins for EM tables");
268  dedxCmd->SetParameterName("binsDEDX",true);
271 
272  lamCmd = new G4UIcmdWithAnInteger("/process/eLoss/binsLambda",this);
273  lamCmd->SetGuidance("Set number of bins for EM tables");
274  lamCmd->SetParameterName("binsL",true);
275  lamCmd->SetDefaultValue(77);
277 
278  amCmd = new G4UIcmdWithAnInteger("/process/eLoss/binsPerDecade",this);
279  amCmd->SetGuidance("Set number of bins per decade for EM tables");
280  amCmd->SetParameterName("bins",true);
283 
284  verCmd = new G4UIcmdWithAnInteger("/process/eLoss/verbose",this);
285  verCmd->SetGuidance("Set verbose level for EM physics");
286  verCmd->SetParameterName("verb",true);
289 
290  ver1Cmd = new G4UIcmdWithAnInteger("/process/em/verbose",this);
291  ver1Cmd->SetGuidance("Set verbose level for EM physics");
292  ver1Cmd->SetParameterName("verb1",true);
295 
296  ver2Cmd = new G4UIcmdWithAnInteger("/process/em/workerVerbose",this);
297  ver2Cmd->SetGuidance("Set worker verbose level for EM physics");
298  ver2Cmd->SetParameterName("verb2",true);
301 
302  mscCmd = new G4UIcmdWithAString("/process/msc/StepLimit",this);
303  mscCmd->SetGuidance("Set msc step limitation type");
304  mscCmd->SetParameterName("StepLim",true);
305  mscCmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
307 
308  msc1Cmd = new G4UIcmdWithAString("/process/msc/StepLimitMuHad",this);
309  msc1Cmd->SetGuidance("Set msc step limitation type for muons/hadrons");
310  msc1Cmd->SetParameterName("StepLim1",true);
311  msc1Cmd->SetCandidates("fMinimal fUseSafety fUseSafetyPlus fUseDistanceToBoundary");
313 
314  pixeXsCmd = new G4UIcmdWithAString("/process/em/pixeXSmodel",this);
315  pixeXsCmd->SetGuidance("The name of PIXE cross section");
316  pixeXsCmd->SetParameterName("pixeXS",true);
317  pixeXsCmd->SetCandidates("ECPSSR_Analytical Empirical ECPSSR_FormFactor");
319 
320  pixeeXsCmd = new G4UIcmdWithAString("/process/em/pixeElecXSmodel",this);
321  pixeeXsCmd->SetGuidance("The name of PIXE cross section for electron");
322  pixeeXsCmd->SetParameterName("pixeEXS",true);
323  pixeeXsCmd->SetCandidates("ECPSSR_Analytical Empirical Livermore Penelope");
325 
326  paiCmd = new G4UIcommand("/process/em/AddPAIRegion",this);
327  paiCmd->SetGuidance("Activate PAI in the G4Region.");
328  paiCmd->SetGuidance(" partName : particle name (default - all)");
329  paiCmd->SetGuidance(" regName : G4Region name");
330  paiCmd->SetGuidance(" paiType : PAI, PAIphoton");
331 
332  G4UIparameter* part = new G4UIparameter("partName",'s',false);
333  paiCmd->SetParameter(part);
334 
335  G4UIparameter* pregName = new G4UIparameter("regName",'s',false);
336  paiCmd->SetParameter(pregName);
337 
338  G4UIparameter* ptype = new G4UIparameter("type",'s',false);
339  paiCmd->SetParameter(ptype);
340 
341  meCmd = new G4UIcmdWithAString("/process/em/AddMicroElecRegion",this);
342  meCmd->SetGuidance("Activate MicroElec model in the G4Region");
343  meCmd->SetParameterName("MicroElec",true);
345 
346  dnaCmd = new G4UIcommand("/process/em/AddDNARegion",this);
347  dnaCmd->SetGuidance("Activate DNA in the G4Region.");
348  dnaCmd->SetGuidance(" regName : G4Region name");
349  dnaCmd->SetGuidance(" dnaType : opt0, opt1, opt2");
350 
351  G4UIparameter* regName = new G4UIparameter("regName",'s',false);
352  dnaCmd->SetParameter(regName);
353 
354  G4UIparameter* type = new G4UIparameter("type",'s',false);
355  dnaCmd->SetParameter(type);
356 
357  dumpCmd = new G4UIcommand("/process/em/printParameters",this);
358  dumpCmd->SetGuidance("Print all EM parameters.");
359 
360 }
361 
362 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
363 
365 {
366  delete eLossDirectory;
367  delete mscDirectory;
368  delete emDirectory;
369 
370  delete flucCmd;
371  delete rangeCmd;
372  delete lpmCmd;
373  delete splCmd;
374  delete rsCmd;
375  delete aplCmd;
376  delete deCmd;
377  delete dirFluoCmd;
378  delete auCmd;
379  delete auCascadeCmd;
380  delete pixeCmd;
381  delete dcutCmd;
382  delete latCmd;
383  delete mulatCmd;
384  delete catCmd;
385  delete delCmd;
386  delete mottCmd;
387 
388  delete minSubSecCmd;
389  delete minEnCmd;
390  delete maxEnCmd;
391  delete cenCmd;
392  delete lowEnCmd;
393  delete lowhEnCmd;
394  delete lllCmd;
395  delete brCmd;
396  delete labCmd;
397  delete mscfCmd;
398  delete angCmd;
399  delete frCmd;
400  delete fr1Cmd;
401  delete fgCmd;
402  delete skinCmd;
403 
404  delete dedxCmd;
405  delete lamCmd;
406  delete amCmd;
407  delete verCmd;
408  delete ver1Cmd;
409  delete ver2Cmd;
410 
411  delete mscCmd;
412  delete msc1Cmd;
413 
414  delete pixeXsCmd;
415  delete pixeeXsCmd;
416 
417  delete paiCmd;
418  delete meCmd;
419  delete dnaCmd;
420  delete dumpCmd;
421 }
422 
423 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
424 
426  G4String newValue)
427 {
428  G4bool physicsModified = false;
429  if (command == flucCmd) {
431  physicsModified = true;
432  } else if (command == rangeCmd) {
434  } else if (command == lpmCmd) {
436  physicsModified = true;
437  } else if (command == splCmd) {
439  physicsModified = true;
440  } else if (command == rsCmd) {
442  physicsModified = true;
443  } else if (command == aplCmd) {
445  physicsModified = true;
446  } else if (command == deCmd) {
448  physicsModified = true;
449  } else if (command == dirFluoCmd) {
451  } else if (command == auCmd) {
453  physicsModified = true;
454  } else if (command == auCascadeCmd) {
456  physicsModified = true;
457  } else if (command == pixeCmd) {
459  physicsModified = true;
460  } else if (command == dcutCmd) {
462  physicsModified = true;
463  } else if (command == latCmd) {
465  physicsModified = true;
466  } else if (command == mulatCmd) {
468  physicsModified = true;
469  } else if (command == catCmd) {
471  physicsModified = true;
472  } else if (command == delCmd) {
474  } else if (command == mottCmd) {
476 
477  } else if (command == minSubSecCmd) {
479  } else if (command == minEnCmd) {
481  physicsModified = true;
482  } else if (command == maxEnCmd) {
484  physicsModified = true;
485  } else if (command == cenCmd) {
487  physicsModified = true;
488  } else if (command == lowEnCmd) {
490  physicsModified = true;
491  } else if (command == lowhEnCmd) {
493  physicsModified = true;
494  } else if (command == lllCmd) {
496  physicsModified = true;
497  } else if (command == brCmd) {
499  physicsModified = true;
500  } else if (command == labCmd) {
502  physicsModified = true;
503  } else if (command == mscfCmd) {
505  physicsModified = true;
506  } else if (command == angCmd) {
508  physicsModified = true;
509  } else if (command == frCmd) {
511  physicsModified = true;
512  } else if (command == fr1Cmd) {
514  physicsModified = true;
515  } else if (command == fgCmd) {
517  physicsModified = true;
518  } else if (command == skinCmd) {
520  physicsModified = true;
521 
522  } else if (command == dedxCmd) {
524  physicsModified = true;
525  } else if (command == lamCmd) {
527  physicsModified = true;
528  } else if (command == amCmd) {
530  physicsModified = true;
531  } else if (command == verCmd) {
533  } else if (command == ver1Cmd) {
535  } else if (command == ver2Cmd) {
537 
538  } else if (command == mscCmd || command == msc1Cmd) {
539  G4MscStepLimitType msctype = fUseSafety;
540  if(newValue == "Minimal") {
541  msctype = fMinimal;
542  } else if(newValue == "UseDistanceToBoundary") {
543  msctype = fUseDistanceToBoundary;
544  } else if(newValue == "UseSafety") {
545  msctype = fUseSafety;
546  } else if(newValue == "UseSafetyPlus") {
547  msctype = fUseSafetyPlus;
548  } else {
549  G4cout << "### G4EmParametersMessenger WARNING: StepLimit type <"
550  << newValue << "> unknown!" << G4endl;
551  return;
552  }
553  if (command == mscCmd) {
555  } else {
557  }
558  physicsModified = true;
559  } else if (command == pixeXsCmd) {
561  physicsModified = true;
562  } else if (command == pixeeXsCmd) {
564  physicsModified = true;
565  } else if (command == paiCmd) {
566  G4String s1(""),s2(""),s3("");
567  std::istringstream is(newValue);
568  is >> s1 >> s2 >> s3;
569  theParameters->AddPAIModel(s1, s2, s3);
570  } else if (command == meCmd) {
571  theParameters->AddMicroElec(newValue);
572  } else if (command == dnaCmd) {
573  G4String s1(""),s2("");
574  std::istringstream is(newValue);
575  is >> s1 >> s2;
576  theParameters->AddDNA(s1, s2);
577  } else if (command == dumpCmd) {
578  theParameters->Dump();
579  }
580  if(physicsModified) {
581  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
582  }
583 }
584 
585 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetLossFluctuations(G4bool val)
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:152
void SetApplyCuts(G4bool val)
void SetVerbose(G4int val)
void SetDeexcitationIgnoreCut(G4bool val)
void SetUseMottCorrection(G4bool val)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetLowestElectronEnergy(G4double val)
void SetLatDisplacementBeyondSafety(G4bool val)
G4UIcmdWithAnInteger * ver2Cmd
static G4int GetNewIntValue(const char *paramString)
void SetMscStepLimitType(G4MscStepLimitType val)
G4EmParametersMessenger(G4EmParameters *)
G4UIcmdWithADoubleAndUnit * brCmd
G4UIcmdWithAnInteger * ver1Cmd
void SetBeardenFluoDir(G4bool val)
void SetLinearLossLimit(G4double val)
void SetAuger(G4bool val)
G4UIcmdWithADoubleAndUnit * lowhEnCmd
void SetNumberOfBins(G4int val)
void SetMinSubRange(G4double val)
void SetMaxEnergyForCSDARange(G4double val)
void SetUnitCategory(const char *unitCategory)
static G4double GetNewDoubleValue(const char *paramString)
void SetPIXEElectronCrossSectionModel(const G4String &)
void SetMaxEnergy(G4double val)
void SetBremsstrahlungTh(G4double val)
static G4bool GetNewBoolValue(const char *paramString)
void SetDefaultValue(G4bool defVal)
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:58
void AddPAIModel(const G4String &particle, const G4String &region, const G4String &type)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetLateralDisplacement(G4bool val)
void SetWorkerVerbose(G4int val)
G4UIcmdWithADoubleAndUnit * angCmd
void SetPIXECrossSectionModel(const G4String &)
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
void SetMscRangeFactor(G4double val)
static G4double GetNewDoubleValue(const char *paramString)
void SetRange(const char *rs)
Definition: G4UIcommand.hh:125
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:161
void SetAugerCascade(G4bool val)
void Dump() const
void SetNumberOfBinsPerDecade(G4int val)
void SetLowestMuHadEnergy(G4double val)
void SetMscGeomFactor(G4double val)
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:239
void SetMscMuHadStepLimitType(G4MscStepLimitType val)
void AddDNA(const G4String &region, const G4String &type)
G4UIcmdWithADouble * minSubSecCmd
G4UIcmdWithADoubleAndUnit * minEnCmd
G4UIcmdWithAnInteger * lamCmd
G4UIcmdWithADoubleAndUnit * maxEnCmd
void SetBuildCSDARange(G4bool val)
G4UIcmdWithAnInteger * verCmd
void AddMicroElec(const G4String &region)
void SetMscMuHadRangeFactor(G4double val)
void SetPixe(G4bool val)
void SetSpline(G4bool val)
void SetMuHadLateralDisplacement(G4bool val)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetMinEnergy(G4double val)
void SetLPM(G4bool val)
G4MscStepLimitType
G4UIcmdWithADoubleAndUnit * lowEnCmd
void SetDefaultValue(G4double defVal)
void ActivateAngularGeneratorForIonisation(G4bool val)
void SetUseCutAsFinalRange(G4bool val)
void SetCandidates(const char *candidateList)
#define G4endl
Definition: G4ios.hh:61
void SetDefaultValue(G4int defVal)
void SetMscThetaLimit(G4double val)
void SetLambdaFactor(G4double val)
void SetFactorForAngleLimit(G4double val)
G4UIcmdWithADoubleAndUnit * cenCmd
virtual void SetNewValue(G4UIcommand *, G4String)
G4UIcmdWithAnInteger * dedxCmd
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
G4UIcmdWithAnInteger * amCmd
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:446
void SetFluo(G4bool val)
void SetMscSkin(G4double val)