Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4OpticalParametersMessenger.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26//
27//----------------------------------------------------------------------------
28//
29// ClassName: G4OpticalParametersMessenger
30//
31// Author: P.Gumplinger 30.09.2009 //
32//
33// Modified: P.Gumplinger 29.09.2011
34// (based on code from I. Hrivnacova)
35//
36//----------------------------------------------------------------------------
37//
38
41
42#include "G4UIcommand.hh"
43#include "G4UIdirectory.hh"
44#include "G4UIcmdWithABool.hh"
45#include "G4UIcmdWithAString.hh"
46#include "G4UIcmdWithADouble.hh"
49#include "G4UImanager.hh"
50#include "G4UIparameter.hh"
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
55 G4OpticalParameters* opticalParameters)
56 : params(opticalParameters)
57
58{
59 G4bool toBeBroadcasted = false;
60 fDir = new G4UIdirectory("/process/optical/", toBeBroadcasted);
61 fDir->SetGuidance(
62 "Commands related to the optical physics simulation engine.");
63
64 fCerenkovDir =
65 new G4UIdirectory("/process/optical/cerenkov/", toBeBroadcasted);
66 fCerenkovDir->SetGuidance("Cerenkov process commands");
67 fScintDir =
68 new G4UIdirectory("/process/optical/scintillation/", toBeBroadcasted);
69 fScintDir->SetGuidance("Scintillation process commands");
70 fWlsDir = new G4UIdirectory("/process/optical/wls/", toBeBroadcasted);
71 fWlsDir->SetGuidance("Wave length shifting process commands");
72 fWls2Dir = new G4UIdirectory("/process/optical/wls2/", toBeBroadcasted);
73 fWls2Dir->SetGuidance("Second Wave length shifting process commands");
74 fBoundaryDir =
75 new G4UIdirectory("/process/optical/boundary/", toBeBroadcasted);
76 fBoundaryDir->SetGuidance("Boundary scattering commands");
77 fMieDir = new G4UIdirectory("/process/optical/mie/", toBeBroadcasted);
78 fMieDir->SetGuidance("Mie scattering process commands");
79 fAbsDir = new G4UIdirectory("/process/optical/absorption/", toBeBroadcasted);
80 fAbsDir->SetGuidance("absorption process commands");
81 fRaylDir = new G4UIdirectory("/process/optical/rayleigh/", toBeBroadcasted);
82 fRaylDir->SetGuidance("Rayleigh scattering commands");
83
84 // general commands
85 fActivateProcessCmd =
86 new G4UIcommand("/process/optical/processActivation", this);
87 fActivateProcessCmd->SetGuidance(
88 "Activate/deactivate the specified optical process");
89 auto par = new G4UIparameter("proc_name", 's', false);
90 G4String candidates;
91 for(G4int i = 0; i < kNoProcess; ++i)
92 {
93 candidates += G4OpticalProcessName(i);
94 candidates += G4String(" ");
95 }
96 par->SetParameterCandidates(candidates);
97 par->SetGuidance("the process name");
98 fActivateProcessCmd->SetParameter(par);
99 par = new G4UIparameter("flag", 'b', true);
100 par->SetDefaultValue(true);
101 par->SetGuidance("activation flag");
102 fActivateProcessCmd->SetParameter(par);
103 fActivateProcessCmd->AvailableForStates(G4State_PreInit);
104
105 fVerboseCmd = new G4UIcmdWithAnInteger("/process/optical/verbose", this);
106 fVerboseCmd->SetGuidance("Set default verbose level for optical processes");
107 fVerboseCmd->SetParameterName("ver", true);
108 fVerboseCmd->SetDefaultValue(1);
109 fVerboseCmd->SetRange("ver>=0");
111
112 fDumpCmd = new G4UIcommand("/process/optical/printParameters", this);
113 fDumpCmd->SetGuidance("Print all optical parameters.");
114
115 // Cerenkov ////////////////////
116 fCerenkovMaxPhotonsCmd =
117 new G4UIcmdWithAnInteger("/process/optical/cerenkov/setMaxPhotons", this);
118 fCerenkovMaxPhotonsCmd->SetGuidance("Set maximum number of photons per step");
119 fCerenkovMaxPhotonsCmd->SetParameterName("CerenkovMaxPhotons", false);
120 fCerenkovMaxPhotonsCmd->SetRange("CerenkovMaxPhotons>=0");
121 fCerenkovMaxPhotonsCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
122
123 fCerenkovMaxBetaChangeCmd =
124 new G4UIcmdWithADouble("/process/optical/cerenkov/setMaxBetaChange", this);
125 fCerenkovMaxBetaChangeCmd->SetGuidance(
126 "Set maximum change of beta of parent particle per step (in percent)");
127 fCerenkovMaxBetaChangeCmd->SetParameterName("CerenkovMaxBetaChange", false);
128 fCerenkovMaxBetaChangeCmd->SetRange("CerenkovMaxBetaChange>=0");
129 fCerenkovMaxBetaChangeCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
130
131 fCerenkovStackPhotonsCmd =
132 new G4UIcmdWithABool("/process/optical/cerenkov/setStackPhotons", this);
133 fCerenkovStackPhotonsCmd->SetGuidance(
134 "Set whether or not to stack secondary Cerenkov photons");
135 fCerenkovStackPhotonsCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
136
137 fCerenkovTrackSecondariesFirstCmd = new G4UIcmdWithABool(
138 "/process/optical/cerenkov/setTrackSecondariesFirst", this);
139 fCerenkovTrackSecondariesFirstCmd->SetGuidance(
140 "Whether to track secondary Cerenkov photons before the primary.");
141 fCerenkovTrackSecondariesFirstCmd->AvailableForStates(G4State_PreInit,
143
144 fCerenkovVerboseLevelCmd =
145 new G4UIcmdWithAnInteger("/process/optical/cerenkov/verbose", this);
146 fCerenkovVerboseLevelCmd->SetGuidance("Verbose level for Cerenkov process.");
147 fCerenkovVerboseLevelCmd->SetParameterName("verbose", true);
148 fCerenkovVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
149 fCerenkovVerboseLevelCmd->SetDefaultValue(2);
150 fCerenkovVerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
151
152 // Scintillation //////////////////////////
153 fScintByParticleTypeCmd = new G4UIcmdWithABool(
154 "/process/optical/scintillation/setByParticleType", this);
155 fScintByParticleTypeCmd->SetGuidance(
156 "Activate/Inactivate scintillation process by particle type");
157 fScintByParticleTypeCmd->SetParameterName(
158 "ScintillationByParticleTypeActivation", false);
159 fScintByParticleTypeCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
160
161 fScintTrackInfoCmd =
162 new G4UIcmdWithABool("/process/optical/scintillation/setTrackInfo", this);
163 fScintTrackInfoCmd->SetGuidance(
164 "Activate/Inactivate scintillation TrackInformation");
165 fScintTrackInfoCmd->SetParameterName("ScintillationTrackInfo", false);
166 fScintTrackInfoCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
167
168 fScintFiniteRiseTimeCmd = new G4UIcmdWithABool(
169 "/process/optical/scintillation/setFiniteRiseTime", this);
170 fScintFiniteRiseTimeCmd->SetGuidance(
171 "Set option of a finite rise-time for G4Scintillation");
172 fScintFiniteRiseTimeCmd->SetGuidance(
173 "If set, the G4Scintillation process expects the user to have set the");
174 fScintFiniteRiseTimeCmd->SetGuidance(
175 "constant material property SCINTILLATIONRISETIME{1,2,3}");
176 fScintFiniteRiseTimeCmd->SetParameterName("FiniteRiseTime", false);
177 fScintFiniteRiseTimeCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
178
179 fScintStackPhotonsCmd = new G4UIcmdWithABool(
180 "/process/optical/scintillation/setStackPhotons", this);
181 fScintStackPhotonsCmd->SetGuidance(
182 "Set whether or not to stack secondary Scintillation photons");
183 fScintStackPhotonsCmd->SetParameterName("ScintillationStackPhotons", true);
184 fScintStackPhotonsCmd->SetDefaultValue(true);
185 fScintStackPhotonsCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
186
187 fScintTrackSecondariesFirstCmd = new G4UIcmdWithABool(
188 "/process/optical/scintillation/setTrackSecondariesFirst", this);
189 fScintTrackSecondariesFirstCmd->SetGuidance(
190 "Whether to track scintillation secondaries before primary.");
191 fScintTrackSecondariesFirstCmd->AvailableForStates(G4State_PreInit,
193
194 fScintVerboseLevelCmd =
195 new G4UIcmdWithAnInteger("/process/optical/scintillation/verbose", this);
196 fScintVerboseLevelCmd->SetGuidance(
197 "Verbose level for scintillation process.");
198 fScintVerboseLevelCmd->SetParameterName("verbose", true);
199 fScintVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
200 fScintVerboseLevelCmd->AvailableForStates(G4State_Idle, G4State_PreInit);
201
202 // WLS //////////////////////////////////
203 fWLSTimeProfileCmd =
204 new G4UIcmdWithAString("/process/optical/wls/setTimeProfile", this);
205 fWLSTimeProfileCmd->SetGuidance(
206 "Set the WLS time profile (delta or exponential)");
207 fWLSTimeProfileCmd->SetParameterName("WLSTimeProfile", false);
208 fWLSTimeProfileCmd->SetCandidates("delta exponential");
209 fWLSTimeProfileCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
210
211 fWLSVerboseLevelCmd =
212 new G4UIcmdWithAnInteger("/process/optical/wls/verbose", this);
213 fWLSVerboseLevelCmd->SetGuidance("Verbose level for WLS process.");
214 fWLSVerboseLevelCmd->SetParameterName("verbose", true);
215 fWLSVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
216 fWLSVerboseLevelCmd->SetDefaultValue(1);
217 fWLSVerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
218
219 // WLS2 //////////////////////////////////
220 fWLS2TimeProfileCmd =
221 new G4UIcmdWithAString("/process/optical/wls2/setTimeProfile", this);
222 fWLS2TimeProfileCmd->SetGuidance(
223 "Set the WLS2 time profile (delta or exponential)");
224 fWLS2TimeProfileCmd->SetParameterName("WLS2TimeProfile", false);
225 fWLS2TimeProfileCmd->SetCandidates("delta exponential");
226 fWLS2TimeProfileCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
227
228 fWLS2VerboseLevelCmd =
229 new G4UIcmdWithAnInteger("/process/optical/wls2/verbose", this);
230 fWLS2VerboseLevelCmd->SetGuidance("Verbose level for WLS2 process.");
231 fWLS2VerboseLevelCmd->SetParameterName("verbose", true);
232 fWLS2VerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
233 fWLS2VerboseLevelCmd->SetDefaultValue(1);
234 fWLS2VerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
235
236 // boundary //////////////////////////////////////
237 fBoundaryInvokeSDCmd =
238 new G4UIcmdWithABool("/process/optical/boundary/setInvokeSD", this);
239 fBoundaryInvokeSDCmd->SetGuidance(
240 "Set option for calling InvokeSD in G4OpBoundaryProcess");
241 fBoundaryInvokeSDCmd->SetParameterName("InvokeSD", false);
242 fBoundaryInvokeSDCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
243
244 fBoundaryVerboseLevelCmd =
245 new G4UIcmdWithAnInteger("/process/optical/boundary/verbose", this);
246 fBoundaryVerboseLevelCmd->SetGuidance("Verbose level for boundary process.");
247 fBoundaryVerboseLevelCmd->SetParameterName("verbose", true);
248 fBoundaryVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
249 fBoundaryVerboseLevelCmd->SetDefaultValue(1);
250 fBoundaryVerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
251
252 // absorption //////////////////////////////////////
253 fAbsorptionVerboseLevelCmd =
254 new G4UIcmdWithAnInteger("/process/optical/absorption/verbose", this);
255 fAbsorptionVerboseLevelCmd->SetGuidance(
256 "Verbose level for absorption process.");
257 fAbsorptionVerboseLevelCmd->SetParameterName("verbose", true);
258 fAbsorptionVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
259 fAbsorptionVerboseLevelCmd->SetDefaultValue(1);
260 fAbsorptionVerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
261
262 // rayleigh //////////////////////////////////////
263 fRayleighVerboseLevelCmd =
264 new G4UIcmdWithAnInteger("/process/optical/rayleigh/verbose", this);
265 fRayleighVerboseLevelCmd->SetGuidance("Verbose level for Rayleigh process.");
266 fRayleighVerboseLevelCmd->SetParameterName("verbose", true);
267 fRayleighVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
268 fRayleighVerboseLevelCmd->SetDefaultValue(1);
269 fRayleighVerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
270
271 // mie //////////////////////////////////////
272 fMieVerboseLevelCmd =
273 new G4UIcmdWithAnInteger("/process/optical/mie/verbose", this);
274 fMieVerboseLevelCmd->SetGuidance("Verbose level for Mie process.");
275 fMieVerboseLevelCmd->SetParameterName("verbose", true);
276 fMieVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
277 fMieVerboseLevelCmd->SetDefaultValue(1);
278 fMieVerboseLevelCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
279}
280
282{
283 delete fDir;
284 delete fCerenkovDir;
285 delete fScintDir;
286 delete fWlsDir;
287 delete fBoundaryDir;
288 delete fMieDir;
289 delete fAbsDir;
290 delete fRaylDir;
291 delete fActivateProcessCmd;
292 delete fVerboseCmd;
293 delete fDumpCmd;
294 delete fCerenkovMaxPhotonsCmd;
295 delete fCerenkovMaxBetaChangeCmd;
296 delete fCerenkovStackPhotonsCmd;
297 delete fCerenkovTrackSecondariesFirstCmd;
298 delete fCerenkovVerboseLevelCmd;
299 delete fScintByParticleTypeCmd;
300 delete fScintTrackInfoCmd;
301 delete fScintStackPhotonsCmd;
302 delete fScintVerboseLevelCmd;
303 delete fScintFiniteRiseTimeCmd;
304 delete fScintTrackSecondariesFirstCmd;
305 delete fWLSTimeProfileCmd;
306 delete fWLSVerboseLevelCmd;
307 delete fWLS2TimeProfileCmd;
308 delete fWLS2VerboseLevelCmd;
309 delete fAbsorptionVerboseLevelCmd;
310 delete fRayleighVerboseLevelCmd;
311 delete fMieVerboseLevelCmd;
312 delete fBoundaryVerboseLevelCmd;
313 delete fBoundaryInvokeSDCmd;
314}
315
317 G4String newValue)
318{
319 // physics needs to be rebuilt for all commands
320 G4bool physicsModified = true;
321
322 /// Apply command to the associated object.
323 if(command == fActivateProcessCmd)
324 {
325 std::istringstream is(newValue.data());
326 G4String pn;
327 G4String flag;
328 is >> pn >> flag;
330 params->SetProcessActivation(pn, value);
331 }
332 else if(command == fVerboseCmd)
333 {
334 params->SetVerboseLevel(fVerboseCmd->GetNewIntValue(newValue));
335 }
336 else if(command == fDumpCmd)
337 {
338 params->Dump();
339 }
340 else if(command == fCerenkovMaxPhotonsCmd)
341 {
343 fCerenkovMaxPhotonsCmd->GetNewIntValue(newValue));
344 G4cout << "Cerenkov max photons: " << params->GetCerenkovMaxPhotonsPerStep()
345 << G4endl;
346 }
347 else if(command == fCerenkovMaxBetaChangeCmd)
348 {
350 fCerenkovMaxBetaChangeCmd->GetNewDoubleValue(newValue));
351 }
352 else if(command == fCerenkovStackPhotonsCmd)
353 {
355 fCerenkovStackPhotonsCmd->GetNewBoolValue(newValue));
356 }
357 else if(command == fCerenkovTrackSecondariesFirstCmd)
358 {
360 fCerenkovTrackSecondariesFirstCmd->GetNewBoolValue(newValue));
361 }
362 else if(command == fCerenkovVerboseLevelCmd)
363 {
365 fCerenkovVerboseLevelCmd->GetNewIntValue(newValue));
366 }
367 else if(command == fScintByParticleTypeCmd)
368 {
370 fScintByParticleTypeCmd->GetNewBoolValue(newValue));
371 }
372 else if(command == fScintTrackInfoCmd)
373 {
374 params->SetScintTrackInfo(fScintTrackInfoCmd->GetNewBoolValue(newValue));
375 }
376 else if(command == fScintFiniteRiseTimeCmd)
377 {
379 fScintFiniteRiseTimeCmd->GetNewBoolValue(newValue));
380 }
381 else if(command == fScintStackPhotonsCmd)
382 {
383 params->SetScintStackPhotons(
384 fScintStackPhotonsCmd->GetNewBoolValue(newValue));
385 }
386 else if(command == fScintTrackSecondariesFirstCmd)
387 {
389 fScintTrackSecondariesFirstCmd->GetNewBoolValue(newValue));
390 }
391 else if(command == fScintVerboseLevelCmd)
392 {
393 params->SetScintVerboseLevel(
394 fScintVerboseLevelCmd->GetNewIntValue(newValue));
395 }
396 else if(command == fWLSTimeProfileCmd)
397 {
398 params->SetWLSTimeProfile(newValue);
399 }
400 else if(command == fWLSVerboseLevelCmd)
401 {
402 params->SetWLSVerboseLevel(fWLSVerboseLevelCmd->GetNewIntValue(newValue));
403 }
404 else if(command == fWLS2TimeProfileCmd)
405 {
406 params->SetWLS2TimeProfile(newValue);
407 }
408 else if(command == fWLS2VerboseLevelCmd)
409 {
410 params->SetWLS2VerboseLevel(fWLS2VerboseLevelCmd->GetNewIntValue(newValue));
411 }
412 else if(command == fAbsorptionVerboseLevelCmd)
413 {
415 fAbsorptionVerboseLevelCmd->GetNewIntValue(newValue));
416 }
417 else if(command == fRayleighVerboseLevelCmd)
418 {
420 fRayleighVerboseLevelCmd->GetNewIntValue(newValue));
421 }
422 else if(command == fMieVerboseLevelCmd)
423 {
424 params->SetMieVerboseLevel(fMieVerboseLevelCmd->GetNewIntValue(newValue));
425 }
426 else if(command == fBoundaryVerboseLevelCmd)
427 {
429 fBoundaryVerboseLevelCmd->GetNewIntValue(newValue));
430 }
431 else if(command == fBoundaryInvokeSDCmd)
432 {
433 params->SetBoundaryInvokeSD(
434 fBoundaryInvokeSDCmd->GetNewBoolValue(newValue));
435 }
436 if(physicsModified)
437 {
438 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
439 }
440}
@ 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
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual void SetNewValue(G4UIcommand *, G4String)
G4OpticalParametersMessenger(G4OpticalParameters *)
void SetScintByParticleType(G4bool)
void SetCerenkovMaxBetaChange(G4double)
void SetRayleighVerboseLevel(G4int)
void SetCerenkovMaxPhotonsPerStep(G4int)
void SetBoundaryInvokeSD(G4bool)
void SetBoundaryVerboseLevel(G4int)
void SetScintTrackSecondariesFirst(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 SetProcessActivation(const G4String &, G4bool)
static G4bool GetNewBoolValue(const char *paramString)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetDefaultValue(G4bool defVal)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
static G4double GetNewDoubleValue(const char *paramString)
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)
static G4int GetNewIntValue(const char *paramString)
void SetDefaultValue(G4int defVal)
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:147
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:157
static G4bool ConvertToBool(const char *st)
Definition: G4UIcommand.cc:549
void SetRange(const char *rs)
Definition: G4UIcommand.hh:121
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:287
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:495
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77