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
G4OpWLS.hh
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// $Id$
28//
29////////////////////////////////////////////////////////////////////////
30// Optical Photon WaveLength Shifting (WLS) Class Definition
31////////////////////////////////////////////////////////////////////////
32//
33// File: G4OpWLS.hh
34// Description: Discrete Process -- Wavelength Shifting of Optical Photons
35// Version: 1.0
36// Created: 2003-05-13
37// Author: John Paul Archambault
38// (Adaptation of G4Scintillation and G4OpAbsorption)
39// Updated: 2005-07-28 add G4ProcessType to constructor
40// 2006-05-07 - add G4VWLSTimeGeneratorProfile
41// mail: gum@triumf.ca
42// jparcham@phys.ualberta.ca
43//
44////////////////////////////////////////////////////////////////////////
45
46#ifndef G4OpWLS_h
47#define G4OpWLS_h 1
48
49/////////////
50// Includes
51/////////////
52
53#include "globals.hh"
54#include "templates.hh"
55#include "Randomize.hh"
56#include "G4Poisson.hh"
57#include "G4ThreeVector.hh"
58#include "G4ParticleMomentum.hh"
59#include "G4Step.hh"
60#include "G4VDiscreteProcess.hh"
61#include "G4DynamicParticle.hh"
62#include "G4Material.hh"
63#include "G4OpticalPhoton.hh"
64#include "G4PhysicsTable.hh"
68
69// Class Description:
70// Discrete Process -- Bulk absorption of Optical Photons.
71// Class inherits publicly from G4VDiscreteProcess
72// Class Description - End:
73
74/////////////////////
75// Class Definition
76/////////////////////
77
79
81{
82
83public:
84
85 ////////////////////////////////
86 // Constructors and Destructor
87 ////////////////////////////////
88
89 G4OpWLS(const G4String& processName = "OpWLS",
90 G4ProcessType type = fOptical);
91 ~G4OpWLS();
92
93private:
94
95 G4OpWLS(const G4OpWLS &right);
96
97 //////////////
98 // Operators
99 //////////////
100
101 G4OpWLS& operator=(const G4OpWLS &right);
102
103public:
104
105 ////////////
106 // Methods
107 ////////////
108
109 G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
110 // Returns true -> 'is applicable' only for an optical photon.
111
112 G4double GetMeanFreePath(const G4Track& aTrack,
113 G4double ,
115 // Returns the absorption length for bulk absorption of optical
116 // photons in media with a specified attenuation length.
117
119 const G4Step& aStep);
120 // This is the method implementing bulk absorption of optical
121 // photons.
122
124 // Returns the address of the WLS integral table.
125
126 void DumpPhysicsTable() const;
127 // Prints the WLS integral table.
128
129 void UseTimeProfile(const G4String name);
130 // Selects the time profile generator
131
132private:
133
134 void BuildThePhysicsTable();
135 // Is the WLS integral table;
136
137protected:
138
141
142};
143
144////////////////////
145// Inline methods
146////////////////////
147
148inline
150{
151 return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
152}
153
154inline
156{
157 return theIntegralTable;
158}
159
160inline
162{
163 G4int PhysicsTableSize = theIntegralTable->entries();
165
166 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
167 {
169 v->DumpValues();
170 }
171}
172
173#endif /* G4OpWLS_h */
G4ForceCondition
G4ProcessType
@ fOptical
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4OpWLS.cc:351
G4PhysicsTable * GetIntegralTable() const
Definition: G4OpWLS.hh:155
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
Definition: G4OpWLS.hh:149
void DumpPhysicsTable() const
Definition: G4OpWLS.hh:161
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4OpWLS.cc:102
G4PhysicsTable * theIntegralTable
Definition: G4OpWLS.hh:140
~G4OpWLS()
Definition: G4OpWLS.cc:85
G4VWLSTimeGeneratorProfile * WLSTimeGeneratorProfile
Definition: G4OpWLS.hh:139
void UseTimeProfile(const G4String name)
Definition: G4OpWLS.cc:385
static G4OpticalPhoton * OpticalPhoton()
size_t entries() const
Definition: G4Step.hh:78