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 /// \file DetectorConstructionMessenger.cc 27 /// \brief Implementation of the DetectorConstruction messenger class 28 // 29 // 30 // 31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 32 33 #include "DetectorConstructionMessenger.hh" 34 #include "DetectorConstruction.hh" 35 36 #include "G4UIdirectory.hh" 37 #include "G4UIcmdWithADoubleAndUnit.hh" 38 #include "G4UIcmdWithADouble.hh" 39 #include "G4UIcmdWithAnInteger.hh" 40 #include "G4UIcmdWith3VectorAndUnit.hh" 41 #include "G4UIcmdWithABool.hh" 42 #include "G4UIcmdWithAString.hh" 43 44 #include "G4RunManager.hh" 45 #include "G4ios.hh" 46 47 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 48 49 DetectorConstructionMessenger::DetectorConstructionMessenger(DetectorConstruction* det): 50 fDetector(det) 51 { 52 fCmdDir = new G4UIdirectory("/crystal/"); 53 fCmdDir->SetGuidance("crystal Control"); 54 55 fCrystalMaterialCmd = new G4UIcmdWithAString("/crystal/setCrystalMaterial",this); 56 fCrystalMaterialCmd->SetGuidance("Set Crystal Material"); 57 fCrystalMaterialCmd->SetParameterName("matname",false); 58 fCrystalMaterialCmd->SetDefaultValue("G4_Si"); 59 60 fCrystalSizeCmd = new G4UIcmdWith3VectorAndUnit("/crystal/setCrystalSize",this); 61 fCrystalSizeCmd->SetGuidance("Set Crystal size"); 62 fCrystalSizeCmd->SetParameterName("dimCrX","dimCrY","dimCrZ",false); 63 fCrystalSizeCmd->SetUnitCategory("Length"); 64 fCrystalSizeCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 65 66 fCrystalLatticeCmd = new G4UIcmdWithAString("/crystal/setCrystalLattice",this); 67 fCrystalLatticeCmd-> 68 SetGuidance("Set Crystal Lattice, use brackets (...) for planes and <...> for axes"); 69 fCrystalLatticeCmd->SetParameterName("lattice",false); 70 fCrystalLatticeCmd->SetDefaultValue("(111)"); 71 72 fCrystalAngleXCmd = new G4UIcmdWithADoubleAndUnit("/crystal/setCrystalAngleX",this); 73 fCrystalAngleXCmd->SetGuidance("Set crystal orientation with respect to the beam"); 74 fCrystalAngleXCmd->SetUnitCategory("Angle"); 75 fCrystalAngleXCmd->SetParameterName("angX",false); 76 fCrystalAngleXCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 77 78 fCrystalAngleYCmd = new G4UIcmdWithADoubleAndUnit("/crystal/setCrystalAngleY",this); 79 fCrystalAngleYCmd->SetGuidance("Set crystal orientation with respect to the beam"); 80 fCrystalAngleYCmd->SetUnitCategory("Angle"); 81 fCrystalAngleYCmd->SetParameterName("angY",false); 82 fCrystalAngleYCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 83 84 fCrystalBendingAngleCmd = 85 new G4UIcmdWithADoubleAndUnit("/crystal/setCrystalBendingAngle",this); 86 fCrystalBendingAngleCmd->SetGuidance("Set crystal bending angle"); 87 fCrystalBendingAngleCmd->SetParameterName("bendingAngle",false); 88 fCrystalBendingAngleCmd->SetUnitCategory("Angle"); 89 fCrystalBendingAngleCmd->SetRange("bendingAngle >= 0"); 90 fCrystalBendingAngleCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 91 92 fCrystallineUndulatorAmplitudeCmd = 93 new G4UIcmdWithADoubleAndUnit 94 ("/crystal/setCrystallineUndulatorAmplitude",this); 95 fCrystallineUndulatorAmplitudeCmd-> 96 SetGuidance("Set crystalline undulator amplitude"); 97 fCrystallineUndulatorAmplitudeCmd->SetUnitCategory("Length"); 98 fCrystallineUndulatorAmplitudeCmd-> 99 SetParameterName("CrystallineUndulatorAmplitude",false); 100 fCrystallineUndulatorAmplitudeCmd-> 101 SetRange("CrystallineUndulatorAmplitude > 0"); 102 fCrystallineUndulatorAmplitudeCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 103 104 fCrystallineUndulatorPeriodCmd = 105 new G4UIcmdWithADoubleAndUnit("/crystal/setCrystallineUndulatorPeriod",this); 106 fCrystallineUndulatorPeriodCmd-> 107 SetGuidance("Set crystalline undulator Period"); 108 fCrystallineUndulatorPeriodCmd->SetUnitCategory("Length"); 109 fCrystallineUndulatorPeriodCmd-> 110 SetParameterName("CrystallineUndulatorPeriod",false); 111 fCrystallineUndulatorPeriodCmd-> 112 SetRange("CrystallineUndulatorPeriod > 0"); 113 fCrystallineUndulatorPeriodCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 114 115 fCrystallineUndulatorPhaseCmd = 116 new G4UIcmdWithADouble("/crystal/setCrystallineUndulatorPhase",this); 117 fCrystallineUndulatorPhaseCmd-> 118 SetGuidance("Set crystalline undulator phase"); 119 fCrystallineUndulatorPhaseCmd-> 120 SetParameterName("CrystallineUndulatorPhase",false); 121 fCrystallineUndulatorPhaseCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 122 123 fDetectorSizeCmd = new G4UIcmdWith3VectorAndUnit("/crystal/setDetectorSize",this); 124 fDetectorSizeCmd->SetGuidance("Set detector size"); 125 fDetectorSizeCmd->SetParameterName("dimDetX","dimDetY","dimDetZ",false); 126 fDetectorSizeCmd->SetUnitCategory("Length"); 127 fDetectorSizeCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 128 129 fDetectorFrontPosZCmd = 130 new G4UIcmdWithADoubleAndUnit("/crystal/setFrontPositionZ",this); 131 fDetectorFrontPosZCmd->SetGuidance("Set detector front position Z"); 132 fDetectorFrontPosZCmd->SetParameterName("frontPosDetZ",false); 133 fDetectorFrontPosZCmd->SetUnitCategory("Length"); 134 fDetectorFrontPosZCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 135 136 fPotentialPathCmd = new G4UIcmdWithAString("/crystal/setChannelingDataPath",this); 137 fPotentialPathCmd-> 138 SetGuidance("Set the path where to find the available data " 139 "for the G4ChannelingFastSimModel " 140 "if different from G4CHANNELINGDATA"); 141 fPotentialPathCmd->SetParameterName("channelingDataPath",false); 142 fPotentialPathCmd->SetDefaultValue(""); 143 144 fChannelingModelCmd = new G4UIcmdWithABool("/crystal/setChannelingModel", this); 145 fChannelingModelCmd->SetGuidance("Activate/deactivate G4ChannelingFastSimModel"); 146 fChannelingModelCmd->SetParameterName("ChannelingModel",true); 147 fChannelingModelCmd->SetDefaultValue(false); 148 149 fRadModelCmd = new G4UIcmdWithABool("/crystal/setRadiationModel", this); 150 fRadModelCmd->SetGuidance("Activate/deactivate G4BaierKatkov"); 151 fRadModelCmd->SetParameterName("ActivateRadiationModel",true); 152 fRadModelCmd->SetDefaultValue(false); 153 154 fMinPhotonEnergyCmd = 155 new G4UIcmdWithADoubleAndUnit("/crystal/setMinPhotonEnergy",this); 156 fMinPhotonEnergyCmd-> 157 SetGuidance("Set the low energy threshold for " 158 "the photons emitted in G4BaierKatkov"); 159 fMinPhotonEnergyCmd->SetParameterName("MinPhotonEnergy",false); 160 fMinPhotonEnergyCmd->SetUnitCategory("Energy"); 161 fMinPhotonEnergyCmd->SetRange("MinPhotonEnergy > 0"); 162 fMinPhotonEnergyCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 163 164 fSamplingPhotonsNumberCmd = 165 new G4UIcmdWithAnInteger("/crystal/setSamplingPhotonsNumber",this); 166 fSamplingPhotonsNumberCmd-> 167 SetGuidance("Set SamplingPhotonsNumber in G4BaierKatkov"); 168 fSamplingPhotonsNumberCmd->SetParameterName("SamplingPhotonsNumber",false); 169 fSamplingPhotonsNumberCmd->SetRange("SamplingPhotonsNumber>1"); 170 fSamplingPhotonsNumberCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 171 172 fNSmallTrajectoryStepsCmd = 173 new G4UIcmdWithAnInteger("/crystal/setNSmallTrajectorySteps",this); 174 fNSmallTrajectoryStepsCmd-> 175 SetGuidance("Set NSmallTrajectorySteps in G4BaierKatkov"); 176 fNSmallTrajectoryStepsCmd->SetParameterName("NSmallTrajectorySteps",false); 177 fNSmallTrajectoryStepsCmd->SetRange("NSmallTrajectorySteps>1"); 178 fNSmallTrajectoryStepsCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 179 180 fRadiationAngleFactorCmd = 181 new G4UIcmdWithADouble("/crystal/setRadiationAngleFactor",this); 182 fRadiationAngleFactorCmd->SetGuidance("Set Radiation Angle Factor"); 183 fRadiationAngleFactorCmd->SetParameterName("RadiationAngleFactor",false); 184 fRadiationAngleFactorCmd->SetRange("RadiationAngleFactor > 0"); 185 fRadiationAngleFactorCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 186 187 fMinPhotonEnergyAddStatCmd = 188 new G4UIcmdWithADoubleAndUnit("/crystal/AddPhotonStatistics/setMinPhotonEnergy",this); 189 fMinPhotonEnergyAddStatCmd-> 190 SetGuidance("Set the min energy in the range to increase " 191 "the sampling photon statistics in G4BaierKatkov"); 192 fMinPhotonEnergyAddStatCmd->SetParameterName("addStatMinPhotonEnergy",false); 193 fMinPhotonEnergyAddStatCmd->SetUnitCategory("Energy"); 194 fMinPhotonEnergyAddStatCmd->SetRange("addStatMinPhotonEnergy > 0"); 195 fMinPhotonEnergyAddStatCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 196 197 fMaxPhotonEnergyAddStatCmd = 198 new G4UIcmdWithADoubleAndUnit("/crystal/AddPhotonStatistics/setMaxPhotonEnergy",this); 199 fMaxPhotonEnergyAddStatCmd-> 200 SetGuidance("Set the max energy in the range to increase " 201 "the sampling photon statistics in G4BaierKatkov"); 202 fMaxPhotonEnergyAddStatCmd->SetParameterName("addStatMaxPhotonEnergy",false); 203 fMaxPhotonEnergyAddStatCmd->SetUnitCategory("Energy"); 204 fMaxPhotonEnergyAddStatCmd->SetRange("addStatMaxPhotonEnergy > 0"); 205 fMaxPhotonEnergyAddStatCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 206 207 fTimesPhotonStatisticsCmd = 208 new G4UIcmdWithAnInteger("/crystal/AddPhotonStatistics/setMultiplePhotonStatistics", 209 this); 210 fTimesPhotonStatisticsCmd-> 211 SetGuidance("Set multiple of the sampling photon statistics in G4BaierKatkov"); 212 fTimesPhotonStatisticsCmd->SetParameterName("timesPhotonStatistics",false); 213 fTimesPhotonStatisticsCmd->SetRange("timesPhotonStatistics > 1"); 214 fTimesPhotonStatisticsCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 215 216 fParticleMinKinEnergyCmd = 217 new G4UIcmdWithADoubleAndUnit("/crystal/setParticleMinKinEnergy",this); 218 fParticleMinKinEnergyCmd-> 219 SetGuidance("Set the low energy threshold for particle " 220 "to enter the G4ChannelingFastSimModel"); 221 fParticleMinKinEnergyCmd->SetParameterName("partLEth",false); 222 fParticleMinKinEnergyCmd->SetUnitCategory("Energy"); 223 fParticleMinKinEnergyCmd->SetRange("partLEth > 0"); 224 fParticleMinKinEnergyCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 225 226 fProtonMinKinEnergyCmd = 227 new G4UIcmdWithADoubleAndUnit("/crystal/setParticleMinKinEnergy/proton",this); 228 fProtonMinKinEnergyCmd-> 229 SetGuidance("Set the low energy threshold for proton " 230 "to enter the G4ChannelingFastSimModel"); 231 fProtonMinKinEnergyCmd->SetParameterName("protonLEth",false); 232 fProtonMinKinEnergyCmd->SetUnitCategory("Energy"); 233 fProtonMinKinEnergyCmd->SetRange("protonLEth > 0"); 234 fProtonMinKinEnergyCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 235 236 fAntiprotonMinKinEnergyCmd = 237 new G4UIcmdWithADoubleAndUnit("/crystal/setParticleMinKinEnergy/anti_proton",this); 238 fAntiprotonMinKinEnergyCmd->SetGuidance("Set the low energy threshold for anti_proton" 239 " to enter the G4ChannelingFastSimModel"); 240 fAntiprotonMinKinEnergyCmd->SetParameterName("anti_protonLEth",false); 241 fAntiprotonMinKinEnergyCmd->SetUnitCategory("Energy"); 242 fAntiprotonMinKinEnergyCmd->SetRange("anti_protonLEth > 0"); 243 fAntiprotonMinKinEnergyCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 244 245 fPiPlusMinKinEnergyCmd = 246 new G4UIcmdWithADoubleAndUnit("/crystal/setParticleMinKinEnergy/pi+",this); 247 fPiPlusMinKinEnergyCmd-> 248 SetGuidance("Set the low energy threshold for pi+ " 249 "to enter the G4ChannelingFastSimModel"); 250 fPiPlusMinKinEnergyCmd->SetParameterName("piPlusLEth",false); 251 fPiPlusMinKinEnergyCmd->SetUnitCategory("Energy"); 252 fPiPlusMinKinEnergyCmd->SetRange("piPlusLEth > 0"); 253 fPiPlusMinKinEnergyCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 254 255 fPiMinusMinKinEnergyCmd = 256 new G4UIcmdWithADoubleAndUnit("/crystal/setParticleMinKinEnergy/pi-",this); 257 fPiMinusMinKinEnergyCmd->SetGuidance("Set the low energy threshold for pi- " 258 "to enter the G4ChannelingFastSimModel"); 259 fPiMinusMinKinEnergyCmd->SetParameterName("piMinusLEth",false); 260 fPiMinusMinKinEnergyCmd->SetUnitCategory("Energy"); 261 fPiMinusMinKinEnergyCmd->SetRange("piMinusLEth > 0"); 262 fPiMinusMinKinEnergyCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 263 264 fPositronMinKinEnergyCmd = 265 new G4UIcmdWithADoubleAndUnit("/crystal/setParticleMinKinEnergy/e+",this); 266 fPositronMinKinEnergyCmd->SetGuidance("Set the low energy threshold for e+ " 267 "to enter the G4ChannelingFastSimModel"); 268 fPositronMinKinEnergyCmd->SetParameterName("ePlusLEth",false); 269 fPositronMinKinEnergyCmd->SetUnitCategory("Energy"); 270 fPositronMinKinEnergyCmd->SetRange("ePlusLEth > 0"); 271 fPositronMinKinEnergyCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 272 273 fElectronMinKinEnergyCmd = 274 new G4UIcmdWithADoubleAndUnit("/crystal/setParticleMinKinEnergy/e-",this); 275 fElectronMinKinEnergyCmd->SetGuidance("Set the low energy threshold for e- " 276 "to enter the G4ChannelingFastSimModel"); 277 fElectronMinKinEnergyCmd->SetParameterName("eMinusLEth",false); 278 fElectronMinKinEnergyCmd->SetUnitCategory("Energy"); 279 fElectronMinKinEnergyCmd->SetRange("eMinusLEth > 0"); 280 fElectronMinKinEnergyCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 281 282 fMuPlusMinKinEnergyCmd = 283 new G4UIcmdWithADoubleAndUnit("/crystal/setParticleMinKinEnergy/mu+",this); 284 fMuPlusMinKinEnergyCmd->SetGuidance("Set the low energy threshold for mu+ " 285 "to enter the G4ChannelingFastSimModel"); 286 fMuPlusMinKinEnergyCmd->SetParameterName("muPlusLEth",false); 287 fMuPlusMinKinEnergyCmd->SetUnitCategory("Energy"); 288 fMuPlusMinKinEnergyCmd->SetRange("muPlusLEth > 0"); 289 fMuPlusMinKinEnergyCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 290 291 fMuMinusMinKinEnergyCmd = 292 new G4UIcmdWithADoubleAndUnit("/crystal/setParticleMinKinEnergy/mu-",this); 293 fMuMinusMinKinEnergyCmd->SetGuidance("Set the low energy threshold for mu- " 294 "to enter the G4ChannelingFastSimModel"); 295 fMuMinusMinKinEnergyCmd->SetParameterName("muMinusLEth",false); 296 fMuMinusMinKinEnergyCmd->SetUnitCategory("Energy"); 297 fMuMinusMinKinEnergyCmd->SetRange("muMinusLEth > 0"); 298 fMuMinusMinKinEnergyCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 299 300 fLindhardAnglesCmd = new G4UIcmdWithADouble("/crystal/setLindhardAngles",this); 301 fLindhardAnglesCmd-> 302 SetGuidance("Set high angular threshold for particle to enter " 303 "the G4ChannelingFastSimModel expressed in Lindhard angles"); 304 fLindhardAnglesCmd->SetParameterName("LindhardAngles",false); 305 fLindhardAnglesCmd->SetRange("LindhardAngles >= 0"); 306 fLindhardAnglesCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 307 308 fLindhardAnglesProtonCmd = 309 new G4UIcmdWithADouble("/crystal/setLindhardAngles/proton",this); 310 fLindhardAnglesProtonCmd-> 311 SetGuidance("Set high angular threshold for proton to enter " 312 "the G4ChannelingFastSimModel expressed in Lindhard angles"); 313 fLindhardAnglesProtonCmd->SetParameterName("LindhardAnglesProton",false); 314 fLindhardAnglesProtonCmd->SetRange("LindhardAnglesProton >= 0"); 315 fLindhardAnglesProtonCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 316 317 fLindhardAnglesAntiprotonCmd = 318 new G4UIcmdWithADouble("/crystal/setLindhardAngles/anti_proton",this); 319 fLindhardAnglesAntiprotonCmd-> 320 SetGuidance("Set high angular threshold for anti_proton to enter " 321 "the G4ChannelingFastSimModel expressed in Lindhard angles"); 322 fLindhardAnglesAntiprotonCmd->SetParameterName("LindhardAnglesAnti_proton",false); 323 fLindhardAnglesAntiprotonCmd->SetRange("LindhardAnglesAnti_proton >= 0"); 324 fLindhardAnglesAntiprotonCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 325 326 fLindhardAnglesPiPlusCmd = 327 new G4UIcmdWithADouble("/crystal/setLindhardAngles/pi+",this); 328 fLindhardAnglesPiPlusCmd-> 329 SetGuidance("Set high angular threshold for pi+ to enter " 330 "the G4ChannelingFastSimModel expressed in Lindhard angles"); 331 fLindhardAnglesPiPlusCmd->SetParameterName("LindhardAnglesPiPlus",false); 332 fLindhardAnglesPiPlusCmd->SetRange("LindhardAnglesPiPlus >= 0"); 333 fLindhardAnglesPiPlusCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 334 335 fLindhardAnglesPiMinusCmd = 336 new G4UIcmdWithADouble("/crystal/setLindhardAngles/pi-",this); 337 fLindhardAnglesPiMinusCmd-> 338 SetGuidance("Set high angular threshold for pi- to enter " 339 "the G4ChannelingFastSimModel expressed in Lindhard angles"); 340 fLindhardAnglesPiMinusCmd->SetParameterName("LindhardAnglesPiMinus",false); 341 fLindhardAnglesPiMinusCmd->SetRange("LindhardAnglesPiMinus >= 0"); 342 fLindhardAnglesPiMinusCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 343 344 fLindhardAnglesPositronCmd = 345 new G4UIcmdWithADouble("/crystal/setLindhardAngles/e+",this); 346 fLindhardAnglesPositronCmd-> 347 SetGuidance("Set high angular threshold for e+ to enter " 348 "the G4ChannelingFastSimModel expressed in Lindhard angles"); 349 fLindhardAnglesPositronCmd->SetParameterName("LindhardAnglesPositron",false); 350 fLindhardAnglesPositronCmd->SetRange("LindhardAnglesPositron >= 0"); 351 fLindhardAnglesPositronCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 352 353 fLindhardAnglesElectronCmd = 354 new G4UIcmdWithADouble("/crystal/setLindhardAngles/e-",this); 355 fLindhardAnglesElectronCmd-> 356 SetGuidance("Set high angular threshold for e- to enter " 357 "the G4ChannelingFastSimModel expressed in Lindhard angles"); 358 fLindhardAnglesElectronCmd->SetParameterName("LindhardAnglesElectron",false); 359 fLindhardAnglesElectronCmd->SetRange("LindhardAnglesElectron >= 0"); 360 fLindhardAnglesElectronCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 361 362 fLindhardAnglesMuPlusCmd = 363 new G4UIcmdWithADouble("/crystal/setLindhardAngles/mu+",this); 364 fLindhardAnglesMuPlusCmd-> 365 SetGuidance("Set high angular threshold for mu+ to enter " 366 "the G4ChannelingFastSimModel expressed in Lindhard angles"); 367 fLindhardAnglesMuPlusCmd->SetParameterName("LindhardAnglesMuPlus",false); 368 fLindhardAnglesMuPlusCmd->SetRange("LindhardAnglesMuPlus >= 0"); 369 fLindhardAnglesMuPlusCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 370 371 fLindhardAnglesMuMinusCmd = 372 new G4UIcmdWithADouble("/crystal/setLindhardAngles/mu-",this); 373 fLindhardAnglesMuMinusCmd-> 374 SetGuidance("Set high angular threshold for mu- to enter " 375 "the G4ChannelingFastSimModel expressed in Lindhard angles"); 376 fLindhardAnglesMuMinusCmd->SetParameterName("LindhardAnglesMuMinus",false); 377 fLindhardAnglesMuMinusCmd->SetRange("LindhardAnglesMuMinus >= 0"); 378 fLindhardAnglesMuMinusCmd->AvailableForStates(G4State_PreInit,G4State_Idle); 379 380 } 381 382 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 383 384 DetectorConstructionMessenger::~DetectorConstructionMessenger() 385 { 386 delete fCmdDir; 387 delete fCrystalMaterialCmd; 388 delete fCrystalSizeCmd; 389 delete fCrystalLatticeCmd; 390 delete fCrystalAngleXCmd; 391 delete fCrystalAngleYCmd; 392 delete fCrystalBendingAngleCmd; 393 delete fRadModelCmd; 394 delete fChannelingModelCmd; 395 396 delete fCrystallineUndulatorAmplitudeCmd; 397 delete fCrystallineUndulatorPeriodCmd; 398 delete fCrystallineUndulatorPhaseCmd; 399 400 delete fDetectorSizeCmd; 401 delete fDetectorFrontPosZCmd; 402 403 delete fPotentialPathCmd; 404 405 delete fMinPhotonEnergyCmd; 406 delete fSamplingPhotonsNumberCmd; 407 delete fNSmallTrajectoryStepsCmd; 408 delete fRadiationAngleFactorCmd; 409 delete fMinPhotonEnergyAddStatCmd; 410 delete fMaxPhotonEnergyAddStatCmd; 411 delete fTimesPhotonStatisticsCmd; 412 413 delete fParticleMinKinEnergyCmd; 414 delete fProtonMinKinEnergyCmd; 415 delete fAntiprotonMinKinEnergyCmd; 416 delete fPiPlusMinKinEnergyCmd; 417 delete fPiMinusMinKinEnergyCmd; 418 delete fElectronMinKinEnergyCmd; 419 delete fPositronMinKinEnergyCmd; 420 delete fMuPlusMinKinEnergyCmd; 421 delete fMuMinusMinKinEnergyCmd; 422 423 delete fLindhardAnglesCmd; 424 delete fLindhardAnglesProtonCmd; 425 delete fLindhardAnglesAntiprotonCmd; 426 delete fLindhardAnglesPiPlusCmd; 427 delete fLindhardAnglesPiMinusCmd; 428 delete fLindhardAnglesElectronCmd; 429 delete fLindhardAnglesPositronCmd; 430 delete fLindhardAnglesMuPlusCmd; 431 delete fLindhardAnglesMuMinusCmd; 432 } 433 434 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 435 436 void DetectorConstructionMessenger::SetNewValue(G4UIcommand* command, G4String newValue) 437 { 438 if (command == fCrystalMaterialCmd) 439 {fDetector->SetCrystalMaterial(newValue);} 440 if (command == fCrystalSizeCmd) 441 {fDetector->SetCrystalSize(fCrystalSizeCmd->GetNew3VectorValue(newValue));} 442 if (command == fCrystalLatticeCmd) 443 {fDetector->SetCrystalLattice(newValue);} 444 if (command == fCrystalAngleXCmd) 445 {fDetector->SetCrystalAngleX(fCrystalAngleXCmd->GetNewDoubleValue(newValue));} 446 if (command == fCrystalAngleYCmd) 447 {fDetector->SetCrystalAngleY(fCrystalAngleYCmd->GetNewDoubleValue(newValue));} 448 if (command == fCrystalBendingAngleCmd) 449 {fDetector->SetCrystalBendingAngle( 450 fCrystalBendingAngleCmd->GetNewDoubleValue(newValue));} 451 if (command == fRadModelCmd) 452 {fDetector->SetRadiationModel(fRadModelCmd->GetNewBoolValue(newValue));} 453 if (command == fChannelingModelCmd) 454 {fDetector->SetChannelingModel(fChannelingModelCmd->GetNewBoolValue(newValue));} 455 456 if (command == fCrystallineUndulatorAmplitudeCmd) 457 {fDetector->SetCrystallineUndulatorAmplitude 458 (fCrystallineUndulatorAmplitudeCmd->GetNewDoubleValue(newValue));} 459 if (command == fCrystallineUndulatorPeriodCmd) 460 {fDetector->SetCrystallineUndulatorPeriod 461 (fCrystallineUndulatorPeriodCmd->GetNewDoubleValue(newValue));} 462 if (command == fCrystallineUndulatorPhaseCmd) 463 {fDetector->SetCrystallineUndulatorPhase 464 (fCrystallineUndulatorPhaseCmd->GetNewDoubleValue(newValue));} 465 466 if (command == fDetectorSizeCmd) 467 {fDetector->SetDetectorSize(fDetectorSizeCmd->GetNew3VectorValue(newValue));} 468 if (command == fDetectorFrontPosZCmd) 469 {fDetector->SetDetectorFrontPositionZ( 470 fDetectorFrontPosZCmd->GetNewDoubleValue(newValue));} 471 472 if (command == fPotentialPathCmd) 473 {fDetector->SetPotentialPath(newValue);} 474 475 if (command == fMinPhotonEnergyCmd) 476 {fDetector->SetMinPhotonEnergy( 477 fMinPhotonEnergyCmd->GetNewDoubleValue(newValue));} 478 if (command == fSamplingPhotonsNumberCmd) 479 {fDetector->SetSamplingPhotonsNumber( 480 fSamplingPhotonsNumberCmd->GetNewIntValue(newValue));} 481 if (command == fNSmallTrajectoryStepsCmd) 482 {fDetector->SetNSmallTrajectorySteps( 483 fNSmallTrajectoryStepsCmd->GetNewIntValue(newValue));} 484 if (command == fRadiationAngleFactorCmd) 485 {fDetector->SetRadiationAngleFactor( 486 fRadiationAngleFactorCmd->GetNewDoubleValue(newValue));} 487 488 if (command == fMinPhotonEnergyAddStatCmd) 489 {fDetector->SetMinPhotonEnergyAddStat( 490 fMinPhotonEnergyAddStatCmd->GetNewDoubleValue(newValue));} 491 if (command == fMaxPhotonEnergyAddStatCmd) 492 {fDetector->SetMaxPhotonEnergyAddStat( 493 fMaxPhotonEnergyAddStatCmd->GetNewDoubleValue(newValue));} 494 if (command == fTimesPhotonStatisticsCmd) 495 {fDetector->SetMultiplePhotonStatistics( 496 fTimesPhotonStatisticsCmd->GetNewIntValue(newValue));} 497 498 if (command == fParticleMinKinEnergyCmd) 499 {fDetector->SetParticleMinKinEnergy( 500 fParticleMinKinEnergyCmd->GetNewDoubleValue(newValue));} 501 if (command == fProtonMinKinEnergyCmd) 502 {fDetector->SetProtonMinKinEnergy( 503 fProtonMinKinEnergyCmd->GetNewDoubleValue(newValue));} 504 if (command == fAntiprotonMinKinEnergyCmd) 505 {fDetector->SetAntiprotonMinKinEnergy( 506 fAntiprotonMinKinEnergyCmd->GetNewDoubleValue(newValue));} 507 if (command == fPiPlusMinKinEnergyCmd) 508 {fDetector->SetPiPlusMinKinEnergy( 509 fPiPlusMinKinEnergyCmd->GetNewDoubleValue(newValue));} 510 if (command == fPiMinusMinKinEnergyCmd) 511 {fDetector->SetPiMinusMinKinEnergy( 512 fPiMinusMinKinEnergyCmd->GetNewDoubleValue(newValue));} 513 if (command == fElectronMinKinEnergyCmd) 514 {fDetector->SetElectronMinKinEnergy( 515 fElectronMinKinEnergyCmd->GetNewDoubleValue(newValue));} 516 if (command == fPositronMinKinEnergyCmd) 517 {fDetector->SetPositronMinKinEnergy( 518 fPositronMinKinEnergyCmd->GetNewDoubleValue(newValue));} 519 if (command == fMuPlusMinKinEnergyCmd) 520 {fDetector->SetMuPlusMinKinEnergy( 521 fMuPlusMinKinEnergyCmd->GetNewDoubleValue(newValue));} 522 if (command == fMuMinusMinKinEnergyCmd) 523 {fDetector->SetMuMinusMinKinEnergy( 524 fMuMinusMinKinEnergyCmd->GetNewDoubleValue(newValue));} 525 526 if (command == fLindhardAnglesCmd) 527 {fDetector->SetLindhardAngles( 528 fLindhardAnglesCmd->GetNewDoubleValue(newValue));} 529 if (command == fLindhardAnglesProtonCmd) 530 {fDetector->SetLindhardAnglesProton( 531 fLindhardAnglesProtonCmd->GetNewDoubleValue(newValue));} 532 if (command == fLindhardAnglesAntiprotonCmd) 533 {fDetector->SetLindhardAnglesAntiproton( 534 fLindhardAnglesAntiprotonCmd->GetNewDoubleValue(newValue));} 535 if (command == fLindhardAnglesPiPlusCmd) 536 {fDetector->SetLindhardAnglesPiPlus( 537 fLindhardAnglesPiPlusCmd->GetNewDoubleValue(newValue));} 538 if (command == fLindhardAnglesPiMinusCmd) 539 {fDetector->SetLindhardAnglesPiMinus( 540 fLindhardAnglesPiMinusCmd->GetNewDoubleValue(newValue));} 541 if (command == fLindhardAnglesElectronCmd) 542 {fDetector->SetLindhardAnglesElectron( 543 fLindhardAnglesElectronCmd->GetNewDoubleValue(newValue));} 544 if (command == fLindhardAnglesPositronCmd) 545 {fDetector->SetLindhardAnglesPositron( 546 fLindhardAnglesPositronCmd->GetNewDoubleValue(newValue));} 547 if (command == fLindhardAnglesMuPlusCmd) 548 {fDetector->SetLindhardAnglesMuPlus( 549 fLindhardAnglesMuPlusCmd->GetNewDoubleValue(newValue));} 550 if (command == fLindhardAnglesMuMinusCmd) 551 {fDetector->SetLindhardAnglesMuMinus( 552 fLindhardAnglesMuMinusCmd->GetNewDoubleValue(newValue));} 553 } 554 555 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... 556