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.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// $Id$
28//
29////////////////////////////////////////////////////////////////////////
30// Optical Photon WaveLength Shifting (WLS) Class Implementation
31////////////////////////////////////////////////////////////////////////
32//
33// File: G4OpWLS.cc
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#include "G4OpWLS.hh"
47
48#include "G4ios.hh"
50#include "G4SystemOfUnits.hh"
51#include "G4OpProcessSubType.hh"
52
55
56/////////////////////////
57// Class Implementation
58/////////////////////////
59
60/////////////////
61// Constructors
62/////////////////
63
64G4OpWLS::G4OpWLS(const G4String& processName, G4ProcessType type)
65 : G4VDiscreteProcess(processName, type)
66{
68
70
71 if (verboseLevel>0) {
72 G4cout << GetProcessName() << " is created " << G4endl;
73 }
74
76 new G4WLSTimeGeneratorProfileDelta("WLSTimeGeneratorProfileDelta");
77
78 BuildThePhysicsTable();
79}
80
81////////////////
82// Destructors
83////////////////
84
86{
87 if (theIntegralTable != 0) {
89 delete theIntegralTable;
90 }
92}
93
94////////////
95// Methods
96////////////
97
98// PostStepDoIt
99// -------------
100//
102G4OpWLS::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
103{
105
107
108 if (verboseLevel>0) {
109 G4cout << "\n** Photon absorbed! **" << G4endl;
110 }
111
112 const G4Material* aMaterial = aTrack.GetMaterial();
113
114 G4StepPoint* pPostStepPoint = aStep.GetPostStepPoint();
115
116 G4MaterialPropertiesTable* aMaterialPropertiesTable =
117 aMaterial->GetMaterialPropertiesTable();
118 if (!aMaterialPropertiesTable)
119 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
120
121 const G4MaterialPropertyVector* WLS_Intensity =
122 aMaterialPropertiesTable->GetProperty("WLSCOMPONENT");
123
124 if (!WLS_Intensity)
125 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
126
127 G4int NumPhotons = 1;
128
129 if (aMaterialPropertiesTable->ConstPropertyExists("WLSMEANNUMBERPHOTONS")) {
130
131 G4double MeanNumberOfPhotons = aMaterialPropertiesTable->
132 GetConstProperty("WLSMEANNUMBERPHOTONS");
133
134 NumPhotons = G4int(G4Poisson(MeanNumberOfPhotons));
135
136 if (NumPhotons <= 0) {
137
138 // return unchanged particle and no secondaries
139
141
142 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
143
144 }
145
146 }
147
149
150 G4int materialIndex = aMaterial->GetIndex();
151
152 // Retrieve the WLS Integral for this material
153 // new G4PhysicsOrderedFreeVector allocated to hold CII's
154
155 G4double WLSTime = 0.*ns;
156 G4PhysicsOrderedFreeVector* WLSIntegral = 0;
157
158 WLSTime = aMaterialPropertiesTable->
159 GetConstProperty("WLSTIMECONSTANT");
160 WLSIntegral =
161 (G4PhysicsOrderedFreeVector*)((*theIntegralTable)(materialIndex));
162
163 // Max WLS Integral
164
165 G4double CIImax = WLSIntegral->GetMaxValue();
166
167 for (G4int i = 0; i < NumPhotons; i++) {
168
169 // Determine photon energy
170
171 G4double CIIvalue = G4UniformRand()*CIImax;
172 G4double sampledEnergy =
173 WLSIntegral->GetEnergy(CIIvalue);
174
175 if (verboseLevel>1) {
176 G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
177 G4cout << "CIIvalue = " << CIIvalue << G4endl;
178 }
179
180 // Generate random photon direction
181
182 G4double cost = 1. - 2.*G4UniformRand();
183 G4double sint = std::sqrt((1.-cost)*(1.+cost));
184
185 G4double phi = twopi*G4UniformRand();
186 G4double sinp = std::sin(phi);
187 G4double cosp = std::cos(phi);
188
189 G4double px = sint*cosp;
190 G4double py = sint*sinp;
191 G4double pz = cost;
192
193 // Create photon momentum direction vector
194
195 G4ParticleMomentum photonMomentum(px, py, pz);
196
197 // Determine polarization of new photon
198
199 G4double sx = cost*cosp;
200 G4double sy = cost*sinp;
201 G4double sz = -sint;
202
203 G4ThreeVector photonPolarization(sx, sy, sz);
204
205 G4ThreeVector perp = photonMomentum.cross(photonPolarization);
206
207 phi = twopi*G4UniformRand();
208 sinp = std::sin(phi);
209 cosp = std::cos(phi);
210
211 photonPolarization = cosp * photonPolarization + sinp * perp;
212
213 photonPolarization = photonPolarization.unit();
214
215 // Generate a new photon:
216
217 G4DynamicParticle* aWLSPhoton =
219 photonMomentum);
220 aWLSPhoton->SetPolarization
221 (photonPolarization.x(),
222 photonPolarization.y(),
223 photonPolarization.z());
224
225 aWLSPhoton->SetKineticEnergy(sampledEnergy);
226
227 // Generate new G4Track object:
228
229 // Must give position of WLS optical photon
230
231 G4double TimeDelay = WLSTimeGeneratorProfile->GenerateTime(WLSTime);
232 G4double aSecondaryTime = (pPostStepPoint->GetGlobalTime()) + TimeDelay;
233
234 G4ThreeVector aSecondaryPosition = pPostStepPoint->GetPosition();
235
236 G4Track* aSecondaryTrack =
237 new G4Track(aWLSPhoton,aSecondaryTime,aSecondaryPosition);
238
239 aSecondaryTrack->SetTouchableHandle(aTrack.GetTouchableHandle());
240 // aSecondaryTrack->SetTouchableHandle((G4VTouchable*)0);
241
242 aSecondaryTrack->SetParentID(aTrack.GetTrackID());
243
244 aParticleChange.AddSecondary(aSecondaryTrack);
245 }
246
247 if (verboseLevel>0) {
248 G4cout << "\n Exiting from G4OpWLS::DoIt -- NumberOfSecondaries = "
250 }
251
252 return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
253}
254
255// BuildThePhysicsTable for the wavelength shifting process
256// --------------------------------------------------
257//
258
259void G4OpWLS::BuildThePhysicsTable()
260{
261 if (theIntegralTable) return;
262
263 const G4MaterialTable* theMaterialTable =
265 G4int numOfMaterials = G4Material::GetNumberOfMaterials();
266
267 // create new physics table
268
269 if(!theIntegralTable)theIntegralTable = new G4PhysicsTable(numOfMaterials);
270
271 // loop for materials
272
273 for (G4int i=0 ; i < numOfMaterials; i++)
274 {
275 G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector =
277
278 // Retrieve vector of WLS wavelength intensity for
279 // the material from the material's optical properties table.
280
281 G4Material* aMaterial = (*theMaterialTable)[i];
282
283 G4MaterialPropertiesTable* aMaterialPropertiesTable =
284 aMaterial->GetMaterialPropertiesTable();
285
286 if (aMaterialPropertiesTable) {
287
288 G4MaterialPropertyVector* theWLSVector =
289 aMaterialPropertiesTable->GetProperty("WLSCOMPONENT");
290
291 if (theWLSVector) {
292
293 // Retrieve the first intensity point in vector
294 // of (photon energy, intensity) pairs
295
296 G4double currentIN = (*theWLSVector)[0];
297
298 if (currentIN >= 0.0) {
299
300 // Create first (photon energy)
301
302 G4double currentPM = theWLSVector->Energy(0);
303
304 G4double currentCII = 0.0;
305
306 aPhysicsOrderedFreeVector->
307 InsertValues(currentPM , currentCII);
308
309 // Set previous values to current ones prior to loop
310
311 G4double prevPM = currentPM;
312 G4double prevCII = currentCII;
313 G4double prevIN = currentIN;
314
315 // loop over all (photon energy, intensity)
316 // pairs stored for this material
317
318 for (size_t j = 1;
319 j < theWLSVector->GetVectorLength();
320 j++)
321 {
322 currentPM = theWLSVector->Energy(j);
323 currentIN = (*theWLSVector)[j];
324
325 currentCII = 0.5 * (prevIN + currentIN);
326
327 currentCII = prevCII +
328 (currentPM - prevPM) * currentCII;
329
330 aPhysicsOrderedFreeVector->
331 InsertValues(currentPM, currentCII);
332
333 prevPM = currentPM;
334 prevCII = currentCII;
335 prevIN = currentIN;
336 }
337 }
338 }
339 }
340 // The WLS integral for a given material
341 // will be inserted in the table according to the
342 // position of the material in the material table.
343
344 theIntegralTable->insertAt(i,aPhysicsOrderedFreeVector);
345 }
346}
347
348// GetMeanFreePath
349// ---------------
350//
352 G4double ,
354{
355 const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
356 const G4Material* aMaterial = aTrack.GetMaterial();
357
358 G4double thePhotonEnergy = aParticle->GetTotalEnergy();
359
360 G4MaterialPropertiesTable* aMaterialPropertyTable;
361 G4MaterialPropertyVector* AttenuationLengthVector;
362
363 G4double AttenuationLength = DBL_MAX;
364
365 aMaterialPropertyTable = aMaterial->GetMaterialPropertiesTable();
366
367 if ( aMaterialPropertyTable ) {
368 AttenuationLengthVector = aMaterialPropertyTable->
369 GetProperty("WLSABSLENGTH");
370 if ( AttenuationLengthVector ){
371 AttenuationLength = AttenuationLengthVector->
372 Value(thePhotonEnergy);
373 }
374 else {
375 // G4cout << "No WLS absorption length specified" << G4endl;
376 }
377 }
378 else {
379 // G4cout << "No WLS absortion length specified" << G4endl;
380 }
381
382 return AttenuationLength;
383}
384
386{
387 if (name == "delta")
388 {
392 }
393 else if (name == "exponential")
394 {
397 new G4WLSTimeGeneratorProfileExponential("exponential");
398 }
399 else
400 {
401 G4Exception("G4OpWLS::UseTimeProfile", "em0202",
403 "generator does not exist");
404 }
405}
@ FatalException
G4ForceCondition
std::vector< G4Material * > G4MaterialTable
@ fOpWLS
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:50
G4ProcessType
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
double z() const
Hep3Vector unit() const
double x() const
double y() const
Hep3Vector cross(const Hep3Vector &) const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4double GetTotalEnergy() const
void SetKineticEnergy(G4double aEnergy)
G4MaterialPropertyVector * GetProperty(const char *key)
G4bool ConstPropertyExists(const char *key)
static const G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:562
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:251
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:569
size_t GetIndex() const
Definition: G4Material.hh:261
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4OpWLS.cc:351
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
G4OpWLS(const G4String &processName="OpWLS", G4ProcessType type=fOptical)
Definition: G4OpWLS.cc:64
void UseTimeProfile(const G4String name)
Definition: G4OpWLS.cc:385
static G4OpticalPhoton * OpticalPhoton()
void AddSecondary(G4Track *aSecondary)
virtual void Initialize(const G4Track &)
G4double GetEnergy(G4double aValue)
void insertAt(size_t, G4PhysicsVector *)
void clearAndDestroy()
size_t GetVectorLength() const
G4double Energy(size_t index) const
G4double GetGlobalTime() const
const G4ThreeVector & GetPosition() const
Definition: G4Step.hh:78
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 &)
void ProposeTrackStatus(G4TrackStatus status)
G4int GetNumberOfSecondaries() const
void SetNumberOfSecondaries(G4int totSecondaries)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4int verboseLevel
Definition: G4VProcess.hh:368
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
virtual G4double GenerateTime(const G4double time_constant)=0
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define DBL_MAX
Definition: templates.hh:83
#define ns
Definition: xmlparse.cc:597