Geant4  10.00.p03
G4EnergyLossMessenger.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: G4EnergyLossMessenger.cc 81864 2014-06-06 11:30:54Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 // File name: G4EnergyLossMessenger
33 //
34 // Author: Michel Maire
35 //
36 // Creation date: 17-03-2011 (original version of 22-06-2000)
37 //
38 // Modifications:
39 // 10-01-06 SetStepLimits -> SetStepFunction (V.Ivanchenko)
40 // 10-01-06 PreciseRange -> CSDARange (V.Ivanchenko)
41 // 10-05-06 Add command MscStepLimit (V.Ivanchenko)
42 // 10-10-06 Add DEDXBinning command (V.Ivanchenko)
43 // 07-02-07 Add MscLateralDisplacement command (V.Ivanchenko)
44 // 12-02-07 Add SetSkin, SetLinearLossLimit (V.Ivanchenko)
45 // 15-03-07 Send a message "/run/physicsModified" if reinitialisation
46 // is needed after the command (V.Ivanchenko)
47 // 16-03-07 modify /process/eLoss/minsubsec command (V.Ivanchenko)
48 // 18-05-07 add /process/msc directory and commands (V.Ivanchenko)
49 // 11-03-08 add /process/em directory and commands (V.Ivanchenko)
50 //
51 // -------------------------------------------------------------------
52 //
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
56 
57 #include "G4EnergyLossMessenger.hh"
58 
59 #include "G4UIdirectory.hh"
60 #include "G4UIcommand.hh"
61 #include "G4UIparameter.hh"
62 #include "G4UIcmdWithABool.hh"
63 #include "G4UIcmdWithAnInteger.hh"
64 #include "G4UIcmdWithADouble.hh"
66 #include "G4UIcmdWithAString.hh"
67 #include "G4EmProcessOptions.hh"
68 #include "G4UImanager.hh"
69 #include "G4MscStepLimitType.hh"
70 #include "G4EmProcessOptions.hh"
71 
72 #include <sstream>
73 
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75 
77 {
78  opt = 0;
79  eLossDirectory = new G4UIdirectory("/process/eLoss/");
80  eLossDirectory->SetGuidance("Commands for EM processes.");
81  mscDirectory = new G4UIdirectory("/process/msc/");
82  mscDirectory->SetGuidance("Commands for EM scattering processes.");
83  emDirectory = new G4UIdirectory("/process/em/");
84  emDirectory->SetGuidance("General commands for EM processes.");
85 
86  RndmStepCmd = new G4UIcmdWithABool("/process/eLoss/useCutAsFinalRange",this);
87  RndmStepCmd->SetGuidance("Use cut in range as a final range");
88  RndmStepCmd->SetParameterName("choice",true);
91 
92  EnlossFlucCmd = new G4UIcmdWithABool("/process/eLoss/fluct",this);
93  EnlossFlucCmd->SetGuidance("Switch true/false the energy loss fluctuations.");
94  EnlossFlucCmd->SetParameterName("choice",true);
97 
98  SubSecCmd = new G4UIcmdWithABool("/process/eLoss/subsec",this);
99  SubSecCmd->SetGuidance("Switch true/false the subcutoff generation.");
100  SubSecCmd->SetParameterName("choice",true);
101  SubSecCmd->SetDefaultValue(true);
103 
104  MinSubSecCmd = new G4UIcmdWithADouble("/process/eLoss/minsubsec",this);
105  MinSubSecCmd->SetGuidance("Set the ratio subcut/cut ");
106  MinSubSecCmd->SetParameterName("rcmin",true);
108 
109  StepFuncCmd = new G4UIcommand("/process/eLoss/StepFunction",this);
110  StepFuncCmd->SetGuidance("Set the energy loss step limitation parameters.");
111  StepFuncCmd->SetGuidance(" dRoverR : max Range variation per step");
112  StepFuncCmd->SetGuidance(" finalRange: range for final step");
113 
114  G4UIparameter* dRoverRPrm = new G4UIparameter("dRoverR",'d',false);
115  dRoverRPrm->SetGuidance("max Range variation per step (fractional number)");
116  dRoverRPrm->SetParameterRange("dRoverR>0. && dRoverR<=1.");
117  StepFuncCmd->SetParameter(dRoverRPrm);
118 
119  G4UIparameter* finalRangePrm = new G4UIparameter("finalRange",'d',false);
120  finalRangePrm->SetGuidance("range for final step");
121  finalRangePrm->SetParameterRange("finalRange>0.");
122  StepFuncCmd->SetParameter(finalRangePrm);
123 
124  G4UIparameter* unitPrm = new G4UIparameter("unit",'s',true);
125  unitPrm->SetGuidance("unit of finalRange");
126  unitPrm->SetDefaultValue("mm");
127  G4String unitCandidates =
129  unitPrm->SetParameterCandidates(unitCandidates);
130 
131  StepFuncCmd->SetParameter(unitPrm);
133 
134  MinEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/minKinEnergy",this);
135  MinEnCmd->SetGuidance("Set the min kinetic energy");
136  MinEnCmd->SetParameterName("emin",true);
137  MinEnCmd->SetUnitCategory("Energy");
139 
140  MaxEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergy",this);
141  MaxEnCmd->SetGuidance("Set the max kinetic energy");
142  MaxEnCmd->SetParameterName("emax",true);
143  MaxEnCmd->SetUnitCategory("Energy");
145 
146  IntegCmd = new G4UIcmdWithABool("/process/eLoss/integral",this);
147  IntegCmd->SetGuidance("Switch true/false the integral option");
148  IntegCmd->SetParameterName("integ",true);
149  IntegCmd->SetDefaultValue(true);
151 
152  rangeCmd = new G4UIcmdWithABool("/process/eLoss/CSDARange",this);
153  rangeCmd->SetGuidance("Switch true/false the CSDA range calculation");
154  rangeCmd->SetParameterName("range",true);
155  rangeCmd->SetDefaultValue(true);
157 
158  lpmCmd = new G4UIcmdWithABool("/process/eLoss/LPM",this);
159  lpmCmd->SetGuidance("The flag of the LPM effect calculation");
160  lpmCmd->SetParameterName("lpm",true);
161  lpmCmd->SetDefaultValue(true);
163 
164  splCmd = new G4UIcmdWithABool("/process/em/spline",this);
165  splCmd->SetGuidance("The flag of usage spline for Physics Vectors");
166  splCmd->SetParameterName("spl",true);
167  splCmd->SetDefaultValue(false);
169 
170  aplCmd = new G4UIcmdWithABool("/process/em/applyCuts",this);
171  aplCmd->SetGuidance("The flag to Apply Cuts for gamma processes");
172  aplCmd->SetParameterName("apl",true);
173  aplCmd->SetDefaultValue(false);
175 
176  deCmd = new G4UIcmdWithABool("/process/em/fluo",this);
177  deCmd->SetGuidance("The flag to enable/disable deexcitation");
178  deCmd->SetParameterName("fluoFlag",true);
179  deCmd->SetDefaultValue(false);
181 
182  auCmd = new G4UIcmdWithABool("/process/em/auger",this);
183  auCmd->SetGuidance("The flag to enable/disable Auger electrons");
184  auCmd->SetParameterName("augerFlag",true);
185  auCmd->SetDefaultValue(false);
187 
188  pixeCmd = new G4UIcmdWithABool("/process/em/pixe",this);
189  pixeCmd->SetGuidance("The flag to enable/disable PIXE");
190  pixeCmd->SetParameterName("pixeFlag",true);
191  pixeCmd->SetDefaultValue(false);
193 
194  pixeXsCmd = new G4UIcmdWithAString("/process/em/pixeXSmodel",this);
195  pixeXsCmd->SetGuidance("The name of PIXE cross section");
196  pixeXsCmd->SetParameterName("pixeXS",true);
198 
199  pixeeXsCmd = new G4UIcmdWithAString("/process/em/pixeElecXSmodel",this);
200  pixeeXsCmd->SetGuidance("The name of PIXE cross section for electron");
201  pixeeXsCmd->SetParameterName("pixeEXS",true);
203 
204  deexCmd = new G4UIcommand("/process/em/deexcitation",this);
205  deexCmd->SetGuidance("Set deexcitation flags per G4Region.");
206  deexCmd->SetGuidance(" regName : G4Region name");
207  deexCmd->SetGuidance(" flagFluo : Fluorescence");
208  deexCmd->SetGuidance(" flagAuger : Auger");
209  deexCmd->SetGuidance(" flagPIXE : PIXE");
210 
211  G4UIparameter* regName = new G4UIparameter("regName",'s',false);
212  deexCmd->SetParameter(regName);
213 
214  G4UIparameter* flagFluo = new G4UIparameter("flagFluo",'s',false);
215  deexCmd->SetParameter(flagFluo);
216 
217  G4UIparameter* flagAuger = new G4UIparameter("flagAuger",'s',false);
218  deexCmd->SetParameter(flagAuger);
219 
220  G4UIparameter* flagPIXE = new G4UIparameter("flagPIXE",'s',false);
221  deexCmd->SetParameter(flagPIXE);
222 
223  dedxCmd = new G4UIcmdWithAnInteger("/process/eLoss/binsDEDX",this);
224  dedxCmd->SetGuidance("Set number of bins for DEDX tables");
225  dedxCmd->SetParameterName("binsDEDX",true);
228 
229  lamCmd = new G4UIcmdWithAnInteger("/process/eLoss/binsLambda",this);
230  lamCmd->SetGuidance("Set number of bins for Lambda tables");
231  lamCmd->SetParameterName("binsL",true);
232  lamCmd->SetDefaultValue(77);
234 
235  verCmd = new G4UIcmdWithAnInteger("/process/eLoss/verbose",this);
236  verCmd->SetGuidance("Set verbose level for EM physics");
237  verCmd->SetParameterName("verb",true);
240 
241  ver1Cmd = new G4UIcmdWithAnInteger("/process/em/verbose",this);
242  ver1Cmd->SetGuidance("Set verbose level for EM physics");
243  ver1Cmd->SetParameterName("verb1",true);
246 
247  ver2Cmd = new G4UIcmdWithAnInteger("/process/em/workerVerbose",this);
248  ver2Cmd->SetGuidance("Set worker verbose level for EM physics");
249  ver2Cmd->SetParameterName("verb2",true);
252 
253  lllCmd = new G4UIcmdWithADouble("/process/eLoss/linLossLimit",this);
254  lllCmd->SetGuidance("Set linearLossLimit parameter");
255  lllCmd->SetParameterName("linlim",true);
257 
258  labCmd = new G4UIcmdWithADouble("/process/eLoss/LambdaFactor",this);
259  labCmd->SetGuidance("Set lambdaFactor parameter for integral option");
260  labCmd->SetParameterName("Fl",true);
262 
263  mscCmd = new G4UIcmdWithAString("/process/msc/StepLimit",this);
264  mscCmd->SetGuidance("Set msc step limitation type");
265  mscCmd->SetParameterName("StepLim",true);
267 
268  latCmd = new G4UIcmdWithABool("/process/msc/LateralDisplacement",this);
269  latCmd->SetGuidance("Set flag of sampling of lateral displacement");
270  latCmd->SetParameterName("lat",true);
271  latCmd->SetDefaultValue(true);
273 
274  frCmd = new G4UIcmdWithADouble("/process/msc/RangeFactor",this);
275  frCmd->SetGuidance("Set RangeFactor parameter for msc processes");
276  frCmd->SetParameterName("Fr",true);
277  frCmd->SetRange("Fr>0");
278  frCmd->SetDefaultValue(0.04);
280 
281  fgCmd = new G4UIcmdWithADouble("/process/msc/GeomFactor",this);
282  fgCmd->SetGuidance("Set GeomFactor parameter for msc processes");
283  fgCmd->SetParameterName("Fg",true);
284  fgCmd->SetRange("Fg>0");
285  fgCmd->SetDefaultValue(3.5);
287 
288  mscfCmd = new G4UIcmdWithADouble("/process/msc/FactorForAngleLimit",this);
289  mscfCmd->SetGuidance("Set factor for computation of a limit for -t (invariant trasfer)");
290  mscfCmd->SetParameterName("Fact",true);
291  mscfCmd->SetRange("Fact>0");
294 
295  skinCmd = new G4UIcmdWithADouble("/process/msc/Skin",this);
296  skinCmd->SetGuidance("Set skin parameter for msc processes");
297  skinCmd->SetParameterName("skin",true);
299 
300  angCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/ThetaLimit",this);
301  angCmd->SetGuidance("Set the limit on the polar angle for msc and single scattering");
302  angCmd->SetParameterName("theta",true);
303  angCmd->SetUnitCategory("Angle");
305 
306  bfCmd = new G4UIcommand("/process/em/setBiasingFactor",this);
307  bfCmd->SetGuidance("Set factor for the process cross section.");
308  bfCmd->SetGuidance(" procName : process name");
309  bfCmd->SetGuidance(" procFact : factor");
310  bfCmd->SetGuidance(" flagFact : flag to change weight");
311 
312  G4UIparameter* procName = new G4UIparameter("procName",'s',false);
313  bfCmd->SetParameter(procName);
314 
315  G4UIparameter* procFact = new G4UIparameter("procFact",'d',false);
316  bfCmd->SetParameter(procFact);
317 
318  G4UIparameter* flagFact = new G4UIparameter("flagFact",'s',false);
319  bfCmd->SetParameter(flagFact);
321 
322  fiCmd = new G4UIcommand("/process/em/setForcedInteraction",this);
323  fiCmd->SetGuidance("Set factor for the process cross section.");
324  fiCmd->SetGuidance(" procNam : process name");
325  fiCmd->SetGuidance(" regNam : region name");
326  fiCmd->SetGuidance(" tlength : fixed target length");
327  fiCmd->SetGuidance(" tflag : flag to change weight");
328 
329  G4UIparameter* procNam = new G4UIparameter("procNam",'s',false);
330  fiCmd->SetParameter(procNam);
331 
332  G4UIparameter* regNam = new G4UIparameter("regNam",'s',false);
333  fiCmd->SetParameter(regNam);
334 
335  G4UIparameter* tlength = new G4UIparameter("tlength",'d',false);
336  fiCmd->SetParameter(tlength);
337 
338  G4UIparameter* unitT = new G4UIparameter("unitT",'s',true);
339  fiCmd->SetParameter(unitT);
340  unitT->SetGuidance("unit of tlength");
341 
342  G4UIparameter* flagT = new G4UIparameter("tflag",'s',true);
343  fiCmd->SetParameter(flagT);
345 
346  brCmd = new G4UIcommand("/process/em/setSecBiasing",this);
347  brCmd->SetGuidance("Set bremsstrahlung or delta-electron splitting/Russian roullette per region.");
348  brCmd->SetGuidance(" bProcNam : process name");
349  brCmd->SetGuidance(" bRegNam : region name");
350  brCmd->SetGuidance(" bFactor : number of splitted gamma or probability of Russian roulette");
351  brCmd->SetGuidance(" bEnergy : max energy of a secondary for this biasing method");
352 
353  G4UIparameter* bProcNam = new G4UIparameter("bProcNam",'s',false);
354  brCmd->SetParameter(bProcNam);
355 
356  G4UIparameter* bRegNam = new G4UIparameter("bRegNam",'s',false);
357  brCmd->SetParameter(bRegNam);
358 
359  G4UIparameter* bFactor = new G4UIparameter("bFactor",'d',false);
360  brCmd->SetParameter(bFactor);
361 
362  G4UIparameter* bEnergy = new G4UIparameter("bEnergy",'d',false);
363  brCmd->SetParameter(bEnergy);
364 
365  G4UIparameter* bUnit = new G4UIparameter("bUnit",'s',true);
366  brCmd->SetParameter(bUnit);
367  brCmd->SetGuidance("unit of energy");
368 
370 }
371 
372 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
373 
375 {
376  delete opt;
377  delete RndmStepCmd;
378  delete EnlossFlucCmd;
379  delete SubSecCmd;
380  delete MinSubSecCmd;
381  delete StepFuncCmd;
382  delete deexCmd;
383  delete eLossDirectory;
384  delete mscDirectory;
385  delete emDirectory;
386  delete MinEnCmd;
387  delete MaxEnCmd;
388  delete IntegCmd;
389  delete rangeCmd;
390  delete lpmCmd;
391  delete splCmd;
392  delete aplCmd;
393  delete latCmd;
394  delete verCmd;
395  delete ver1Cmd;
396  delete ver2Cmd;
397  delete mscCmd;
398  delete dedxCmd;
399  delete deCmd;
400  delete auCmd;
401  delete pixeCmd;
402  delete pixeXsCmd;
403  delete pixeeXsCmd;
404  delete frCmd;
405  delete fgCmd;
406  delete lllCmd;
407  delete lamCmd;
408  delete labCmd;
409  delete skinCmd;
410  delete angCmd;
411  delete mscfCmd;
412  delete bfCmd;
413  delete fiCmd;
414  delete brCmd;
415 }
416 
417 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
418 
420 {
421  if(!opt) { opt = new G4EmProcessOptions(); }
422 
423  if (command == RndmStepCmd) {
425  } else if (command == EnlossFlucCmd) {
427  } else if(command == SubSecCmd) {
429  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
430  } else if (command == MinSubSecCmd) {
432  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
433  } else if (command == StepFuncCmd) {
434  G4double v1,v2;
435  G4String unt;
436  std::istringstream is(newValue);
437  is >> v1 >> v2 >> unt;
438  v2 *= G4UIcommand::ValueOf(unt);
439  opt->SetStepFunction(v1,v2);
440  } else if (command == deexCmd) {
441  G4String s1 (""), s2(""), s3(""), s4("");
442  G4bool b2(false), b3(false), b4(false);
443  std::istringstream is(newValue);
444  is >> s1 >> s2 >> s3 >> s4;
445  if(s2 == "true") { b2 = true; }
446  if(s3 == "true") { b3 = true; }
447  if(s4 == "true") { b4 = true; }
449  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
450  } else if (command == deCmd) {
451  opt->SetFluo(deCmd->GetNewBoolValue(newValue));
452  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
453  } else if (command == auCmd) {
454  opt->SetAuger(auCmd->GetNewBoolValue(newValue));
455  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
456  } else if (command == pixeCmd) {
457  opt->SetPIXE(pixeCmd->GetNewBoolValue(newValue));
458  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
459  } else if (command == pixeXsCmd) {
460  G4String name;
461  if (newValue == "ecpssr_analytical")
462  {name = "ECPSSR_Analytical";}
463  else if (newValue == "ecpssr_interpolated")
464  {name = "ECPSSR_FormFactor";}
465  else
466  {name = newValue;}
468  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
469  } else if (command == pixeeXsCmd) {
471  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
472  } else if (command == mscCmd) {
473  if(newValue == "Minimal")
475 
476  else if(newValue == "UseDistanceToBoundary")
478 
479  else if(newValue == "UseSafety")
481 
482  else {
483  G4cout << "### G4EnergyLossMessenger WARNING: StepLimit type <"
484  << newValue << "> unknown!" << G4endl;
485  return;
486  }
487  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
488  } else if (command == MinEnCmd) {
490  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
491  } else if (command == MaxEnCmd) {
493  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
494  } else if (command == IntegCmd) {
496  } else if (command == rangeCmd) {
498  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
499  } else if (command == lpmCmd) {
500  opt->SetLPMFlag(lpmCmd->GetNewBoolValue(newValue));
501  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
502  } else if (command == splCmd) {
504  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
505  } else if (command == aplCmd) {
506  opt->SetApplyCuts(aplCmd->GetNewBoolValue(newValue));
507  } else if (command == latCmd) {
509  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
510  } else if (command == verCmd) {
511  opt->SetVerbose(verCmd->GetNewIntValue(newValue),"all",false);
512  } else if (command == ver1Cmd) {
513  opt->SetVerbose(ver1Cmd->GetNewIntValue(newValue),"all",false);
514  } else if (command == ver2Cmd) {
515  opt->SetVerbose(ver2Cmd->GetNewIntValue(newValue),"all",true);
516  } else if (command == lllCmd) {
518  } else if (command == labCmd) {
520  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
521  } else if (command == skinCmd) {
522  opt->SetSkin(skinCmd->GetNewDoubleValue(newValue));
523  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
524  } else if (command == dedxCmd) {
526  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
527  } else if (command == lamCmd) {
529  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
530  } else if (command == frCmd) {
532  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
533  } else if (command == fgCmd) {
535  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
536  } else if (command == mscfCmd) {
538  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
539  } else if (command == angCmd) {
541  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
542  } else if (command == bfCmd) {
543  G4double v1(1.0);
544  G4String s0(""),s1("");
545  std::istringstream is(newValue);
546  is >> s0 >> v1 >> s1;
547  G4bool yes = false;
548  if(s1 == "true") { yes = true; }
549  opt->SetProcessBiasingFactor(s0,v1,yes);
550  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
551  } else if (command == fiCmd) {
552  G4double v1(0.0);
553  G4String s1(""),s2(""),s3(""),unt("mm");
554  std::istringstream is(newValue);
555  is >> s1 >> s2 >> v1 >> unt >> s3;
556  G4bool yes = false;
557  if(s3 == "true") { yes = true; }
558  v1 *= G4UIcommand::ValueOf(unt);
559  opt->ActivateForcedInteraction(s1,v1,s2,yes);
560  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
561  } else if (command == brCmd) {
562  G4double fb(1.0),en(1.e+30);
563  G4String s1(""),s2(""),unt("MeV");
564  std::istringstream is(newValue);
565  is >> s1 >> s2 >> fb >> en >> unt;
566  en *= G4UIcommand::ValueOf(unt);
567  if (s1=="phot"||s1=="compt"||s1=="conv")
569  else opt->ActivateSecondaryBiasing(s1,s2,fb,en);
570  G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
571  }
572 }
573 
574 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:152
void SetSkin(G4double val)
G4UIcmdWithADouble * mscfCmd
void SetLambdaFactor(G4double val)
G4UIcmdWithAnInteger * verCmd
void SetPIXE(G4bool val)
G4UIcmdWithADouble * skinCmd
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
G4UIcmdWithADouble * labCmd
static G4int GetNewIntValue(const char *paramString)
G4UIcmdWithAnInteger * lamCmd
void SetParameterRange(const char *theRange)
void SetSplineFlag(G4bool val)
void SetParameterCandidates(const char *theString)
G4String name
Definition: TRTMaterials.hh:40
G4UIcmdWithADoubleAndUnit * MaxEnCmd
void SetMinEnergy(G4double val)
void SetNewValue(G4UIcommand *, G4String)
void SetDefaultValue(const char *theDefaultValue)
void SetMscGeomFactor(G4double val)
void SetStepFunction(G4double v1, G4double v2)
G4UIcmdWithAnInteger * ver1Cmd
G4UIcmdWithAString * pixeeXsCmd
void ActivateSecondaryBiasingForGamma(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
G4UIcmdWithABool * RndmStepCmd
void SetUnitCategory(const char *unitCategory)
static G4double GetNewDoubleValue(const char *paramString)
void SetMscLateralDisplacement(G4bool val)
static G4bool GetNewBoolValue(const char *paramString)
void SetFluo(G4bool val)
G4UIcmdWithAString * pixeXsCmd
G4UIcmdWithABool * EnlossFlucCmd
void SetMinSubRange(G4double val)
void SetDefaultValue(G4bool defVal)
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:58
void SetAuger(G4bool val)
void SetDEDXBinning(G4int val)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
G4UIcmdWithADouble * lllCmd
G4GLOB_DLL std::ostream G4cout
void SetLambdaBinning(G4int val)
static G4String UnitsList(const char *unitCategory)
Definition: G4UIcommand.cc:306
static const G4double b3
G4UIcmdWithADoubleAndUnit * MinEnCmd
bool G4bool
Definition: G4Types.hh:79
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 SetLossFluctuations(G4bool val)
static const G4double b2
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:225
void SetLinearLossLimit(G4double val)
G4UIcmdWithADouble * frCmd
void SetMaxEnergy(G4double val)
G4UIcmdWithAString * mscCmd
G4UIcmdWithADoubleAndUnit * angCmd
G4UIcmdWithADouble * fgCmd
G4UIcmdWithAnInteger * dedxCmd
void SetDeexcitationActiveRegion(const G4String &rname="", G4bool valDeexcitation=true, G4bool valAuger=true, G4bool valPIXE=true)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetFactorForAngleLimit(G4double val)
static G4double ValueOf(const char *unitName)
Definition: G4UIcommand.cc:294
G4EmProcessOptions * opt
G4UIcmdWithABool * SubSecCmd
G4UIcmdWithAnInteger * ver2Cmd
void ActivateSecondaryBiasing(const G4String &name, const G4String &region, G4double factor, G4double energyLimit)
void ActivateForcedInteraction(const G4String &name, G4double length=0.0, const G4String &region="", G4bool flag=true)
void SetDefaultValue(G4double defVal)
void SetPIXECrossSectionModel(const G4String &val)
void SetBuildCSDARange(G4bool val)
void SetIntegral(G4bool val)
#define G4endl
Definition: G4ios.hh:61
void SetMscStepLimitation(G4MscStepLimitType val)
void SetDefaultValue(G4int defVal)
void SetRandomStep(G4bool val)
double G4double
Definition: G4Types.hh:76
void SetGuidance(const char *theGuidance)
void SetMscRangeFactor(G4double val)
void SetSubCutoff(G4bool val, const G4Region *r=0)
static const G4double b4
G4UIcmdWithADouble * MinSubSecCmd
static G4String CategoryOf(const char *unitName)
Definition: G4UIcommand.cc:301
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetLPMFlag(G4bool val)
void SetProcessBiasingFactor(const G4String &name, G4double val, G4bool flag=true)
void SetPIXEElectronCrossSectionModel(const G4String &val)
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:432
void SetApplyCuts(G4bool val)
void SetVerbose(G4int val, const G4String &name="all", G4bool worker=false)
void SetPolarAngleLimit(G4double val)