source: trunk/source/processes/optical/src/G4OpMieHG.cc @ 1350

Last change on this file since 1350 was 1350, checked in by garnier, 13 years ago

update to last version 4.9.4

File size: 7.0 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////////////////////////////////////////////////////////////////////////
28//
29// File G4OpMieHG.hh
30// Description: Discrete Process -- Mie Scattering of Optical Photons
31// Created: 2010-07-03
32// Author: Xin Qian
33// Based on work from Vlasios Vasileiou
34//
35// This subroutine will mimic the Mie scattering based on
36// Henyey-Greenstein phase function
37// Forward and backward angles are treated separately.
38//
39// mail:  gum@triumf.ca
40//
41////////////////////////////////////////////////////////////////////////
42
43#include "G4OpProcessSubType.hh"
44
45#include "G4OpMieHG.hh"
46
47G4OpMieHG::G4OpMieHG(const G4String& processName, G4ProcessType type)
48           : G4VDiscreteProcess(processName, type)
49{
50        if (verboseLevel>0) {
51           G4cout << GetProcessName() << " is created " << G4endl;
52        }
53
54        SetProcessSubType(fOpMieHG);
55}
56
57G4OpMieHG::~G4OpMieHG(){}
58
59        ////////////
60        // Methods
61        ////////////
62
63// PostStepDoIt
64// -------------
65//
66G4VParticleChange* 
67G4OpMieHG::PostStepDoIt(const G4Track& aTrack, const G4Step& aStep)
68{
69        aParticleChange.Initialize(aTrack);
70
71        const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
72        const G4Material* aMaterial = aTrack.GetMaterial();
73        G4MaterialPropertiesTable* aMaterialPropertyTable =
74          aMaterial->GetMaterialPropertiesTable();
75
76        G4double forward_g =
77              aMaterialPropertyTable->GetConstProperty("MIEHG_FORWARD");
78        G4double backward_g =
79              aMaterialPropertyTable->GetConstProperty("MIEHG_BACKWARD");
80        G4double ForwardRatio =
81              aMaterialPropertyTable->GetConstProperty("MIEHG_FORWARD_RATIO");
82
83        if (verboseLevel>0) {
84                G4cout << "MIE Scattering Photon!" << G4endl;
85                G4cout << "MIE Old Momentum Direction: "
86                     << aParticle->GetMomentumDirection() << G4endl;
87                G4cout << "MIE Old Polarization: "
88                     << aParticle->GetPolarization() << G4endl;
89        }
90
91        G4double g;
92        G4int direction;
93        if (G4UniformRand()<=ForwardRatio){
94           g = forward_g;
95           direction = 1;
96        } else {
97           g = backward_g;
98           direction = -1;
99        }
100
101        G4double r = G4UniformRand();
102
103        G4double Theta;
104        //sample the direction
105        if (g!=0) {
106          Theta = std::acos(2*r*(1+g)*(1+g)*(1-g+g*r)/((1-g+2*g*r)*(1-g+2*g*r)) -1);
107        } else {
108          Theta = std::acos(2*r-1.);
109        }
110        G4double Phi = G4UniformRand()*2*pi;
111
112        if (direction==-1) Theta = pi - Theta; //backward scattering
113
114        G4ThreeVector NewMomentumDirection, OldMomentumDirection;
115        G4ThreeVector OldPolarization, NewPolarization;
116
117        NewMomentumDirection.set
118                       (std::sin(Theta)*std::cos(Phi), std::sin(Theta)*std::sin(Phi), std::cos(Theta));
119        OldMomentumDirection = aParticle->GetMomentumDirection();
120        NewMomentumDirection.rotateUz(OldMomentumDirection);
121        NewMomentumDirection = NewMomentumDirection.unit();
122
123        OldPolarization = aParticle->GetPolarization();
124        G4double constant = -1./NewMomentumDirection.dot(OldPolarization);
125
126        NewPolarization = NewMomentumDirection + constant*OldPolarization;
127        NewPolarization = NewPolarization.unit();
128
129        if (NewPolarization.mag()==0) {
130           r = G4UniformRand()*twopi;
131           NewPolarization.set(std::cos(r),std::sin(r),0.);
132           NewPolarization.rotateUz(NewMomentumDirection);
133        } else {
134           // There are two directions which perpendicular
135           // new momentum direction
136           if (G4UniformRand() < 0.5) NewPolarization = -NewPolarization;
137        }
138
139        aParticleChange.ProposePolarization(NewPolarization);
140        aParticleChange.ProposeMomentumDirection(NewMomentumDirection);
141
142        if (verboseLevel>0) {
143              G4cout << "MIE New Polarization: " 
144                     << NewPolarization << G4endl;
145              G4cout << "MIE Polarization Change: "
146                     << *(aParticleChange.GetPolarization()) << G4endl; 
147              G4cout << "MIE New Momentum Direction: " 
148                     << NewMomentumDirection << G4endl;
149              G4cout << "MIE Momentum Change: "
150                     << *(aParticleChange.GetMomentumDirection()) << G4endl; 
151        }
152
153        return G4VDiscreteProcess::PostStepDoIt(aTrack, aStep);
154}
155
156// GetMeanFreePath()
157// -----------------
158//
159G4double G4OpMieHG::GetMeanFreePath(const G4Track& aTrack,
160                                    G4double ,
161                                    G4ForceCondition* )
162{
163        const G4DynamicParticle* aParticle = aTrack.GetDynamicParticle();
164        const G4Material* aMaterial = aTrack.GetMaterial();
165
166        G4double thePhotonEnergy = aParticle->GetTotalEnergy();
167
168        G4double AttenuationLength = DBL_MAX;
169
170        G4MaterialPropertiesTable* aMaterialPropertyTable =
171          aMaterial->GetMaterialPropertiesTable();
172
173        if (aMaterialPropertyTable) {
174           G4MaterialPropertyVector* AttenuationLengthVector =
175                                 aMaterialPropertyTable->GetProperty("MIEHG");
176           if (AttenuationLengthVector) {
177              AttenuationLength = AttenuationLengthVector ->
178                                    GetProperty(thePhotonEnergy);
179           } else {
180//              G4cout << "No Mie scattering length specified" << G4endl;
181           }
182        } else {
183//           G4cout << "No Mie scattering length specified" << G4endl;
184        }
185
186//      G4cout << thePhotonEnergy/GeV << " \t" << AttenuationLength/m << G4endl;
187
188        return AttenuationLength;
189}
Note: See TracBrowser for help on using the repository browser.