Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4OpWLS2.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// Optical Photon WaveLength Shifting (WLS) Class Implementation
30////////////////////////////////////////////////////////////////////////
31//
32// File: G4OpWLS2.cc
33// Description: Discrete Process -- Wavelength Shifting of Optical Photons
34// Version: 1.0
35// Created: 2003-05-13
36// Author: John Paul Archambault
37// (Adaptation of G4Scintillation and G4OpAbsorption)
38// Updated: 2005-07-28 - add G4ProcessType to constructor
39// 2006-05-07 - add G4VWLSTimeGeneratorProfile
40//
41////////////////////////////////////////////////////////////////////////
42
43#include "G4OpWLS2.hh"
44#include "G4ios.hh"
46#include "G4SystemOfUnits.hh"
47#include "G4OpProcessSubType.hh"
48#include "G4Poisson.hh"
52
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
54G4OpWLS2::G4OpWLS2(const G4String& processName, G4ProcessType type)
55 : G4VDiscreteProcess(processName, type)
56{
58 Initialise();
60 theIntegralTable = nullptr;
61
62 if(verboseLevel > 0)
63 G4cout << GetProcessName() << " is created " << G4endl;
64}
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
68{
70 {
71 theIntegralTable->clearAndDestroy();
72 delete theIntegralTable;
73 }
75}
76
77//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82
83//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
90
91//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93 const G4Step& aStep)
94{
95 aParticleChange.Initialize(aTrack);
96 aParticleChange.ProposeTrackStatus(fStopAndKill);
97
98 if(verboseLevel > 1)
99 {
100 G4cout << "\n** G4OpWLS2: Photon absorbed! **" << G4endl;
101 }
102
105 if(!MPT)
106 {
107 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
108 }
109
110 G4PhysicsFreeVector* WLSIntegral = nullptr;
111 if(!MPT->GetProperty(kWLSCOMPONENT2))
112 {
113 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
114 }
115 else
116 {
117 // Retrieve the WLS Integral for this material
118 WLSIntegral = (G4PhysicsFreeVector*) ((*theIntegralTable)
119 (aTrack.GetMaterial()->GetIndex()));
120 }
121
122 G4double primaryEnergy = aTrack.GetDynamicParticle()->GetKineticEnergy();
123 // No WLS photons are produced if the primary photon's energy is below
124 // the lower bound of the WLS integral range
125 if(primaryEnergy < WLSIntegral->GetMinValue())
126 {
127 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
128 }
129
130 G4int NumPhotons = 1;
132 {
133 G4double MeanNumberOfPhotons =
135 NumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
136 if(NumPhotons <= 0)
137 {
138 // return unchanged particle and no secondaries
139 aParticleChange.SetNumberOfSecondaries(0);
140 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
141 }
142 }
143
145 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
146 std::vector<G4Track*> proposedSecondaries;
147
148 // Max WLS Integral - force sampling of a WLS photon with energy below
149 // the primary photon's energy
150 G4double CIImax = (primaryEnergy > WLSIntegral->GetMaxEnergy()) ?
151 WLSIntegral->GetMaxValue() : WLSIntegral->Value(primaryEnergy);
152
153 for(G4int i = 0; i < NumPhotons; ++i)
154 {
155 // Determine photon energy
156 G4double sampledEnergy = WLSIntegral->GetEnergy(G4UniformRand() * CIImax);
157
158 // Make sure the energy of the secondary is less than that of the primary
159 if(sampledEnergy > primaryEnergy)
160 {
162 ed << "Sampled photon energy " << sampledEnergy << " is greater than "
163 << "the primary photon energy " << primaryEnergy << G4endl;
164 G4Exception("G4OpWLS2::PostStepDoIt", "WSL02", FatalException, ed);
165 }
166 else if(verboseLevel > 1)
167 {
168 G4cout << "G4OpWLS2: Created photon with energy: " << sampledEnergy
169 << G4endl;
170 }
171
172 // Generate random photon direction
173 G4double cost = 1. - 2. * G4UniformRand();
174 G4double sint = std::sqrt((1. - cost) * (1. + cost));
175 G4double phi = twopi * G4UniformRand();
176 G4double sinp = std::sin(phi);
177 G4double cosp = std::cos(phi);
178 G4ParticleMomentum photonMomentum(sint * cosp, sint * sinp, cost);
179
180 G4ThreeVector photonPolarization(cost * cosp, cost * sinp, -sint);
181 G4ThreeVector perp = photonMomentum.cross(photonPolarization);
182
183 phi = twopi * G4UniformRand();
184 sinp = std::sin(phi);
185 cosp = std::cos(phi);
186 photonPolarization = (cosp * photonPolarization + sinp * perp).unit();
187
188 // Generate a new photon:
189 auto sec_dp =
191 sec_dp->SetPolarization(photonPolarization);
192 sec_dp->SetKineticEnergy(sampledEnergy);
193
194 G4double secTime = pPostStepPoint->GetGlobalTime() +
195 WLSTimeGeneratorProfile->GenerateTime(WLSTime);
196 G4ThreeVector secPos = pPostStepPoint->GetPosition();
197 G4Track* secTrack = new G4Track(sec_dp, secTime, secPos);
198
199 secTrack->SetTouchableHandle(aTrack.GetTouchableHandle());
200 secTrack->SetParentID(aTrack.GetTrackID());
201
202 proposedSecondaries.push_back(secTrack);
203 }
204
205 aParticleChange.SetNumberOfSecondaries((G4int)proposedSecondaries.size());
206 for(auto sec : proposedSecondaries)
207 {
208 aParticleChange.AddSecondary(sec);
209 }
210 if(verboseLevel > 1)
211 {
212 G4cout << "\n Exiting from G4OpWLS2::DoIt -- NumberOfSecondaries = "
213 << aParticleChange.GetNumberOfSecondaries() << G4endl;
214 }
215
216 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
217}
218
219//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
221{
223 {
224 theIntegralTable->clearAndDestroy();
225 delete theIntegralTable;
226 theIntegralTable = nullptr;
227 }
228
229 const G4MaterialTable* materialTable = G4Material::GetMaterialTable();
230 std::size_t numOfMaterials = G4Material::GetNumberOfMaterials();
231 theIntegralTable = new G4PhysicsTable(numOfMaterials);
232
233 // loop for materials
234 for(std::size_t i = 0; i < numOfMaterials; ++i)
235 {
236 auto physVector = new G4PhysicsFreeVector();
237
238 // Retrieve vector of WLS2 wavelength intensity for
239 // the material from the material's optical properties table.
241 (*materialTable)[i]->GetMaterialPropertiesTable();
242 if(MPT)
243 {
245 if(wlsVector)
246 {
247 // Retrieve the first intensity point in vector
248 // of (photon energy, intensity) pairs
249 G4double currentIN = (*wlsVector)[0];
250 if(currentIN >= 0.0)
251 {
252 // Create first (photon energy)
253 G4double currentPM = wlsVector->Energy(0);
254 G4double currentCII = 0.0;
255 physVector->InsertValues(currentPM, currentCII);
256
257 // Set previous values to current ones prior to loop
258 G4double prevPM = currentPM;
259 G4double prevCII = currentCII;
260 G4double prevIN = currentIN;
261
262 // loop over all (photon energy, intensity)
263 // pairs stored for this material
264 for(std::size_t j = 1; j < wlsVector->GetVectorLength(); ++j)
265 {
266 currentPM = wlsVector->Energy(j);
267 currentIN = (*wlsVector)[j];
268 currentCII =
269 prevCII + 0.5 * (currentPM - prevPM) * (prevIN + currentIN);
270
271 physVector->InsertValues(currentPM, currentCII);
272
273 prevPM = currentPM;
274 prevCII = currentCII;
275 prevIN = currentIN;
276 }
277 }
278 }
279 }
280 theIntegralTable->insertAt(i, physVector);
281 }
282}
283
284//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
287{
288 G4double thePhotonEnergy = aTrack.GetDynamicParticle()->GetTotalEnergy();
289 G4double attLength = DBL_MAX;
292
293 if(MPT)
294 {
296 if(attVector)
297 {
298 attLength = attVector->Value(thePhotonEnergy, idx_wls2);
299 }
300 }
301 return attLength;
302}
303
304//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
306{
308 {
310 WLSTimeGeneratorProfile = nullptr;
311 }
312 if(name == "delta")
313 {
315 }
316 else if(name == "exponential")
317 {
319 new G4WLSTimeGeneratorProfileExponential("exponential");
320 }
321 else
322 {
323 G4Exception("G4OpWLS::UseTimeProfile", "em0202", FatalException,
324 "generator does not exist");
325 }
327}
328
329//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
G4ForceCondition
G4PhysicsFreeVector G4MaterialPropertyVector
std::vector< G4Material * > G4MaterialTable
G4ThreeVector G4ParticleMomentum
G4long G4Poisson(G4double mean)
Definition G4Poisson.hh:50
G4ProcessType
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
Hep3Vector cross(const Hep3Vector &) const
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
G4bool ConstPropertyExists(const G4String &key) const
G4double GetConstProperty(const G4String &key) const
G4MaterialPropertyVector * GetProperty(const char *key) const
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
static std::size_t GetNumberOfMaterials()
std::size_t GetIndex() const
static G4MaterialTable * GetMaterialTable()
void SetVerboseLevel(G4int)
Definition G4OpWLS2.cc:330
virtual void Initialise()
Definition G4OpWLS2.cc:84
virtual void UseTimeProfile(const G4String name)
Definition G4OpWLS2.cc:305
virtual ~G4OpWLS2()
Definition G4OpWLS2.cc:67
virtual void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
Definition G4OpWLS2.cc:220
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
Definition G4OpWLS2.cc:78
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
Definition G4OpWLS2.cc:92
G4PhysicsTable * theIntegralTable
Definition G4OpWLS2.hh:91
G4VWLSTimeGeneratorProfile * WLSTimeGeneratorProfile
Definition G4OpWLS2.hh:90
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
Definition G4OpWLS2.cc:285
G4OpWLS2(const G4String &processName="OpWLS2", G4ProcessType type=fOptical)
Definition G4OpWLS2.cc:54
G4int GetWLS2VerboseLevel() const
G4String GetWLS2TimeProfile() const
void SetWLS2TimeProfile(const G4String &)
static G4OpticalParameters * Instance()
static G4OpticalPhoton * OpticalPhoton()
G4double GetEnergy(const G4double value) const
G4double GetMaxEnergy() const
G4double GetMaxValue() const
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
G4double GetGlobalTime() const
const G4ThreeVector & GetPosition() const
G4StepPoint * GetPostStepPoint() const
G4int GetTrackID() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
void SetParentID(const G4int aValue)
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4VDiscreteProcess(const G4String &aName, G4ProcessType aType=fNotDefined)
G4ParticleChange aParticleChange
G4int verboseLevel
void SetProcessSubType(G4int)
const G4String & GetProcessName() const
#define DBL_MAX
Definition templates.hh:62