Geant4  10.02.p01
G4OpticalPhysics.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: G4OpticalPhysics
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 
39 #include "G4OpticalPhysics.hh"
40 
41 #include "G4OpAbsorption.hh"
42 #include "G4OpRayleigh.hh"
43 #include "G4OpMieHG.hh"
44 
45 #include "G4OpBoundaryProcess.hh"
46 
47 #include "G4OpWLS.hh"
48 #include "G4Scintillation.hh"
49 #include "G4Cerenkov.hh"
50 
51 #include "G4LossTableManager.hh"
52 #include "G4EmSaturation.hh"
53 
54 #include "G4GenericMessenger.hh"
55 
56 // factory
58 //
60 
61 
62 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
63 
65  : G4VPhysicsConstructor(name),
66 
67  fMaxNumPhotons(100),
68  fMaxBetaChange(10.0),
69  fYieldFactor(1.),
70  fExcitationRatio(0.0),
71  fProfile("delta"),
72  fFiniteRiseTime(false),
73  fScintillationByParticleType(false)
74 {
75  verboseLevel = verbose;
77 
78  for ( G4int i=0; i<kNoProcess; i++ ) {
79  fProcessUse.push_back(true);
80  fProcessTrackSecondariesFirst.push_back(true);
81  }
82 }
83 
84 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85 
87 {
88  delete fMessenger;
89 }
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92 
94 {
95 // Print all processes activation and their parameters
96 
97  for ( G4int i=0; i<kNoProcess; i++ ) {
98  G4cout << " " << G4OpticalProcessName(i) << " process: ";
99  if ( ! fProcessUse[i] ) {
100  G4cout << "not used" << G4endl;
101  }
102  else {
103  G4cout << "used" << G4endl;
104  if ( i == kCerenkov ) {
105  G4cout << " Max number of photons per step: " << fMaxNumPhotons << G4endl;
106  G4cout << " Max beta change per step: " << fMaxBetaChange << G4endl;
108  G4cout << " Track secondaries first: activated" << G4endl;
109  }
110  else {
111  G4cout << " Track secondaries first: inactivated" << G4endl;
112  }
113  }
114  if ( i == kScintillation ) {
116  G4cout << " Scintillation by Particle Type: activated " << G4endl;
117  G4cout << " Yield factor: " << fYieldFactor << G4endl;
118  G4cout << " ExcitationRatio: " << fExcitationRatio << G4endl;
120  G4cout << " Track secondaries first: activated" << G4endl;
121  }
122  else {
123  G4cout << " Track secondaries first: inactivated" << G4endl;
124  }
125  }
126  if ( i == kWLS ) {
127  G4cout << " WLS process time profile: " << fProfile << G4endl;
128  }
129  }
130  }
131 }
132 
133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134 
135 #include "G4OpticalPhoton.hh"
136 
138 {
140 
141  // optical photon
143 }
144 
145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
146 #define DIR_CMDS "/process/optical"
147 #define GUIDANCE "Commands related to the optical physics simulation engine for "
148 
149 namespace UIhelpers {
150  //Helper functions to create UI commands for processes
151 
152  template<class T>
153  void buildCommands(T* proc,const char* dir, const char* guidance)
154  {
155  //Generic function to add a "verbose" command for a process
156  G4GenericMessenger* mess = new G4GenericMessenger(proc,dir,guidance);
158  G4GenericMessenger::Command& wlscmd4 = mess->DeclareMethod("verbose",&T::SetVerboseLevel,
159  "Set the verbose level");
160  wlscmd4.SetParameterName("ver",true);
161  wlscmd4.SetDefaultValue("1");
162  wlscmd4.SetStates(G4State_Idle);
163  }
164 
165  void buildCommands( G4OpWLS* op )
166  {
167  //Build UI commands for WLS
168  G4GenericMessenger* mess = new G4GenericMessenger(op,DIR_CMDS"/wls/",GUIDANCE" for WLS process.");
170  //Here, more correctly DeclareProperty should be used, but to do that, I would need public/friend
171  //access of G4GenericMessenger to private members of G4Scintillation. Thus I use the approach
172  //of DeclareMethod, that has a draw back: range checking does not work
173  G4GenericMessenger::Command& wlscmd1 = mess->DeclareMethod("setTimeProfile",&G4OpWLS::UseTimeProfile,
174  "Set the WLS time profile (delta or exponential)");
175  wlscmd1.SetParameterName("profile",false);
176  wlscmd1.SetCandidates("delta exponential");
177  wlscmd1.SetStates(G4State_Idle);
178  buildCommands(op,DIR_CMDS"/wls/",GUIDANCE" for WLS process.");
179  }
180 
181  void buildCommands(G4Scintillation* ScintillationProcess)
182  {
183  //Build UI commands for scintillation
184  G4GenericMessenger* mess = new G4GenericMessenger(ScintillationProcess,DIR_CMDS"/scintillation/",GUIDANCE" for scintillation process.");
186  G4GenericMessenger::Command& sccmd1 = mess->DeclareMethod("setFiniteRiseTime",
188  "Set option of a finite rise-time for G4Scintillation - If set, the G4Scintillation process expects the user to have set the constant material property FAST/SLOWSCINTILLATIONRISETIME");
189  sccmd1.SetParameterName("time",false);
190  sccmd1.SetStates(G4State_Idle);
191 
192  G4GenericMessenger::Command& sccmd2 = mess->DeclareMethod("setYieldFactor",
194  "Set scintillation yield factor");
195  sccmd2.SetParameterName("factor",false);
196  //sccmd2.SetRange("factor>=0."); //LIMITATION: w/ DeclareMethod range checking does not work
197  sccmd2.SetStates(G4State_Idle);
198 
199  G4GenericMessenger::Command& sccmd3 = mess->DeclareMethod("setExcitationRatio",
201  "Set scintillation excitation ratio");
202  sccmd3.SetParameterName("ratio",false);
203  //sccmd3.SetRange("ratio>=0.&&ratio<=1.");//LIMITATION: w/ DeclareMethod range checking does not work
204  sccmd3.SetStates(G4State_Idle);
205  G4GenericMessenger::Command& sccmd4 = mess->DeclareMethod("setByParticleType",
207  "Activate/Inactivate scintillation process by particle type");
208  sccmd4.SetParameterName("flag", false);
209  sccmd4.SetStates(G4State_Idle);
210 
211  G4GenericMessenger::Command& sccmd5 = mess->DeclareMethod("setTrackSecondariesFirst",
213  "Set option to track secondaries before finishing their parent track");
214  sccmd5.SetParameterName("flag",false);
215  sccmd5.SetStates(G4State_Idle);
216 
217  buildCommands(ScintillationProcess,DIR_CMDS"/scintillation/",GUIDANCE" for scintillation process.");
218  }
219 
220  void buildCommands(G4Cerenkov* CerenkovProcess)
221  {
222  //BUild UI commands for cerenkov
223  G4GenericMessenger* mess = new G4GenericMessenger(CerenkovProcess,DIR_CMDS"/cerenkov/",GUIDANCE" for Cerenkov process.");
225  G4GenericMessenger::Command& cecmd1 = mess->DeclareMethod("setMaxPhotons",
227  "Set maximum number of photons per step");
228  cecmd1.SetParameterName("max",false);
229  //cecmd1.SetRange("max>=0");//LIMITATION: w/ DeclareMethod range checking does not work
230  cecmd1.SetStates(G4State_Idle);
231 
232  G4GenericMessenger::Command& cecmd2 = mess->DeclareMethod("setMaxBetaChange",
234  "Set maximum change of beta of parent particle per step");
235  cecmd2.SetParameterName("max",false);
236  //cecmd2.SetRange("max>=0.");//LIMITATION: w/ DeclareMethod range checking does not work
237  cecmd2.SetStates(G4State_Idle);
238 
239  G4GenericMessenger::Command& cecmd3 = mess->DeclareMethod("setTrackSecondariesFirst",
241  "Set option to track secondaries before finishing their parent track");
242  cecmd3.SetParameterName("flag",false);
243  cecmd3.SetStates(G4State_Idle);
244 
245  buildCommands(CerenkovProcess,DIR_CMDS"/cerenkov/",GUIDANCE" for Cerenkov process.");
246  }
247 }
248 
249 
250 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
251 
252 #include "G4Threading.hh"
253 #include "G4ParticleDefinition.hh"
254 #include "G4ProcessManager.hh"
255 #include "G4AutoDelete.hh"
256 
258 {
259 // Construct optical processes.
260 
261  if(verboseLevel>0)
262  G4cout <<"G4OpticalPhysics:: Add Optical Physics Processes"<< G4endl;
263 
264  // A vector of optical processes
265  std::vector<G4VProcess*> OpProcesses;
266 
267  for ( G4int i=0; i<kNoProcess; i++ ) OpProcesses.push_back(NULL);
268 
269  // Add Optical Processes
270 
271  G4OpAbsorption* OpAbsorptionProcess = new G4OpAbsorption();
272  UIhelpers::buildCommands(OpAbsorptionProcess,DIR_CMDS"/absorption/",GUIDANCE" for absorption process");
273  OpProcesses[kAbsorption] = OpAbsorptionProcess;
274 
275  G4OpRayleigh* OpRayleighScatteringProcess = new G4OpRayleigh();
276  UIhelpers::buildCommands(OpRayleighScatteringProcess,DIR_CMDS"/rayleigh/",GUIDANCE" for Reyleigh scattering process");
277  OpProcesses[kRayleigh] = OpRayleighScatteringProcess;
278 
279  G4OpMieHG* OpMieHGScatteringProcess = new G4OpMieHG();
280  UIhelpers::buildCommands(OpMieHGScatteringProcess,DIR_CMDS"/mie/",GUIDANCE" for Mie cattering process");
281  OpProcesses[kMieHG] = OpMieHGScatteringProcess;
282 
283  G4OpBoundaryProcess* OpBoundaryProcess = new G4OpBoundaryProcess();
284  UIhelpers::buildCommands(OpBoundaryProcess,DIR_CMDS"/boundary/",GUIDANCE" for boundary process");
285  OpProcesses[kBoundary] = OpBoundaryProcess;
286 
287  G4OpWLS* OpWLSProcess = new G4OpWLS();
288  OpWLSProcess->UseTimeProfile(fProfile);
289  UIhelpers::buildCommands(OpWLSProcess);
290  OpProcesses[kWLS] = OpWLSProcess;
291 
292  G4ProcessManager * pManager = 0;
294 
295  if (!pManager) {
296  std::ostringstream o;
297  o << "Optical Photon without a Process Manager";
298  G4Exception("G4OpticalPhysics::ConstructProcess()","",
299  FatalException,o.str().c_str());
300  return;
301  }
302 
303  for ( G4int i=kAbsorption; i<=kWLS; i++ ) {
304  if ( fProcessUse[i] ) {
305  pManager->AddDiscreteProcess(OpProcesses[i]);
306  }
307  }
308 
309  G4Scintillation* ScintillationProcess = new G4Scintillation();
310  ScintillationProcess->SetScintillationYieldFactor(fYieldFactor);
311  ScintillationProcess->SetScintillationExcitationRatio(fExcitationRatio);
312  ScintillationProcess->SetFiniteRiseTime(fFiniteRiseTime);
316  ScintillationProcess->AddSaturation(emSaturation);
317  UIhelpers::buildCommands(ScintillationProcess);
318  OpProcesses[kScintillation] = ScintillationProcess;
319 
320  G4Cerenkov* CerenkovProcess = new G4Cerenkov();
321  CerenkovProcess->SetMaxNumPhotonsPerStep(fMaxNumPhotons);
322  CerenkovProcess->SetMaxBetaChangePerStep(fMaxBetaChange);
324  UIhelpers::buildCommands(CerenkovProcess);
325  OpProcesses[kCerenkov] = CerenkovProcess;
326 
327  aParticleIterator->reset();
328 
329  while( (*aParticleIterator)() ){
330 
331  G4ParticleDefinition* particle = aParticleIterator->value();
332  G4String particleName = particle->GetParticleName();
333 
334  pManager = particle->GetProcessManager();
335  if (!pManager) {
336  std::ostringstream o;
337  o << "Particle " << particleName << "without a Process Manager";
338  G4Exception("G4OpticalPhysics::ConstructProcess()","",
339  FatalException,o.str().c_str());
340  return; // else coverity complains for pManager use below
341  }
342 
343  if( CerenkovProcess->IsApplicable(*particle) &&
345  pManager->AddProcess(CerenkovProcess);
346  pManager->SetProcessOrdering(CerenkovProcess,idxPostStep);
347  }
348  if( ScintillationProcess->IsApplicable(*particle) &&
350  pManager->AddProcess(ScintillationProcess);
351  pManager->SetProcessOrderingToLast(ScintillationProcess,idxAtRest);
352  pManager->SetProcessOrderingToLast(ScintillationProcess,idxPostStep);
353  }
354 
355  }
356 
357  // Add verbose
358  for ( G4int i=0; i<kNoProcess; i++ ) {
359  if ( fProcessUse[i] ) OpProcesses[i]->SetVerboseLevel(verboseLevel);
360  }
361 
362  if (verboseLevel > 1) PrintStatistics();
363  if (verboseLevel > 0)
364  G4cout << "### " << namePhysics << " physics constructed." << G4endl;
365 }
366 
368 {
370 
371  fYieldFactor = yieldFactor;
372  //G4Scintillation::SetScintillationYieldFactor(yieldFactor);
373 }
374 
376 {
378 
379  fExcitationRatio = excitationRatio;
380 // G4Scintillation::SetScintillationExcitationRatio(excitationRatio);
381 }
382 
384 {
386 
387  fMaxNumPhotons = maxNumPhotons;
388 // G4Cerenkov::SetMaxNumPhotonsPerStep(maxNumPhotons);
389 }
390 
392 {
394 
395  fMaxBetaChange = maxBetaChange;
396 // G4Cerenkov::SetMaxBetaChangePerStep(maxBetaChange);
397 }
398 
400 {
402  fProfile = profile;
403 }
404 
405 //void G4OpticalPhysics::AddScintillationSaturation(G4EmSaturation* saturation)
406 //{
408 // G4Scintillation::AddSaturation(saturation);
409 //}
410 
412  SetScintillationByParticleType(G4bool scintillationByParticleType)
413 {
414  fScintillationByParticleType = scintillationByParticleType;
415  //G4Scintillation::SetScintillationByParticleType(scintillationByParticleType);
416 }
417 
419  G4bool trackSecondariesFirst)
420 {
421  if ( index >= kNoProcess ) return;
422  if ( fProcessTrackSecondariesFirst[index] == trackSecondariesFirst ) return;
423  fProcessTrackSecondariesFirst[index] = trackSecondariesFirst;
424 
425 // if ( index == kCerenkov )
426 // G4Cerenkov::SetTrackSecondariesFirst(trackSecondariesFirst);
427 // if ( index == kScintillation)
428 // G4Scintillation::SetTrackSecondariesFirst(trackSecondariesFirst);
429 }
430 
432 {
433  fFiniteRiseTime = finiteRiseTime;
434  //G4Scintillation::SetFiniteRiseTime(finiteRiseTime);
435 }
436 
438 {
439  // Configure the physics constructor to use/not use a selected process.
440  // This method can only be called in PreInit> phase (before execution of
441  // ConstructProcess). The process is not added to particle's process manager
442  // and so it cannot be re-activated later in Idle> phase with the command
443  // /process/activate.
444 
445  if ( index >= kNoProcess ) return;
446  if ( fProcessUse[index] == isUse ) return;
447  fProcessUse[index] = isUse;
448 }
449 
450 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void SetScintillationByParticleType(const G4bool)
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.cc:151
static G4LossTableManager * Instance()
void SetFiniteRiseTime(const G4bool state)
Number of processes, no selected process.
virtual void ConstructParticle()
G4OpticalProcessIndex
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4Cerenkov.cc:146
void SetMaxBetaChangePerStep(G4double)
This class is generic messenger.
G4String name
Definition: TRTMaterials.hh:40
G4OpticalPhysicsMessenger * fMessenger
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
Command & DeclareMethod(const G4String &name, const G4AnyMethod &fun, const G4String &doc="")
Scintillation process index.
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
Mie scattering process index.
void SetFiniteRiseTime(G4bool)
Absorption process index.
G4ProcessManager * GetProcessManager() const
int G4int
Definition: G4Types.hh:78
#define DIR_CMDS
const G4String & GetParticleName() const
void SetScintillationByParticleType(G4bool)
std::vector< G4bool > fProcessTrackSecondariesFirst
G4double fMaxBetaChange
max change of beta per step
void UseTimeProfile(const G4String name)
Definition: G4OpWLS.cc:417
virtual void ConstructProcess()
virtual ~G4OpticalPhysics()
Command & SetDefaultValue(const G4String &)
void SetScintillationYieldFactor(const G4double yieldfactor)
void Register(T *inst)
Definition: G4AutoDelete.hh:65
std::vector< G4bool > fProcessUse
G4double fYieldFactor
scintillation yield factor
G4GLOB_DLL std::ostream G4cout
void SetMaxNumPhotonsPerStep(G4int)
void AddSaturation(G4EmSaturation *)
G4bool fFiniteRiseTime
option to set a finite rise-time; Note: the G4Scintillation process expects the user to have set the ...
bool G4bool
Definition: G4Types.hh:79
#define aParticleIterator
void Configure(G4OpticalProcessIndex, G4bool)
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
G4EmSaturation * EmSaturation()
G4int fMaxNumPhotons
max number of Cerenkov photons per step
G4_DECLARE_PHYSCONSTR_FACTORY(G4OpticalPhysics)
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.cc:156
G4double fExcitationRatio
scintillation excitation ratio
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
Command & SetStates(G4ApplicationState s0)
Wave Length Shifting process index.
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
void SetScintillationExcitationRatio(const G4double ratio)
static G4OpticalPhoton * OpticalPhoton()
void SetScintillationExcitationRatio(G4double)
G4String G4OpticalProcessName(G4int)
Return the name for a given optical process index.
G4OpticalPhysics(G4int verbose=0, const G4String &name="Optical")
Boundary process index.
void SetScintillationYieldFactor(G4double)
void SetTrackSecondariesFirst(const G4bool state)
void SetProcessOrderingToLast(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
Command & SetCandidates(const G4String &)
void buildCommands(T *proc, const char *dir, const char *guidance)
Cerenkov process index.
void SetWLSTimeProfile(G4String)
void SetTrackSecondariesFirst(G4OpticalProcessIndex, G4bool)
Command & SetParameterName(const G4String &, G4bool, G4bool=false)
#define G4endl
Definition: G4ios.hh:61
static G4OpticalPhoton * OpticalPhotonDefinition()
Rayleigh scattering process index.
#define GUIDANCE
double G4double
Definition: G4Types.hh:76
G4String fProfile
the WLS process time profile
void PrintStatistics() const
G4bool fScintillationByParticleType
option to allow for the light yield to be a function of particle type and deposited energy in case of...
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
Definition: G4Cerenkov.cc:135