Geant4  10.02.p01
G4OpRayleigh.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: G4OpRayleigh.hh 84596 2014-10-17 07:34:41Z gcosmo $
28 //
29 //
31 // Optical Photon Rayleigh Scattering Class Definition
33 //
34 // File: G4OpRayleigh.hh
35 // Description: Discrete Process -- Rayleigh scattering of optical photons
36 // Version: 1.0
37 // Created: 1996-05-31
38 // Author: Juliet Armstrong
39 // Updated: 2014-08-20 allow for more material types
40 // 2005-07-28 add G4ProcessType to constructor
41 // 1999-10-29 add method and class descriptors
42 // 1997-04-09 by Peter Gumplinger
43 // > new physics/tracking scheme
44 // mail: gum@triumf.ca
45 //
47 
48 #ifndef G4OpRayleigh_h
49 #define G4OpRayleigh_h 1
50 
52 // Includes
54 
55 #include "globals.hh"
56 #include "templates.hh"
57 #include "Randomize.hh"
58 #include "G4ThreeVector.hh"
59 #include "G4ParticleMomentum.hh"
60 #include "G4Step.hh"
61 #include "G4VDiscreteProcess.hh"
62 #include "G4DynamicParticle.hh"
63 #include "G4Material.hh"
64 #include "G4OpticalPhoton.hh"
65 #include "G4PhysicsTable.hh"
67 
68 // Class Description:
69 // Discrete Process -- Rayleigh scattering of optical photons.
70 // Class inherits publicly from G4VDiscreteProcess.
71 // Class Description - End:
72 
74 // Class Definition
76 
78 {
79 
80 public:
81 
83  // Constructors and Destructor
85 
86  G4OpRayleigh(const G4String& processName = "OpRayleigh",
87  G4ProcessType type = fOptical);
88  ~G4OpRayleigh();
89 
90 private:
91 
93 
95  // Operators
97 
98  G4OpRayleigh& operator=(const G4OpRayleigh &right);
99 
100 public:
101 
103  // Methods
105 
106  G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
107  // Returns true -> 'is applicable' only for an optical photon.
108 
109  void BuildPhysicsTable(const G4ParticleDefinition& aParticleType);
110  // Build thePhysicsTable at a right time
111 
112  G4double GetMeanFreePath(const G4Track& aTrack,
113  G4double ,
114  G4ForceCondition* );
115  // Returns the mean free path for Rayleigh scattering in water.
116  // --- Not yet implemented for other materials! ---
117 
118  G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
119  const G4Step& aStep);
120  // This is the method implementing Rayleigh scattering.
121 
123  // Returns the address of the physics table.
124 
125  void DumpPhysicsTable() const;
126  // Prints the physics table.
127 
128 private:
129 
131  // Helper Functions
133 
140  CalculateRayleighMeanFreePaths( const G4Material* material ) const;
141 
143  // Class Data Members
145 
146 protected:
147 
149  // A Physics Table can be either a cross-sections table or
150  // an energy table (or can be used for other specific
151  // purposes).
152 
153 private:
154 };
155 
157 // Inline methods
159 
160 inline
162 {
163  return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
164 }
165 
166 inline
168 
169 {
170  G4int PhysicsTableSize = thePhysicsTable->entries();
172 
173  for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
174  {
176  v->DumpValues();
177  }
178 }
179 
181 {
182  return thePhysicsTable;
183 }
184 
185 
186 #endif /* G4OpRayleigh_h */
G4OpRayleigh & operator=(const G4OpRayleigh &right)
void DumpPhysicsTable() const
int G4int
Definition: G4Types.hh:78
G4PhysicsTable * thePhysicsTable
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
bool G4bool
Definition: G4Types.hh:79
Definition: G4Step.hh:76
static G4OpticalPhoton * OpticalPhoton()
G4OpRayleigh(const G4String &processName="OpRayleigh", G4ProcessType type=fOptical)
Definition: G4OpRayleigh.cc:88
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
double G4double
Definition: G4Types.hh:76
G4PhysicsOrderedFreeVector * CalculateRayleighMeanFreePaths(const G4Material *material) const
Calculates the mean free paths for a material as a function of photon energy.
G4ForceCondition
size_t entries() const
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
G4PhysicsTable * GetPhysicsTable() const
G4ProcessType