Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4OpticalParametersMessenger Class Reference

#include <G4OpticalParametersMessenger.hh>

Inheritance diagram for G4OpticalParametersMessenger:
G4UImessenger

Public Member Functions

G4bool CommandsShouldBeInMaster () const
 
 G4OpticalParametersMessenger (G4OpticalParameters *)
 
virtual G4String GetCurrentValue (G4UIcommand *command)
 
G4bool operator!= (const G4UImessenger &messenger) const
 
G4bool operator== (const G4UImessenger &messenger) const
 
virtual void SetNewValue (G4UIcommand *, G4String)
 
virtual ~G4OpticalParametersMessenger ()
 

Protected Member Functions

void AddUIcommand (G4UIcommand *newCommand)
 
G4String BtoS (G4bool b)
 
template<typename T >
T * CreateCommand (const G4String &cname, const G4String &dsc)
 
void CreateDirectory (const G4String &path, const G4String &dsc, G4bool commandsToBeBroadcasted=true)
 
G4String DtoS (G4double a)
 
G4String ItoS (G4int i)
 
G4bool StoB (G4String s)
 
G4double StoD (G4String s)
 
G4int StoI (G4String s)
 
G4long StoL (G4String s)
 

Protected Attributes

G4UIdirectorybaseDir = nullptr
 
G4String baseDirName = ""
 
G4bool commandsShouldBeInMaster = false
 

Private Member Functions

 G4OpticalParametersMessenger ()=delete
 
 G4OpticalParametersMessenger (const G4OpticalParametersMessenger &right)=delete
 
G4OpticalParametersMessengeroperator= (const G4OpticalParametersMessenger &right)=delete
 

Private Attributes

G4UIdirectoryfAbsDir
 
G4UIcmdWithAnIntegerfAbsorptionVerboseLevelCmd
 
G4UIcommandfActivateProcessCmd
 selectOpProcess command More...
 
G4UIdirectoryfBoundaryDir
 
G4UIcmdWithABoolfBoundaryInvokeSDCmd
 setInvokeSD command More...
 
G4UIcmdWithAnIntegerfBoundaryVerboseLevelCmd
 
G4UIdirectoryfCerenkovDir
 
G4UIcmdWithADoublefCerenkovMaxBetaChangeCmd
 setCerenkovMaxBetaChange command More...
 
G4UIcmdWithAnIntegerfCerenkovMaxPhotonsCmd
 
G4UIcmdWithABoolfCerenkovStackPhotonsCmd
 setStackPhotons command More...
 
G4UIcmdWithABoolfCerenkovTrackSecondariesFirstCmd
 
G4UIcmdWithAnIntegerfCerenkovVerboseLevelCmd
 
G4UIdirectoryfDir
 command directory More...
 
G4UIcommandfDumpCmd
 
G4UIdirectoryfMieDir
 
G4UIcmdWithAnIntegerfMieVerboseLevelCmd
 
G4UIdirectoryfRaylDir
 
G4UIcmdWithAnIntegerfRayleighVerboseLevelCmd
 
G4UIcmdWithABoolfScintByParticleTypeCmd
 setScintillationByParticleType command More...
 
G4UIdirectoryfScintDir
 
G4UIcmdWithABoolfScintFiniteRiseTimeCmd
 setFiniteRiseTime command More...
 
G4UIcmdWithABoolfScintStackPhotonsCmd
 setStackPhotons command More...
 
G4UIcmdWithABoolfScintTrackInfoCmd
 setScintillationTrackInfo command More...
 
G4UIcmdWithABoolfScintTrackSecondariesFirstCmd
 
G4UIcmdWithAnIntegerfScintVerboseLevelCmd
 
G4UIcmdWithAnIntegerfVerboseCmd
 setProcessVerbose command More...
 
G4UIdirectoryfWls2Dir
 
G4UIcmdWithAStringfWLS2TimeProfileCmd
 setWLS2TimeProfile command More...
 
G4UIcmdWithAnIntegerfWLS2VerboseLevelCmd
 
G4UIdirectoryfWlsDir
 
G4UIcmdWithAStringfWLSTimeProfileCmd
 setWLSTimeProfile command More...
 
G4UIcmdWithAnIntegerfWLSVerboseLevelCmd
 
G4OpticalParametersparams
 associated class More...
 

Detailed Description

Definition at line 61 of file G4OpticalParametersMessenger.hh.

Constructor & Destructor Documentation

◆ G4OpticalParametersMessenger() [1/3]

G4OpticalParametersMessenger::G4OpticalParametersMessenger ( G4OpticalParameters opticalParameters)

Definition at line 54 of file G4OpticalParametersMessenger.cc.

56 : params(opticalParameters)
57
58{
59 G4bool toBeBroadcasted = false;
60 fDir = new G4UIdirectory("/process/optical/", toBeBroadcasted);
62 "Commands related to the optical physics simulation engine.");
63
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");
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
86 new G4UIcommand("/process/optical/processActivation", this);
88 "Activate/deactivate the specified optical process");
89 G4UIparameter* 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");
99 par = new G4UIparameter("flag", 'b', true);
100 par->SetDefaultValue(true);
101 par->SetGuidance("activation flag");
104
105 fVerboseCmd = new G4UIcmdWithAnInteger("/process/optical/verbose", this);
106 fVerboseCmd->SetGuidance("Set default verbose level for optical processes");
107 fVerboseCmd->SetParameterName("ver", true);
109 fVerboseCmd->SetRange("ver>=0");
111
112 fDumpCmd = new G4UIcommand("/process/optical/printParameters", this);
113 fDumpCmd->SetGuidance("Print all optical parameters.");
114
115 // Cerenkov ////////////////////
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");
122
124 new G4UIcmdWithADouble("/process/optical/cerenkov/setMaxBetaChange", this);
126 "Set maximum change of beta of parent particle per step (in percent)");
127 fCerenkovMaxBetaChangeCmd->SetParameterName("CerenkovMaxBetaChange", false);
128 fCerenkovMaxBetaChangeCmd->SetRange("CerenkovMaxBetaChange>=0");
130
132 new G4UIcmdWithABool("/process/optical/cerenkov/setStackPhotons", this);
134 "Set whether or not to stack secondary Cerenkov photons");
136
138 "/process/optical/cerenkov/setTrackSecondariesFirst", this);
140 "Whether to track secondary Cerenkov photons before the primary.");
143
145 new G4UIcmdWithAnInteger("/process/optical/cerenkov/verbose", this);
146 fCerenkovVerboseLevelCmd->SetGuidance("Verbose level for Cerenkov process.");
148 fCerenkovVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
151
152 // Scintillation //////////////////////////
154 "/process/optical/scintillation/setByParticleType", this);
156 "Activate/Inactivate scintillation process by particle type");
158 "ScintillationByParticleTypeActivation", false);
160
162 new G4UIcmdWithABool("/process/optical/scintillation/setTrackInfo", this);
164 "Activate/Inactivate scintillation TrackInformation");
165 fScintTrackInfoCmd->SetParameterName("ScintillationTrackInfo", false);
167
169 "/process/optical/scintillation/setFiniteRiseTime", this);
171 "Set option of a finite rise-time for G4Scintillation");
173 "If set, the G4Scintillation process expects the user to have set the");
175 "constant material property SCINTILLATIONRISETIME{1,2,3}");
176 fScintFiniteRiseTimeCmd->SetParameterName("FiniteRiseTime", false);
178
180 "/process/optical/scintillation/setStackPhotons", this);
182 "Set whether or not to stack secondary Scintillation photons");
183 fScintStackPhotonsCmd->SetParameterName("ScintillationStackPhotons", true);
186
188 "/process/optical/scintillation/setTrackSecondariesFirst", this);
190 "Whether to track scintillation secondaries before primary.");
193
195 new G4UIcmdWithAnInteger("/process/optical/scintillation/verbose", this);
197 "Verbose level for scintillation process.");
198 fScintVerboseLevelCmd->SetParameterName("verbose", true);
199 fScintVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
201
202 // WLS //////////////////////////////////
204 new G4UIcmdWithAString("/process/optical/wls/setTimeProfile", this);
206 "Set the WLS time profile (delta or exponential)");
207 fWLSTimeProfileCmd->SetParameterName("WLSTimeProfile", false);
208 fWLSTimeProfileCmd->SetCandidates("delta exponential");
210
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");
218
219 // WLS2 //////////////////////////////////
221 new G4UIcmdWithAString("/process/optical/wls2/setTimeProfile", this);
223 "Set the WLS2 time profile (delta or exponential)");
224 fWLS2TimeProfileCmd->SetParameterName("WLS2TimeProfile", false);
225 fWLS2TimeProfileCmd->SetCandidates("delta exponential");
227
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");
235
236 // boundary //////////////////////////////////////
238 new G4UIcmdWithABool("/process/optical/boundary/setInvokeSD", this);
240 "Set option for calling InvokeSD in G4OpBoundaryProcess");
241 fBoundaryInvokeSDCmd->SetParameterName("InvokeSD", false);
243
245 new G4UIcmdWithAnInteger("/process/optical/boundary/verbose", this);
246 fBoundaryVerboseLevelCmd->SetGuidance("Verbose level for boundary process.");
248 fBoundaryVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
251
252 // absorption //////////////////////////////////////
254 new G4UIcmdWithAnInteger("/process/optical/absorption/verbose", this);
256 "Verbose level for absorption process.");
258 fAbsorptionVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
261
262 // rayleigh //////////////////////////////////////
264 new G4UIcmdWithAnInteger("/process/optical/rayleigh/verbose", this);
265 fRayleighVerboseLevelCmd->SetGuidance("Verbose level for Rayleigh process.");
267 fRayleighVerboseLevelCmd->SetRange("verbose >= 0 && verbose <= 2");
270
271 // mie //////////////////////////////////////
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");
279}
@ 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
G4UIcmdWithABool * fBoundaryInvokeSDCmd
setInvokeSD command
G4UIcmdWithABool * fScintStackPhotonsCmd
setStackPhotons command
G4UIcmdWithABool * fScintFiniteRiseTimeCmd
setFiniteRiseTime command
G4UIcmdWithADouble * fCerenkovMaxBetaChangeCmd
setCerenkovMaxBetaChange command
G4UIcommand * fActivateProcessCmd
selectOpProcess command
G4UIcmdWithAnInteger * fVerboseCmd
setProcessVerbose command
G4UIcmdWithABool * fScintByParticleTypeCmd
setScintillationByParticleType command
G4UIcmdWithAnInteger * fAbsorptionVerboseLevelCmd
G4UIdirectory * fDir
command directory
G4UIcmdWithABool * fScintTrackInfoCmd
setScintillationTrackInfo command
G4UIcmdWithAString * fWLSTimeProfileCmd
setWLSTimeProfile command
G4UIcmdWithAString * fWLS2TimeProfileCmd
setWLS2TimeProfile command
G4UIcmdWithABool * fCerenkovStackPhotonsCmd
setStackPhotons command
G4OpticalParameters * params
associated class
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:288
void SetDefaultValue(const char *theDefaultValue)
void SetGuidance(const char *theGuidance)
void SetParameterCandidates(const char *theString)

References G4UIcommand::AvailableForStates(), fAbsDir, fAbsorptionVerboseLevelCmd, fActivateProcessCmd, fBoundaryDir, fBoundaryInvokeSDCmd, fBoundaryVerboseLevelCmd, fCerenkovDir, fCerenkovMaxBetaChangeCmd, fCerenkovMaxPhotonsCmd, fCerenkovStackPhotonsCmd, fCerenkovTrackSecondariesFirstCmd, fCerenkovVerboseLevelCmd, fDir, fDumpCmd, fMieDir, fMieVerboseLevelCmd, fRaylDir, fRayleighVerboseLevelCmd, fScintByParticleTypeCmd, fScintDir, fScintFiniteRiseTimeCmd, fScintStackPhotonsCmd, fScintTrackInfoCmd, fScintTrackSecondariesFirstCmd, fScintVerboseLevelCmd, fVerboseCmd, fWls2Dir, fWLS2TimeProfileCmd, fWLS2VerboseLevelCmd, fWlsDir, fWLSTimeProfileCmd, fWLSVerboseLevelCmd, G4OpticalProcessName(), G4State_Idle, G4State_PreInit, kNoProcess, G4UIcmdWithAString::SetCandidates(), G4UIparameter::SetDefaultValue(), G4UIcmdWithABool::SetDefaultValue(), G4UIcmdWithAnInteger::SetDefaultValue(), G4UIcommand::SetGuidance(), G4UIparameter::SetGuidance(), G4UIcommand::SetParameter(), G4UIparameter::SetParameterCandidates(), G4UIcmdWithABool::SetParameterName(), G4UIcmdWithADouble::SetParameterName(), G4UIcmdWithAnInteger::SetParameterName(), G4UIcmdWithAString::SetParameterName(), and G4UIcommand::SetRange().

◆ ~G4OpticalParametersMessenger()

G4OpticalParametersMessenger::~G4OpticalParametersMessenger ( )
virtual

Definition at line 281 of file G4OpticalParametersMessenger.cc.

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;
300 delete fScintTrackInfoCmd;
305 delete fWLSTimeProfileCmd;
306 delete fWLSVerboseLevelCmd;
307 delete fWLS2TimeProfileCmd;
311 delete fMieVerboseLevelCmd;
314}

References fAbsDir, fAbsorptionVerboseLevelCmd, fActivateProcessCmd, fBoundaryDir, fBoundaryInvokeSDCmd, fBoundaryVerboseLevelCmd, fCerenkovDir, fCerenkovMaxBetaChangeCmd, fCerenkovMaxPhotonsCmd, fCerenkovStackPhotonsCmd, fCerenkovTrackSecondariesFirstCmd, fCerenkovVerboseLevelCmd, fDir, fDumpCmd, fMieDir, fMieVerboseLevelCmd, fRaylDir, fRayleighVerboseLevelCmd, fScintByParticleTypeCmd, fScintDir, fScintFiniteRiseTimeCmd, fScintStackPhotonsCmd, fScintTrackInfoCmd, fScintTrackSecondariesFirstCmd, fScintVerboseLevelCmd, fVerboseCmd, fWLS2TimeProfileCmd, fWLS2VerboseLevelCmd, fWlsDir, fWLSTimeProfileCmd, and fWLSVerboseLevelCmd.

◆ G4OpticalParametersMessenger() [2/3]

G4OpticalParametersMessenger::G4OpticalParametersMessenger ( )
privatedelete

◆ G4OpticalParametersMessenger() [3/3]

G4OpticalParametersMessenger::G4OpticalParametersMessenger ( const G4OpticalParametersMessenger right)
privatedelete

Member Function Documentation

◆ AddUIcommand()

void G4UImessenger::AddUIcommand ( G4UIcommand newCommand)
protectedinherited

Definition at line 149 of file G4UImessenger.cc.

150{
151 G4cerr << "Warning : Old style definition of G4UIcommand <"
152 << newCommand->GetCommandPath() << ">." << G4endl;
153}
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
const G4String & GetCommandPath() const
Definition: G4UIcommand.hh:136

References G4cerr, G4endl, and G4UIcommand::GetCommandPath().

◆ BtoS()

G4String G4UImessenger::BtoS ( G4bool  b)
protectedinherited

Definition at line 98 of file G4UImessenger.cc.

99{
100 G4String vl = "0";
101 if(b)
102 vl = "true";
103 return vl;
104}

◆ CommandsShouldBeInMaster()

G4bool G4UImessenger::CommandsShouldBeInMaster ( ) const
inlineinherited

Definition at line 77 of file G4UImessenger.hh.

78 {
80 }
G4bool commandsShouldBeInMaster

References G4UImessenger::commandsShouldBeInMaster.

Referenced by G4UIcommand::G4UIcommandCommonConstructorCode().

◆ CreateCommand()

template<typename T >
T * G4UImessenger::CreateCommand ( const G4String cname,
const G4String dsc 
)
protectedinherited

Definition at line 110 of file G4UImessenger.hh.

111{
112 G4String path;
113 if(cname[0] != '/')
114 {
115 path = baseDirName + cname;
116 if(path[0] != '/')
117 path = "/" + path;
118 }
119
120 T* command = new T(path.c_str(), this);
121 command->SetGuidance(dsc.c_str());
122
123 return command;
124}
G4String baseDirName

References G4UImessenger::baseDirName.

◆ CreateDirectory()

void G4UImessenger::CreateDirectory ( const G4String path,
const G4String dsc,
G4bool  commandsToBeBroadcasted = true 
)
protectedinherited

Definition at line 156 of file G4UImessenger.cc.

158{
160
161 G4String fullpath = path;
162 if(fullpath.back() != '/')
163 fullpath.append("/");
164
165 G4UIcommandTree* tree = ui->GetTree()->FindCommandTree(fullpath.c_str());
166 if(tree != nullptr)
167 {
168 baseDirName = tree->GetPathName();
169 }
170 else
171 {
172 baseDir = new G4UIdirectory(fullpath.c_str(), commandsToBeBroadcasted);
173 baseDirName = fullpath;
174 baseDir->SetGuidance(dsc.c_str());
175 }
176}
const G4String & GetPathName() const
G4UIcommandTree * FindCommandTree(const char *commandPath)
G4UIcommandTree * GetTree() const
Definition: G4UImanager.hh:186
static G4UImanager * GetUIpointer()
Definition: G4UImanager.cc:77
G4UIdirectory * baseDir

References G4UImessenger::baseDir, G4UImessenger::baseDirName, G4UIcommandTree::FindCommandTree(), G4UIcommandTree::GetPathName(), G4UImanager::GetTree(), G4UImanager::GetUIpointer(), and G4UIcommand::SetGuidance().

Referenced by G4MoleculeShootMessenger::G4MoleculeShootMessenger(), and G4UImessenger::G4UImessenger().

◆ DtoS()

G4String G4UImessenger::DtoS ( G4double  a)
protectedinherited

Definition at line 90 of file G4UImessenger.cc.

91{
92 std::ostringstream os;
93 os << a;
94 return G4String(os.str());
95}

Referenced by G4ScoreQuantityMessenger::FilterCommands(), and G4UIcontrolMessenger::SetNewValue().

◆ GetCurrentValue()

G4String G4UImessenger::GetCurrentValue ( G4UIcommand command)
virtualinherited

Reimplemented in G4ScoreQuantityMessenger, G4VisCommandModelCreate< Factory >, G4VisCommandListManagerList< Manager >, G4VisCommandListManagerSelect< Manager >, G4VisCommandManagerMode< Manager >, G4ToolsAnalysisMessenger, G4ScoringMessenger, G4EvManMessenger, G4GeneralParticleSourceMessenger, G4ParticleGunMessenger, G4GeometryMessenger, G4GenericMessenger, G4UIcontrolMessenger, GFlashShowerModelMessenger, G4DecayTableMessenger, G4ParticleMessenger, G4ParticlePropertyMessenger, G4tgrMessenger, G4PersistencyCenterMessenger, G4ProductionCutsTableMessenger, G4SchedulerMessenger, G4VITSteppingVerbose, G4MoleculeShootMessenger, G4MoleculeGunMessenger, G4ProcessManagerMessenger, G4ProcessTableMessenger, G4MatScanMessenger, G4RunMessenger, G4UserPhysicsListMessenger, G4TrackingMessenger, G4GMocrenMessenger, G4HepRepMessenger, G4VisCommandAbortReviewKeptEvents, G4VisCommandDrawOnlyToBeKeptEvents, G4VisCommandEnable, G4VisCommandList, G4VisCommandReviewKeptEvents, G4VisCommandVerbose, G4VisCommandGeometryList, G4VisCommandGeometryRestore, G4VisCommandGeometrySetColour, G4VisCommandGeometrySetDaughtersInvisible, G4VisCommandGeometrySetForceAuxEdgeVisible, G4VisCommandGeometrySetForceCloud, G4VisCommandGeometrySetForceSolid, G4VisCommandGeometrySetForceLineSegmentsPerCircle, G4VisCommandGeometrySetForceWireframe, G4VisCommandGeometrySetLineStyle, G4VisCommandGeometrySetLineWidth, G4VisCommandGeometrySetVisibility, G4VisCommandSceneActivateModel, G4VisCommandSceneCreate, G4VisCommandSceneEndOfEventAction, G4VisCommandSceneEndOfRunAction, G4VisCommandSceneList, G4VisCommandSceneNotifyHandlers, G4VisCommandSceneRemoveModel, G4VisCommandSceneSelect, G4VisCommandSceneShowExtents, G4VisCommandSceneAddArrow, G4VisCommandSceneAddArrow2D, G4VisCommandSceneAddAxes, G4VisCommandSceneAddDate, G4VisCommandSceneAddDigis, G4VisCommandSceneAddEventID, G4VisCommandSceneAddExtent, G4VisCommandSceneAddElectricField, G4VisCommandSceneAddFrame, G4VisCommandSceneAddGPS, G4VisCommandSceneAddGhosts, G4VisCommandSceneAddHits, G4VisCommandSceneAddLine, G4VisCommandSceneAddLine2D, G4VisCommandSceneAddLocalAxes, G4VisCommandSceneAddLogicalVolume, G4VisCommandSceneAddLogo, G4VisCommandSceneAddLogo2D, G4VisCommandSceneAddMagneticField, G4VisCommandSceneAddPSHits, G4VisCommandSceneAddScale, G4VisCommandSceneAddText, G4VisCommandSceneAddText2D, G4VisCommandSceneAddTrajectories, G4VisCommandSceneAddUserAction, G4VisCommandSceneAddVolume, G4VisCommandSceneAddPlotter, G4VisCommandSceneHandlerAttach, G4VisCommandSceneHandlerCreate, G4VisCommandSceneHandlerList, G4VisCommandSceneHandlerSelect, G4VisCommandSetArrow3DLineSegmentsPerCircle, G4VisCommandSetColour, G4VisCommandSetExtentForField, G4VisCommandSetLineWidth, G4VisCommandSetTextColour, G4VisCommandSetTextLayout, G4VisCommandSetTextSize, G4VisCommandSetTouchable, G4VisCommandSetVolumeForField, G4VisCommandsTouchable, G4VisCommandsTouchableSet, G4VisCommandViewerAddCutawayPlane, G4VisCommandViewerCentreOn, G4VisCommandViewerChangeCutawayPlane, G4VisCommandViewerClear, G4VisCommandViewerClearCutawayPlanes, G4VisCommandViewerClearTransients, G4VisCommandViewerClearVisAttributesModifiers, G4VisCommandViewerClone, G4VisCommandViewerColourByDensity, G4VisCommandViewerCopyViewFrom, G4VisCommandViewerCreate, G4VisCommandViewerDolly, G4VisCommandViewerFlush, G4VisCommandViewerInterpolate, G4VisCommandViewerList, G4VisCommandViewerPan, G4VisCommandViewerReset, G4VisCommandViewerRefresh, G4VisCommandViewerRebuild, G4VisCommandViewerSave, G4VisCommandViewerScale, G4VisCommandViewerSelect, G4VisCommandViewerUpdate, G4VisCommandViewerZoom, G4VisCommandViewerDefaultHiddenEdge, G4VisCommandViewerDefaultStyle, G4VisCommandsViewerSet, G4VModelCommand< T >, G4VModelCommand< M >, G4RTMessenger, G4ASCIITreeMessenger, G4VtkMessenger, G4PolarizationMessenger, and G4DNAChemistryManager.

Definition at line 58 of file G4UImessenger.cc.

59{
60 G4String nullString;
61 return nullString;
62}

Referenced by G4UIcommand::DoIt(), and G4UIcommand::GetCurrentValue().

◆ ItoS()

G4String G4UImessenger::ItoS ( G4int  i)
protectedinherited

Definition at line 82 of file G4UImessenger.cc.

83{
84 std::ostringstream os;
85 os << i;
86 return G4String(os.str());
87}

Referenced by G4GenericMessenger::DeclareMethod(), and G4ParticleGunMessenger::GetCurrentValue().

◆ operator!=()

G4bool G4UImessenger::operator!= ( const G4UImessenger messenger) const
inherited

Definition at line 76 of file G4UImessenger.cc.

77{
78 return this != &messenger;
79}

◆ operator=()

G4OpticalParametersMessenger & G4OpticalParametersMessenger::operator= ( const G4OpticalParametersMessenger right)
privatedelete

◆ operator==()

G4bool G4UImessenger::operator== ( const G4UImessenger messenger) const
inherited

Definition at line 70 of file G4UImessenger.cc.

71{
72 return this == &messenger;
73}

◆ SetNewValue()

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

Apply command to the associated object.

Reimplemented from G4UImessenger.

Definition at line 316 of file G4OpticalParametersMessenger.cc.

318{
319 // physics needs to be rebuilt for all commands
320 G4bool physicsModified = true;
321
323 if(command == fActivateProcessCmd)
324 {
325 std::istringstream is(newValue.data());
326 G4String pn;
327 G4String flag;
328 is >> pn >> flag;
331 }
332 else if(command == fVerboseCmd)
333 {
335 }
336 else if(command == fDumpCmd)
337 {
338 params->Dump();
339 }
340 else if(command == fCerenkovMaxPhotonsCmd)
341 {
344 G4cout << "Cerenkov max photons: " << params->GetCerenkovMaxPhotonsPerStep()
345 << G4endl;
346 }
347 else if(command == fCerenkovMaxBetaChangeCmd)
348 {
351 }
352 else if(command == fCerenkovStackPhotonsCmd)
353 {
356 }
357 else if(command == fCerenkovTrackSecondariesFirstCmd)
358 {
361 }
362 else if(command == fCerenkovVerboseLevelCmd)
363 {
366 }
367 else if(command == fScintByParticleTypeCmd)
368 {
371 }
372 else if(command == fScintTrackInfoCmd)
373 {
375 }
376 else if(command == fScintFiniteRiseTimeCmd)
377 {
380 }
381 else if(command == fScintStackPhotonsCmd)
382 {
385 }
386 else if(command == fScintTrackSecondariesFirstCmd)
387 {
390 }
391 else if(command == fScintVerboseLevelCmd)
392 {
395 }
396 else if(command == fWLSTimeProfileCmd)
397 {
398 params->SetWLSTimeProfile(newValue);
399 }
400 else if(command == fWLSVerboseLevelCmd)
401 {
403 }
404 else if(command == fWLS2TimeProfileCmd)
405 {
406 params->SetWLS2TimeProfile(newValue);
407 }
408 else if(command == fWLS2VerboseLevelCmd)
409 {
411 }
412 else if(command == fAbsorptionVerboseLevelCmd)
413 {
416 }
417 else if(command == fRayleighVerboseLevelCmd)
418 {
421 }
422 else if(command == fMieVerboseLevelCmd)
423 {
425 }
426 else if(command == fBoundaryVerboseLevelCmd)
427 {
430 }
431 else if(command == fBoundaryInvokeSDCmd)
432 {
435 }
436 if(physicsModified)
437 {
438 G4UImanager::GetUIpointer()->ApplyCommand("/run/physicsModified");
439 }
440}
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 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)
static G4double GetNewDoubleValue(const char *paramString)
static G4int GetNewIntValue(const char *paramString)
static G4bool ConvertToBool(const char *st)
Definition: G4UIcommand.cc:545
G4int ApplyCommand(const char *aCommand)
Definition: G4UImanager.cc:485

References G4UImanager::ApplyCommand(), G4UIcommand::ConvertToBool(), G4OpticalParameters::Dump(), fAbsorptionVerboseLevelCmd, fActivateProcessCmd, fBoundaryInvokeSDCmd, fBoundaryVerboseLevelCmd, fCerenkovMaxBetaChangeCmd, fCerenkovMaxPhotonsCmd, fCerenkovStackPhotonsCmd, fCerenkovTrackSecondariesFirstCmd, fCerenkovVerboseLevelCmd, fDumpCmd, fMieVerboseLevelCmd, fRayleighVerboseLevelCmd, fScintByParticleTypeCmd, fScintFiniteRiseTimeCmd, fScintStackPhotonsCmd, fScintTrackInfoCmd, fScintTrackSecondariesFirstCmd, fScintVerboseLevelCmd, fVerboseCmd, fWLS2TimeProfileCmd, fWLS2VerboseLevelCmd, fWLSTimeProfileCmd, fWLSVerboseLevelCmd, G4cout, G4endl, G4OpticalParameters::GetCerenkovMaxPhotonsPerStep(), G4UIcmdWithABool::GetNewBoolValue(), G4UIcmdWithADouble::GetNewDoubleValue(), G4UIcmdWithAnInteger::GetNewIntValue(), G4UImanager::GetUIpointer(), params, G4InuclParticleNames::pn, G4OpticalParameters::SetAbsorptionVerboseLevel(), G4OpticalParameters::SetBoundaryInvokeSD(), G4OpticalParameters::SetBoundaryVerboseLevel(), G4OpticalParameters::SetCerenkovMaxBetaChange(), G4OpticalParameters::SetCerenkovMaxPhotonsPerStep(), G4OpticalParameters::SetCerenkovStackPhotons(), G4OpticalParameters::SetCerenkovTrackSecondariesFirst(), G4OpticalParameters::SetCerenkovVerboseLevel(), G4OpticalParameters::SetMieVerboseLevel(), G4OpticalParameters::SetProcessActivation(), G4OpticalParameters::SetRayleighVerboseLevel(), G4OpticalParameters::SetScintByParticleType(), G4OpticalParameters::SetScintFiniteRiseTime(), G4OpticalParameters::SetScintStackPhotons(), G4OpticalParameters::SetScintTrackInfo(), G4OpticalParameters::SetScintTrackSecondariesFirst(), G4OpticalParameters::SetScintVerboseLevel(), G4OpticalParameters::SetVerboseLevel(), G4OpticalParameters::SetWLS2TimeProfile(), G4OpticalParameters::SetWLS2VerboseLevel(), G4OpticalParameters::SetWLSTimeProfile(), and G4OpticalParameters::SetWLSVerboseLevel().

◆ StoB()

G4bool G4UImessenger::StoB ( G4String  s)
protectedinherited

Definition at line 137 of file G4UImessenger.cc.

138{
140 G4bool vl = false;
141 if(v == "Y" || v == "YES" || v == "1" || v == "T" || v == "TRUE")
142 {
143 vl = true;
144 }
145 return vl;
146}
G4String to_upper_copy(G4String str)
Return uppercase copy of string.

References G4StrUtil::to_upper_copy().

Referenced by G4LocalThreadCoutMessenger::SetNewValue(), G4CascadeParamMessenger::SetNewValue(), G4ScoreQuantityMessenger::SetNewValue(), and G4ScoringMessenger::SetNewValue().

◆ StoD()

G4double G4UImessenger::StoD ( G4String  s)
protectedinherited

◆ StoI()

G4int G4UImessenger::StoI ( G4String  s)
protectedinherited

◆ StoL()

G4long G4UImessenger::StoL ( G4String  s)
protectedinherited

Definition at line 117 of file G4UImessenger.cc.

118{
119 G4long vl;
120 const char* t = str;
121 std::istringstream is(t);
122 is >> vl;
123 return vl;
124}
long G4long
Definition: G4Types.hh:87

Referenced by G4RunMessenger::SetNewValue().

Field Documentation

◆ baseDir

G4UIdirectory* G4UImessenger::baseDir = nullptr
protectedinherited

◆ baseDirName

G4String G4UImessenger::baseDirName = ""
protectedinherited

◆ commandsShouldBeInMaster

G4bool G4UImessenger::commandsShouldBeInMaster = false
protectedinherited

◆ fAbsDir

G4UIdirectory* G4OpticalParametersMessenger::fAbsDir
private

◆ fAbsorptionVerboseLevelCmd

G4UIcmdWithAnInteger* G4OpticalParametersMessenger::fAbsorptionVerboseLevelCmd
private

◆ fActivateProcessCmd

G4UIcommand* G4OpticalParametersMessenger::fActivateProcessCmd
private

selectOpProcess command

Definition at line 94 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger(), SetNewValue(), and ~G4OpticalParametersMessenger().

◆ fBoundaryDir

G4UIdirectory* G4OpticalParametersMessenger::fBoundaryDir
private

◆ fBoundaryInvokeSDCmd

G4UIcmdWithABool* G4OpticalParametersMessenger::fBoundaryInvokeSDCmd
private

setInvokeSD command

Definition at line 144 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger(), SetNewValue(), and ~G4OpticalParametersMessenger().

◆ fBoundaryVerboseLevelCmd

G4UIcmdWithAnInteger* G4OpticalParametersMessenger::fBoundaryVerboseLevelCmd
private

◆ fCerenkovDir

G4UIdirectory* G4OpticalParametersMessenger::fCerenkovDir
private

◆ fCerenkovMaxBetaChangeCmd

G4UIcmdWithADouble* G4OpticalParametersMessenger::fCerenkovMaxBetaChangeCmd
private

setCerenkovMaxBetaChange command

Definition at line 105 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger(), SetNewValue(), and ~G4OpticalParametersMessenger().

◆ fCerenkovMaxPhotonsCmd

G4UIcmdWithAnInteger* G4OpticalParametersMessenger::fCerenkovMaxPhotonsCmd
private

◆ fCerenkovStackPhotonsCmd

G4UIcmdWithABool* G4OpticalParametersMessenger::fCerenkovStackPhotonsCmd
private

setStackPhotons command

Definition at line 108 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger(), SetNewValue(), and ~G4OpticalParametersMessenger().

◆ fCerenkovTrackSecondariesFirstCmd

G4UIcmdWithABool* G4OpticalParametersMessenger::fCerenkovTrackSecondariesFirstCmd
private

◆ fCerenkovVerboseLevelCmd

G4UIcmdWithAnInteger* G4OpticalParametersMessenger::fCerenkovVerboseLevelCmd
private

◆ fDir

G4UIdirectory* G4OpticalParametersMessenger::fDir
private

command directory

Definition at line 83 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger(), and ~G4OpticalParametersMessenger().

◆ fDumpCmd

G4UIcommand* G4OpticalParametersMessenger::fDumpCmd
private

◆ fMieDir

G4UIdirectory* G4OpticalParametersMessenger::fMieDir
private

◆ fMieVerboseLevelCmd

G4UIcmdWithAnInteger* G4OpticalParametersMessenger::fMieVerboseLevelCmd
private

◆ fRaylDir

G4UIdirectory* G4OpticalParametersMessenger::fRaylDir
private

◆ fRayleighVerboseLevelCmd

G4UIcmdWithAnInteger* G4OpticalParametersMessenger::fRayleighVerboseLevelCmd
private

◆ fScintByParticleTypeCmd

G4UIcmdWithABool* G4OpticalParametersMessenger::fScintByParticleTypeCmd
private

setScintillationByParticleType command

Definition at line 116 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger(), SetNewValue(), and ~G4OpticalParametersMessenger().

◆ fScintDir

G4UIdirectory* G4OpticalParametersMessenger::fScintDir
private

◆ fScintFiniteRiseTimeCmd

G4UIcmdWithABool* G4OpticalParametersMessenger::fScintFiniteRiseTimeCmd
private

setFiniteRiseTime command

Definition at line 127 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger(), SetNewValue(), and ~G4OpticalParametersMessenger().

◆ fScintStackPhotonsCmd

G4UIcmdWithABool* G4OpticalParametersMessenger::fScintStackPhotonsCmd
private

setStackPhotons command

Definition at line 122 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger(), SetNewValue(), and ~G4OpticalParametersMessenger().

◆ fScintTrackInfoCmd

G4UIcmdWithABool* G4OpticalParametersMessenger::fScintTrackInfoCmd
private

setScintillationTrackInfo command

Definition at line 119 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger(), SetNewValue(), and ~G4OpticalParametersMessenger().

◆ fScintTrackSecondariesFirstCmd

G4UIcmdWithABool* G4OpticalParametersMessenger::fScintTrackSecondariesFirstCmd
private

◆ fScintVerboseLevelCmd

G4UIcmdWithAnInteger* G4OpticalParametersMessenger::fScintVerboseLevelCmd
private

◆ fVerboseCmd

G4UIcmdWithAnInteger* G4OpticalParametersMessenger::fVerboseCmd
private

setProcessVerbose command

Definition at line 97 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger(), SetNewValue(), and ~G4OpticalParametersMessenger().

◆ fWls2Dir

G4UIdirectory* G4OpticalParametersMessenger::fWls2Dir
private

Definition at line 87 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger().

◆ fWLS2TimeProfileCmd

G4UIcmdWithAString* G4OpticalParametersMessenger::fWLS2TimeProfileCmd
private

setWLS2TimeProfile command

Definition at line 140 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger(), SetNewValue(), and ~G4OpticalParametersMessenger().

◆ fWLS2VerboseLevelCmd

G4UIcmdWithAnInteger* G4OpticalParametersMessenger::fWLS2VerboseLevelCmd
private

◆ fWlsDir

G4UIdirectory* G4OpticalParametersMessenger::fWlsDir
private

◆ fWLSTimeProfileCmd

G4UIcmdWithAString* G4OpticalParametersMessenger::fWLSTimeProfileCmd
private

setWLSTimeProfile command

Definition at line 134 of file G4OpticalParametersMessenger.hh.

Referenced by G4OpticalParametersMessenger(), SetNewValue(), and ~G4OpticalParametersMessenger().

◆ fWLSVerboseLevelCmd

G4UIcmdWithAnInteger* G4OpticalParametersMessenger::fWLSVerboseLevelCmd
private

◆ params

G4OpticalParameters* G4OpticalParametersMessenger::params
private

associated class

Definition at line 80 of file G4OpticalParametersMessenger.hh.

Referenced by SetNewValue().


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