source: trunk/source/processes/electromagnetic/lowenergy/src/G4LowEnergyPolarizedRayleigh.cc @ 961

Last change on this file since 961 was 961, checked in by garnier, 15 years ago

update processes

File size: 7.8 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// $Id: G4LowEnergyPolarizedRayleigh.cc,v 1.7 2006/06/29 19:40:27 gunter Exp $
27// GEANT4 tag $Name: geant4-09-02-ref-02 $
28//
29// --------------------------------------------------------------
30//
31// File name:     G4LowEnergyPolarizedRayleigh.cc
32//
33// Author:        Capra Riccardo
34//
35// Creation date: May 2005
36//
37// History:
38// -----------
39// 02 May 2005  R. Capra         1st implementation
40// 20 May 2005  MGP              Changed name of a local variable hiding
41//                               a data member of the base class
42//
43//----------------------------------------------------------------
44     
45
46#include "G4LowEnergyPolarizedRayleigh.hh"
47
48#include "G4LogLogInterpolation.hh"
49#include "G4VCrossSectionHandler.hh"
50#include "G4VEMDataSet.hh"
51#include "globals.hh" //G4Exception
52
53
54G4LowEnergyPolarizedRayleigh::G4LowEnergyPolarizedRayleigh(const G4String& processName)
55  :
56  G4VLowEnergyDiscretePhotonProcess(processName, "rayl/re-cs-", "rayl/re-ff-", new G4LogLogInterpolation, 250*eV, 100*GeV),
57  intrinsicLowEnergyLimit(10*eV),
58  intrinsicHighEnergyLimit(100*GeV)
59{
60  if (GetLowEnergyLimit() < intrinsicLowEnergyLimit || 
61      GetHighEnergyLimit() > intrinsicHighEnergyLimit)
62    G4Exception("G4LowEnergyPolarizedRayleigh::G4LowEnergyPolarizedRayleigh - Energy limit outside intrinsic process validity range");
63 
64 
65}
66
67
68G4VParticleChange* G4LowEnergyPolarizedRayleigh::PostStepDoIt(const G4Track&  aTrack, const G4Step&  aStep)
69{
70  // aParticleChange comes from G4VProcess
71  aParticleChange.Initialize(aTrack);
72 
73  const G4DynamicParticle* incomingPhoton = aTrack.GetDynamicParticle();
74  G4double incomingPhotonEnergy = incomingPhoton->GetKineticEnergy();
75 
76  if (incomingPhotonEnergy <= GetLowEnergyLimit())
77    {
78      aParticleChange.ProposeTrackStatus(fStopAndKill);
79      aParticleChange.ProposeEnergy(0.);
80      aParticleChange.ProposeLocalEnergyDeposit(incomingPhotonEnergy);
81 
82      return G4VLowEnergyDiscretePhotonProcess::PostStepDoIt(aTrack, aStep);
83    }
84
85  const G4VCrossSectionHandler* crossSectionHandle = GetCrossSectionHandler();
86  const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
87  G4int zAtom = crossSectionHandle->SelectRandomAtom(couple, incomingPhotonEnergy);
88
89  G4double outcomingPhotonCosTheta = GenerateCosTheta(incomingPhotonEnergy, zAtom);
90  G4double outcomingPhotonPhi = GeneratePhi(outcomingPhotonCosTheta);
91  G4double beta=GeneratePolarizationAngle();
92 
93  // incomingPhoton reference frame:
94  // z = versor parallel to the incomingPhotonDirection
95  // x = versor parallel to the incomingPhotonPolarization
96  // y = defined as z^x
97 
98  // outgoingPhoton reference frame:
99  // z' = versor parallel to the outgoingPhotonDirection
100  // x' = defined as x-x*z'z' normalized
101  // y' = defined as z'^x'
102 
103  G4ThreeVector z(incomingPhoton->GetMomentumDirection().unit()); 
104  G4ThreeVector x(GetPhotonPolarization(*incomingPhoton));
105  G4ThreeVector y(z.cross(x));
106 
107  // z' = std::cos(phi)*std::sin(theta) x + std::sin(phi)*std::sin(theta) y + std::cos(theta) z
108  G4double xDir;
109  G4double yDir;
110  G4double zDir;
111  zDir=outcomingPhotonCosTheta;
112  xDir=std::sqrt(1-outcomingPhotonCosTheta*outcomingPhotonCosTheta);
113  yDir=xDir;
114  xDir*=std::cos(outcomingPhotonPhi);
115  yDir*=std::sin(outcomingPhotonPhi);
116 
117  G4ThreeVector zPrime((xDir*x + yDir*y + zDir*z).unit());
118  G4ThreeVector xPrime(x.perpPart(zPrime).unit());
119  G4ThreeVector yPrime(zPrime.cross(xPrime));
120 
121  // outgoingPhotonPolarization is directed as x' std::cos(beta) + y' std::sin(beta)
122  G4ThreeVector outcomingPhotonPolarization(xPrime*std::cos(beta) + yPrime*std::sin(beta));
123 
124  aParticleChange.ProposeEnergy(incomingPhotonEnergy);
125  aParticleChange.ProposeMomentumDirection(zPrime);
126  aParticleChange.ProposePolarization(outcomingPhotonPolarization);
127  aParticleChange.SetNumberOfSecondaries(0);
128
129  // returns aParticleChange though pParticleChange and G4VProcess::PostStepDoIt
130  return G4VLowEnergyDiscretePhotonProcess::PostStepDoIt(aTrack, aStep);
131}
132
133
134
135
136G4double G4LowEnergyPolarizedRayleigh::GenerateCosTheta(G4double incomingPhotonEnergy, G4int zAtom) const
137{
138  //  d sigma                                                                    k0
139  // --------- =  r0^2 * pi * F^2(x, Z) * ( 2 - sin^2 theta) * std::sin (theta), x = ---- std::sin(theta/2)
140  //  d theta                                                                    hc
141 
142  //  d sigma                                             k0          1 - y
143  // --------- = r0^2 * pi * F^2(x, Z) * ( 1 + y^2), x = ---- std::sqrt ( ------- ), y = std::cos(theta)
144  //    d y                                               hc            2
145
146  //              Z
147  // F(x, Z) ~ --------
148  //            a + bx
149  //
150  // The time to exit from the outer loop grows as ~ k0
151  // On pcgeant2 the time is ~ 1 s for k0 ~ 1 MeV on the oxygen element. A 100 GeV
152  // event will take ~ 10 hours.
153  //
154  // On the avarage the inner loop does 1.5 iterations before exiting
155 
156  const G4double xFactor = (incomingPhotonEnergy*cm)/(h_Planck*c_light);
157  const G4VEMDataSet * formFactorData = GetScatterFunctionData();
158
159  G4double cosTheta;
160  G4double fCosTheta;
161  G4double x;
162  G4double fValue;
163
164  do
165    {
166      do
167        {
168          cosTheta = 2.*G4UniformRand()-1.;
169          fCosTheta = (1.+cosTheta*cosTheta)/2.;
170        }
171      while (fCosTheta < G4UniformRand());
172 
173      x = xFactor*std::sqrt((1.-cosTheta)/2.);
174 
175      if (x > 1.e+005)
176        fValue = formFactorData->FindValue(x, zAtom-1);
177      else
178        fValue = formFactorData->FindValue(0., zAtom-1);
179   
180      fValue/=zAtom;
181      fValue*=fValue;
182    }
183  while(fValue < G4UniformRand());
184
185  return cosTheta;
186}
187
188
189
190G4double G4LowEnergyPolarizedRayleigh::GeneratePhi(G4double cosTheta) const
191{
192  //  d sigma
193  // --------- = alpha * ( 1 - sin^2 (theta) * cos^2 (phi) )
194  //   d phi
195 
196  // On the average the loop takes no more than 2 iterations before exiting
197
198  G4double phi;
199  G4double cosPhi;
200  G4double phiProbability;
201  G4double sin2Theta;
202 
203  sin2Theta=1.-cosTheta*cosTheta;
204 
205  do
206    {
207      phi = twopi * G4UniformRand();
208      cosPhi = std::cos(phi);
209      phiProbability= 1. - sin2Theta*cosPhi*cosPhi;
210    }
211  while (phiProbability < G4UniformRand());
212 
213  return phi;
214}
215
216
217
218
219
220G4double G4LowEnergyPolarizedRayleigh::GeneratePolarizationAngle(void) const
221{
222  // Rayleigh polarization is always on the x' direction
223
224  return 0;
225}
Note: See TracBrowser for help on using the repository browser.