Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 #include "G4UIdirectory.hh"
48 #include "G4UIcommand.hh"
49 #include "G4UIparameter.hh"
50 #include "G4UIcmdWithABool.hh"
51 #include "G4UIcmdWithAnInteger.hh"
52 #include "G4UIcmdWithADouble.hh"
54 #include "G4UIcmdWithAString.hh"
55 #include "G4UImanager.hh"
56 #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);
120  dirFluoCmd->AvailableForStates(G4State_PreInit);
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);
131  auCascadeCmd->SetDefaultValue(false);
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  IntegCmd = new G4UIcmdWithABool("/process/eLoss/integral",this);
171  IntegCmd->SetGuidance("Switch true/false the integral option");
172  IntegCmd->SetParameterName("integ",true);
173  IntegCmd->SetDefaultValue(true);
175 
176  mottCmd = new G4UIcmdWithABool("/process/msc/UseMottCorrection",this);
177  mottCmd->SetGuidance("Enable usage of Mott corrections for e- elastic scattering");
178  mottCmd->SetParameterName("mott",true);
179  mottCmd->SetDefaultValue(false);
181 
182  birksCmd = new G4UIcmdWithABool("/process/msc/UseG4EmSaturation",this);
183  birksCmd->SetGuidance("Enable usage of built-in Birks saturation");
184  birksCmd->SetParameterName("birks",true);
185  birksCmd->SetDefaultValue(false);
187 
188  minSubSecCmd = new G4UIcmdWithADouble("/process/eLoss/minsubsec",this);
189  minSubSecCmd->SetGuidance("Set the ratio subcut/cut ");
190  minSubSecCmd->SetParameterName("rcmin",true);
191  minSubSecCmd->AvailableForStates(G4State_PreInit);
192 
193  minEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/minKinEnergy",this);
194  minEnCmd->SetGuidance("Set the min kinetic energy for EM tables");
195  minEnCmd->SetParameterName("emin",true);
196  minEnCmd->SetUnitCategory("Energy");
198 
199  maxEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergy",this);
200  maxEnCmd->SetGuidance("Set the max kinetic energy for EM tables");
201  maxEnCmd->SetParameterName("emax",true);
202  maxEnCmd->SetUnitCategory("Energy");
204 
205  cenCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergyCSDA",this);
206  cenCmd->SetGuidance("Set the max kinetic energy for CSDA table");
207  cenCmd->SetParameterName("emaxCSDA",true);
208  cenCmd->SetUnitCategory("Energy");
210 
211  lowEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestElectronEnergy",this);
212  lowEnCmd->SetGuidance("Set the lowest kinetic energy for e+-");
213  lowEnCmd->SetParameterName("elow",true);
214  lowEnCmd->SetUnitCategory("Energy");
216 
217  lowhEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestMuHadEnergy",this);
218  lowhEnCmd->SetGuidance("Set the lowest kinetic energy for muons and hadrons");
219  lowhEnCmd->SetParameterName("elowh",true);
220  lowhEnCmd->SetUnitCategory("Energy");
222 
223  lllCmd = new G4UIcmdWithADouble("/process/eLoss/linLossLimit",this);
224  lllCmd->SetGuidance("Set linearLossLimit parameter");
225  lllCmd->SetParameterName("linlim",true);
227 
228  brCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/bremThreshold",this);
229  brCmd->SetGuidance("Set bremsstrahlung energy threshold");
230  brCmd->SetParameterName("emaxBrem",true);
231  brCmd->SetUnitCategory("Energy");
233 
234  labCmd = new G4UIcmdWithADouble("/process/eLoss/LambdaFactor",this);
235  labCmd->SetGuidance("Set lambdaFactor parameter for integral option");
236  labCmd->SetParameterName("Fl",true);
238 
239  mscfCmd = new G4UIcmdWithADouble("/process/msc/FactorForAngleLimit",this);
240  mscfCmd->SetGuidance("Set factor for computation of a limit for -t (invariant trasfer)");
241  mscfCmd->SetParameterName("Fact",true);
242  mscfCmd->SetRange("Fact>0");
243  mscfCmd->SetDefaultValue(1.);
245 
246  angCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/ThetaLimit",this);
247  angCmd->SetGuidance("Set the limit on the polar angle for msc and single scattering");
248  angCmd->SetParameterName("theta",true);
249  angCmd->SetUnitCategory("Angle");
251 
252  frCmd = new G4UIcmdWithADouble("/process/msc/RangeFactor",this);
253  frCmd->SetGuidance("Set RangeFactor for msc processes of e+-");
254  frCmd->SetParameterName("Fr",true);
255  frCmd->SetRange("Fr>0");
256  frCmd->SetDefaultValue(0.04);
258 
259  fr1Cmd = new G4UIcmdWithADouble("/process/msc/RangeFactorMuHad",this);
260  fr1Cmd->SetGuidance("Set RangeFactor for msc processes of muons/hadrons");
261  fr1Cmd->SetParameterName("Fr1",true);
262  fr1Cmd->SetRange("Fr1>0");
263  fr1Cmd->SetDefaultValue(0.2);
265 
266  fgCmd = new G4UIcmdWithADouble("/process/msc/GeomFactor",this);
267  fgCmd->SetGuidance("Set GeomFactor parameter for msc processes");
268  fgCmd->SetParameterName("Fg",true);
269  fgCmd->SetRange("Fg>0");
270  fgCmd->SetDefaultValue(3.5);
272 
273  skinCmd = new G4UIcmdWithADouble("/process/msc/Skin",this);
274  skinCmd->SetGuidance("Set skin parameter for msc processes");
275  skinCmd->SetParameterName("skin",true);
277 
278  dedxCmd = new G4UIcmdWithAnInteger("/process/eLoss/binsDEDX",this);
279  dedxCmd->SetGuidance("Set number of bins for EM tables");
280  dedxCmd->SetParameterName("binsDEDX",true);
281  dedxCmd->SetDefaultValue(77);
283 
284  lamCmd = new G4UIcmdWithAnInteger("/process/eLoss/binsLambda",this);
285  lamCmd->SetGuidance("Set number of bins for EM tables");
286  lamCmd->SetParameterName("binsL",true);
287  lamCmd->SetDefaultValue(77);
289 
290  amCmd = new G4UIcmdWithAnInteger("/process/eLoss/binsPerDecade",this);
291  amCmd->SetGuidance("Set number of bins per decade for EM tables");
292  amCmd->SetParameterName("bins",true);
293  amCmd->SetDefaultValue(7);
295 
296  verCmd = new G4UIcmdWithAnInteger("/process/eLoss/verbose",this);
297  verCmd->SetGuidance("Set verbose level for EM physics");
298  verCmd->SetParameterName("verb",true);
299  verCmd->SetDefaultValue(1);
301 
302  ver1Cmd = new G4UIcmdWithAnInteger("/process/em/verbose",this);
303  ver1Cmd->SetGuidance("Set verbose level for EM physics");
304  ver1Cmd->SetParameterName("verb1",true);
305  ver1Cmd->SetDefaultValue(1);
307 
308  ver2Cmd = new G4UIcmdWithAnInteger("/process/em/workerVerbose",this);
309  ver2Cmd->SetGuidance("Set worker verbose level for EM physics");
310  ver2Cmd->SetParameterName("verb2",true);
311  ver2Cmd->SetDefaultValue(1);
313 
314  mscCmd = new G4UIcmdWithAString("/process/msc/StepLimit",this);
315  mscCmd->SetGuidance("Set msc step limitation type");
316  mscCmd->SetParameterName("StepLim",true);
317  mscCmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
319 
320  msc1Cmd = new G4UIcmdWithAString("/process/msc/StepLimitMuHad",this);
321  msc1Cmd->SetGuidance("Set msc step limitation type for muons/hadrons");
322  msc1Cmd->SetParameterName("StepLim1",true);
323  msc1Cmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary");
325 
326  pixeXsCmd = new G4UIcmdWithAString("/process/em/pixeXSmodel",this);
327  pixeXsCmd->SetGuidance("The name of PIXE cross section");
328  pixeXsCmd->SetParameterName("pixeXS",true);
329  pixeXsCmd->SetCandidates("ECPSSR_Analytical Empirical ECPSSR_FormFactor");
331 
332  pixeeXsCmd = new G4UIcmdWithAString("/process/em/pixeElecXSmodel",this);
333  pixeeXsCmd->SetGuidance("The name of PIXE cross section for electron");
334  pixeeXsCmd->SetParameterName("pixeEXS",true);
335  pixeeXsCmd->SetCandidates("ECPSSR_Analytical Empirical Livermore Penelope");
337 
338  paiCmd = new G4UIcommand("/process/em/AddPAIRegion",this);
339  paiCmd->SetGuidance("Activate PAI in the G4Region.");
340  paiCmd->SetGuidance(" partName : particle name (default - all)");
341  paiCmd->SetGuidance(" regName : G4Region name");
342  paiCmd->SetGuidance(" paiType : PAI, PAIphoton");
344 
345  G4UIparameter* part = new G4UIparameter("partName",'s',false);
346  paiCmd->SetParameter(part);
347 
348  G4UIparameter* pregName = new G4UIparameter("regName",'s',false);
349  paiCmd->SetParameter(pregName);
350 
351  G4UIparameter* ptype = new G4UIparameter("type",'s',false);
352  paiCmd->SetParameter(ptype);
353 
354  meCmd = new G4UIcmdWithAString("/process/em/AddMicroElecRegion",this);
355  meCmd->SetGuidance("Activate MicroElec model in the G4Region");
356  meCmd->SetParameterName("MicroElec",true);
358 
359  dnaCmd = new G4UIcommand("/process/em/AddDNARegion",this);
360  dnaCmd->SetGuidance("Activate DNA in a G4Region.");
361  dnaCmd->SetGuidance(" regName : G4Region name");
362  dnaCmd->SetGuidance(" dnaType : DNA_opt0, DNA_opt1, DNA_opt2");
364 
365  G4UIparameter* regName = new G4UIparameter("regName",'s',false);
366  dnaCmd->SetParameter(regName);
367 
368  G4UIparameter* type = new G4UIparameter("dnaType",'s',false);
369  dnaCmd->SetParameter(type);
370 
371  mscoCmd = new G4UIcommand("/process/em/AddEmRegion",this);
372  mscoCmd->SetGuidance("Add optional EM configuration for a G4Region.");
373  mscoCmd->SetGuidance(" regName : G4Region name");
374  mscoCmd->SetGuidance(" mscType : G4EmStandard, G4EmStandard_opt1, ...");
376 
377  G4UIparameter* mregName = new G4UIparameter("regName",'s',false);
378  mscoCmd->SetParameter(mregName);
379 
380  G4UIparameter* mtype = new G4UIparameter("mscType",'s',false);
381  mscoCmd->SetParameter(mtype);
382 
383  dumpCmd = new G4UIcommand("/process/em/printParameters",this);
384  dumpCmd->SetGuidance("Print all EM parameters.");
385 
386  SubSecCmd = new G4UIcommand("/process/eLoss/subsec",this);
387  SubSecCmd->SetGuidance("Switch true/false the subcutoff generation per region.");
388  SubSecCmd->SetGuidance(" subSec : true/false");
389  SubSecCmd->SetGuidance(" Region : region name");
391 
392  G4UIparameter* subSec = new G4UIparameter("subSec",'s',false);
393  SubSecCmd->SetParameter(subSec);
394 
395  G4UIparameter* subSecReg = new G4UIparameter("Region",'s',false);
396  SubSecCmd->SetParameter(subSecReg);
397 
398  StepFuncCmd = new G4UIcommand("/process/eLoss/StepFunction",this);
399  StepFuncCmd->SetGuidance("Set the energy loss step limitation parameters for e+-.");
400  StepFuncCmd->SetGuidance(" dRoverR : max Range variation per step");
401  StepFuncCmd->SetGuidance(" finalRange: range for final step");
403 
404  G4UIparameter* dRoverRPrm = new G4UIparameter("dRoverR",'d',false);
405  dRoverRPrm->SetParameterRange("dRoverR>0. && dRoverR<=1.");
406  StepFuncCmd->SetParameter(dRoverRPrm);
407 
408  G4UIparameter* finalRangePrm = new G4UIparameter("finalRange",'d',false);
409  finalRangePrm->SetParameterRange("finalRange>0.");
410  StepFuncCmd->SetParameter(finalRangePrm);
411 
412  G4UIparameter* unitPrm = new G4UIparameter("unit",'s',true);
413  unitPrm->SetDefaultValue("mm");
414  StepFuncCmd->SetParameter(unitPrm);
415 
416  StepFuncCmd1 = new G4UIcommand("/process/eLoss/StepFunctionMuHad",this);
417  StepFuncCmd1->SetGuidance("Set the energy loss step limitation parameters for muon/hadron.");
418  StepFuncCmd1->SetGuidance(" dRoverR : max Range variation per step");
419  StepFuncCmd1->SetGuidance(" finalRange: range for final step");
421 
422  G4UIparameter* dRoverRPrm1 = new G4UIparameter("dRoverRMuHad",'d',false);
423  dRoverRPrm1->SetParameterRange("dRoverRMuHad>0. && dRoverRMuHad<=1.");
424  StepFuncCmd1->SetParameter(dRoverRPrm1);
425 
426  G4UIparameter* finalRangePrm1 = new G4UIparameter("finalRangeMuHad",'d',false);
427  finalRangePrm1->SetParameterRange("finalRangeMuHad>0.");
428  StepFuncCmd1->SetParameter(finalRangePrm1);
429 
430  G4UIparameter* unitPrm1 = new G4UIparameter("unit",'s',true);
431  unitPrm1->SetDefaultValue("mm");
432  StepFuncCmd1->SetParameter(unitPrm1);
433 
434  deexCmd = new G4UIcommand("/process/em/deexcitation",this);
435  deexCmd->SetGuidance("Set deexcitation flags per G4Region.");
436  deexCmd->SetGuidance(" regName : G4Region name");
437  deexCmd->SetGuidance(" flagFluo : Fluorescence");
438  deexCmd->SetGuidance(" flagAuger : Auger");
439  deexCmd->SetGuidance(" flagPIXE : PIXE");
441 
442  G4UIparameter* regNameD = new G4UIparameter("regName",'s',false);
443  deexCmd->SetParameter(regNameD);
444 
445  G4UIparameter* flagFluo = new G4UIparameter("flagFluo",'s',false);
446  deexCmd->SetParameter(flagFluo);
447 
448  G4UIparameter* flagAuger = new G4UIparameter("flagAuger",'s',false);
449  deexCmd->SetParameter(flagAuger);
450 
451  G4UIparameter* flagPIXE = new G4UIparameter("flagPIXE",'s',false);
452  deexCmd->SetParameter(flagPIXE);
453 
454  bfCmd = new G4UIcommand("/process/em/setBiasingFactor",this);
455  bfCmd->SetGuidance("Set factor for the process cross section.");
456  bfCmd->SetGuidance(" procName : process name");
457  bfCmd->SetGuidance(" procFact : factor");
458  bfCmd->SetGuidance(" flagFact : flag to change weight");
460 
461  G4UIparameter* procName = new G4UIparameter("procName",'s',false);
462  bfCmd->SetParameter(procName);
463 
464  G4UIparameter* procFact = new G4UIparameter("procFact",'d',false);
465  bfCmd->SetParameter(procFact);
466 
467  G4UIparameter* flagFact = new G4UIparameter("flagFact",'s',false);
468  bfCmd->SetParameter(flagFact);
469 
470  fiCmd = new G4UIcommand("/process/em/setForcedInteraction",this);
471  fiCmd->SetGuidance("Set factor for the process cross section.");
472  fiCmd->SetGuidance(" procNam : process name");
473  fiCmd->SetGuidance(" regNam : region name");
474  fiCmd->SetGuidance(" tlength : fixed target length");
475  fiCmd->SetGuidance(" unitT : length unit");
476  fiCmd->SetGuidance(" tflag : flag to change weight");
478 
479  G4UIparameter* procNam = new G4UIparameter("procNam",'s',false);
480  fiCmd->SetParameter(procNam);
481 
482  G4UIparameter* regNam = new G4UIparameter("regNam",'s',false);
483  fiCmd->SetParameter(regNam);
484 
485  G4UIparameter* tlength = new G4UIparameter("tlength",'d',false);
486  fiCmd->SetParameter(tlength);
487 
488  G4UIparameter* unitT = new G4UIparameter("unitT",'s',true);
489  fiCmd->SetParameter(unitT);
490 
491  G4UIparameter* flagT = new G4UIparameter("tflag",'s',true);
492  fiCmd->SetParameter(flagT);
493 
494  bsCmd = new G4UIcommand("/process/em/setSecBiasing",this);
495  bsCmd->SetGuidance("Set bremsstrahlung or delta-e- splitting/Russian roullette per region.");
496  bsCmd->SetGuidance(" bProcNam : process name");
497  bsCmd->SetGuidance(" bRegNam : region name");
498  bsCmd->SetGuidance(" bFactor : number of splitted gamma or probability of Russian roulette");
499  bsCmd->SetGuidance(" bEnergy : max energy of a secondary for this biasing method");
500  bsCmd->SetGuidance(" bUnit : energy unit");
502 
503  G4UIparameter* bProcNam = new G4UIparameter("bProcNam",'s',false);
504  bsCmd->SetParameter(bProcNam);
505 
506  G4UIparameter* bRegNam = new G4UIparameter("bRegNam",'s',false);
507  bsCmd->SetParameter(bRegNam);
508 
509  G4UIparameter* bFactor = new G4UIparameter("bFactor",'d',false);
510  bsCmd->SetParameter(bFactor);
511 
512  G4UIparameter* bEnergy = new G4UIparameter("bEnergy",'d',false);
513  bsCmd->SetParameter(bEnergy);
514 
515  G4UIparameter* bUnit = new G4UIparameter("bUnit",'s',true);
516  bsCmd->SetParameter(bUnit);
517 
518  nffCmd = new G4UIcmdWithAString("/process/em/setNuclearFormFactor",this);
519  nffCmd->SetGuidance("Define typy of nuclear form-factor");
520  nffCmd->SetParameterName("NucFF",true);
521  nffCmd->SetCandidates("None Exponential Gaussian Flat");
523 }
524 
525 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
526 
528 {
529  delete eLossDirectory;
530  delete mscDirectory;
531  delete emDirectory;
532 
533  delete flucCmd;
534  delete rangeCmd;
535  delete lpmCmd;
536  delete splCmd;
537  delete rsCmd;
538  delete aplCmd;
539  delete deCmd;
540  delete dirFluoCmd;
541  delete auCmd;
542  delete auCascadeCmd;
543  delete pixeCmd;
544  delete dcutCmd;
545  delete latCmd;
546  delete mulatCmd;
547  delete catCmd;
548  delete delCmd;
549  delete IntegCmd;
550  delete mottCmd;
551  delete birksCmd;
552 
553  delete minSubSecCmd;
554  delete minEnCmd;
555  delete maxEnCmd;
556  delete cenCmd;
557  delete lowEnCmd;
558  delete lowhEnCmd;
559  delete lllCmd;
560  delete brCmd;
561  delete labCmd;
562  delete mscfCmd;
563  delete angCmd;
564  delete frCmd;
565  delete fr1Cmd;
566  delete fgCmd;
567  delete skinCmd;
568 
569  delete dedxCmd;
570  delete lamCmd;
571  delete amCmd;
572  delete verCmd;
573  delete ver1Cmd;
574  delete ver2Cmd;
575 
576  delete mscCmd;
577  delete msc1Cmd;
578 
579  delete pixeXsCmd;
580  delete pixeeXsCmd;
581 
582  delete paiCmd;
583  delete meCmd;
584  delete dnaCmd;
585  delete mscoCmd;
586  delete dumpCmd;
587 
588  delete SubSecCmd;
589  delete StepFuncCmd;
590  delete StepFuncCmd1;
591  delete deexCmd;
592  delete bfCmd;
593  delete fiCmd;
594  delete bsCmd;
595  delete nffCmd;
596 }
597 
598 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
599 
601  G4String newValue)
602 {
603  G4bool physicsModified = false;
604  if (command == flucCmd) {
605  theParameters->SetLossFluctuations(flucCmd->GetNewBoolValue(newValue));
606  physicsModified = true;
607  } else if (command == rangeCmd) {
608  theParameters->SetBuildCSDARange(rangeCmd->GetNewBoolValue(newValue));
609  } else if (command == lpmCmd) {
610  theParameters->SetLPM(lpmCmd->GetNewBoolValue(newValue));
611  physicsModified = true;
612  } else if (command == splCmd) {
613  theParameters->SetSpline(splCmd->GetNewBoolValue(newValue));
614  physicsModified = true;
615  } else if (command == rsCmd) {
616  theParameters->SetUseCutAsFinalRange(rsCmd->GetNewBoolValue(newValue));
617  physicsModified = true;
618  } else if (command == aplCmd) {
619  theParameters->SetApplyCuts(aplCmd->GetNewBoolValue(newValue));
620  physicsModified = true;
621  } else if (command == deCmd) {
622  theParameters->SetFluo(deCmd->GetNewBoolValue(newValue));
623  physicsModified = true;
624  } else if (command == dirFluoCmd) {
625  theParameters->SetBeardenFluoDir(dirFluoCmd->GetNewBoolValue(newValue));
626  } else if (command == auCmd) {
627  theParameters->SetAuger(auCmd->GetNewBoolValue(newValue));
628  physicsModified = true;
629  } else if (command == auCascadeCmd) {
630  theParameters->SetAugerCascade(auCascadeCmd->GetNewBoolValue(newValue));
631  physicsModified = true;
632  } else if (command == pixeCmd) {
633  theParameters->SetPixe(pixeCmd->GetNewBoolValue(newValue));
634  physicsModified = true;
635  } else if (command == dcutCmd) {
636  theParameters->SetDeexcitationIgnoreCut(dcutCmd->GetNewBoolValue(newValue));
637  physicsModified = true;
638  } else if (command == latCmd) {
639  theParameters->SetLateralDisplacement(latCmd->GetNewBoolValue(newValue));
640  physicsModified = true;
641  } else if (command == mulatCmd) {
642  theParameters->SetMuHadLateralDisplacement(mulatCmd->GetNewBoolValue(newValue));
643  physicsModified = true;
644  } else if (command == catCmd) {
645  theParameters->SetLatDisplacementBeyondSafety(catCmd->GetNewBoolValue(newValue));
646  physicsModified = true;
647  } else if (command == delCmd) {
648  theParameters->ActivateAngularGeneratorForIonisation(delCmd->GetNewBoolValue(newValue));
649  } else if (command == IntegCmd) {
650  theParameters->SetIntegral(IntegCmd->GetNewBoolValue(newValue));
651  physicsModified = true;
652  } else if (command == mottCmd) {
653  theParameters->SetUseMottCorrection(mottCmd->GetNewBoolValue(newValue));
654  } else if (command == birksCmd) {
655  theParameters->SetBirksActive(birksCmd->GetNewBoolValue(newValue));
656 
657  } else if (command == minSubSecCmd) {
658  theParameters->SetMinSubRange(minSubSecCmd->GetNewDoubleValue(newValue));
659  } else if (command == minEnCmd) {
660  theParameters->SetMinEnergy(minEnCmd->GetNewDoubleValue(newValue));
661  physicsModified = true;
662  } else if (command == maxEnCmd) {
663  theParameters->SetMaxEnergy(maxEnCmd->GetNewDoubleValue(newValue));
664  physicsModified = true;
665  } else if (command == cenCmd) {
666  theParameters->SetMaxEnergyForCSDARange(cenCmd->GetNewDoubleValue(newValue));
667  physicsModified = true;
668  } else if (command == lowEnCmd) {
669  theParameters->SetLowestElectronEnergy(lowEnCmd->GetNewDoubleValue(newValue));
670  physicsModified = true;
671  } else if (command == lowhEnCmd) {
672  theParameters->SetLowestMuHadEnergy(lowhEnCmd->GetNewDoubleValue(newValue));
673  physicsModified = true;
674  } else if (command == lllCmd) {
675  theParameters->SetLinearLossLimit(lllCmd->GetNewDoubleValue(newValue));
676  physicsModified = true;
677  } else if (command == brCmd) {
678  theParameters->SetBremsstrahlungTh(brCmd->GetNewDoubleValue(newValue));
679  physicsModified = true;
680  } else if (command == labCmd) {
681  theParameters->SetLambdaFactor(labCmd->GetNewDoubleValue(newValue));
682  physicsModified = true;
683  } else if (command == mscfCmd) {
684  theParameters->SetFactorForAngleLimit(mscfCmd->GetNewDoubleValue(newValue));
685  physicsModified = true;
686  } else if (command == angCmd) {
687  theParameters->SetMscThetaLimit(angCmd->GetNewDoubleValue(newValue));
688  physicsModified = true;
689  } else if (command == frCmd) {
690  theParameters->SetMscRangeFactor(frCmd->GetNewDoubleValue(newValue));
691  physicsModified = true;
692  } else if (command == fr1Cmd) {
693  theParameters->SetMscMuHadRangeFactor(fr1Cmd->GetNewDoubleValue(newValue));
694  physicsModified = true;
695  } else if (command == fgCmd) {
696  theParameters->SetMscGeomFactor(fgCmd->GetNewDoubleValue(newValue));
697  physicsModified = true;
698  } else if (command == skinCmd) {
699  theParameters->SetMscSkin(skinCmd->GetNewDoubleValue(newValue));
700  physicsModified = true;
701 
702  } else if (command == dedxCmd) {
703  theParameters->SetNumberOfBins(dedxCmd->GetNewIntValue(newValue));
704  physicsModified = true;
705  } else if (command == lamCmd) {
706  theParameters->SetNumberOfBins(lamCmd->GetNewIntValue(newValue));
707  physicsModified = true;
708  } else if (command == amCmd) {
709  theParameters->SetNumberOfBinsPerDecade(amCmd->GetNewIntValue(newValue));
710  physicsModified = true;
711  } else if (command == verCmd) {
712  theParameters->SetVerbose(verCmd->GetNewIntValue(newValue));
713  } else if (command == ver1Cmd) {
714  theParameters->SetVerbose(ver1Cmd->GetNewIntValue(newValue));
715  } else if (command == ver2Cmd) {
716  theParameters->SetWorkerVerbose(ver2Cmd->GetNewIntValue(newValue));
717 
718  } else if (command == mscCmd || command == msc1Cmd) {
719  G4MscStepLimitType msctype = fUseSafety;
720  if(newValue == "Minimal") {
721  msctype = fMinimal;
722  } else if(newValue == "UseDistanceToBoundary") {
723  msctype = fUseDistanceToBoundary;
724  } else if(newValue == "UseSafety") {
725  msctype = fUseSafety;
726  } else if(newValue == "UseSafetyPlus") {
727  msctype = fUseSafetyPlus;
728  } else {
729  G4cout << "### G4EmParametersMessenger WARNING: StepLimit type <"
730  << newValue << "> unknown!" << G4endl;
731  return;
732  }
733  if (command == mscCmd) {
734  theParameters->SetMscStepLimitType(msctype);
735  } else {
736  theParameters->SetMscMuHadStepLimitType(msctype);
737  }
738  physicsModified = true;
739  } else if (command == pixeXsCmd) {
740  theParameters->SetPIXECrossSectionModel(newValue);
741  physicsModified = true;
742  } else if (command == pixeeXsCmd) {
743  theParameters->SetPIXEElectronCrossSectionModel(newValue);
744  physicsModified = true;
745  } else if (command == paiCmd) {
746  G4String s1(""),s2(""),s3("");
747  std::istringstream is(newValue);
748  is >> s1 >> s2 >> s3;
749  theParameters->AddPAIModel(s1, s2, s3);
750  } else if (command == meCmd) {
751  theParameters->AddMicroElec(newValue);
752  } else if (command == dnaCmd) {
753  G4String s1(""),s2("");
754  std::istringstream is(newValue);
755  is >> s1 >> s2;
756  theParameters->AddDNA(s1, s2);
757  } else if (command == mscoCmd) {
758  G4String s1(""),s2("");
759  std::istringstream is(newValue);
760  is >> s1 >> s2;
761  theParameters->AddMsc(s1, s2);
762  } else if (command == dumpCmd) {
763  theParameters->Dump();
764  } else if (command == SubSecCmd) {
765  G4String s1, s2;
766  std::istringstream is(newValue);
767  is >> s1 >> s2;
768  G4bool yes = false;
769  if(s1 == "true") { yes = true; }
770  theParameters->SetSubCutoff(yes,s2);
771  } else if (command == StepFuncCmd || command == StepFuncCmd1) {
772  G4double v1,v2;
773  G4String unt;
774  std::istringstream is(newValue);
775  is >> v1 >> v2 >> unt;
776  v2 *= G4UIcommand::ValueOf(unt);
777  if(command == StepFuncCmd) {
778  theParameters->SetStepFunction(v1,v2);
779  } else {
780  theParameters->SetStepFunctionMuHad(v1,v2);
781  }
782  physicsModified = true;
783  } else if (command == deexCmd) {
784  G4String s1 (""), s2(""), s3(""), s4("");
785  G4bool b2(false), b3(false), b4(false);
786  std::istringstream is(newValue);
787  is >> s1 >> s2 >> s3 >> s4;
788  if(s2 == "true") { b2 = true; }
789  if(s3 == "true") { b3 = true; }
790  if(s4 == "true") { b4 = true; }
791  theParameters->SetDeexActiveRegion(s1,b2,b3,b4);
792  physicsModified = true;
793  } else if (command == bfCmd) {
794  G4double v1(1.0);
795  G4String s0(""),s1("");
796  std::istringstream is(newValue);
797  is >> s0 >> v1 >> s1;
798  G4bool yes = false;
799  if(s1 == "true") { yes = true; }
800  theParameters->SetProcessBiasingFactor(s0,v1,yes);
801  physicsModified = true;
802  } else if (command == fiCmd) {
803  G4double v1(0.0);
804  G4String s1(""),s2(""),s3(""),unt("mm");
805  std::istringstream is(newValue);
806  is >> s1 >> s2 >> v1 >> unt >> s3;
807  G4bool yes = false;
808  if(s3 == "true") { yes = true; }
809  v1 *= G4UIcommand::ValueOf(unt);
810  theParameters->ActivateForcedInteraction(s1,s2,v1,yes);
811  physicsModified = true;
812  } else if (command == bsCmd) {
813  G4double fb(1.0),en(1.e+30);
814  G4String s1(""),s2(""),unt("MeV");
815  std::istringstream is(newValue);
816  is >> s1 >> s2 >> fb >> en >> unt;
817  en *= G4UIcommand::ValueOf(unt);
818  theParameters->ActivateSecondaryBiasing(s1,s2,fb,en);
819  physicsModified = true;
820  } else if (command == nffCmd) {
822  if(newValue == "Exponential") { x = fExponentialNF; }
823  else if(newValue == "Gaussian") { x = fGaussianNF; }
824  else if(newValue == "Flat") { x = fFlatNF; }
825  theParameters->SetNuclearFormfactorType(x);
826  physicsModified = true;
827  }
828  if(physicsModified) {
829  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
830  }
831 }
832 
833 //....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
virtual void SetNewValue(G4UIcommand *, G4String) override
void SetApplyCuts(G4bool val)
void SetProcessBiasingFactor(const G4String &procname, G4double val, G4bool wflag)
void SetVerbose(G4int val)
void ActivateSecondaryBiasing(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
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)
static G4int GetNewIntValue(const char *paramString)
void SetMscStepLimitType(G4MscStepLimitType val)
G4EmParametersMessenger(G4EmParameters *)
void SetParameterRange(const char *theRange)
void SetBeardenFluoDir(G4bool val)
void SetLinearLossLimit(G4double val)
void SetAuger(G4bool val)
void SetDefaultValue(const char *theDefaultValue)
void SetNumberOfBins(G4int val)
void SetMinSubRange(G4double val)
void SetMaxEnergyForCSDARange(G4double val)
void SetUnitCategory(const char *unitCategory)
void SetStepFunctionMuHad(G4double v1, G4double v2)
static G4double GetNewDoubleValue(const char *paramString)
void SetPIXEElectronCrossSectionModel(const G4String &)
void SetDeexActiveRegion(const G4String &region, G4bool fdeex, G4bool fauger, G4bool fpixe)
void SetMaxEnergy(G4double val)
void SetBremsstrahlungTh(G4double val)
static G4bool GetNewBoolValue(const char *paramString)
void SetDefaultValue(G4bool defVal)
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:59
void SetBirksActive(G4bool val)
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)
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 SetNuclearFormfactorType(G4NuclearFormfactorType val)
void SetLowestMuHadEnergy(G4double val)
void SetMscGeomFactor(G4double val)
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:240
void SetMscMuHadStepLimitType(G4MscStepLimitType val)
void AddDNA(const G4String &region, const G4String &type)
void SetSubCutoff(G4bool val, const G4String &region="")
void SetBuildCSDARange(G4bool val)
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)
G4NuclearFormfactorType
G4MscStepLimitType
static G4double ValueOf(const char *unitName)
Definition: G4UIcommand.cc:309
void SetDefaultValue(G4double defVal)
void SetIntegral(G4bool val)
void AddMsc(const G4String &region, const G4String &type)
void ActivateAngularGeneratorForIonisation(G4bool val)
void ActivateForcedInteraction(const G4String &procname, const G4String &region, G4double length, G4bool wflag)
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)
double G4double
Definition: G4Types.hh:76
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetStepFunction(G4double v1, G4double v2)
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:447
void SetFluo(G4bool val)
void SetMscSkin(G4double val)