Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
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: G4OpWLS.cc 86052 2014-11-07 08:31:04Z gcosmo $
28 //
30 // Optical Photon WaveLength Shifting (WLS) Class Implementation
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 //
45 
46 #include "G4OpWLS.hh"
47 
48 #include "G4ios.hh"
49 #include "G4PhysicalConstants.hh"
50 #include "G4SystemOfUnits.hh"
51 #include "G4OpProcessSubType.hh"
52 
55 
57 // Class Implementation
59 
61  // static data members
63 
65 // Constructors
67 
68 G4OpWLS::G4OpWLS(const G4String& processName, G4ProcessType type)
69  : G4VDiscreteProcess(processName, type)
70 {
72 
73  theIntegralTable = NULL;
74 
76  new G4WLSTimeGeneratorProfileDelta("WLSTimeGeneratorProfileDelta");
77 
78  if (verboseLevel>0) G4cout << GetProcessName() << " is created " << G4endl;
79 }
80 
82 // Destructors
84 
86 {
87  if (theIntegralTable) {
89  delete theIntegralTable;
90  }
92 }
93 
95 // Methods
97 
98 // PostStepDoIt
99 // -------------
100 //
102 G4OpWLS::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
103 {
104  aParticleChange.Initialize(aTrack);
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  G4double primaryEnergy = aTrack.GetDynamicParticle()->GetKineticEnergy();
151 
152  G4int materialIndex = aMaterial->GetIndex();
153 
154  // Retrieve the WLS Integral for this material
155  // new G4PhysicsOrderedFreeVector allocated to hold CII's
156 
157  G4double WLSTime = 0.*ns;
158  G4PhysicsOrderedFreeVector* WLSIntegral = 0;
159 
160  WLSTime = aMaterialPropertiesTable->
161  GetConstProperty("WLSTIMECONSTANT");
162  WLSIntegral =
163  (G4PhysicsOrderedFreeVector*)((*theIntegralTable)(materialIndex));
164 
165  // Max WLS Integral
166 
167  G4double CIImax = WLSIntegral->GetMaxValue();
168 
169  G4int NumberOfPhotons = NumPhotons;
170 
171  for (G4int i = 0; i < NumPhotons; i++) {
172 
173  G4double sampledEnergy;
174 
175  // Make sure the energy of the secondary is less than that of the primary
176 
177  for (G4int j = 1; j <= 100; j++) {
178 
179  // Determine photon energy
180 
181  G4double CIIvalue = G4UniformRand()*CIImax;
182  sampledEnergy = WLSIntegral->GetEnergy(CIIvalue);
183 
184  if (verboseLevel>1) {
185  G4cout << "sampledEnergy = " << sampledEnergy << G4endl;
186  G4cout << "CIIvalue = " << CIIvalue << G4endl;
187  }
188 
189  if (sampledEnergy <= primaryEnergy) break;
190  }
191 
192  // If no such energy can be sampled, return one less secondary, or none
193 
194  if (sampledEnergy > primaryEnergy) {
195  if (verboseLevel>1)
196  G4cout << " *** One less WLS photon will be returned ***" << G4endl;
197  NumberOfPhotons--;
198  aParticleChange.SetNumberOfSecondaries(NumberOfPhotons);
199  if (NumberOfPhotons == 0) {
200  if (verboseLevel>1)
201  G4cout << " *** No WLS photon can be sampled for this primary ***"
202  << G4endl;
203  // return unchanged particle and no secondaries
204  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
205  }
206  continue;
207  }
208 
209  // Generate random photon direction
210 
211  G4double cost = 1. - 2.*G4UniformRand();
212  G4double sint = std::sqrt((1.-cost)*(1.+cost));
213 
214  G4double phi = twopi*G4UniformRand();
215  G4double sinp = std::sin(phi);
216  G4double cosp = std::cos(phi);
217 
218  G4double px = sint*cosp;
219  G4double py = sint*sinp;
220  G4double pz = cost;
221 
222  // Create photon momentum direction vector
223 
224  G4ParticleMomentum photonMomentum(px, py, pz);
225 
226  // Determine polarization of new photon
227 
228  G4double sx = cost*cosp;
229  G4double sy = cost*sinp;
230  G4double sz = -sint;
231 
232  G4ThreeVector photonPolarization(sx, sy, sz);
233 
234  G4ThreeVector perp = photonMomentum.cross(photonPolarization);
235 
236  phi = twopi*G4UniformRand();
237  sinp = std::sin(phi);
238  cosp = std::cos(phi);
239 
240  photonPolarization = cosp * photonPolarization + sinp * perp;
241 
242  photonPolarization = photonPolarization.unit();
243 
244  // Generate a new photon:
245 
246  G4DynamicParticle* aWLSPhoton =
248  photonMomentum);
249  aWLSPhoton->SetPolarization
250  (photonPolarization.x(),
251  photonPolarization.y(),
252  photonPolarization.z());
253 
254  aWLSPhoton->SetKineticEnergy(sampledEnergy);
255 
256  // Generate new G4Track object:
257 
258  // Must give position of WLS optical photon
259 
260  G4double TimeDelay = WLSTimeGeneratorProfile->GenerateTime(WLSTime);
261  G4double aSecondaryTime = (pPostStepPoint->GetGlobalTime()) + TimeDelay;
262 
263  G4ThreeVector aSecondaryPosition = pPostStepPoint->GetPosition();
264 
265  G4Track* aSecondaryTrack =
266  new G4Track(aWLSPhoton,aSecondaryTime,aSecondaryPosition);
267 
268  aSecondaryTrack->SetTouchableHandle(aTrack.GetTouchableHandle());
269  // aSecondaryTrack->SetTouchableHandle((G4VTouchable*)0);
270 
271  aSecondaryTrack->SetParentID(aTrack.GetTrackID());
272 
273  aParticleChange.AddSecondary(aSecondaryTrack);
274  }
275 
276  if (verboseLevel>0) {
277  G4cout << "\n Exiting from G4OpWLS::DoIt -- NumberOfSecondaries = "
279  }
280 
281  return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
282 }
283 
284 // BuildPhysicsTable for the wavelength shifting process
285 // --------------------------------------------------
286 
288 {
289  if (theIntegralTable) {
291  delete theIntegralTable;
292  theIntegralTable = NULL;
293  }
294 
295  const G4MaterialTable* theMaterialTable =
297  G4int numOfMaterials = G4Material::GetNumberOfMaterials();
298 
299  // create new physics table
300 
301  theIntegralTable = new G4PhysicsTable(numOfMaterials);
302 
303  // loop for materials
304 
305  for (G4int i=0 ; i < numOfMaterials; i++)
306  {
307  G4PhysicsOrderedFreeVector* aPhysicsOrderedFreeVector =
309 
310  // Retrieve vector of WLS wavelength intensity for
311  // the material from the material's optical properties table.
312 
313  G4Material* aMaterial = (*theMaterialTable)[i];
314 
315  G4MaterialPropertiesTable* aMaterialPropertiesTable =
316  aMaterial->GetMaterialPropertiesTable();
317 
318  if (aMaterialPropertiesTable) {
319 
320  G4MaterialPropertyVector* theWLSVector =
321  aMaterialPropertiesTable->GetProperty("WLSCOMPONENT");
322 
323  if (theWLSVector) {
324 
325  // Retrieve the first intensity point in vector
326  // of (photon energy, intensity) pairs
327 
328  G4double currentIN = (*theWLSVector)[0];
329 
330  if (currentIN >= 0.0) {
331 
332  // Create first (photon energy)
333 
334  G4double currentPM = theWLSVector->Energy(0);
335 
336  G4double currentCII = 0.0;
337 
338  aPhysicsOrderedFreeVector->
339  InsertValues(currentPM , currentCII);
340 
341  // Set previous values to current ones prior to loop
342 
343  G4double prevPM = currentPM;
344  G4double prevCII = currentCII;
345  G4double prevIN = currentIN;
346 
347  // loop over all (photon energy, intensity)
348  // pairs stored for this material
349 
350  for (size_t j = 1;
351  j < theWLSVector->GetVectorLength();
352  j++)
353  {
354  currentPM = theWLSVector->Energy(j);
355  currentIN = (*theWLSVector)[j];
356 
357  currentCII = 0.5 * (prevIN + currentIN);
358 
359  currentCII = prevCII +
360  (currentPM - prevPM) * currentCII;
361 
362  aPhysicsOrderedFreeVector->
363  InsertValues(currentPM, currentCII);
364 
365  prevPM = currentPM;
366  prevCII = currentCII;
367  prevIN = currentIN;
368  }
369  }
370  }
371  }
372  // The WLS integral for a given material
373  // will be inserted in the table according to the
374  // position of the material in the material table.
375 
376  theIntegralTable->insertAt(i,aPhysicsOrderedFreeVector);
377  }
378 }
379 
380 // GetMeanFreePath
381 // ---------------
382 //
384  G4double ,
386 {
387  const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
388  const G4Material* aMaterial = aTrack.GetMaterial();
389 
390  G4double thePhotonEnergy = aParticle->GetTotalEnergy();
391 
392  G4MaterialPropertiesTable* aMaterialPropertyTable;
393  G4MaterialPropertyVector* AttenuationLengthVector;
394 
395  G4double AttenuationLength = DBL_MAX;
396 
397  aMaterialPropertyTable = aMaterial->GetMaterialPropertiesTable();
398 
399  if ( aMaterialPropertyTable ) {
400  AttenuationLengthVector = aMaterialPropertyTable->
401  GetProperty("WLSABSLENGTH");
402  if ( AttenuationLengthVector ){
403  AttenuationLength = AttenuationLengthVector->
404  Value(thePhotonEnergy);
405  }
406  else {
407  // G4cout << "No WLS absorption length specified" << G4endl;
408  }
409  }
410  else {
411  // G4cout << "No WLS absortion length specified" << G4endl;
412  }
413 
414  return AttenuationLength;
415 }
416 
418 {
419  if (name == "delta")
420  {
423  new G4WLSTimeGeneratorProfileDelta("delta");
424  }
425  else if (name == "exponential")
426  {
429  new G4WLSTimeGeneratorProfileExponential("exponential");
430  }
431  else
432  {
433  G4Exception("G4OpWLS::UseTimeProfile", "em0202",
435  "generator does not exist");
436  }
437 }
const XML_Char * name
Definition: expat.h:151
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
G4int GetNumberOfSecondaries() const
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double GetKineticEnergy() const
G4double GetTotalEnergy() const
double x() const
const G4DynamicParticle * GetDynamicParticle() const
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
Definition: G4OpWLS.cc:102
size_t GetIndex() const
Definition: G4Material.hh:262
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:587
std::vector< G4Material * > G4MaterialTable
void SetTouchableHandle(const G4TouchableHandle &apValue)
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
double z() const
static constexpr double twopi
Definition: G4SIunits.hh:76
void UseTimeProfile(const G4String name)
Definition: G4OpWLS.cc:417
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
virtual G4double GenerateTime(const G4double time_constant)=0
const G4ThreeVector & GetPosition() const
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
Definition: G4OpWLS.cc:287
G4VWLSTimeGeneratorProfile * WLSTimeGeneratorProfile
Definition: G4OpWLS.hh:137
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
~G4OpWLS()
Definition: G4OpWLS.cc:85
Definition: G4Step.hh:76
G4int GetTrackID() const
void SetPolarization(G4double polX, G4double polY, G4double polZ)
G4double Energy(size_t index) const
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:594
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
void SetKineticEnergy(G4double aEnergy)
const G4TouchableHandle & GetTouchableHandle() const
G4double GetEnergy(G4double aValue)
G4Material * GetMaterial() const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4bool ConstPropertyExists(const char *key) const
static G4OpticalPhoton * OpticalPhoton()
virtual void Initialize(const G4Track &)
G4MaterialPropertiesTable * GetMaterialPropertiesTable() const
Definition: G4Material.hh:252
void SetNumberOfSecondaries(G4int totSecondaries)
Hep3Vector unit() const
void SetParentID(const G4int aValue)
G4StepPoint * GetPostStepPoint() const
double y() const
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
G4OpWLS(const G4String &processName="OpWLS", G4ProcessType type=fOptical)
Definition: G4OpWLS.cc:68
void AddSecondary(G4Track *aSecondary)
void insertAt(size_t, G4PhysicsVector *)
#define G4endl
Definition: G4ios.hh:61
G4double GetGlobalTime() const
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ForceCondition
G4PhysicsTable * theIntegralTable
Definition: G4OpWLS.hh:138
#define DBL_MAX
Definition: templates.hh:83
#define ns
Definition: xmlparse.cc:614
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
G4MaterialPropertyVector * GetProperty(const char *key)
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
Definition: G4OpWLS.cc:383
void clearAndDestroy()
G4ProcessType