Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4OpticalParametersMessenger Class Reference

#include <G4OpticalParametersMessenger.hh>

+ Inheritance diagram for G4OpticalParametersMessenger:

Public Member Functions

 G4OpticalParametersMessenger (G4OpticalParameters *)
 
virtual ~G4OpticalParametersMessenger ()
 
virtual void SetNewValue (G4UIcommand *, G4String)
 
- Public Member Functions inherited from G4UImessenger
 G4UImessenger ()
 
 G4UImessenger (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
virtual ~G4UImessenger ()
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
virtual void SetNewValue (G4UIcommand *command, G4String newValue)
 
G4bool operator== (const G4UImessenger &messenger) const
 
G4bool operator!= (const G4UImessenger &messenger) const
 
G4bool CommandsShouldBeInMaster () const
 

Additional Inherited Members

- Protected Member Functions inherited from G4UImessenger
G4String ItoS (G4int i)
 
G4String DtoS (G4double a)
 
G4String BtoS (G4bool b)
 
G4int StoI (G4String s)
 
G4long StoL (G4String s)
 
G4double StoD (G4String s)
 
G4bool StoB (G4String s)
 
void AddUIcommand (G4UIcommand *newCommand)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 
- Protected Attributes inherited from G4UImessenger
G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Detailed Description

Definition at line 62 of file G4OpticalParametersMessenger.hh.

Constructor & Destructor Documentation

◆ G4OpticalParametersMessenger()

G4OpticalParametersMessenger::G4OpticalParametersMessenger ( G4OpticalParameters opticalParameters)

Definition at line 59 of file G4OpticalParametersMessenger.cc.

61 : params(opticalParameters)
62
63{
64 G4bool toBeBroadcasted = false;
65 fDir = new G4UIdirectory("/process/optical/defaults/",toBeBroadcasted);
66 fDir->SetGuidance("DEPRECATED Commands related to the optical physics simulation engine.");
67 fDir2 = new G4UIdirectory("/process/optical/",toBeBroadcasted);
68 fDir2->SetGuidance("Commands related to the optical physics simulation engine.");
69
70 CreateDirectory("/process/optical/defaults/cerenkov/", "DEPRECATED Cerenkov process commands");
71 CreateDirectory("/process/optical/defaults/scintillation/", "DEPRECATED Scintillation process commands");
72 CreateDirectory("/process/optical/defaults/wls/", "DEPRECATED Wave length shifting process commands");
73 CreateDirectory("/process/optical/defaults/boundary/", "DEPRECATED Boundary scattering commands");
74
75 CreateDirectory("/process/optical/cerenkov/", "Cerenkov process commands");
76 CreateDirectory("/process/optical/scintillation/", "Scintillation process commands");
77 CreateDirectory("/process/optical/wls/", "Wave length shifting process commands");
78 CreateDirectory("/process/optical/wls2/", "Second Wave length shifting process commands");
79 CreateDirectory("/process/optical/boundary/", "Boundary scattering commands");
80 CreateDirectory("/process/optical/mie/", "Mie scattering process commands");
81 CreateDirectory("/process/optical/absorption/", "absorption process commands");
82 CreateDirectory("/process/optical/rayleigh/", "Rayleigh scattering commands");
83
84 // general commands
85 fActivateProcessCmd= new G4UIcommand("/process/optical/processActivation", this);
86 fActivateProcessCmd->SetGuidance("Activate/deactivate the specified optical process");
87 G4UIparameter* par = new G4UIparameter("proc_name",'s',false);
88 G4String candidates;
89 for ( G4int i=0; i<kNoProcess; i++ ) {
90 candidates += G4OpticalProcessName(i);
91 candidates += G4String(" ");
92 }
93 par->SetParameterCandidates(candidates);
94 par->SetGuidance("the process name");
95 fActivateProcessCmd->SetParameter(par);
96 par = new G4UIparameter("flag",'b',true);
97 par->SetDefaultValue(true);
98 par->SetGuidance("activation flag");
99 fActivateProcessCmd->SetParameter(par);
100 fActivateProcessCmd->AvailableForStates(G4State_PreInit);
101
102 // DEPRECATED
103 fTrackSecondariesFirstCmd = new G4UIcommand("/process/optical/setTrackSecondariesFirst", this);
104 fTrackSecondariesFirstCmd->SetGuidance("Activate/deactivate tracking of secondaries before finishing their parent track");
105 fTrackSecondariesFirstCmd->SetGuidance("DEPRECATED: Use /process/optical/cerenkov/setTrackSecondariesFirst and");
106 fTrackSecondariesFirstCmd->SetGuidance("/process/optical/scintillation/setTrackSecondariesFirst and instead.");
107 par = new G4UIparameter("proc_name",'s',false);
108 par->SetParameterCandidates("Cerenkov Scintillation");
109 fTrackSecondariesFirstCmd->SetParameter(par);
110 par = new G4UIparameter("flag",'b',false);
111 par->SetDefaultValue(true);
112 fTrackSecondariesFirstCmd->SetParameter(par);
113 fTrackSecondariesFirstCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
114
115 fVerboseCmd = new G4UIcmdWithAnInteger("/process/optical/verbose", this);
116 fVerboseCmd->SetGuidance("Set default verbose level for optical processes");
117 fVerboseCmd->SetParameterName("ver", true);
118 fVerboseCmd->SetDefaultValue(1);
119 fVerboseCmd->SetRange("ver>=0");
121
122 fDumpCmd = new G4UIcommand("/process/optical/printParameters", this);
123 fDumpCmd->SetGuidance("Print all optical parameters.");
124
125 // Cerenkov ////////////////////
126 fCerenkovMaxPhotons1Cmd = new G4UIcmdWithAnInteger("/process/optical/defaults/cerenkov/setMaxPhotons", this);
127 fCerenkovMaxPhotons1Cmd->SetGuidance("Set maximum number of photons per step");
128 fCerenkovMaxPhotons1Cmd->SetGuidance("DEPRECATED: use /process/optical/cerenkov/setMaxPhotons instead.");
129 fCerenkovMaxPhotons1Cmd->SetParameterName("CerenkovMaxPhotons", false);
130 fCerenkovMaxPhotons1Cmd->SetRange("CerenkovMaxPhotons>=0");
131 fCerenkovMaxPhotons1Cmd->AvailableForStates(G4State_PreInit, G4State_Idle);
132
133 fCerenkovMaxPhotonsCmd = new G4UIcmdWithAnInteger("/process/optical/cerenkov/setMaxPhotons", this);
134 fCerenkovMaxPhotonsCmd->SetGuidance("Set maximum number of photons per step");
135 fCerenkovMaxPhotonsCmd->SetParameterName("CerenkovMaxPhotons", false);
136 fCerenkovMaxPhotonsCmd->SetRange("CerenkovMaxPhotons>=0");
137 fCerenkovMaxPhotonsCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
138
139 fCerenkovMaxBetaChange1Cmd = new G4UIcmdWithADouble("/process/optical/defaults/cerenkov/setMaxBetaChange", this);
140 fCerenkovMaxBetaChange1Cmd->SetGuidance("Set maximum change of beta of parent particle per step");
141 fCerenkovMaxBetaChange1Cmd->SetGuidance("DEPRECATED: use /process/optical/cerenkov/setMaxBetaChange instead.");
142 fCerenkovMaxBetaChange1Cmd->SetParameterName("CerenkovMaxBetaChange", false);
143 fCerenkovMaxBetaChange1Cmd->SetRange("CerenkovMaxBetaChange>=0");
144 fCerenkovMaxBetaChange1Cmd->AvailableForStates(G4State_PreInit, G4State_Idle);
145
146 fCerenkovMaxBetaChangeCmd = new G4UIcmdWithADouble("/process/optical/cerenkov/setMaxBetaChange", this);
147 fCerenkovMaxBetaChangeCmd->SetGuidance("Set maximum change of beta of parent particle per step");
148 fCerenkovMaxBetaChangeCmd->SetParameterName("CerenkovMaxBetaChange", false);
149 fCerenkovMaxBetaChangeCmd->SetRange("CerenkovMaxBetaChange>=0");
150 fCerenkovMaxBetaChangeCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
151
152 fCerenkovStackPhotons1Cmd = new G4UIcmdWithABool("/process/optical/defaults/cerenkov/setStackPhotons", this);
153 fCerenkovStackPhotons1Cmd->SetGuidance("Set whether or not to stack secondary Cerenkov photons");
154 fCerenkovStackPhotons1Cmd->SetGuidance("DEPRECATED: use /process/optical/cerenkov/setStackPhotons instead.");
155 fCerenkovStackPhotons1Cmd->SetParameterName("CerenkovStackPhotons", true);
156 fCerenkovStackPhotons1Cmd->SetDefaultValue(true);
157 fCerenkovStackPhotons1Cmd->AvailableForStates(G4State_PreInit, G4State_Idle);
158
159 fCerenkovStackPhotonsCmd = new G4UIcmdWithABool("/process/optical/cerenkov/setStackPhotons", this);
160 fCerenkovStackPhotonsCmd->SetGuidance("Set whether or not to stack secondary Cerenkov photons");
161 fCerenkovStackPhotonsCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
162
163 fCerenkovTrackSecondariesFirstCmd = new G4UIcmdWithABool("/process/optical/cerenkov/setTrackSecondariesFirst", this);
164 fCerenkovTrackSecondariesFirstCmd->SetGuidance("Whether to track secondary Cerenkov photons before the primary.");
165 fCerenkovTrackSecondariesFirstCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
166
167 fCerenkovVerboseLevelCmd = new G4UIcmdWithAnInteger("/process/optical/cerenkov/verbose", this);
168 fCerenkovVerboseLevelCmd->SetGuidance("Verbose level for Cerenkov process.");
169 fCerenkovVerboseLevelCmd->SetParameterName("verbose", true);
170 fCerenkovVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
171 fCerenkovVerboseLevelCmd->SetDefaultValue(2);
172 fCerenkovVerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
173
174 // Scintillation //////////////////////////
175 fScintYieldFactor1Cmd = new G4UIcmdWithADouble("/process/optical/defaults/scintillation/setYieldFactor", this);
176 fScintYieldFactor1Cmd->SetGuidance("Set scintillation yield factor");
177 fScintYieldFactor1Cmd->SetGuidance("DEPRECATED: use /process/optical/scintillation/setYieldFactorinstead.");
178 fScintYieldFactor1Cmd->SetParameterName("ScintillationYieldFactor", false);
179 fScintYieldFactor1Cmd->SetRange("ScintillationYieldFactor>=0");
180 fScintYieldFactor1Cmd->AvailableForStates(G4State_PreInit, G4State_Idle);
181
182 fScintYieldFactorCmd = new G4UIcmdWithADouble("/process/optical/scintillation/setYieldFactor", this);
183 fScintYieldFactorCmd->SetGuidance("Set scintillation yield factor");
184 fScintYieldFactorCmd->SetParameterName("ScintillationYieldFactor", false);
185 fScintYieldFactorCmd->SetRange("ScintillationYieldFactor>=0");
186 fScintYieldFactorCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
187
188 fScintExcitationRatioCmd = new G4UIcmdWithADouble("/process/optical/scintillation/setExcitationRatio", this);
189 fScintExcitationRatioCmd->SetGuidance("Set scintillation excitation ratio");
190 fScintExcitationRatioCmd->SetParameterName("ExcitationRatio", false);
191 fScintExcitationRatioCmd->SetRange("ExcitationRatio >= 0 && ExcitationRatio <=1");
192 fScintExcitationRatioCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
193
194 fScintByParticleType1Cmd = new G4UIcmdWithABool("/process/optical/defaults/scintillation/setByParticleType", this);
195 fScintByParticleType1Cmd->SetGuidance("Activate/Inactivate scintillation process by particle type");
196 fScintByParticleType1Cmd->SetGuidance("DEPRECATED: use /process/optical/scintillation/setByParticleType instead.");
197 fScintByParticleType1Cmd->SetParameterName("ScintillationByParticleTypeActivation", false);
198 fScintByParticleType1Cmd->AvailableForStates(G4State_PreInit, G4State_Idle);
199
200 fScintByParticleTypeCmd = new G4UIcmdWithABool("/process/optical/scintillation/setByParticleType", this);
201 fScintByParticleTypeCmd->SetGuidance("Activate/Inactivate scintillation process by particle type");
202 fScintByParticleTypeCmd->SetParameterName("ScintillationByParticleTypeActivation", false);
203 fScintByParticleTypeCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
204
205 fScintEnhancedTimeConstantsCmd = new G4UIcmdWithABool("/process/optical/scintillation/setEnhancedTimeConstants", this);
206 fScintEnhancedTimeConstantsCmd->SetGuidance("Activate/Inactivate enhanced time constants for scintillation.");
207 fScintEnhancedTimeConstantsCmd->SetGuidance("This will be the default in the next major release.");
208 fScintEnhancedTimeConstantsCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
209
210 fScintTrackInfo1Cmd = new G4UIcmdWithABool("/process/optical/defaults/scintillation/setTrackInfo", this);
211 fScintTrackInfo1Cmd->SetGuidance("Activate/Inactivate scintillation TrackInformation");
212 fScintTrackInfo1Cmd->SetGuidance("DEPRECATED: use /process/optical/scintillation/setTrackInfo instead.");
213 fScintTrackInfo1Cmd->SetParameterName("ScintillationTrackInfo", false);
214 fScintTrackInfo1Cmd->AvailableForStates(G4State_PreInit, G4State_Idle);
215
216 fScintTrackInfoCmd = new G4UIcmdWithABool("/process/optical/scintillation/setTrackInfo", this);
217 fScintTrackInfoCmd->SetGuidance("Activate/Inactivate scintillation TrackInformation");
218 fScintTrackInfoCmd->SetParameterName("ScintillationTrackInfo", false);
219 fScintTrackInfoCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
220
221 fScintFiniteRiseTime1Cmd = new G4UIcmdWithABool("/process/optical/defaults/scintillation/setFiniteRiseTime", this);
222 fScintFiniteRiseTime1Cmd->SetGuidance("Set option of a finite rise-time for G4Scintillation");
223 fScintFiniteRiseTime1Cmd->SetGuidance("If set, the G4Scintillation process expects the user to have set the");
224 fScintFiniteRiseTime1Cmd->SetGuidance("constant material property FAST/SLOWSCINTILLATIONRISETIME");
225 fScintFiniteRiseTime1Cmd->SetGuidance("DEPRECATED: use /process/optical/scintillation/setFiniteRiseTime instead.");
226 fScintFiniteRiseTime1Cmd->SetParameterName("FiniteRiseTime", false);
227 fScintFiniteRiseTime1Cmd->AvailableForStates(G4State_PreInit, G4State_Idle);
228
229 fScintFiniteRiseTimeCmd = new G4UIcmdWithABool("/process/optical/scintillation/setFiniteRiseTime", this);
230 fScintFiniteRiseTimeCmd->SetGuidance("Set option of a finite rise-time for G4Scintillation");
231 fScintFiniteRiseTimeCmd->SetGuidance("If set, the G4Scintillation process expects the user to have set the");
232 fScintFiniteRiseTimeCmd->SetGuidance("constant material property FAST/SLOWSCINTILLATIONRISETIME");
233 fScintFiniteRiseTimeCmd->SetParameterName("FiniteRiseTime", false);
234 fScintFiniteRiseTimeCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
235
236 fScintStackPhotons1Cmd = new G4UIcmdWithABool("/process/optical/defaults/scintillation/setStackPhotons", this);
237 fScintStackPhotons1Cmd->SetGuidance("Set whether or not to stack secondary Scintillation photons");
238 fScintStackPhotons1Cmd->SetGuidance("DEPRECATED: use /process/optical/scintillation/setStackPhotons instead.");
239 fScintStackPhotons1Cmd->SetParameterName("ScintillationStackPhotons", true);
240 fScintStackPhotons1Cmd->SetDefaultValue(true);
241 fScintStackPhotons1Cmd->AvailableForStates(G4State_PreInit, G4State_Idle);
242
243 fScintStackPhotonsCmd = new G4UIcmdWithABool("/process/optical/scintillation/setStackPhotons", this);
244 fScintStackPhotonsCmd->SetGuidance("Set whether or not to stack secondary Scintillation photons");
245 fScintStackPhotonsCmd->SetParameterName("ScintillationStackPhotons", true);
246 fScintStackPhotonsCmd->SetDefaultValue(true);
247 fScintStackPhotonsCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
248
249 fScintTrackSecondariesFirstCmd = new G4UIcmdWithABool("/process/optical/scintillation/setTrackSecondariesFirst", this);
250 fScintTrackSecondariesFirstCmd->SetGuidance("Whether to track scintillation secondaries before primary.");
251 fScintTrackSecondariesFirstCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
252
253 fScintVerboseLevelCmd = new G4UIcmdWithAnInteger("/process/optical/scintillation/verbose", this);
254 fScintVerboseLevelCmd->SetGuidance("Verbose level for scintillation process.");
255 fScintVerboseLevelCmd->SetParameterName("verbose", true);
256 fScintVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
257 fScintVerboseLevelCmd->AvailableForStates(G4State_Idle, G4State_PreInit);
258
259 // WLS //////////////////////////////////
260 fWLSTimeProfile1Cmd = new G4UIcmdWithAString("/process/optical/defaults/wls/setTimeProfile", this);
261 fWLSTimeProfile1Cmd->SetGuidance("Set the WLS time profile (delta or exponential)");
262 fWLSTimeProfile1Cmd->SetGuidance("DEPRECATED: use /process/optical/wls/setTimeProfile instead.");
263 fWLSTimeProfile1Cmd->SetParameterName("WLSTimeProfile", false);
264 fWLSTimeProfile1Cmd->SetCandidates("delta exponential");
265 fWLSTimeProfile1Cmd->AvailableForStates(G4State_PreInit, G4State_Idle);
266
267 fWLSTimeProfileCmd = new G4UIcmdWithAString("/process/optical/wls/setTimeProfile", this);
268 fWLSTimeProfileCmd->SetGuidance("Set the WLS time profile (delta or exponential)");
269 fWLSTimeProfileCmd->SetParameterName("WLSTimeProfile", false);
270 fWLSTimeProfileCmd->SetCandidates("delta exponential");
271 fWLSTimeProfileCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
272
273 fWLSVerboseLevelCmd = new G4UIcmdWithAnInteger("/process/optical/wls/verbose", this);
274 fWLSVerboseLevelCmd->SetGuidance("Verbose level for WLS process.");
275 fWLSVerboseLevelCmd->SetParameterName("verbose", true);
276 fWLSVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
277 fWLSVerboseLevelCmd->SetDefaultValue(1);
278 fWLSVerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
279
280 // WLS2 //////////////////////////////////
281 fWLS2TimeProfileCmd = new G4UIcmdWithAString("/process/optical/wls2/setTimeProfile", this);
282 fWLS2TimeProfileCmd->SetGuidance("Set the WLS2 time profile (delta or exponential)");
283 fWLS2TimeProfileCmd->SetParameterName("WLS2TimeProfile", false);
284 fWLS2TimeProfileCmd->SetCandidates("delta exponential");
285 fWLS2TimeProfileCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
286
287 fWLS2VerboseLevelCmd = new G4UIcmdWithAnInteger("/process/optical/wls2/verbose", this);
288 fWLS2VerboseLevelCmd->SetGuidance("Verbose level for WLS2 process.");
289 fWLS2VerboseLevelCmd->SetParameterName("verbose", true);
290 fWLS2VerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
291 fWLS2VerboseLevelCmd->SetDefaultValue(1);
292 fWLS2VerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
293
294 // boundary //////////////////////////////////////
295 fBoundaryInvokeSD1Cmd = new G4UIcmdWithABool("/process/optical/defaults/boundary/setInvokeSD", this);
296 fBoundaryInvokeSD1Cmd->SetGuidance("Set option for calling InvokeSD in G4OpBoundaryProcess");
297 fBoundaryInvokeSD1Cmd->SetGuidance("DEPRECATED: use /process/optical/boundary/setInvokeSD instead.");
298 fBoundaryInvokeSD1Cmd->SetParameterName("InvokeSD", false);
299 fBoundaryInvokeSD1Cmd->AvailableForStates(G4State_PreInit, G4State_Idle);
300
301 fBoundaryInvokeSDCmd = new G4UIcmdWithABool("/process/optical/boundary/setInvokeSD", this);
302 fBoundaryInvokeSDCmd->SetGuidance("Set option for calling InvokeSD in G4OpBoundaryProcess");
303 fBoundaryInvokeSDCmd->SetParameterName("InvokeSD", false);
304 fBoundaryInvokeSDCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
305
306 fBoundaryVerboseLevelCmd = new G4UIcmdWithAnInteger("/process/optical/boundary/verbose", this);
307 fBoundaryVerboseLevelCmd->SetGuidance("Verbose level for boundary process.");
308 fBoundaryVerboseLevelCmd->SetParameterName("verbose", true);
309 fBoundaryVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
310 fBoundaryVerboseLevelCmd->SetDefaultValue(1);
311 fBoundaryVerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
312
313 // absorption //////////////////////////////////////
314 fBoundaryInvokeSD1Cmd = new G4UIcmdWithABool("/process/optical/defaults/boundary/setInvokeSD", this);
315 fAbsorptionVerboseLevelCmd = new G4UIcmdWithAnInteger("/process/optical/absorption/verbose", this);
316 fAbsorptionVerboseLevelCmd->SetGuidance("Verbose level for absorption process.");
317 fAbsorptionVerboseLevelCmd->SetParameterName("verbose", true);
318 fAbsorptionVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
319 fAbsorptionVerboseLevelCmd->SetDefaultValue(1);
320 fAbsorptionVerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
321
322 // rayleigh //////////////////////////////////////
323 fBoundaryInvokeSD1Cmd = new G4UIcmdWithABool("/process/optical/defaults/boundary/setInvokeSD", this);
324 fRayleighVerboseLevelCmd = new G4UIcmdWithAnInteger("/process/optical/rayleigh/verbose", this);
325 fRayleighVerboseLevelCmd->SetGuidance("Verbose level for Rayleigh process.");
326 fRayleighVerboseLevelCmd->SetParameterName("verbose", true);
327 fRayleighVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
328 fRayleighVerboseLevelCmd->SetDefaultValue(1);
329 fRayleighVerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
330
331 // mie //////////////////////////////////////
332 fMieVerboseLevelCmd = new G4UIcmdWithAnInteger("/process/optical/mie/verbose", this);
333 fMieVerboseLevelCmd->SetGuidance("Verbose level for Mie process.");
334 fMieVerboseLevelCmd->SetParameterName("verbose", true);
335 fMieVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
336 fMieVerboseLevelCmd->SetDefaultValue(1);
337 fMieVerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
338}
@ G4State_Idle
@ G4State_PreInit
@ kNoProcess
Number of processes, no selected process.
G4String G4OpticalProcessName(G4int)
Return the name for a given optical process index.
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(G4bool defVal)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetCandidates(const char *candidateList)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(G4int defVal)
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:146
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:156
void SetRange(const char *rs)
Definition: G4UIcommand.hh:120
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:273
void CreateDirectory(const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
void SetDefaultValue(const char *theDefaultValue)
void SetGuidance(const char *theGuidance)
void SetParameterCandidates(const char *theString)

Referenced by G4OpticalParametersMessenger().

◆ ~G4OpticalParametersMessenger()

G4OpticalParametersMessenger::~G4OpticalParametersMessenger ( )
virtual

Definition at line 340 of file G4OpticalParametersMessenger.cc.

341{
342 delete fDir;
343 delete fDir2;
344 delete fActivateProcessCmd;
345 delete fVerboseCmd;
346 delete fDumpCmd;
347 delete fCerenkovMaxPhotonsCmd;
348 delete fCerenkovMaxPhotons1Cmd;
349 delete fCerenkovMaxBetaChangeCmd;
350 delete fCerenkovMaxBetaChange1Cmd;
351 delete fCerenkovStackPhotonsCmd;
352 delete fCerenkovStackPhotons1Cmd;
353 delete fCerenkovTrackSecondariesFirstCmd;
354 delete fCerenkovVerboseLevelCmd;
355 delete fScintYieldFactorCmd;
356 delete fScintYieldFactor1Cmd;
357 delete fScintByParticleTypeCmd;
358 delete fScintByParticleType1Cmd;
359 delete fScintEnhancedTimeConstantsCmd;
360 delete fScintTrackInfoCmd;
361 delete fScintTrackInfo1Cmd;
362 delete fScintStackPhotonsCmd;
363 delete fScintStackPhotons1Cmd;
364 delete fScintExcitationRatioCmd;
365 delete fScintVerboseLevelCmd;
366 delete fScintFiniteRiseTimeCmd;
367 delete fScintFiniteRiseTime1Cmd;
368 delete fScintTrackSecondariesFirstCmd;
369 delete fWLSTimeProfileCmd;
370 delete fWLSTimeProfile1Cmd;
371 delete fWLSVerboseLevelCmd;
372 delete fWLS2TimeProfileCmd;
373 delete fWLS2VerboseLevelCmd;
374 delete fAbsorptionVerboseLevelCmd;
375 delete fRayleighVerboseLevelCmd;
376 delete fMieVerboseLevelCmd;
377 delete fBoundaryVerboseLevelCmd;
378 delete fTrackSecondariesFirstCmd;
379 delete fBoundaryInvokeSDCmd;
380 delete fBoundaryInvokeSD1Cmd;
381}

Member Function Documentation

◆ SetNewValue()

void G4OpticalParametersMessenger::SetNewValue ( G4UIcommand command,
G4String  newValue 
)
virtual

Apply command to the associated object.

Reimplemented from G4UImessenger.

Definition at line 383 of file G4OpticalParametersMessenger.cc.

385{
386 // physics needs to be rebuilt for all commands
387 G4bool physicsModified = true;
388
389 /// Apply command to the associated object.
390 if (command == fActivateProcessCmd) {
391 std::istringstream is(newValue.data());
392 G4String pn;
393 G4String flag;
394 is >> pn >> flag;
396 params->SetProcessActivation(pn, value);
397 }
398 else if (command == fTrackSecondariesFirstCmd) {
399 std::istringstream is(newValue.data());
400 G4String pn;
401 G4String flag;
402 is >> pn >> flag;
404 if (pn == "Cerenkov") params->SetCerenkovStackPhotons(value);
405 else if (pn == "Scintillation") params->SetScintStackPhotons(value);
406 else {
408 msg << "Process name not allowed: "<<pn<<" (UI: "<<newValue<<")";
409 G4Exception("G4OpticalParametersMessenger::SetNewValue(...)","Optical001",
410 FatalException,msg);
411 }
412 }
413 else if (command == fVerboseCmd) {
414 params->SetVerboseLevel(fVerboseCmd->GetNewIntValue(newValue));
415 }
416 else if (command == fDumpCmd) {
417 params->Dump();
418 }
419 else if (command == fCerenkovMaxPhotons1Cmd) {
421 fCerenkovMaxPhotons1Cmd->GetNewIntValue(newValue));
422 Deprecated();
423 }
424 else if (command == fCerenkovMaxPhotonsCmd) {
426 fCerenkovMaxPhotonsCmd->GetNewIntValue(newValue));
427 G4cout << "Cerenkov max photons: " << params->GetCerenkovMaxPhotonsPerStep() << G4endl;
428 }
429 else if (command == fCerenkovMaxBetaChange1Cmd) {
431 fCerenkovMaxBetaChange1Cmd->GetNewDoubleValue(newValue));
432 Deprecated();
433 }
434 else if (command == fCerenkovMaxBetaChangeCmd) {
436 fCerenkovMaxBetaChangeCmd->GetNewDoubleValue(newValue));
437 }
438 else if (command == fCerenkovStackPhotons1Cmd) {
440 fCerenkovStackPhotons1Cmd->GetNewBoolValue(newValue));
441 Deprecated();
442 }
443 else if (command == fCerenkovStackPhotonsCmd) {
445 fCerenkovStackPhotonsCmd->GetNewBoolValue(newValue));
446 }
447 else if (command == fCerenkovTrackSecondariesFirstCmd) {
449 fCerenkovTrackSecondariesFirstCmd->GetNewBoolValue(newValue));
450 }
451 else if (command == fCerenkovVerboseLevelCmd) {
453 fCerenkovVerboseLevelCmd->GetNewIntValue(newValue));
454 }
455 else if (command == fScintYieldFactor1Cmd) {
456 params->SetScintYieldFactor(
457 fScintYieldFactor1Cmd->GetNewDoubleValue(newValue));
458 Deprecated();
459 }
460 else if (command == fScintYieldFactorCmd) {
461 params->SetScintYieldFactor(
462 fScintYieldFactorCmd->GetNewDoubleValue(newValue));
463 }
464 else if (command == fScintByParticleType1Cmd) {
466 fScintByParticleType1Cmd->GetNewBoolValue(newValue));
467 Deprecated();
468 }
469 else if (command == fScintByParticleTypeCmd) {
471 fScintByParticleTypeCmd->GetNewBoolValue(newValue));
472 }
473 else if (command == fScintEnhancedTimeConstantsCmd) {
475 fScintEnhancedTimeConstantsCmd->GetNewBoolValue(newValue));
476 }
477 else if (command == fScintTrackInfo1Cmd) {
478 params->SetScintTrackInfo(
479 fScintTrackInfo1Cmd->GetNewBoolValue(newValue));
480 Deprecated();
481 }
482 else if (command == fScintTrackInfoCmd) {
483 params->SetScintTrackInfo(
484 fScintTrackInfoCmd->GetNewBoolValue(newValue));
485 }
486 else if (command == fScintFiniteRiseTime1Cmd) {
488 fScintFiniteRiseTime1Cmd->GetNewBoolValue(newValue));
489 Deprecated();
490 }
491 else if (command == fScintFiniteRiseTimeCmd) {
493 fScintFiniteRiseTimeCmd->GetNewBoolValue(newValue));
494 }
495 else if (command == fScintStackPhotons1Cmd) {
496 params->SetScintStackPhotons(
497 fScintStackPhotons1Cmd->GetNewBoolValue(newValue));
498 Deprecated();
499 }
500 else if (command == fScintStackPhotonsCmd) {
501 params->SetScintStackPhotons(
502 fScintStackPhotonsCmd->GetNewBoolValue(newValue));
503 }
504 else if (command == fScintExcitationRatioCmd) {
506 fScintExcitationRatioCmd->GetNewDoubleValue(newValue));
507 }
508 else if (command == fScintTrackSecondariesFirstCmd) {
510 fScintTrackSecondariesFirstCmd->GetNewBoolValue(newValue));
511 }
512 else if (command == fScintVerboseLevelCmd) {
513 params->SetScintVerboseLevel(
514 fScintVerboseLevelCmd->GetNewIntValue(newValue));
515 }
516 else if (command == fWLSTimeProfile1Cmd) {
517 params->SetWLSTimeProfile(newValue);
518 Deprecated();
519 }
520 else if (command == fWLSTimeProfileCmd) {
521 params->SetWLSTimeProfile(newValue);
522 }
523 else if (command == fWLSVerboseLevelCmd) {
524 params->SetWLSVerboseLevel(fWLSVerboseLevelCmd->GetNewIntValue(newValue));
525 }
526 else if (command == fWLS2TimeProfileCmd) {
527 params->SetWLS2TimeProfile(newValue);
528 }
529 else if (command == fWLS2VerboseLevelCmd) {
530 params->SetWLS2VerboseLevel(fWLS2VerboseLevelCmd->GetNewIntValue(newValue));
531 }
532 else if (command == fAbsorptionVerboseLevelCmd) {
533 params->SetAbsorptionVerboseLevel(fAbsorptionVerboseLevelCmd->GetNewIntValue(newValue));
534 }
535 else if (command == fRayleighVerboseLevelCmd) {
536 params->SetRayleighVerboseLevel(fRayleighVerboseLevelCmd->GetNewIntValue(newValue));
537 }
538 else if (command == fMieVerboseLevelCmd) {
539 params->SetMieVerboseLevel(fMieVerboseLevelCmd->GetNewIntValue(newValue));
540 }
541 else if (command == fBoundaryVerboseLevelCmd) {
542 params->SetBoundaryVerboseLevel(fBoundaryVerboseLevelCmd->GetNewIntValue(newValue));
543 }
544 else if (command == fBoundaryInvokeSD1Cmd) {
545 params->SetBoundaryInvokeSD(fBoundaryInvokeSD1Cmd->GetNewBoolValue(newValue));
546 Deprecated();
547 }
548 else if (command == fBoundaryInvokeSDCmd) {
549 params->SetBoundaryInvokeSD(fBoundaryInvokeSDCmd->GetNewBoolValue(newValue));
550 }
551 if (physicsModified) {
552 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
553 }
554}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetScintByParticleType(G4bool)
void SetCerenkovMaxBetaChange(G4double)
void SetRayleighVerboseLevel(G4int)
void SetCerenkovMaxPhotonsPerStep(G4int)
void SetBoundaryInvokeSD(G4bool)
void SetBoundaryVerboseLevel(G4int)
void SetScintTrackSecondariesFirst(G4bool)
void SetScintEnhancedTimeConstants(G4bool)
void SetScintStackPhotons(G4bool)
void SetWLS2TimeProfile(const G4String &)
G4int GetCerenkovMaxPhotonsPerStep() const
void SetAbsorptionVerboseLevel(G4int)
void SetCerenkovStackPhotons(G4bool)
void SetCerenkovTrackSecondariesFirst(G4bool)
void SetScintFiniteRiseTime(G4bool)
void SetWLSTimeProfile(const G4String &)
void SetCerenkovVerboseLevel(G4int)
void SetScintExcitationRatio(G4double)
void SetScintYieldFactor(G4double)
void SetProcessActivation(const G4String &, G4bool)
const char * data() const
static G4bool GetNewBoolValue(const char *paramString)
static G4double GetNewDoubleValue(const char *paramString)
static G4int GetNewIntValue(const char *paramString)
static G4bool ConvertToBool(const char *st)
Definition: G4UIcommand.cc:530
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:485
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77

The documentation for this class was generated from the following files: