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
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$
28//
29//
30////////////////////////////////////////////////////////////////////////
31// Optical Photon Rayleigh Scattering Class Definition
32////////////////////////////////////////////////////////////////////////
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: 2005-07-28 add G4ProcessType to constructor
40// 1999-10-29 add method and class descriptors
41// 1997-04-09 by Peter Gumplinger
42// > new physics/tracking scheme
43// mail: gum@triumf.ca
44//
45////////////////////////////////////////////////////////////////////////
46
47#ifndef G4OpRayleigh_h
48#define G4OpRayleigh_h 1
49
50/////////////
51// Includes
52/////////////
53
54#include "globals.hh"
55#include "templates.hh"
56#include "Randomize.hh"
57#include "G4ThreeVector.hh"
58#include "G4ParticleMomentum.hh"
59#include "G4Step.hh"
60#include "G4VDiscreteProcess.hh"
61#include "G4DynamicParticle.hh"
62#include "G4Material.hh"
63#include "G4OpticalPhoton.hh"
64#include "G4PhysicsTable.hh"
66
67// Class Description:
68// Discrete Process -- Rayleigh scattering of optical photons.
69// Class inherits publicly from G4VDiscreteProcess.
70// Class Description - End:
71
72/////////////////////
73// Class Definition
74/////////////////////
75
77{
78
79public:
80
81 ////////////////////////////////
82 // Constructors and Destructor
83 ////////////////////////////////
84
85 G4OpRayleigh(const G4String& processName = "OpRayleigh",
86 G4ProcessType type = fOptical);
88
89private:
90
91 G4OpRayleigh(const G4OpRayleigh &right);
92
93 //////////////
94 // Operators
95 //////////////
96
97 G4OpRayleigh& operator=(const G4OpRayleigh &right);
98
99public:
100
101 ////////////
102 // Methods
103 ////////////
104
105 G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
106 // Returns true -> 'is applicable' only for an optical photon.
107
108 G4double GetMeanFreePath(const G4Track& aTrack,
109 G4double ,
111 // Returns the mean free path for Rayleigh scattering in water.
112 // --- Not yet implemented for other materials! ---
113
115 const G4Step& aStep);
116 // This is the method implementing Rayleigh scattering.
117
119 // Returns the address of the physics table.
120
121 void DumpPhysicsTable() const;
122 // Prints the physics table.
123
124private:
125
126 void BuildThePhysicsTable();
127
128 /////////////////////
129 // Helper Functions
130 /////////////////////
131
132 G4PhysicsOrderedFreeVector* RayleighAttenuationLengthGenerator(
134
135 ///////////////////////
136 // Class Data Members
137 ///////////////////////
138
139protected:
140
142 // A Physics Table can be either a cross-sections table or
143 // an energy table (or can be used for other specific
144 // purposes).
145
146private:
147
148 G4bool DefaultWater;
149
150};
151
152////////////////////
153// Inline methods
154////////////////////
155
156inline
158{
159 return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
160}
161
162inline
164
165{
166 G4int PhysicsTableSize = thePhysicsTable->entries();
168
169 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
170 {
172 v->DumpValues();
173 }
174}
175
177{
178 return thePhysicsTable;
179}
180
181
182#endif /* G4OpRayleigh_h */
G4ForceCondition
G4ProcessType
@ fOptical
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
void DumpPhysicsTable() const
G4PhysicsTable * GetPhysicsTable() const
G4PhysicsTable * thePhysicsTable
G4double GetMeanFreePath(const G4Track &aTrack, G4double, G4ForceCondition *)
G4bool IsApplicable(const G4ParticleDefinition &aParticleType)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep)
static G4OpticalPhoton * OpticalPhoton()
size_t entries() const
Definition: G4Step.hh:78