Geant4 Cross Reference |
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 // ------------------------------------------------------------------- 27 // 28 // GEANT4 Class file 29 // 30 // File name: G4EmParametersMessenger 31 // 32 // Author: Vladimir Ivanchenko created from G4EnergyLossMessenger 33 // 34 // Creation date: 22-05-2013 35 // 36 // ------------------------------------------------------------------- 37 // 38 39 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 40 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 41 42 #include "G4EmParametersMessenger.hh" 43 #include "G4UIdirectory.hh" 44 #include "G4UIcommand.hh" 45 #include "G4UIparameter.hh" 46 #include "G4UIcmdWithABool.hh" 47 #include "G4UIcmdWithAnInteger.hh" 48 #include "G4UIcmdWithADouble.hh" 49 #include "G4UIcmdWithADoubleAndUnit.hh" 50 #include "G4UIcmdWithAString.hh" 51 #include "G4UIcmdWith3VectorAndUnit.hh" 52 #include "G4UImanager.hh" 53 #include "G4MscStepLimitType.hh" 54 #include "G4NuclearFormfactorType.hh" 55 #include "G4EmParameters.hh" 56 57 #include <sstream> 58 59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 60 61 G4EmParametersMessenger::G4EmParametersMessenger(G4EmParameters* ptr) 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); 79 flucCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 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); 86 rangeCmd->AvailableForStates(G4State_PreInit); 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); 93 lpmCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 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); 100 rsCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 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); 107 aplCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 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); 114 intCmd->AvailableForStates(G4State_PreInit); 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); 121 latCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 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); 128 lat96Cmd->AvailableForStates(G4State_PreInit,G4State_Idle); 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); 135 mulatCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 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); 142 delCmd->AvailableForStates(G4State_PreInit); 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); 149 mottCmd->AvailableForStates(G4State_PreInit); 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); 156 birksCmd->AvailableForStates(G4State_PreInit,G4State_Init); 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); 163 sharkCmd->AvailableForStates(G4State_PreInit); 164 sharkCmd->SetToBeBroadcasted(false); 165 166 poCmd = new G4UIcmdWithABool("/process/em/Polarisation",this); 167 poCmd->SetGuidance("Enable polarisation"); 168 poCmd->AvailableForStates(G4State_PreInit); 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); 175 sampleTCmd->AvailableForStates(G4State_PreInit); 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); 182 icru90Cmd->AvailableForStates(G4State_PreInit); 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); 189 mudatCmd->AvailableForStates(G4State_PreInit); 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); 196 peKCmd->AvailableForStates(G4State_PreInit); 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); 203 mscPCmd->AvailableForStates(G4State_PreInit, G4State_Idle); 204 mscPCmd->SetToBeBroadcasted(false); 205 206 pepicsCmd = new G4UIcmdWithABool("/process/em/UseEPICS2017XS",this); 207 pepicsCmd->SetGuidance("Use EPICS2017 data for gamma x-ections"); 208 pepicsCmd->SetParameterName("pepics",true); 209 pepicsCmd->SetDefaultValue(false); 210 pepicsCmd->AvailableForStates(G4State_PreInit); 211 pepicsCmd->SetToBeBroadcasted(false); 212 213 f3gCmd = new G4UIcmdWithABool("/process/em/3GammaAnnihilationOnFly",this); 214 f3gCmd->SetGuidance("Enable/disable 3 gamma annihilation on fly"); 215 f3gCmd->SetParameterName("f3gamma",true); 216 f3gCmd->SetDefaultValue(false); 217 f3gCmd->AvailableForStates(G4State_PreInit); 218 f3gCmd->SetToBeBroadcasted(false); 219 220 minEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/minKinEnergy",this); 221 minEnCmd->SetGuidance("Set the min kinetic energy for EM tables"); 222 minEnCmd->SetParameterName("emin",true); 223 minEnCmd->SetUnitCategory("Energy"); 224 minEnCmd->AvailableForStates(G4State_PreInit); 225 minEnCmd->SetToBeBroadcasted(false); 226 227 maxEnCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergy",this); 228 maxEnCmd->SetGuidance("Set the max kinetic energy for EM tables"); 229 maxEnCmd->SetParameterName("emax",true); 230 maxEnCmd->SetUnitCategory("Energy"); 231 maxEnCmd->AvailableForStates(G4State_PreInit); 232 maxEnCmd->SetToBeBroadcasted(false); 233 234 cenCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/maxKinEnergyCSDA",this); 235 cenCmd->SetGuidance("Set the max kinetic energy for CSDA table"); 236 cenCmd->SetParameterName("emaxCSDA",true); 237 cenCmd->SetUnitCategory("Energy"); 238 cenCmd->AvailableForStates(G4State_PreInit); 239 cenCmd->SetToBeBroadcasted(false); 240 241 max5DCmd = new G4UIcmdWithADoubleAndUnit("/process/em/max5DMuPairEnergy",this); 242 max5DCmd->SetGuidance("Set the max kinetic energy for 5D muon pair production"); 243 max5DCmd->SetParameterName("emax5D",true); 244 max5DCmd->SetUnitCategory("Energy"); 245 max5DCmd->AvailableForStates(G4State_PreInit); 246 max5DCmd->SetToBeBroadcasted(false); 247 248 lowEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestElectronEnergy",this); 249 lowEnCmd->SetGuidance("Set the lowest kinetic energy for e+-"); 250 lowEnCmd->SetParameterName("elow",true); 251 lowEnCmd->SetUnitCategory("Energy"); 252 lowEnCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 253 lowEnCmd->SetToBeBroadcasted(false); 254 255 lowhEnCmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestMuHadEnergy",this); 256 lowhEnCmd->SetGuidance("Set the lowest kinetic energy for muons and hadrons"); 257 lowhEnCmd->SetParameterName("elowh",true); 258 lowhEnCmd->SetUnitCategory("Energy"); 259 lowhEnCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 260 lowhEnCmd->SetToBeBroadcasted(false); 261 262 lowEn3Cmd = new G4UIcmdWithADoubleAndUnit("/process/em/lowestTripletEnergy",this); 263 lowEn3Cmd->SetGuidance("Set the lowest kinetic energy for triplet production"); 264 lowEn3Cmd->SetParameterName("elow3",true); 265 lowEn3Cmd->SetUnitCategory("Energy"); 266 lowEn3Cmd->AvailableForStates(G4State_PreInit,G4State_Idle); 267 lowEn3Cmd->SetToBeBroadcasted(false); 268 269 lllCmd = new G4UIcmdWithADouble("/process/eLoss/linLossLimit",this); 270 lllCmd->SetGuidance("Set linearLossLimit parameter"); 271 lllCmd->SetParameterName("linlim",true); 272 lllCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 273 lllCmd->SetToBeBroadcasted(false); 274 275 brCmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/bremThreshold",this); 276 brCmd->SetGuidance("Set e+- bremsstrahlung energy threshold"); 277 brCmd->SetParameterName("emaxBrem",true); 278 brCmd->SetUnitCategory("Energy"); 279 brCmd->AvailableForStates(G4State_PreInit); 280 brCmd->SetToBeBroadcasted(false); 281 282 br1Cmd = new G4UIcmdWithADoubleAndUnit("/process/eLoss/bremMuHadThreshold",this); 283 br1Cmd->SetGuidance("Set muon/hadron bremsstrahlung energy threshold"); 284 br1Cmd->SetParameterName("emaxMuHadBrem",true); 285 br1Cmd->SetUnitCategory("Energy"); 286 br1Cmd->AvailableForStates(G4State_PreInit); 287 br1Cmd->SetToBeBroadcasted(false); 288 289 labCmd = new G4UIcmdWithADouble("/process/eLoss/LambdaFactor",this); 290 labCmd->SetGuidance("Set lambdaFactor parameter for integral option"); 291 labCmd->SetParameterName("Fl",true); 292 labCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 293 labCmd->SetToBeBroadcasted(false); 294 295 mscfCmd = new G4UIcmdWithADouble("/process/msc/FactorForAngleLimit",this); 296 mscfCmd->SetGuidance("Set factor for computation of a limit for -t (invariant transfer)"); 297 mscfCmd->SetParameterName("Fact",true); 298 mscfCmd->SetRange("Fact>0"); 299 mscfCmd->SetDefaultValue(1.); 300 mscfCmd->AvailableForStates(G4State_PreInit); 301 mscfCmd->SetToBeBroadcasted(false); 302 303 angCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/ThetaLimit",this); 304 angCmd->SetGuidance("Set the limit on the polar angle for msc and single scattering"); 305 angCmd->SetParameterName("theta",true); 306 angCmd->SetUnitCategory("Angle"); 307 angCmd->AvailableForStates(G4State_PreInit); 308 angCmd->SetToBeBroadcasted(false); 309 310 msceCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/EnergyLimit",this); 311 msceCmd->SetGuidance("Set the upper energy limit for msc"); 312 msceCmd->SetParameterName("mscE",true); 313 msceCmd->SetUnitCategory("Energy"); 314 msceCmd->AvailableForStates(G4State_PreInit); 315 msceCmd->SetToBeBroadcasted(false); 316 317 nielCmd = new G4UIcmdWithADoubleAndUnit("/process/em/MaxEnergyNIEL",this); 318 nielCmd->SetGuidance("Set the upper energy limit for NIEL"); 319 nielCmd->SetParameterName("niel",true); 320 nielCmd->SetUnitCategory("Energy"); 321 nielCmd->AvailableForStates(G4State_PreInit); 322 nielCmd->SetToBeBroadcasted(false); 323 324 frCmd = new G4UIcmdWithADouble("/process/msc/RangeFactor",this); 325 frCmd->SetGuidance("Set RangeFactor for msc processes of e+-"); 326 frCmd->SetParameterName("Fr",true); 327 frCmd->SetRange("Fr>0"); 328 frCmd->SetDefaultValue(0.04); 329 frCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 330 frCmd->SetToBeBroadcasted(false); 331 332 fr1Cmd = new G4UIcmdWithADouble("/process/msc/RangeFactorMuHad",this); 333 fr1Cmd->SetGuidance("Set RangeFactor for msc processes of muons/hadrons"); 334 fr1Cmd->SetParameterName("Fr1",true); 335 fr1Cmd->SetRange("Fr1>0"); 336 fr1Cmd->SetDefaultValue(0.2); 337 fr1Cmd->AvailableForStates(G4State_PreInit,G4State_Idle); 338 fr1Cmd->SetToBeBroadcasted(false); 339 340 fgCmd = new G4UIcmdWithADouble("/process/msc/GeomFactor",this); 341 fgCmd->SetGuidance("Set GeomFactor parameter for msc processes"); 342 fgCmd->SetParameterName("Fg",true); 343 fgCmd->SetRange("Fg>0"); 344 fgCmd->SetDefaultValue(2.5); 345 fgCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 346 fgCmd->SetToBeBroadcasted(false); 347 348 skinCmd = new G4UIcmdWithADouble("/process/msc/Skin",this); 349 skinCmd->SetGuidance("Set skin parameter for msc processes"); 350 skinCmd->SetParameterName("skin",true); 351 skinCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 352 skinCmd->SetToBeBroadcasted(false); 353 354 screCmd = new G4UIcmdWithADouble("/process/msc/ScreeningFactor",this); 355 screCmd->SetGuidance("Set screening factor"); 356 screCmd->SetParameterName("screen",true); 357 screCmd->AvailableForStates(G4State_PreInit); 358 screCmd->SetToBeBroadcasted(false); 359 360 safCmd = new G4UIcmdWithADouble("/process/msc/SafetyFactor",this); 361 safCmd->SetGuidance("Set safety factor"); 362 safCmd->SetParameterName("fsafe",true); 363 safCmd->AvailableForStates(G4State_PreInit); 364 safCmd->SetToBeBroadcasted(false); 365 366 llimCmd = new G4UIcmdWithADoubleAndUnit("/process/msc/LambdaLimit",this); 367 llimCmd->SetGuidance("Set the upper energy limit for NIEL"); 368 llimCmd->SetParameterName("ll",true); 369 llimCmd->SetUnitCategory("Length"); 370 llimCmd->AvailableForStates(G4State_PreInit); 371 llimCmd->SetToBeBroadcasted(false); 372 373 amCmd = new G4UIcmdWithAnInteger("/process/em/binsPerDecade",this); 374 amCmd->SetGuidance("Set number of bins per decade for EM tables"); 375 amCmd->SetParameterName("bins",true); 376 amCmd->SetDefaultValue(7); 377 amCmd->AvailableForStates(G4State_PreInit); 378 amCmd->SetToBeBroadcasted(false); 379 380 verCmd = new G4UIcmdWithAnInteger("/process/eLoss/verbose",this); 381 verCmd->SetGuidance("Set verbose level for EM physics"); 382 verCmd->SetParameterName("verb",true); 383 verCmd->SetDefaultValue(1); 384 verCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 385 verCmd->SetToBeBroadcasted(false); 386 387 ver1Cmd = new G4UIcmdWithAnInteger("/process/em/verbose",this); 388 ver1Cmd->SetGuidance("Set verbose level for EM physics"); 389 ver1Cmd->SetParameterName("verb1",true); 390 ver1Cmd->SetDefaultValue(1); 391 ver1Cmd->AvailableForStates(G4State_PreInit,G4State_Idle); 392 ver1Cmd->SetToBeBroadcasted(false); 393 394 ver2Cmd = new G4UIcmdWithAnInteger("/process/em/workerVerbose",this); 395 ver2Cmd->SetGuidance("Set worker verbose level for EM physics"); 396 ver2Cmd->SetParameterName("verb2",true); 397 ver2Cmd->SetDefaultValue(0); 398 ver2Cmd->AvailableForStates(G4State_PreInit,G4State_Idle); 399 ver2Cmd->SetToBeBroadcasted(false); 400 401 nFreeCmd = new G4UIcmdWithAnInteger("/process/em/nForFreeVector",this); 402 nFreeCmd->SetGuidance("Set number for logarithmic bin search algorithm"); 403 nFreeCmd->SetParameterName("nFree",true); 404 nFreeCmd->SetDefaultValue(2); 405 nFreeCmd->AvailableForStates(G4State_PreInit); 406 nFreeCmd->SetToBeBroadcasted(false); 407 408 transWithMscCmd = new G4UIcmdWithAString("/process/em/transportationWithMsc",this); 409 transWithMscCmd->SetGuidance("Enable/disable the G4TransportationWithMsc process"); 410 transWithMscCmd->SetParameterName("trans",true); 411 transWithMscCmd->SetCandidates("Disabled Enabled MultipleSteps"); 412 transWithMscCmd->AvailableForStates(G4State_PreInit); 413 transWithMscCmd->SetToBeBroadcasted(false); 414 415 mscCmd = new G4UIcmdWithAString("/process/msc/StepLimit",this); 416 mscCmd->SetGuidance("Set msc step limitation type"); 417 mscCmd->SetParameterName("StepLim",true); 418 mscCmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary"); 419 mscCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 420 mscCmd->SetToBeBroadcasted(false); 421 422 msc1Cmd = new G4UIcmdWithAString("/process/msc/StepLimitMuHad",this); 423 msc1Cmd->SetGuidance("Set msc step limitation type for muons/hadrons"); 424 msc1Cmd->SetParameterName("StepLim1",true); 425 msc1Cmd->SetCandidates("Minimal UseSafety UseSafetyPlus UseDistanceToBoundary"); 426 msc1Cmd->AvailableForStates(G4State_PreInit,G4State_Idle); 427 msc1Cmd->SetToBeBroadcasted(false); 428 429 dumpCmd = new G4UIcommand("/process/em/printParameters",this); 430 dumpCmd->SetGuidance("Print all EM parameters."); 431 dumpCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 432 dumpCmd->SetToBeBroadcasted(false); 433 434 nffCmd = new G4UIcmdWithAString("/process/em/setNuclearFormFactor",this); 435 nffCmd->SetGuidance("Define type of nuclear form-factor"); 436 nffCmd->SetParameterName("NucFF",true); 437 nffCmd->SetCandidates("None Exponential Gaussian Flat"); 438 nffCmd->AvailableForStates(G4State_PreInit); 439 nffCmd->SetToBeBroadcasted(false); 440 441 ssCmd = new G4UIcmdWithAString("/process/em/setSingleScattering",this); 442 ssCmd->SetGuidance("Define type of e+- single scattering model"); 443 ssCmd->SetParameterName("SS",true); 444 ssCmd->SetCandidates("WVI Mott DPWA"); 445 ssCmd->AvailableForStates(G4State_PreInit); 446 ssCmd->SetToBeBroadcasted(false); 447 448 fluc1Cmd = new G4UIcmdWithAString("/process/eLoss/setFluctModel",this); 449 fluc1Cmd->SetGuidance("Define type of energy loss fluctuation model"); 450 fluc1Cmd->SetParameterName("Fluc1",true); 451 fluc1Cmd->SetCandidates("Dummy Universal Urban"); 452 fluc1Cmd->AvailableForStates(G4State_PreInit); 453 fluc1Cmd->SetToBeBroadcasted(false); 454 455 posiCmd = new G4UIcmdWithAString("/process/em/setPositronAtRestModel",this); 456 posiCmd->SetGuidance("Define model of positron annihilation at rest"); 457 posiCmd->SetParameterName("Posi",true); 458 posiCmd->SetCandidates("Simple Allison OrePawell OrePowellPolar"); 459 posiCmd->AvailableForStates(G4State_PreInit); 460 posiCmd->SetToBeBroadcasted(false); 461 462 tripletCmd = new G4UIcmdWithAnInteger("/process/gconv/conversionType",this); 463 tripletCmd->SetGuidance("gamma conversion triplet/nuclear generation type:"); 464 tripletCmd->SetGuidance("0 - (default) both triplet and nuclear"); 465 tripletCmd->SetGuidance("1 - force nuclear"); 466 tripletCmd->SetGuidance("2 - force triplet"); 467 tripletCmd->SetParameterName("type",false); 468 tripletCmd->SetRange("type >= 0 && type <= 2"); 469 tripletCmd->SetDefaultValue(0); 470 tripletCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 471 tripletCmd->SetToBeBroadcasted(false); 472 473 onIsolatedCmd = new G4UIcmdWithABool("/process/gconv/onIsolated",this); 474 onIsolatedCmd->SetGuidance("Conversion on isolated charged particles"); 475 onIsolatedCmd->SetGuidance("false (default) : atomic electron screening"); 476 onIsolatedCmd->SetGuidance("true : conversion on isolated particles."); 477 onIsolatedCmd->SetParameterName("flag",false); 478 onIsolatedCmd->SetDefaultValue(false); 479 onIsolatedCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 480 onIsolatedCmd->SetToBeBroadcasted(false); 481 } 482 483 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 484 485 G4EmParametersMessenger::~G4EmParametersMessenger() 486 { 487 delete gconvDirectory; 488 delete eLossDirectory; 489 delete mscDirectory; 490 delete emDirectory; 491 delete dnaDirectory; 492 493 delete flucCmd; 494 delete rangeCmd; 495 delete lpmCmd; 496 delete rsCmd; 497 delete aplCmd; 498 delete intCmd; 499 delete latCmd; 500 delete lat96Cmd; 501 delete mulatCmd; 502 delete delCmd; 503 delete mottCmd; 504 delete birksCmd; 505 delete sharkCmd; 506 delete onIsolatedCmd; 507 delete sampleTCmd; 508 delete poCmd; 509 delete icru90Cmd; 510 delete mudatCmd; 511 delete peKCmd; 512 delete f3gCmd; 513 delete mscPCmd; 514 delete pepicsCmd; 515 516 delete minEnCmd; 517 delete maxEnCmd; 518 delete max5DCmd; 519 delete cenCmd; 520 delete lowEnCmd; 521 delete lowhEnCmd; 522 delete lowEn3Cmd; 523 delete lllCmd; 524 delete brCmd; 525 delete br1Cmd; 526 delete labCmd; 527 delete mscfCmd; 528 delete angCmd; 529 delete msceCmd; 530 delete nielCmd; 531 delete frCmd; 532 delete fr1Cmd; 533 delete fgCmd; 534 delete skinCmd; 535 delete safCmd; 536 delete llimCmd; 537 delete screCmd; 538 539 delete amCmd; 540 delete verCmd; 541 delete ver1Cmd; 542 delete ver2Cmd; 543 delete transWithMscCmd; 544 delete nFreeCmd; 545 delete tripletCmd; 546 547 delete mscCmd; 548 delete msc1Cmd; 549 delete nffCmd; 550 delete ssCmd; 551 delete fluc1Cmd; 552 delete posiCmd; 553 554 delete dumpCmd; 555 } 556 557 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 558 559 void G4EmParametersMessenger::SetNewValue(G4UIcommand* command, 560 G4String newValue) 561 { 562 G4bool physicsModified = false; 563 if (command == flucCmd) { 564 theParameters->SetLossFluctuations(flucCmd->GetNewBoolValue(newValue)); 565 physicsModified = true; 566 } else if (command == rangeCmd) { 567 theParameters->SetBuildCSDARange(rangeCmd->GetNewBoolValue(newValue)); 568 } else if (command == lpmCmd) { 569 theParameters->SetLPM(lpmCmd->GetNewBoolValue(newValue)); 570 physicsModified = true; 571 } else if (command == rsCmd) { 572 theParameters->SetUseCutAsFinalRange(rsCmd->GetNewBoolValue(newValue)); 573 physicsModified = true; 574 } else if (command == aplCmd) { 575 theParameters->SetApplyCuts(aplCmd->GetNewBoolValue(newValue)); 576 physicsModified = true; 577 } else if (command == intCmd) { 578 theParameters->SetIntegral(intCmd->GetNewBoolValue(newValue)); 579 } else if (command == latCmd) { 580 theParameters->SetLateralDisplacement(latCmd->GetNewBoolValue(newValue)); 581 physicsModified = true; 582 } else if (command == lat96Cmd) { 583 theParameters->SetLateralDisplacementAlg96(lat96Cmd->GetNewBoolValue(newValue)); 584 physicsModified = true; 585 } else if (command == mulatCmd) { 586 theParameters->SetMuHadLateralDisplacement(mulatCmd->GetNewBoolValue(newValue)); 587 physicsModified = true; 588 } else if (command == delCmd) { 589 theParameters->ActivateAngularGeneratorForIonisation(delCmd->GetNewBoolValue(newValue)); 590 } else if (command == mottCmd) { 591 theParameters->SetUseMottCorrection(mottCmd->GetNewBoolValue(newValue)); 592 } else if (command == birksCmd) { 593 theParameters->SetBirksActive(birksCmd->GetNewBoolValue(newValue)); 594 } else if (command == icru90Cmd) { 595 theParameters->SetUseICRU90Data(icru90Cmd->GetNewBoolValue(newValue)); 596 } else if (command == sharkCmd) { 597 theParameters->SetGeneralProcessActive(sharkCmd->GetNewBoolValue(newValue)); 598 } else if (command == poCmd) { 599 theParameters->SetEnablePolarisation(poCmd->GetNewBoolValue(newValue)); 600 } else if (command == sampleTCmd) { 601 theParameters->SetEnableSamplingTable(sampleTCmd->GetNewBoolValue(newValue)); 602 } else if (command == mudatCmd) { 603 theParameters->SetRetrieveMuDataFromFile(mudatCmd->GetNewBoolValue(newValue)); 604 } else if (command == peKCmd) { 605 theParameters->SetPhotoeffectBelowKShell(peKCmd->GetNewBoolValue(newValue)); 606 } else if (command == f3gCmd) { 607 theParameters->Set3GammaAnnihilationOnFly(f3gCmd->GetNewBoolValue(newValue)); 608 } else if (command == mscPCmd) { 609 theParameters->SetMscPositronCorrection(mscPCmd->GetNewBoolValue(newValue)); 610 } else if (command == pepicsCmd) { 611 theParameters->SetUseEPICS2017XS(pepicsCmd->GetNewBoolValue(newValue)); 612 613 } else if (command == minEnCmd) { 614 theParameters->SetMinEnergy(minEnCmd->GetNewDoubleValue(newValue)); 615 } else if (command == maxEnCmd) { 616 theParameters->SetMaxEnergy(maxEnCmd->GetNewDoubleValue(newValue)); 617 } else if (command == max5DCmd) { 618 theParameters->SetMaxEnergyFor5DMuPair(max5DCmd->GetNewDoubleValue(newValue)); 619 } else if (command == cenCmd) { 620 theParameters->SetMaxEnergyForCSDARange(cenCmd->GetNewDoubleValue(newValue)); 621 physicsModified = true; 622 } else if (command == lowEnCmd) { 623 theParameters->SetLowestElectronEnergy(lowEnCmd->GetNewDoubleValue(newValue)); 624 physicsModified = true; 625 } else if (command == lowEn3Cmd) { 626 theParameters->SetLowestTripletEnergy(lowEn3Cmd->GetNewDoubleValue(newValue)); 627 physicsModified = true; 628 } else if (command == lowhEnCmd) { 629 theParameters->SetLowestMuHadEnergy(lowhEnCmd->GetNewDoubleValue(newValue)); 630 physicsModified = true; 631 } else if (command == lllCmd) { 632 theParameters->SetLinearLossLimit(lllCmd->GetNewDoubleValue(newValue)); 633 physicsModified = true; 634 } else if (command == brCmd) { 635 theParameters->SetBremsstrahlungTh(brCmd->GetNewDoubleValue(newValue)); 636 physicsModified = true; 637 } else if (command == br1Cmd) { 638 theParameters->SetMuHadBremsstrahlungTh(br1Cmd->GetNewDoubleValue(newValue)); 639 physicsModified = true; 640 } else if (command == labCmd) { 641 theParameters->SetLambdaFactor(labCmd->GetNewDoubleValue(newValue)); 642 physicsModified = true; 643 } else if (command == mscfCmd) { 644 theParameters->SetFactorForAngleLimit(mscfCmd->GetNewDoubleValue(newValue)); 645 } else if (command == angCmd) { 646 theParameters->SetMscThetaLimit(angCmd->GetNewDoubleValue(newValue)); 647 } else if (command == msceCmd) { 648 theParameters->SetMscEnergyLimit(msceCmd->GetNewDoubleValue(newValue)); 649 } else if (command == nielCmd) { 650 theParameters->SetMaxNIELEnergy(nielCmd->GetNewDoubleValue(newValue)); 651 } else if (command == frCmd) { 652 theParameters->SetMscRangeFactor(frCmd->GetNewDoubleValue(newValue)); 653 physicsModified = true; 654 } else if (command == fr1Cmd) { 655 theParameters->SetMscMuHadRangeFactor(fr1Cmd->GetNewDoubleValue(newValue)); 656 physicsModified = true; 657 } else if (command == fgCmd) { 658 theParameters->SetMscGeomFactor(fgCmd->GetNewDoubleValue(newValue)); 659 physicsModified = true; 660 } else if (command == skinCmd) { 661 theParameters->SetMscSkin(skinCmd->GetNewDoubleValue(newValue)); 662 physicsModified = true; 663 } else if (command == safCmd) { 664 theParameters->SetMscSafetyFactor(safCmd->GetNewDoubleValue(newValue)); 665 } else if (command == llimCmd) { 666 theParameters->SetMscLambdaLimit(llimCmd->GetNewDoubleValue(newValue)); 667 } else if (command == screCmd) { 668 theParameters->SetScreeningFactor(screCmd->GetNewDoubleValue(newValue)); 669 } else if (command == amCmd) { 670 theParameters->SetNumberOfBinsPerDecade(amCmd->GetNewIntValue(newValue)); 671 } else if (command == verCmd) { 672 theParameters->SetVerbose(verCmd->GetNewIntValue(newValue)); 673 } else if (command == ver1Cmd) { 674 theParameters->SetVerbose(ver1Cmd->GetNewIntValue(newValue)); 675 } else if (command == ver2Cmd) { 676 theParameters->SetWorkerVerbose(ver2Cmd->GetNewIntValue(newValue)); 677 } else if (command == nFreeCmd) { 678 theParameters->SetNumberForFreeVector(nFreeCmd->GetNewIntValue(newValue)); 679 } else if (command == dumpCmd) { 680 theParameters->SetIsPrintedFlag(false); 681 theParameters->Dump(); 682 } else if (command == transWithMscCmd) { 683 G4TransportationWithMscType type = G4TransportationWithMscType::fDisabled; 684 if(newValue == "Disabled") { 685 type = G4TransportationWithMscType::fDisabled; 686 } else if(newValue == "Enabled") { 687 type = G4TransportationWithMscType::fEnabled; 688 } else if(newValue == "MultipleSteps") { 689 type = G4TransportationWithMscType::fMultipleSteps; 690 } else { 691 G4ExceptionDescription ed; 692 ed << " TransportationWithMsc type <" << newValue << "> unknown!"; 693 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed); 694 } 695 theParameters->SetTransportationWithMsc(type); 696 } else if (command == mscCmd || command == msc1Cmd) { 697 G4MscStepLimitType msctype = fUseSafety; 698 if(newValue == "Minimal") { 699 msctype = fMinimal; 700 } else if(newValue == "UseDistanceToBoundary") { 701 msctype = fUseDistanceToBoundary; 702 } else if(newValue == "UseSafety") { 703 msctype = fUseSafety; 704 } else if(newValue == "UseSafetyPlus") { 705 msctype = fUseSafetyPlus; 706 } else { 707 G4ExceptionDescription ed; 708 ed << " StepLimit type <" << newValue << "> unknown!"; 709 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed); 710 return; 711 } 712 if (command == mscCmd) { 713 theParameters->SetMscStepLimitType(msctype); 714 } else { 715 theParameters->SetMscMuHadStepLimitType(msctype); 716 } 717 physicsModified = true; 718 } else if (command == nffCmd) { 719 G4NuclearFormfactorType x = fNoneNF; 720 if(newValue == "Exponential") { x = fExponentialNF; } 721 else if(newValue == "Gaussian") { x = fGaussianNF; } 722 else if(newValue == "Flat") { x = fFlatNF; } 723 else if(newValue != "None") { 724 G4ExceptionDescription ed; 725 ed << " NuclearFormFactor type <" << newValue << "> unknown!"; 726 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed); 727 return; 728 } 729 theParameters->SetNuclearFormfactorType(x); 730 } else if (command == ssCmd) { 731 G4eSingleScatteringType x = fWVI; 732 if(newValue == "DPWA") { x = fDPWA; } 733 else if(newValue == "Mott") { x = fMott; } 734 else if(newValue != "WVI") { 735 G4ExceptionDescription ed; 736 ed << " G4eSingleScatteringType type <" << newValue << "> unknown!"; 737 G4Exception("G4EmParametersMessenger", "em0044", JustWarning, ed); 738 return; 739 } 740 theParameters->SetSingleScatteringType(x); 741 } else if (command == fluc1Cmd) { 742 G4EmFluctuationType x = fUniversalFluctuation; 743 if(newValue == "Dummy") { x = fDummyFluctuation; } 744 else if(newValue == "Urban") { x = fUrbanFluctuation; } 745 theParameters->SetFluctuationType(x); 746 } else if (command == posiCmd) { 747 G4PositronAtRestModelType x = fSimplePositronium; 748 if (newValue == "Allison") { x = fAllisonPositronium; } 749 else if (newValue == "OrePowell") { x = fOrePowell; } 750 else if (newValue == "OrePowellPolar") { x = fOrePowellPolar; } 751 theParameters->SetPositronAtRestModelType(x); 752 } else if ( command==tripletCmd ) { 753 theParameters->SetConversionType(tripletCmd->GetNewIntValue(newValue)); 754 } else if ( command==onIsolatedCmd ) { 755 theParameters->SetOnIsolated(onIsolatedCmd->GetNewBoolValue(newValue)); 756 physicsModified = true; 757 } 758 759 if(physicsModified) { 760 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified"); 761 } 762 } 763 764 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 765