source: trunk/source/processes/optical/src/G4OpRayleigh.cc @ 1337

Last change on this file since 1337 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 11.4 KB
Line 
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.cc,v 1.17 2008/10/24 19:51:12 gum Exp $
28// GEANT4 tag $Name: geant4-09-04-beta-01 $
29//
30//
31////////////////////////////////////////////////////////////////////////
32// Optical Photon Rayleigh Scattering Class Implementation
33////////////////////////////////////////////////////////////////////////
34//
35// File:        G4OpRayleigh.cc
36// Description: Discrete Process -- Rayleigh scattering of optical
37//              photons 
38// Version:     1.0
39// Created:     1996-05-31 
40// Author:      Juliet Armstrong
41// Updated:     2005-07-28 - add G4ProcessType to constructor
42//              2001-10-18 by Peter Gumplinger
43//              eliminate unused variable warning on Linux (gcc-2.95.2)
44//              2001-09-18 by mma
45//              >numOfMaterials=G4Material::GetNumberOfMaterials() in BuildPhy
46//              2001-01-30 by Peter Gumplinger
47//              > allow for positiv and negative CosTheta and force the
48//              > new momentum direction to be in the same plane as the
49//              > new and old polarization vectors
50//              2001-01-29 by Peter Gumplinger
51//              > fix calculation of SinTheta (from CosTheta)
52//              1997-04-09 by Peter Gumplinger
53//              > new physics/tracking scheme
54// mail:        gum@triumf.ca
55//
56////////////////////////////////////////////////////////////////////////
57
58#include "G4ios.hh"
59#include "G4OpProcessSubType.hh"
60
61#include "G4OpRayleigh.hh"
62
63/////////////////////////
64// Class Implementation
65/////////////////////////
66
67        //////////////
68        // Operators
69        //////////////
70
71// G4OpRayleigh::operator=(const G4OpRayleigh &right)
72// {
73// }
74
75        /////////////////
76        // Constructors
77        /////////////////
78
79G4OpRayleigh::G4OpRayleigh(const G4String& processName, G4ProcessType type)
80           : G4VDiscreteProcess(processName, type)
81{
82        SetProcessSubType(fOpRayleigh);
83
84        thePhysicsTable = 0;
85
86        DefaultWater = false;
87
88        if (verboseLevel>0) {
89           G4cout << GetProcessName() << " is created " << G4endl;
90        }
91
92        BuildThePhysicsTable();
93}
94
95// G4OpRayleigh::G4OpRayleigh(const G4OpRayleigh &right)
96// {
97// }
98
99        ////////////////
100        // Destructors
101        ////////////////
102
103G4OpRayleigh::~G4OpRayleigh()
104{
105        if (thePhysicsTable!= 0) {
106           thePhysicsTable->clearAndDestroy();
107           delete thePhysicsTable;
108        }
109}
110
111        ////////////
112        // Methods
113        ////////////
114
115// PostStepDoIt
116// -------------
117//
118G4VParticleChange* 
119G4OpRayleigh::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
120{
121        aParticleChange.Initialize(aTrack);
122
123        const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
124
125        if (verboseLevel>0) {
126                G4cout << "Scattering Photon!" << G4endl;
127                G4cout << "Old Momentum Direction: "
128                     << aParticle->GetMomentumDirection() << G4endl;
129                G4cout << "Old Polarization: "
130                     << aParticle->GetPolarization() << G4endl;
131        }
132
133        // find polar angle w.r.t. old polarization vector
134
135        G4double rand = G4UniformRand();
136
137        G4double CosTheta = std::pow(rand, 1./3.);
138        G4double SinTheta = std::sqrt(1.-CosTheta*CosTheta);
139
140        if(G4UniformRand() < 0.5)CosTheta = -CosTheta;
141
142        // find azimuthal angle w.r.t old polarization vector
143
144        rand = G4UniformRand();
145
146        G4double Phi = twopi*rand;
147        G4double SinPhi = std::sin(Phi); 
148        G4double CosPhi = std::cos(Phi); 
149       
150        G4double unit_x = SinTheta * CosPhi; 
151        G4double unit_y = SinTheta * SinPhi; 
152        G4double unit_z = CosTheta; 
153       
154        G4ThreeVector NewPolarization (unit_x,unit_y,unit_z);
155
156        // Rotate new polarization direction into global reference system
157
158        G4ThreeVector OldPolarization = aParticle->GetPolarization();
159        OldPolarization = OldPolarization.unit();
160
161        NewPolarization.rotateUz(OldPolarization);
162        NewPolarization = NewPolarization.unit();
163       
164        // -- new momentum direction is normal to the new
165        // polarization vector and in the same plane as the
166        // old and new polarization vectors --
167
168        G4ThreeVector NewMomentumDirection = 
169                              OldPolarization - NewPolarization * CosTheta;
170
171        if(G4UniformRand() < 0.5)NewMomentumDirection = -NewMomentumDirection;
172        NewMomentumDirection = NewMomentumDirection.unit();
173
174        aParticleChange.ProposePolarization(NewPolarization);
175
176        aParticleChange.ProposeMomentumDirection(NewMomentumDirection);
177
178        if (verboseLevel>0) {
179                G4cout << "New Polarization: " 
180                     << NewPolarization << G4endl;
181                G4cout << "Polarization Change: "
182                     << *(aParticleChange.GetPolarization()) << G4endl; 
183                G4cout << "New Momentum Direction: " 
184                     << NewMomentumDirection << G4endl;
185                G4cout << "Momentum Change: "
186                     << *(aParticleChange.GetMomentumDirection()) << G4endl; 
187        }
188
189        return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
190}
191
192// BuildThePhysicsTable for the Rayleigh Scattering process
193// --------------------------------------------------------
194//
195void G4OpRayleigh::BuildThePhysicsTable()
196{
197//      Builds a table of scattering lengths for each material
198
199        if (thePhysicsTable) return;
200
201        const G4MaterialTable* theMaterialTable=
202                               G4Material::GetMaterialTable();
203        G4int numOfMaterials = G4Material::GetNumberOfMaterials();
204
205        // create a new physics table
206
207        thePhysicsTable = new G4PhysicsTable(numOfMaterials);
208
209        // loop for materials
210
211        for (G4int i=0 ; i < numOfMaterials; i++)
212        {
213            G4PhysicsOrderedFreeVector* ScatteringLengths =
214                                new G4PhysicsOrderedFreeVector();
215
216            G4MaterialPropertiesTable *aMaterialPropertiesTable =
217                         (*theMaterialTable)[i]->GetMaterialPropertiesTable();
218                                                                               
219            if(aMaterialPropertiesTable){
220
221              G4MaterialPropertyVector* AttenuationLengthVector =
222                            aMaterialPropertiesTable->GetProperty("RAYLEIGH");
223
224              if(!AttenuationLengthVector){
225
226                if ((*theMaterialTable)[i]->GetName() == "Water")
227                {
228                   // Call utility routine to Generate
229                   // Rayleigh Scattering Lengths
230
231                   DefaultWater = true;
232
233                   ScatteringLengths =
234                   RayleighAttenuationLengthGenerator(aMaterialPropertiesTable);
235                }
236              }
237            }
238
239            thePhysicsTable->insertAt(i,ScatteringLengths);
240        } 
241}
242
243// GetMeanFreePath()
244// -----------------
245//
246G4double G4OpRayleigh::GetMeanFreePath(const G4Track& aTrack,
247                                     G4double ,
248                                     G4ForceCondition* )
249{
250        const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
251        const G4Material* aMaterial = aTrack.GetMaterial();
252
253        G4double thePhotonEnergy = aParticle->GetTotalEnergy();
254
255        G4double AttenuationLength = DBL_MAX;
256
257        if (aMaterial->GetName() == "Water" && DefaultWater){
258
259           G4bool isOutRange;
260
261           AttenuationLength =
262                (*thePhysicsTable)(aMaterial->GetIndex())->
263                           GetValue(thePhotonEnergy, isOutRange);
264        }
265        else {
266
267           G4MaterialPropertiesTable* aMaterialPropertyTable =
268                           aMaterial->GetMaterialPropertiesTable();
269
270           if(aMaterialPropertyTable){
271             G4MaterialPropertyVector* AttenuationLengthVector =
272                   aMaterialPropertyTable->GetProperty("RAYLEIGH");
273             if(AttenuationLengthVector){
274               AttenuationLength = AttenuationLengthVector ->
275                                    GetProperty(thePhotonEnergy);
276             }
277             else{
278//               G4cout << "No Rayleigh scattering length specified" << G4endl;
279             }
280           }
281           else{
282//             G4cout << "No Rayleigh scattering length specified" << G4endl;
283           }
284        }
285
286        return AttenuationLength;
287}
288
289// RayleighAttenuationLengthGenerator()
290// ------------------------------------
291// Private method to compute Rayleigh Scattering Lengths (for water)
292//
293G4PhysicsOrderedFreeVector* 
294G4OpRayleigh::RayleighAttenuationLengthGenerator(G4MaterialPropertiesTable *aMPT) 
295{
296        // Physical Constants
297
298        // isothermal compressibility of water
299        G4double betat = 7.658e-23*m3/MeV;
300
301        // K Boltzman
302        G4double kboltz = 8.61739e-11*MeV/kelvin;
303
304        // Temperature of water is 10 degrees celsius
305        // conversion to kelvin:
306        // TCelsius = TKelvin - 273.15 => 273.15 + 10 = 283.15
307        G4double temp = 283.15*kelvin;
308
309        // Retrieve vectors for refraction index
310        // and photon energy from the material properties table
311
312        G4MaterialPropertyVector* Rindex = aMPT->GetProperty("RINDEX");
313
314        G4double refsq;
315        G4double e;
316        G4double xlambda;
317        G4double c1, c2, c3, c4;
318        G4double Dist;
319        G4double refraction_index;
320
321        G4PhysicsOrderedFreeVector *RayleighScatteringLengths = 
322                                new G4PhysicsOrderedFreeVector();
323
324        if (Rindex ) {
325
326           Rindex->ResetIterator();
327
328           while (++(*Rindex)) {
329
330                e = (Rindex->GetPhotonEnergy());
331
332                refraction_index = Rindex->GetProperty();
333                refsq = refraction_index*refraction_index;
334                xlambda = h_Planck*c_light/e;
335
336                if (verboseLevel>0) {
337                        G4cout << Rindex->GetPhotonEnergy() << " MeV\t";
338                        G4cout << xlambda << " mm\t";
339                }
340
341                c1 = 1 / (6.0 * pi);
342                c2 = std::pow((2.0 * pi / xlambda), 4);
343                c3 = std::pow( ( (refsq - 1.0) * (refsq + 2.0) / 3.0 ), 2);
344                c4 = betat * temp * kboltz;
345
346                Dist = 1.0 / (c1*c2*c3*c4);
347
348                if (verboseLevel>0) {
349                        G4cout << Dist << " mm" << G4endl;
350                }
351                RayleighScatteringLengths->
352                        InsertValues(Rindex->GetPhotonEnergy(), Dist);
353           }
354
355        }
356
357        return RayleighScatteringLengths;
358}
Note: See TracBrowser for help on using the repository browser.