Geant4  10.02.p03
G4OpticalPhysicsMessenger Class Reference

#include <G4OpticalPhysicsMessenger.hh>

Inheritance diagram for G4OpticalPhysicsMessenger:
Collaboration diagram for G4OpticalPhysicsMessenger:

Public Member Functions

 G4OpticalPhysicsMessenger (G4OpticalPhysics *)
 
virtual ~G4OpticalPhysicsMessenger ()
 
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)
 
G4bool operator== (const G4UImessenger &messenger) const
 
G4bool CommandsShouldBeInMaster () const
 

Private Member Functions

 G4OpticalPhysicsMessenger ()
 Not implemented. More...
 
 G4OpticalPhysicsMessenger (const G4OpticalPhysicsMessenger &right)
 Not implemented. More...
 
G4OpticalPhysicsMessengeroperator= (const G4OpticalPhysicsMessenger &right)
 Not implemented. More...
 

Private Attributes

G4OpticalPhysicsfOpticalPhysics
 associated class More...
 
G4UIdirectoryfDir
 command directory More...
 
G4UIdirectoryfDir2
 
G4OpticalProcessIndex fSelectedProcessIndex
 selected optical process More...
 
G4UIcommandfActivateProcessCmd
 selectOpProcess command More...
 
G4UIcmdWithAnIntegerfSetOpProcessVerboseCmd
 setProcessVerbose command More...
 
G4UIcmdWithAnIntegerfSetCerenkovMaxPhotonsCmd
 setCerenkovMaxPhotons command More...
 
G4UIcmdWithADoublefSetCerenkovMaxBetaChangeCmd
 setCerenkovMaxBetaChange command More...
 
G4UIcmdWithADoublefSetScintillationYieldFactorCmd
 setScintillationYieldFactor command More...
 
G4UIcmdWithABoolfSetScintillationByParticleTypeCmd
 setScintillationByParticleType command More...
 
G4UIcmdWithAStringfSetWLSTimeProfileCmd
 setWLSTimeProfile command More...
 
G4UIcommandfSetTrackSecondariesFirstCmd
 setTrackSecondariesFirst command More...
 
G4UIcmdWithABoolfSetFiniteRiseTimeCmd
 setFiniteRiseTime command More...
 

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)
 
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
 
G4String baseDirName
 
G4bool commandsShouldBeInMaster
 

Detailed Description

Definition at line 73 of file G4OpticalPhysicsMessenger.hh.

Constructor & Destructor Documentation

◆ G4OpticalPhysicsMessenger() [1/3]

G4OpticalPhysicsMessenger::G4OpticalPhysicsMessenger ( G4OpticalPhysics opticalPhysics)

Definition at line 60 of file G4OpticalPhysicsMessenger.cc.

62  : G4UImessenger(),
63  fOpticalPhysics(opticalPhysics),
74 {
75  G4bool toBeBroadcasted = false;
76  fDir = new G4UIdirectory("/process/optical/defaults/",toBeBroadcasted);
77  fDir->SetGuidance("Commands related to the optical physics simulation engine.");
78  fDir2 = new G4UIdirectory("/process/optical/",toBeBroadcasted);
79  fDir2->SetGuidance("Commands related to the optical physics simulation engine.");
80 
81  fActivateProcessCmd= new G4UIcommand("/process/optical/processActivation", this);
82  fActivateProcessCmd->SetGuidance("Activate/deactivate the specified optical process");
83  G4UIparameter* par = new G4UIparameter("proc_name",'s',false);
84  G4String candidates;
85  for ( G4int i=0; i<kNoProcess; i++ ) {
86  candidates += G4OpticalProcessName(i);
87  candidates += G4String(" ");
88  }
89  par->SetParameterCandidates(candidates);
90  par->SetGuidance("the process name");
92  par = new G4UIparameter("flag",'b',true);
93  par->SetDefaultValue(true);
94  par->SetGuidance("activation flag");
97 
98 
99  fSetOpProcessVerboseCmd = new G4UIcmdWithAnInteger("/process/optical/verbose", this);
100  fSetOpProcessVerboseCmd->SetGuidance("Set default verbosity level for optical processes");
105 
106 
107  fSetTrackSecondariesFirstCmd = new G4UIcommand("/process/optical/setTrackSecondariesFirst", this);
108  fSetTrackSecondariesFirstCmd->SetGuidance("Activate/deactivate tracking of secondaries before finishing their parent track");
109  par = new G4UIparameter("proc_name",'s',false);
110  par->SetParameterCandidates(candidates);
112  par = new G4UIparameter("flag",'b',false);
113  par->SetDefaultValue(true);
116 
117  //This are repetition of process specific ui commands needed by ATLICE and VMC
118  //(UI messages needed to exist in PreInit state before actual processes are instantiated)
119  fSetCerenkovMaxPhotonsCmd = new G4UIcmdWithAnInteger("/process/optical/defaults/cerenkov/setMaxPhotons", this);
120  fSetCerenkovMaxPhotonsCmd->SetGuidance("Set default maximum number of photons per step");
121  fSetCerenkovMaxPhotonsCmd->SetGuidance("Note this command is used to set the default value,");
122  fSetCerenkovMaxPhotonsCmd->SetGuidance("if process is not active command will not have effect.");
123  fSetCerenkovMaxPhotonsCmd->SetParameterName("CerenkovMaxPhotons", false);
124  fSetCerenkovMaxPhotonsCmd->SetRange("CerenkovMaxPhotons>=0");
126 
127  fSetCerenkovMaxBetaChangeCmd = new G4UIcmdWithADouble("/process/optical/defaults/cerenkov/setMaxBetaChange", this);
128  fSetCerenkovMaxBetaChangeCmd->SetGuidance("Set default maximum change of beta of parent particle per step");
129  fSetCerenkovMaxBetaChangeCmd->SetGuidance("Note this command is used to set the default value,");
130  fSetCerenkovMaxBetaChangeCmd->SetGuidance("if process is not active command will not have effect.");
131  fSetCerenkovMaxBetaChangeCmd->SetParameterName("CerenkovMaxBetaChange", false);
132  fSetCerenkovMaxBetaChangeCmd->SetRange("CerenkovMaxBetaChange>=0");
134 
135  fSetScintillationYieldFactorCmd = new G4UIcmdWithADouble("/process/optical/defaults/scintillation/setYieldFactor", this);
136  fSetScintillationYieldFactorCmd->SetGuidance("Set scintillation yield factor");
137  fSetScintillationYieldFactorCmd->SetGuidance("Note this command is used to set the default value,");
138  fSetScintillationYieldFactorCmd->SetGuidance("if process is not active command will not have effect.");
139  fSetScintillationYieldFactorCmd->SetParameterName("ScintillationYieldFactor", false);
140  fSetScintillationYieldFactorCmd->SetRange("ScintillationYieldFactor>=0");
142 
143  fSetScintillationByParticleTypeCmd = new G4UIcmdWithABool("/process/optical/defaults/scintillation/setByParticleType", this);
144  fSetScintillationByParticleTypeCmd->SetGuidance("Activate/Inactivate scintillation process by particle type");
145  fSetScintillationByParticleTypeCmd->SetGuidance("Note this command is used to set the default value,");
146  fSetScintillationByParticleTypeCmd->SetGuidance("if process is not active command will not have effect.");
147  fSetScintillationByParticleTypeCmd->SetParameterName("ScintillationByParticleTypeActivation", false);
149 
150  fSetFiniteRiseTimeCmd = new G4UIcmdWithABool("/process/optical/defaults/scintillation/setFiniteRiseTime", this);
151  fSetFiniteRiseTimeCmd->SetGuidance("Set option of a finite rise-time for G4Scintillation");
152  fSetFiniteRiseTimeCmd->SetGuidance("If set, the G4Scintillation process expects the user to have set the constant material property FAST/SLOWSCINTILLATIONRISETIME");
153  fSetFiniteRiseTimeCmd->SetGuidance("Note this command is used to set the default value,");
154  fSetFiniteRiseTimeCmd->SetGuidance("if process is not active command will not have effect.");
155  fSetFiniteRiseTimeCmd->SetParameterName("FiniteRiseTime", false);
157 
158  fSetWLSTimeProfileCmd = new G4UIcmdWithAString("/process/optical/defaults/wls/setTimeProfile", this);
159  fSetWLSTimeProfileCmd->SetGuidance("Set the WLS time profile (delta or exponential)");
160  fSetWLSTimeProfileCmd->SetParameterName("WLSTimeProfile", false);
161  fSetWLSTimeProfileCmd->SetCandidates("delta exponential");
163 }
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
G4OpticalPhysics * fOpticalPhysics
associated class
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:152
G4UIcmdWithABool * fSetFiniteRiseTimeCmd
setFiniteRiseTime command
Number of processes, no selected process.
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
G4OpticalProcessIndex fSelectedProcessIndex
selected optical process
void SetParameterCandidates(const char *theString)
void SetDefaultValue(const char *theDefaultValue)
G4UIcmdWithAString * fSetWLSTimeProfileCmd
setWLSTimeProfile command
G4UIcommand * fActivateProcessCmd
selectOpProcess command
G4UIcmdWithAnInteger * fSetOpProcessVerboseCmd
setProcessVerbose command
int G4int
Definition: G4Types.hh:78
G4UIcmdWithAnInteger * fSetCerenkovMaxPhotonsCmd
setCerenkovMaxPhotons command
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
bool G4bool
Definition: G4Types.hh:79
void SetRange(const char *rs)
Definition: G4UIcommand.hh:125
void SetGuidance(const char *aGuidance)
Definition: G4UIcommand.hh:161
void AvailableForStates(G4ApplicationState s1)
Definition: G4UIcommand.cc:239
G4UIcmdWithABool * fSetScintillationByParticleTypeCmd
setScintillationByParticleType command
G4UIdirectory * fDir
command directory
G4String G4OpticalProcessName(G4int)
Return the name for a given optical process index.
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
G4UIcmdWithADouble * fSetScintillationYieldFactorCmd
setScintillationYieldFactor command
void SetCandidates(const char *candidateList)
G4UIcmdWithADouble * fSetCerenkovMaxBetaChangeCmd
setCerenkovMaxBetaChange command
void SetDefaultValue(G4int defVal)
G4UIcommand * fSetTrackSecondariesFirstCmd
setTrackSecondariesFirst command
void SetGuidance(const char *theGuidance)
Here is the call graph for this function:

◆ ~G4OpticalPhysicsMessenger()

G4OpticalPhysicsMessenger::~G4OpticalPhysicsMessenger ( )
virtual

Definition at line 165 of file G4OpticalPhysicsMessenger.cc.

166 {
167 // Destructor
168 
169  delete fDir;
170  delete fDir2;
171  delete fActivateProcessCmd;
177  delete fSetWLSTimeProfileCmd;
178  delete fSetFiniteRiseTimeCmd;
179 }
G4UIcmdWithABool * fSetFiniteRiseTimeCmd
setFiniteRiseTime command
G4UIcmdWithAString * fSetWLSTimeProfileCmd
setWLSTimeProfile command
G4UIcommand * fActivateProcessCmd
selectOpProcess command
G4UIcmdWithAnInteger * fSetOpProcessVerboseCmd
setProcessVerbose command
G4UIcmdWithAnInteger * fSetCerenkovMaxPhotonsCmd
setCerenkovMaxPhotons command
G4UIcmdWithABool * fSetScintillationByParticleTypeCmd
setScintillationByParticleType command
G4UIdirectory * fDir
command directory
G4UIcmdWithADouble * fSetScintillationYieldFactorCmd
setScintillationYieldFactor command
G4UIcmdWithADouble * fSetCerenkovMaxBetaChangeCmd
setCerenkovMaxBetaChange command

◆ G4OpticalPhysicsMessenger() [2/3]

G4OpticalPhysicsMessenger::G4OpticalPhysicsMessenger ( )
private

Not implemented.

◆ G4OpticalPhysicsMessenger() [3/3]

G4OpticalPhysicsMessenger::G4OpticalPhysicsMessenger ( const G4OpticalPhysicsMessenger right)
private

Not implemented.

Member Function Documentation

◆ operator=()

G4OpticalPhysicsMessenger& G4OpticalPhysicsMessenger::operator= ( const G4OpticalPhysicsMessenger right)
private

Not implemented.

◆ SetNewValue()

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

Apply command to the associated object.

Reimplemented from G4UImessenger.

Definition at line 182 of file G4OpticalPhysicsMessenger.cc.

184 {
186  if (command == fActivateProcessCmd) {
187  std::istringstream is(newValue.data());
188  G4String pn;
189  G4int flag;
190  is >> pn >> flag;
191  if ( pn == "Cerenkov" ) {
193  } else if ( pn == "Scintillation" ) {
195  } else if ( pn == "OpAbsorption" ) {
197  } else if ( pn == "OpRayleigh" ) {
199  } else if ( pn == "OpMieHG" ) {
201  } else if ( pn == "OpBoundary" ) {
203  } else if ( pn == "OpWLS" ) {
205  } else {
207  msg << "Not allowed process name: "<<pn<<" (UI: "<<newValue<<")";
208  G4Exception("G4OpticalPhysicsMessenger::SetNewValue(...)","Optical001",FatalException,pn);
209  }
211  }
212  else if (command == fSetTrackSecondariesFirstCmd )
213  {
214  std::istringstream is(newValue.data());
215  G4String pn;
216  G4int flag;
217  is >> pn >> flag;
218  if ( pn == "Cerenkov" ) {
220  } else if ( pn == "Scintillation" ) {
222  } else if ( pn == "OpAbsorption" ) {
224  } else if ( pn == "OpRayleigh" ) {
226  } else if ( pn == "OpMieHG" ) {
228  } else if ( pn == "OpBoundary" ) {
230  } else if ( pn == "OpWLS" ) {
232  } else {
234  msg << "Not allowed process name: "<<pn<<" (UI: "<<newValue<<")";
235  G4Exception("G4OpticalPhysicsMessenger::SetNewValue(...)","Optical001",FatalException,pn);
236  }
238  }
239  else if (command == fSetOpProcessVerboseCmd) {
241  }
242  else if (command == fSetCerenkovMaxPhotonsCmd) {
246  }
247  else if (command == fSetCerenkovMaxBetaChangeCmd) {
251  }
252  else if (command == fSetScintillationYieldFactorCmd) {
256  }
257  else if (command == fSetScintillationByParticleTypeCmd) {
261  }
262  else if (command == fSetFiniteRiseTimeCmd) {
266  }
267  else if (command == fSetWLSTimeProfileCmd) {
269  }
270 }
G4OpticalPhysics * fOpticalPhysics
associated class
G4UIcmdWithABool * fSetFiniteRiseTimeCmd
setFiniteRiseTime command
const char * data() const
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4OpticalProcessIndex fSelectedProcessIndex
selected optical process
static G4int GetNewIntValue(const char *paramString)
void SetMaxBetaChangePerStep(G4double)
G4UIcmdWithAString * fSetWLSTimeProfileCmd
setWLSTimeProfile command
Scintillation process index.
G4UIcommand * fActivateProcessCmd
selectOpProcess command
Mie scattering process index.
G4UIcmdWithAnInteger * fSetOpProcessVerboseCmd
setProcessVerbose command
void SetFiniteRiseTime(G4bool)
Absorption process index.
int G4int
Definition: G4Types.hh:78
static G4bool GetNewBoolValue(const char *paramString)
void SetScintillationByParticleType(G4bool)
G4UIcmdWithAnInteger * fSetCerenkovMaxPhotonsCmd
setCerenkovMaxPhotons command
void SetMaxNumPhotonsPerStep(G4int)
void Configure(G4OpticalProcessIndex, G4bool)
static G4double GetNewDoubleValue(const char *paramString)
void SetVerboseLevel(G4int value)
G4UIcmdWithABool * fSetScintillationByParticleTypeCmd
setScintillationByParticleType command
Wave Length Shifting process index.
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
Boundary process index.
void SetScintillationYieldFactor(G4double)
Cerenkov process index.
void SetWLSTimeProfile(G4String)
void SetTrackSecondariesFirst(G4OpticalProcessIndex, G4bool)
G4UIcmdWithADouble * fSetScintillationYieldFactorCmd
setScintillationYieldFactor command
G4UIcmdWithADouble * fSetCerenkovMaxBetaChangeCmd
setCerenkovMaxBetaChange command
G4UIcommand * fSetTrackSecondariesFirstCmd
setTrackSecondariesFirst command
Rayleigh scattering process index.
Here is the call graph for this function:

Member Data Documentation

◆ fActivateProcessCmd

G4UIcommand* G4OpticalPhysicsMessenger::fActivateProcessCmd
private

selectOpProcess command

Definition at line 105 of file G4OpticalPhysicsMessenger.hh.

◆ fDir

G4UIdirectory* G4OpticalPhysicsMessenger::fDir
private

command directory

Definition at line 98 of file G4OpticalPhysicsMessenger.hh.

◆ fDir2

G4UIdirectory* G4OpticalPhysicsMessenger::fDir2
private

Definition at line 99 of file G4OpticalPhysicsMessenger.hh.

◆ fOpticalPhysics

G4OpticalPhysics* G4OpticalPhysicsMessenger::fOpticalPhysics
private

associated class

Definition at line 95 of file G4OpticalPhysicsMessenger.hh.

◆ fSelectedProcessIndex

G4OpticalProcessIndex G4OpticalPhysicsMessenger::fSelectedProcessIndex
private

selected optical process

Definition at line 102 of file G4OpticalPhysicsMessenger.hh.

◆ fSetCerenkovMaxBetaChangeCmd

G4UIcmdWithADouble* G4OpticalPhysicsMessenger::fSetCerenkovMaxBetaChangeCmd
private

setCerenkovMaxBetaChange command

Definition at line 115 of file G4OpticalPhysicsMessenger.hh.

◆ fSetCerenkovMaxPhotonsCmd

G4UIcmdWithAnInteger* G4OpticalPhysicsMessenger::fSetCerenkovMaxPhotonsCmd
private

setCerenkovMaxPhotons command

Definition at line 112 of file G4OpticalPhysicsMessenger.hh.

◆ fSetFiniteRiseTimeCmd

G4UIcmdWithABool* G4OpticalPhysicsMessenger::fSetFiniteRiseTimeCmd
private

setFiniteRiseTime command

Definition at line 133 of file G4OpticalPhysicsMessenger.hh.

◆ fSetOpProcessVerboseCmd

G4UIcmdWithAnInteger* G4OpticalPhysicsMessenger::fSetOpProcessVerboseCmd
private

setProcessVerbose command

Definition at line 109 of file G4OpticalPhysicsMessenger.hh.

◆ fSetScintillationByParticleTypeCmd

G4UIcmdWithABool* G4OpticalPhysicsMessenger::fSetScintillationByParticleTypeCmd
private

setScintillationByParticleType command

Definition at line 121 of file G4OpticalPhysicsMessenger.hh.

◆ fSetScintillationYieldFactorCmd

G4UIcmdWithADouble* G4OpticalPhysicsMessenger::fSetScintillationYieldFactorCmd
private

setScintillationYieldFactor command

Definition at line 118 of file G4OpticalPhysicsMessenger.hh.

◆ fSetTrackSecondariesFirstCmd

G4UIcommand* G4OpticalPhysicsMessenger::fSetTrackSecondariesFirstCmd
private

setTrackSecondariesFirst command

Definition at line 130 of file G4OpticalPhysicsMessenger.hh.

◆ fSetWLSTimeProfileCmd

G4UIcmdWithAString* G4OpticalPhysicsMessenger::fSetWLSTimeProfileCmd
private

setWLSTimeProfile command

Definition at line 127 of file G4OpticalPhysicsMessenger.hh.


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