source: trunk/source/processes/optical/include/G4OpRayleigh.hh@ 1090

Last change on this file since 1090 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 5.7 KB
RevLine 
[819]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,v 1.9 2006/06/29 21:08:40 gunter Exp $
[1007]28// GEANT4 tag $Name: geant4-09-02 $
[819]29//
30//
31////////////////////////////////////////////////////////////////////////
32// Optical Photon Rayleigh Scattering Class Definition
33////////////////////////////////////////////////////////////////////////
34//
35// File: G4OpRayleigh.hh
36// Description: Discrete Process -- Rayleigh scattering of optical photons
37// Version: 1.0
38// Created: 1996-05-31
39// Author: Juliet Armstrong
40// Updated: 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//
46////////////////////////////////////////////////////////////////////////
47
48#ifndef G4OpRayleigh_h
49#define G4OpRayleigh_h 1
50
51/////////////
52// Includes
53/////////////
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"
66#include "G4PhysicsOrderedFreeVector.hh"
67
68// Class Description:
69// Discrete Process -- Rayleigh scattering of optical photons.
70// Class inherits publicly from G4VDiscreteProcess.
71// Class Description - End:
72
73/////////////////////
74// Class Definition
75/////////////////////
76
77class G4OpRayleigh : public G4VDiscreteProcess
78{
79
80private:
81
82 //////////////
83 // Operators
84 //////////////
85
86 // G4OpRayleigh& operator=(const G4OpRayleigh &right);
87
88public: // Without description
89
90 ////////////////////////////////
91 // Constructors and Destructor
92 ////////////////////////////////
93
94 G4OpRayleigh(const G4String& processName = "OpRayleigh",
95 G4ProcessType type = fOptical);
96
97 // G4OpRayleigh(const G4OpRayleigh &right);
98
99 ~G4OpRayleigh();
100
101 ////////////
102 // Methods
103 ////////////
104
105public: // With description
106
107 G4bool IsApplicable(const G4ParticleDefinition& aParticleType);
108 // Returns true -> 'is applicable' only for an optical photon.
109
110 G4double GetMeanFreePath(const G4Track& aTrack,
111 G4double ,
112 G4ForceCondition* );
113 // Returns the mean free path for Rayleigh scattering in water.
114 // --- Not yet implemented for other materials! ---
115
116 G4VParticleChange* PostStepDoIt(const G4Track& aTrack,
117 const G4Step& aStep);
118 // This is the method implementing Rayleigh scattering.
119
120 G4PhysicsTable* GetPhysicsTable() const;
121 // Returns the address of the physics table.
122
123 void DumpPhysicsTable() const;
124 // Prints the physics table.
125
126private:
127
128 void BuildThePhysicsTable();
129
130 /////////////////////
131 // Helper Functions
132 /////////////////////
133
134 G4PhysicsOrderedFreeVector* RayleighAttenuationLengthGenerator(
135 G4MaterialPropertiesTable *aMPT);
136
137 ///////////////////////
138 // Class Data Members
139 ///////////////////////
140
141protected:
142
143 G4PhysicsTable* thePhysicsTable;
144 // A Physics Table can be either a cross-sections table or
145 // an energy table (or can be used for other specific
146 // purposes).
147
148private:
149
150 G4bool DefaultWater;
151
152};
153
154////////////////////
155// Inline methods
156////////////////////
157
158inline
159G4bool G4OpRayleigh::IsApplicable(const G4ParticleDefinition& aParticleType)
160{
161 return ( &aParticleType == G4OpticalPhoton::OpticalPhoton() );
162}
163
164inline
165void G4OpRayleigh::DumpPhysicsTable() const
166
167{
168 G4int PhysicsTableSize = thePhysicsTable->entries();
169 G4PhysicsOrderedFreeVector *v;
170
171 for (G4int i = 0 ; i < PhysicsTableSize ; i++ )
172 {
173 v = (G4PhysicsOrderedFreeVector*)(*thePhysicsTable)[i];
174 v->DumpValues();
175 }
176}
177
178inline G4PhysicsTable* G4OpRayleigh::GetPhysicsTable() const
179{
180 return thePhysicsTable;
181}
182
183
184#endif /* G4OpRayleigh_h */
Note: See TracBrowser for help on using the repository browser.