Geant4  10.02
G4OpticalPhysicsMessenger.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: G4OpticalPhysicsMessenger
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 
40 #include "G4OpticalPhysics.hh"
41 
42 #include "G4UIcommand.hh"
43 #include "G4UIdirectory.hh"
44 
45 #include "G4UIcommand.hh"
46 #include "G4UIdirectory.hh"
47 #include "G4UIcmdWithABool.hh"
48 #include "G4UIcmdWithAString.hh"
49 #include "G4UIcmdWithADouble.hh"
50 #include "G4UIcmdWithAnInteger.hh"
52 #include "G4UIparameter.hh"
53 
54 #include "G4ParticleTable.hh"
55 #include "G4ProcessManager.hh"
56 #include "G4ParticleDefinition.hh"
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59 
61  G4OpticalPhysics* opticalPhysics)
62  : G4UImessenger(),
63  fOpticalPhysics(opticalPhysics),
64  fSelectedProcessIndex(kNoProcess),
65  fActivateProcessCmd(0),
66  fSetOpProcessVerboseCmd(0),
67  fSetCerenkovMaxPhotonsCmd(0),
68  fSetCerenkovMaxBetaChangeCmd(0),
69  fSetScintillationYieldFactorCmd(0),
70  fSetScintillationByParticleTypeCmd(0),
71  fSetWLSTimeProfileCmd(0),
72  fSetTrackSecondariesFirstCmd(0),
73  fSetFiniteRiseTimeCmd(0)
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 }
164 
166 {
167 // Destructor
168 
169  delete fDir;
170  delete fDir2;
171  delete fActivateProcessCmd;
177  delete fSetWLSTimeProfileCmd;
178  delete fSetFiniteRiseTimeCmd;
179 }
180 
181 #include <iostream>
183  G4String newValue)
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 }
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
G4OpticalPhysics * fOpticalPhysics
associated class
virtual void SetNewValue(G4UIcommand *, G4String)
void SetParameter(G4UIparameter *const newParameter)
Definition: G4UIcommand.hh:152
G4UIcmdWithABool * fSetFiniteRiseTimeCmd
setFiniteRiseTime command
Number of processes, no selected process.
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
G4OpticalProcessIndex fSelectedProcessIndex
selected optical process
static G4int GetNewIntValue(const char *paramString)
void SetMaxBetaChangePerStep(G4double)
void SetParameterCandidates(const char *theString)
void SetDefaultValue(const char *theDefaultValue)
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 SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
void SetMaxNumPhotonsPerStep(G4int)
bool G4bool
Definition: G4Types.hh:79
void Configure(G4OpticalProcessIndex, G4bool)
static G4double GetNewDoubleValue(const char *paramString)
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
void SetVerboseLevel(G4int value)
G4UIcmdWithABool * fSetScintillationByParticleTypeCmd
setScintillationByParticleType command
Wave Length Shifting process index.
G4OpticalPhysicsMessenger()
Not implemented.
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4UIdirectory * fDir
command directory
const char * data() const
G4String G4OpticalProcessName(G4int)
Return the name for a given optical process index.
Boundary process index.
void SetScintillationYieldFactor(G4double)
void SetParameterName(const char *theName, G4bool omittable, G4bool currentAsDefault=false)
Cerenkov process index.
void SetWLSTimeProfile(G4String)
void SetTrackSecondariesFirst(G4OpticalProcessIndex, G4bool)
G4UIcmdWithADouble * fSetScintillationYieldFactorCmd
setScintillationYieldFactor command
void SetCandidates(const char *candidateList)
G4UIcmdWithADouble * fSetCerenkovMaxBetaChangeCmd
setCerenkovMaxBetaChange command
void SetDefaultValue(G4int defVal)
G4UIcommand * fSetTrackSecondariesFirstCmd
setTrackSecondariesFirst command
Rayleigh scattering process index.
void SetGuidance(const char *theGuidance)