Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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 "G4LossTableManager.hh"
42#include "G4EmSaturation.hh"
43
44
45// factory
47//
49
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52
55
56 wasActivated(false),
57
58 fScintillationProcess(NULL),
59 fCerenkovProcess(NULL),
60 fOpWLSProcess(NULL),
61 fOpAbsorptionProcess(NULL),
62 fOpRayleighScatteringProcess(NULL),
63 fOpMieHGScatteringProcess(NULL),
64 fOpBoundaryProcess(NULL),
65 fMaxNumPhotons(100),
66 fMaxBetaChange(10.0),
67 fYieldFactor(1.),
68 fExcitationRatio(0.0),
69 fProfile("delta"),
70 fFiniteRiseTime(false),
71 fScintillationByParticleType(false)
72{
73 verboseLevel = verbose;
74 fMessenger = new G4OpticalPhysicsMessenger(this);
75
76 for ( G4int i=0; i<kNoProcess; i++ ) {
77 fProcesses.push_back(NULL);
78 fProcessUse.push_back(true);
79 fProcessVerbose.push_back(verbose);
80 fProcessTrackSecondariesFirst.push_back(true);
81 }
82}
83
84//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
85
87{
88 delete fMessenger;
89
90 if (wasActivated) {
91
92 if (fScintillationProcess) delete fScintillationProcess;
93 if (fCerenkovProcess) delete fCerenkovProcess;
94 if (fOpWLSProcess) delete fOpWLSProcess;
95
96 if (fOpAbsorptionProcess) delete fOpAbsorptionProcess;
97 if (fOpRayleighScatteringProcess) delete fOpRayleighScatteringProcess;
98 if (fOpMieHGScatteringProcess) delete fOpMieHGScatteringProcess;
99 if (fOpBoundaryProcess) delete fOpBoundaryProcess;
100
101 }
102}
103
104//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
105
106void G4OpticalPhysics::PrintStatistics() const
107{
108// Print all processes activation and their parameters
109
110 for ( G4int i=0; i<kNoProcess; i++ ) {
111 G4cout << " " << G4OpticalProcessName(i) << " process: ";
112 if ( ! fProcessUse[i] ) {
113 G4cout << "not used" << G4endl;
114 }
115 else {
116 G4cout << "used" << G4endl;
117 if ( i == kCerenkov ) {
118 G4cout << " Max number of photons per step: " << fMaxNumPhotons << G4endl;
119 G4cout << " Max beta change per step: " << fMaxBetaChange << G4endl;
120 if ( fProcessTrackSecondariesFirst[kCerenkov] ) G4cout << " Track secondaries first: activated" << G4endl;
121 }
122 if ( i == kScintillation ) {
123 if (fScintillationByParticleType)
124 G4cout << " Scintillation by Particle Type: activated " << G4endl;
125 G4cout << " Yield factor: " << fYieldFactor << G4endl;
126 G4cout << " ExcitationRatio: " << fExcitationRatio << G4endl;
127 if ( fProcessTrackSecondariesFirst[kScintillation] ) G4cout << " Track secondaries first: activated" << G4endl;
128 }
129 if ( i == kWLS ) {
130 G4cout << " WLS process time profile: " << fProfile << G4endl;
131 }
132 }
133 }
134}
135
136//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
137
138#include "G4OpticalPhoton.hh"
139
141{
142/// Instantiate particles.
143
144 // optical photon
146}
147
148//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
149
151#include "G4ProcessManager.hh"
152
154{
155// Construct optical processes.
156
157 if (wasActivated) return;
158
159 if(verboseLevel>0)
160 G4cout <<"G4OpticalPhysics:: Add Optical Physics Processes"<< G4endl;
161
162 // Add Optical Processes
163
164 fProcesses[kAbsorption] = fOpAbsorptionProcess = new G4OpAbsorption();
165 fProcesses[kRayleigh] = fOpRayleighScatteringProcess = new G4OpRayleigh();
166 fProcesses[kMieHG] = fOpMieHGScatteringProcess = new G4OpMieHG();
167
168 fProcesses[kBoundary] = fOpBoundaryProcess = new G4OpBoundaryProcess();
169
170 fProcesses[kWLS] = fOpWLSProcess = new G4OpWLS();
171 fOpWLSProcess->UseTimeProfile(fProfile);
172
173 G4ProcessManager * pManager = 0;
175
176 if (!pManager) {
177 std::ostringstream o;
178 o << "Optical Photon without a Process Manager";
179 G4Exception("G4OpticalPhysics::ConstructProcess()","",
180 FatalException,o.str().c_str());
181 return;
182 }
183
184 for ( G4int i=kAbsorption; i<=kWLS; i++ ) {
185 if ( fProcessUse[i] ) {
186 pManager->AddDiscreteProcess(fProcesses[i]);
187 }
188 }
189
190 fProcesses[kScintillation] = fScintillationProcess = new G4Scintillation();
191 fScintillationProcess->SetScintillationYieldFactor(fYieldFactor);
192 fScintillationProcess->SetScintillationExcitationRatio(fExcitationRatio);
193 fScintillationProcess->SetFiniteRiseTime(fFiniteRiseTime);
194 fScintillationProcess->SetScintillationByParticleType(fScintillationByParticleType);
195 fScintillationProcess->
196 SetTrackSecondariesFirst(fProcessTrackSecondariesFirst[kScintillation]);
197
198 // Use Birks Correction in the Scintillation process
199
201 fScintillationProcess->AddSaturation(emSaturation);
202
203 fProcesses[kCerenkov] = fCerenkovProcess = new G4Cerenkov();
204 fCerenkovProcess->SetMaxNumPhotonsPerStep(fMaxNumPhotons);
205 fCerenkovProcess->SetMaxBetaChangePerStep(fMaxBetaChange);
206 fCerenkovProcess->
207 SetTrackSecondariesFirst(fProcessTrackSecondariesFirst[kCerenkov]);
208
210
211 while( (*theParticleIterator)() ){
212
214 G4String particleName = particle->GetParticleName();
215
216 pManager = particle->GetProcessManager();
217 if (!pManager) {
218 std::ostringstream o;
219 o << "Particle " << particleName << "without a Process Manager";
220 G4Exception("G4OpticalPhysics::ConstructProcess()","",
221 FatalException,o.str().c_str());
222 return; // else coverity complains for pManager use below
223 }
224
225 if( fCerenkovProcess->IsApplicable(*particle) &&
226 fProcessUse[kCerenkov] ) {
227 pManager->AddProcess(fCerenkovProcess);
228 pManager->SetProcessOrdering(fCerenkovProcess,idxPostStep);
229 }
230 if( fScintillationProcess->IsApplicable(*particle) &&
231 fProcessUse[kScintillation] ){
232 pManager->AddProcess(fScintillationProcess);
233 pManager->SetProcessOrderingToLast(fScintillationProcess,idxAtRest);
234 pManager->SetProcessOrderingToLast(fScintillationProcess,idxPostStep);
235 }
236
237 }
238
239 // Add verbose
240 for ( G4int i=0; i<kNoProcess; i++ ) {
241 fProcesses[i]->SetVerboseLevel(fProcessVerbose[i]);
242 }
243
244 if (verboseLevel > 1) PrintStatistics();
245 if (verboseLevel > 0)
246 G4cout << "### " << namePhysics << " physics constructed." << G4endl;
247
248 wasActivated = true;
249
250}
251
253{
254/// Set the scintillation yield factor
255
256 fYieldFactor = yieldFactor;
257
258 if(fScintillationProcess)
259 fScintillationProcess->SetScintillationYieldFactor(yieldFactor);
260}
261
263{
264/// Set the scintillation excitation ratio
265
266 fExcitationRatio = excitationRatio;
267
268 if(fScintillationProcess)
269 fScintillationProcess->SetScintillationExcitationRatio(excitationRatio);
270}
271
273{
274/// Limit step to the specified maximum number of Cherenkov photons
275
276 fMaxNumPhotons = maxNumPhotons;
277
278 if(fCerenkovProcess)
279 fCerenkovProcess->SetMaxNumPhotonsPerStep(maxNumPhotons);
280}
281
283{
284/// Limit step to the specified maximum change of beta of the parent particle
285
286 fMaxBetaChange = maxBetaChange;
287
288 if(fCerenkovProcess)
289 fCerenkovProcess->SetMaxBetaChangePerStep(maxBetaChange);
290}
291
293{
294/// Set the WLS time profile (delta or exponential)
295
296 fProfile = profile;
297
298 if(fOpWLSProcess)
299 fOpWLSProcess->UseTimeProfile(fProfile);
300}
301
303{
304/// Adds Birks Saturation to the G4Scintillation Process
305
306 if(fScintillationProcess)
307 fScintillationProcess->AddSaturation(saturation);
308}
309
311{
312 fScintillationByParticleType = scintillationByParticleType;
313
314 if (fScintillationProcess)
315 fScintillationProcess->SetScintillationByParticleType(scintillationByParticleType);
316}
317
319 G4bool trackSecondariesFirst)
320{
321 if ( index >= kNoProcess ) return;
322 if ( fProcessTrackSecondariesFirst[index] == trackSecondariesFirst ) return;
323 fProcessTrackSecondariesFirst[index] = trackSecondariesFirst;
324
325 if(fCerenkovProcess && index == kCerenkov )
326 fCerenkovProcess->SetTrackSecondariesFirst(trackSecondariesFirst);
327
328 if(fScintillationProcess && index == kScintillation)
329 fScintillationProcess->SetTrackSecondariesFirst(trackSecondariesFirst);
330}
331
333{
334 fFiniteRiseTime = finiteRiseTime;
335 if(fScintillationProcess)
336 fScintillationProcess->SetFiniteRiseTime(finiteRiseTime);
337}
338
340{
341 // Configure the physics constructor to use/not use a selected process.
342 // This method can only be called in PreInit> phase (before execution of
343 // ConstructProcess). The process is not added to particle's process manager
344 // and so it cannot be re-activated later in Idle> phase with the command
345 // /process/activate.
346
347 if ( index >= kNoProcess ) return;
348 if ( fProcessUse[index] == isUse ) return;
349 fProcessUse[index] = isUse;
350}
351
353 G4int inputVerboseLevel)
354{
355 // Set new verbose level to a selected process
356
357 if ( index >= kNoProcess ) return;
358 if ( fProcessVerbose[index] == inputVerboseLevel ) return;
359
360 fProcessVerbose[index] = inputVerboseLevel;
361
362 if ( fProcesses[index] ) fProcesses[index]->SetVerboseLevel(inputVerboseLevel);
363}
364
365//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ FatalException
G4OpticalProcessIndex
@ kWLS
Wave Length Shifting process index.
@ kScintillation
Scintillation process index.
@ kRayleigh
Rayleigh scattering process index.
@ kAbsorption
Absorption process index.
@ kBoundary
Boundary process index.
@ kNoProcess
Number of processes, no selected process.
@ kCerenkov
Cerenkov process index.
@ kMieHG
Mie scattering process index.
G4String G4OpticalProcessName(G4int)
Return the name for a given optical process index.
#define G4_DECLARE_PHYSCONSTR_FACTORY(physics_constructor)
@ idxPostStep
@ idxAtRest
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void SetMaxBetaChangePerStep(const G4double d)
Definition: G4Cerenkov.hh:229
void SetTrackSecondariesFirst(const G4bool state)
Definition: G4Cerenkov.hh:223
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
Definition: G4Cerenkov.hh:214
void SetMaxNumPhotonsPerStep(const G4int NumPhotons)
Definition: G4Cerenkov.hh:235
static G4LossTableManager * Instance()
G4EmSaturation * EmSaturation()
void UseTimeProfile(const G4String name)
Definition: G4OpWLS.cc:385
static G4OpticalPhoton * OpticalPhotonDefinition()
static G4OpticalPhoton * OpticalPhoton()
void SetProcessVerbose(G4int, G4int)
void SetMaxBetaChangePerStep(G4double)
void SetScintillationYieldFactor(G4double)
void Configure(G4OpticalProcessIndex, G4bool)
void SetScintillationByParticleType(G4bool)
void SetMaxNumPhotonsPerStep(G4int)
virtual ~G4OpticalPhysics()
void SetTrackSecondariesFirst(G4OpticalProcessIndex, G4bool)
virtual void ConstructProcess()
virtual void ConstructParticle()
void SetWLSTimeProfile(G4String)
void AddScintillationSaturation(G4EmSaturation *)
void SetFiniteRiseTime(G4bool)
void SetScintillationExcitationRatio(G4double)
G4OpticalPhysics(G4int verbose=0, const G4String &name="Optical")
G4ProcessManager * GetProcessManager() const
const G4String & GetParticleName() const
void SetProcessOrdering(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt, G4int ordDoIt=ordDefault)
G4int AddDiscreteProcess(G4VProcess *aProcess, G4int ord=ordDefault)
G4int AddProcess(G4VProcess *aProcess, G4int ordAtRestDoIt=ordInActive, G4int ordAlongSteptDoIt=ordInActive, G4int ordPostStepDoIt=ordInActive)
void SetProcessOrderingToLast(G4VProcess *aProcess, G4ProcessVectorDoItIndex idDoIt)
void SetScintillationYieldFactor(const G4double yieldfactor)
void SetTrackSecondariesFirst(const G4bool state)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
void SetScintillationExcitationRatio(const G4double excitationratio)
void AddSaturation(G4EmSaturation *sat)
void SetFiniteRiseTime(const G4bool state)
void SetScintillationByParticleType(const G4bool)
G4ParticleTable::G4PTblDicIterator * theParticleIterator
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41