Geant4  10.03
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  fScintillationTrackInfo(false),
75  fInvokeSD(true)
76 {
77  verboseLevel = verbose;
79 
80  for ( G4int i=0; i<kNoProcess; i++ ) {
81  fProcessUse.push_back(true);
82  fProcessTrackSecondariesFirst.push_back(true);
83  }
84 }
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
87 
89 {
90  delete fMessenger;
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
94 
96 {
97 // Print all processes activation and their parameters
98 
99  for ( G4int i=0; i<kNoProcess; i++ ) {
100  G4cout << " " << G4OpticalProcessName(i) << " process: ";
101  if ( ! fProcessUse[i] ) {
102  G4cout << "not used" << G4endl;
103  }
104  else {
105  G4cout << "used" << G4endl;
106  if ( i == kCerenkov ) {
107  G4cout << " Max number of photons per step: " << fMaxNumPhotons << G4endl;
108  G4cout << " Max beta change per step: " << fMaxBetaChange << G4endl;
110  G4cout << " Track secondaries first: activated" << G4endl;
111  }
112  else {
113  G4cout << " Track secondaries first: inactivated" << G4endl;
114  }
115  }
116  if ( i == kScintillation ) {
118  G4cout << " Scintillation by Particle Type: activated " << G4endl;
119  G4cout << " Yield factor: " << fYieldFactor << G4endl;
120  G4cout << " ExcitationRatio: " << fExcitationRatio << G4endl;
122  G4cout << " Track secondaries first: activated" << G4endl;
123  }
124  else {
125  G4cout << " Track secondaries first: inactivated" << G4endl;
126  }
127  }
128  if ( i == kWLS ) {
129  G4cout << " WLS process time profile: " << fProfile << G4endl;
130  }
131  }
132  }
133 }
134 
135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136 
137 #include "G4OpticalPhoton.hh"
138 
140 {
142 
143  // optical photon
145 }
146 
147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148 #define DIR_CMDS "/process/optical"
149 #define GUIDANCE "Commands related to the optical physics simulation engine for "
150 
151 namespace UIhelpers {
152  //Helper functions to create UI commands for processes
153 
154  template<class T>
155  void buildCommands(T* proc,const char* dir, const char* guidance)
156  {
157  //Generic function to add a "verbose" command for a process
158  G4GenericMessenger* mess = new G4GenericMessenger(proc,dir,guidance);
160  G4GenericMessenger::Command& wlscmd4 = mess->DeclareMethod("verbose",&T::SetVerboseLevel,
161  "Set the verbose level");
162  wlscmd4.SetParameterName("ver",true);
163  wlscmd4.SetDefaultValue("1");
164  wlscmd4.SetStates(G4State_Idle);
165  }
166 
167  void buildCommands( G4OpWLS* op )
168  {
169  //Build UI commands for WLS
170  G4GenericMessenger* mess = new G4GenericMessenger(op,DIR_CMDS"/wls/",GUIDANCE" for WLS process.");
172  //Here, more correctly DeclareProperty should be used, but to do that, I would need public/friend
173  //access of G4GenericMessenger to private members of G4Scintillation. Thus I use the approach
174  //of DeclareMethod, that has a draw back: range checking does not work
175  G4GenericMessenger::Command& wlscmd1 = mess->DeclareMethod("setTimeProfile",&G4OpWLS::UseTimeProfile,
176  "Set the WLS time profile (delta or exponential)");
177  wlscmd1.SetParameterName("profile",false);
178  wlscmd1.SetCandidates("delta exponential");
179  wlscmd1.SetStates(G4State_Idle);
180  buildCommands(op,DIR_CMDS"/wls/",GUIDANCE" for WLS process.");
181  }
182 
183  void buildCommands(G4Scintillation* ScintillationProcess)
184  {
185  //Build UI commands for scintillation
186  G4GenericMessenger* mess = new G4GenericMessenger(ScintillationProcess,DIR_CMDS"/scintillation/",GUIDANCE" for scintillation process.");
188  G4GenericMessenger::Command& sccmd1 = mess->DeclareMethod("setFiniteRiseTime",
190  "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");
191  sccmd1.SetParameterName("time",false);
192  sccmd1.SetStates(G4State_Idle);
193 
194  G4GenericMessenger::Command& sccmd2 = mess->DeclareMethod("setYieldFactor",
196  "Set scintillation yield factor");
197  sccmd2.SetParameterName("factor",false);
198  //sccmd2.SetRange("factor>=0."); //LIMITATION: w/ DeclareMethod range checking does not work
199  sccmd2.SetStates(G4State_Idle);
200 
201  G4GenericMessenger::Command& sccmd3 = mess->DeclareMethod("setExcitationRatio",
203  "Set scintillation excitation ratio");
204  sccmd3.SetParameterName("ratio",false);
205  //sccmd3.SetRange("ratio>=0.&&ratio<=1.");//LIMITATION: w/ DeclareMethod range checking does not work
206  sccmd3.SetStates(G4State_Idle);
207 
208  G4GenericMessenger::Command& sccmd4 = mess->DeclareMethod("setByParticleType",
210  "Activate/Inactivate scintillation process by particle type");
211  sccmd4.SetParameterName("flag", false);
212  sccmd4.SetStates(G4State_Idle);
213 
214  G4GenericMessenger::Command& sccmd5 = mess->DeclareMethod("setTrackInfo",
216  "Activate/Inactivate scintillation track info");
217  sccmd5.SetParameterName("flag", false);
218  sccmd5.SetStates(G4State_Idle);
219 
220  G4GenericMessenger::Command& sccmd6 = mess->DeclareMethod("setTrackSecondariesFirst",
222  "Set option to track secondaries before finishing their parent track");
223  sccmd6.SetParameterName("flag",false);
224  sccmd6.SetStates(G4State_Idle);
225 
226  buildCommands(ScintillationProcess,DIR_CMDS"/scintillation/",GUIDANCE" for scintillation process.");
227  }
228 
229  void buildCommands(G4Cerenkov* CerenkovProcess)
230  {
231  //BUild UI commands for cerenkov
232  G4GenericMessenger* mess = new G4GenericMessenger(CerenkovProcess,DIR_CMDS"/cerenkov/",GUIDANCE" for Cerenkov process.");
234  G4GenericMessenger::Command& cecmd1 = mess->DeclareMethod("setMaxPhotons",
236  "Set maximum number of photons per step");
237  cecmd1.SetParameterName("max",false);
238  //cecmd1.SetRange("max>=0");//LIMITATION: w/ DeclareMethod range checking does not work
239  cecmd1.SetStates(G4State_Idle);
240 
241  G4GenericMessenger::Command& cecmd2 = mess->DeclareMethod("setMaxBetaChange",
243  "Set maximum change of beta of parent particle per step");
244  cecmd2.SetParameterName("max",false);
245  //cecmd2.SetRange("max>=0.");//LIMITATION: w/ DeclareMethod range checking does not work
246  cecmd2.SetStates(G4State_Idle);
247 
248  G4GenericMessenger::Command& cecmd3 = mess->DeclareMethod("setTrackSecondariesFirst",
250  "Set option to track secondaries before finishing their parent track");
251  cecmd3.SetParameterName("flag",false);
252  cecmd3.SetStates(G4State_Idle);
253 
254  buildCommands(CerenkovProcess,DIR_CMDS"/cerenkov/",GUIDANCE" for Cerenkov process.");
255  }
256 }
257 
258 
259 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
260 
261 #include "G4Threading.hh"
262 #include "G4ParticleDefinition.hh"
263 #include "G4ProcessManager.hh"
264 #include "G4AutoDelete.hh"
265 
267 {
268 // Construct optical processes.
269 
270  if(verboseLevel>0)
271  G4cout <<"G4OpticalPhysics:: Add Optical Physics Processes"<< G4endl;
272 
273  // A vector of optical processes
274  std::vector<G4VProcess*> OpProcesses;
275 
276  for ( G4int i=0; i<kNoProcess; i++ ) OpProcesses.push_back(NULL);
277 
278  // Add Optical Processes
279 
280  G4OpAbsorption* OpAbsorptionProcess = new G4OpAbsorption();
281  UIhelpers::buildCommands(OpAbsorptionProcess,DIR_CMDS"/absorption/",GUIDANCE" for absorption process");
282  OpProcesses[kAbsorption] = OpAbsorptionProcess;
283 
284  G4OpRayleigh* OpRayleighScatteringProcess = new G4OpRayleigh();
285  UIhelpers::buildCommands(OpRayleighScatteringProcess,DIR_CMDS"/rayleigh/",GUIDANCE" for Reyleigh scattering process");
286  OpProcesses[kRayleigh] = OpRayleighScatteringProcess;
287 
288  G4OpMieHG* OpMieHGScatteringProcess = new G4OpMieHG();
289  UIhelpers::buildCommands(OpMieHGScatteringProcess,DIR_CMDS"/mie/",GUIDANCE" for Mie cattering process");
290  OpProcesses[kMieHG] = OpMieHGScatteringProcess;
291 
292  G4OpBoundaryProcess* OpBoundaryProcess = new G4OpBoundaryProcess();
293  UIhelpers::buildCommands(OpBoundaryProcess,DIR_CMDS"/boundary/",GUIDANCE" for boundary process");
294  OpBoundaryProcess->SetInvokeSD(fInvokeSD);
295  OpProcesses[kBoundary] = OpBoundaryProcess;
296 
297  G4OpWLS* OpWLSProcess = new G4OpWLS();
298  OpWLSProcess->UseTimeProfile(fProfile);
299  UIhelpers::buildCommands(OpWLSProcess);
300  OpProcesses[kWLS] = OpWLSProcess;
301 
302  G4ProcessManager * pManager = 0;
304 
305  if (!pManager) {
306  std::ostringstream o;
307  o << "Optical Photon without a Process Manager";
308  G4Exception("G4OpticalPhysics::ConstructProcess()","",
309  FatalException,o.str().c_str());
310  return;
311  }
312 
313  for ( G4int i=kAbsorption; i<=kWLS; i++ ) {
314  if ( fProcessUse[i] ) {
315  pManager->AddDiscreteProcess(OpProcesses[i]);
316  }
317  }
318 
319  G4Scintillation* ScintillationProcess = new G4Scintillation();
320  ScintillationProcess->SetScintillationYieldFactor(fYieldFactor);
321  ScintillationProcess->SetScintillationExcitationRatio(fExcitationRatio);
322  ScintillationProcess->SetFiniteRiseTime(fFiniteRiseTime);
324  ScintillationProcess->SetScintillationTrackInfo(fScintillationTrackInfo);
327  ScintillationProcess->AddSaturation(emSaturation);
328  UIhelpers::buildCommands(ScintillationProcess);
329  OpProcesses[kScintillation] = ScintillationProcess;
330 
331  G4Cerenkov* CerenkovProcess = new G4Cerenkov();
332  CerenkovProcess->SetMaxNumPhotonsPerStep(fMaxNumPhotons);
333  CerenkovProcess->SetMaxBetaChangePerStep(fMaxBetaChange);
335  UIhelpers::buildCommands(CerenkovProcess);
336  OpProcesses[kCerenkov] = CerenkovProcess;
337 
338  auto myParticleIterator=GetParticleIterator();
339  myParticleIterator->reset();
340 
341  while( (*myParticleIterator)() ){
342 
343  G4ParticleDefinition* particle = myParticleIterator->value();
344  G4String particleName = particle->GetParticleName();
345 
346  pManager = particle->GetProcessManager();
347  if (!pManager) {
348  std::ostringstream o;
349  o << "Particle " << particleName << "without a Process Manager";
350  G4Exception("G4OpticalPhysics::ConstructProcess()","",
351  FatalException,o.str().c_str());
352  return; // else coverity complains for pManager use below
353  }
354 
355  if( CerenkovProcess->IsApplicable(*particle) &&
357  pManager->AddProcess(CerenkovProcess);
358  pManager->SetProcessOrdering(CerenkovProcess,idxPostStep);
359  }
360  if( ScintillationProcess->IsApplicable(*particle) &&
362  pManager->AddProcess(ScintillationProcess);
363  pManager->SetProcessOrderingToLast(ScintillationProcess,idxAtRest);
364  pManager->SetProcessOrderingToLast(ScintillationProcess,idxPostStep);
365  }
366 
367  }
368 
369  // Add verbose
370  for ( G4int i=0; i<kNoProcess; i++ ) {
371  if ( fProcessUse[i] ) OpProcesses[i]->SetVerboseLevel(verboseLevel);
372  }
373 
374  if (verboseLevel > 1) PrintStatistics();
375  if (verboseLevel > 0)
376  G4cout << "### " << namePhysics << " physics constructed." << G4endl;
377 }
378 
380 {
382 
383  fYieldFactor = yieldFactor;
384  //G4Scintillation::SetScintillationYieldFactor(yieldFactor);
385 }
386 
388 {
390 
391  fExcitationRatio = excitationRatio;
392 // G4Scintillation::SetScintillationExcitationRatio(excitationRatio);
393 }
394 
396 {
398 
399  fMaxNumPhotons = maxNumPhotons;
400 // G4Cerenkov::SetMaxNumPhotonsPerStep(maxNumPhotons);
401 }
402 
404 {
406 
407  fMaxBetaChange = maxBetaChange;
408 // G4Cerenkov::SetMaxBetaChangePerStep(maxBetaChange);
409 }
410 
412 {
414  fProfile = profile;
415 }
416 
417 //void G4OpticalPhysics::AddScintillationSaturation(G4EmSaturation* saturation)
418 //{
420 // G4Scintillation::AddSaturation(saturation);
421 //}
422 
424  SetScintillationByParticleType(G4bool scintillationByParticleType)
425 {
426  fScintillationByParticleType = scintillationByParticleType;
427  //G4Scintillation::SetScintillationByParticleType(scintillationByParticleType);
428 }
429 
431  SetScintillationTrackInfo(G4bool scintillationTrackInfo)
432 {
433  fScintillationTrackInfo = scintillationTrackInfo;
434  //G4Scintillation::SetScintillationTrackInfo(scintillationTrackInfo);
435 }
436 
438  G4bool trackSecondariesFirst)
439 {
440  if ( index >= kNoProcess ) return;
441  if ( fProcessTrackSecondariesFirst[index] == trackSecondariesFirst ) return;
442  fProcessTrackSecondariesFirst[index] = trackSecondariesFirst;
443 
444 // if ( index == kCerenkov )
445 // G4Cerenkov::SetTrackSecondariesFirst(trackSecondariesFirst);
446 // if ( index == kScintillation)
447 // G4Scintillation::SetTrackSecondariesFirst(trackSecondariesFirst);
448 }
449 
451 {
452  fFiniteRiseTime = finiteRiseTime;
453  //G4Scintillation::SetFiniteRiseTime(finiteRiseTime);
454 }
455 
457 {
458  fInvokeSD = invokeSD;
459 }
460 
462 {
463  // Configure the physics constructor to use/not use a selected process.
464  // This method can only be called in PreInit> phase (before execution of
465  // ConstructProcess). The process is not added to particle's process manager
466  // and so it cannot be re-activated later in Idle> phase with the command
467  // /process/activate.
468 
469  if ( index >= kNoProcess ) return;
470  if ( fProcessUse[index] == isUse ) return;
471  fProcessUse[index] = isUse;
472 }
473 
474 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
void SetScintillationByParticleType(const G4bool)
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.cc:150
void SetScintillationTrackInfo(const G4bool trackType)
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:145
void SetMaxBetaChangePerStep(G4double)
This class is generic messenger.
void SetInvokeSD(G4bool)
G4OpticalPhysicsMessenger * fMessenger
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
Command & DeclareMethod(const G4String &name, const G4AnyMethod &fun, const G4String &doc="")
Scintillation process index.
Mie scattering process index.
const char * name(G4int ptype)
void SetFiniteRiseTime(G4bool)
Absorption process index.
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 &)
G4ParticleTable::G4PTblDicIterator * GetParticleIterator() const
void SetScintillationYieldFactor(const G4double yieldfactor)
void Register(T *inst)
Definition: G4AutoDelete.hh:65
G4bool fInvokeSD
option to allow for G4OpBoundaryProcess to call/not call InvokeSD method
std::vector< G4bool > fProcessUse
G4double fYieldFactor
scintillation yield factor
G4GLOB_DLL std::ostream G4cout
G4bool fScintillationTrackInfo
option to allow for G4ScintillationTrackInformation to be attached to a scintillation photon's track ...
void SetMaxNumPhotonsPerStep(G4int)
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
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:155
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)
G4ProcessManager * GetProcessManager() const
Command & SetCandidates(const G4String &)
void buildCommands(T *proc, const char *dir, const char *guidance)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
Definition: G4Cerenkov.cc:133
Cerenkov process index.
void SetScintillationTrackInfo(G4bool)
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
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
G4String fProfile
the WLS process time profile
void PrintStatistics() const
void AddSaturation(G4EmSaturation *sat)
G4bool fScintillationByParticleType
option to allow for the light yield to be a function of particle type and deposited energy in case of...