Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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//
28////////////////////////////////////////////////////////////////////////
29// Optical Photon WaveLength Shifting (WLS) Class Definition
30////////////////////////////////////////////////////////////////////////
31//
32// File: G4OpWLS.hh
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#ifndef G4OpWLS_h
44#define G4OpWLS_h 1
45
46#include "G4VDiscreteProcess.hh"
47#include "G4OpticalPhoton.hh"
48
50
52{
53 public:
54 explicit G4OpWLS(const G4String& processName = "OpWLS",
55 G4ProcessType type = fOptical);
56 virtual ~G4OpWLS();
57
58 virtual G4bool IsApplicable(
59 const G4ParticleDefinition& aParticleType) override;
60 // Returns true -> 'is applicable' only for an optical photon.
61
62 virtual void BuildPhysicsTable(
63 const G4ParticleDefinition& aParticleType) override;
64 // Build the WLS integral table at the right time
65
66 virtual G4double GetMeanFreePath(const G4Track& aTrack, G4double,
67 G4ForceCondition*) override;
68 // Returns the absorption length for WLS absorption of optical
69 // photons in media with a specified attenuation length.
70
71 virtual G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
72 const G4Step& aStep) override;
73 // This is the method implementing WLS for optical photons.
74
75 virtual G4PhysicsTable* GetIntegralTable() const;
76 // Returns the address of the WLS integral table.
77
78 virtual void DumpPhysicsTable() const;
79 // Prints the WLS integral table.
80
81 virtual void UseTimeProfile(const G4String name);
82 // Selects the time profile generator
83
84 virtual void PreparePhysicsTable(const G4ParticleDefinition&) override;
85 virtual void Initialise();
86
88
89 protected:
92
93 private:
94 G4OpWLS(const G4OpWLS& right) = delete;
95 G4OpWLS& operator=(const G4OpWLS& right) = delete;
96
97 std::size_t idx_wls = 0;
98};
99
100////////////////////
101// Inline methods
102////////////////////
103
105{
106 return (&aParticleType == G4OpticalPhoton::OpticalPhoton());
107}
108
110{
111 return theIntegralTable;
112}
113
114inline void G4OpWLS::DumpPhysicsTable() const
115{
116 std::size_t PhysicsTableSize = theIntegralTable->entries();
118
119 for(std::size_t i = 0; i < PhysicsTableSize; ++i)
120 {
122 v->DumpValues();
123 }
124}
125
126#endif /* G4OpWLS_h */
G4ForceCondition
G4ProcessType
@ fOptical
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual void BuildPhysicsTable(const G4ParticleDefinition &aParticleType) override
Definition: G4OpWLS.cc:229
virtual void Initialise()
Definition: G4OpWLS.cc:81
virtual G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
Definition: G4OpWLS.hh:104
virtual G4PhysicsTable * GetIntegralTable() const
Definition: G4OpWLS.hh:109
void SetVerboseLevel(G4int)
Definition: G4OpWLS.cc:339
virtual void DumpPhysicsTable() const
Definition: G4OpWLS.hh:114
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *) override
Definition: G4OpWLS.cc:294
G4PhysicsTable * theIntegralTable
Definition: G4OpWLS.hh:91
virtual G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
Definition: G4OpWLS.cc:89
virtual ~G4OpWLS()
Definition: G4OpWLS.cc:67
G4VWLSTimeGeneratorProfile * WLSTimeGeneratorProfile
Definition: G4OpWLS.hh:90
virtual void UseTimeProfile(const G4String name)
Definition: G4OpWLS.cc:314
virtual void PreparePhysicsTable(const G4ParticleDefinition &) override
Definition: G4OpWLS.cc:78
static G4OpticalPhoton * OpticalPhoton()
std::size_t entries() const
void DumpValues(G4double unitE=1.0, G4double unitV=1.0) const
Definition: G4Step.hh:62