source: trunk/source/processes/electromagnetic/lowenergy/src/G4VLowEnergyDiscretePhotonProcess.cc @ 1063

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

update to geant4.9.2

File size: 5.7 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: G4VLowEnergyDiscretePhotonProcess.cc,v 1.5 2006/06/29 19:41:44 gunter Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29// --------------------------------------------------------------
30//
31// File name:     G4VLowEnergyDiscretePhotonProcess.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//
41//----------------------------------------------------------------
42
43 
44
45#include "G4VLowEnergyDiscretePhotonProcess.hh"
46
47#include "G4String.hh"
48#include "G4CrossSectionHandler.hh"
49#include "G4CompositeEMDataSet.hh"
50#include "G4Gamma.hh"
51#include "G4Track.hh"
52#include "G4DynamicParticle.hh"
53#include "G4ThreeVector.hh"
54#include "Randomize.hh" // G4UniformRand
55
56G4VLowEnergyDiscretePhotonProcess :: G4VLowEnergyDiscretePhotonProcess(const G4String& processName, 
57                                                                       const G4String& aCrossSectionFileName, 
58                                                                       const G4String& aScatterFileName, 
59                                                                       G4VDataSetAlgorithm* aScatterInterpolation, 
60                                                                       G4double aLowEnergyLimit, 
61                                                                       G4double aHighEnergyLimit)
62  :
63  G4VLowEnergyTestableDiscreteProcess(processName),
64  lowEnergyLimit(aLowEnergyLimit),
65  highEnergyLimit(aHighEnergyLimit),
66  crossSectionFileName(aCrossSectionFileName),
67  meanFreePathTable(0)
68{
69  crossSectionHandler = new G4CrossSectionHandler();
70  scatterFunctionData = new G4CompositeEMDataSet(aScatterInterpolation, 1., 1.);
71  scatterFunctionData->LoadData(aScatterFileName);
72 
73  if (verboseLevel > 0)
74    {
75      G4cout << GetProcessName() << " is created " << G4endl
76             << "Energy range: "
77             << lowEnergyLimit / keV << " keV - "
78             << highEnergyLimit / GeV << " GeV"
79             << G4endl;
80    }
81}
82
83
84
85G4VLowEnergyDiscretePhotonProcess::~G4VLowEnergyDiscretePhotonProcess(void)
86{
87  if (meanFreePathTable)
88    delete meanFreePathTable;
89 
90  delete crossSectionHandler;
91  delete scatterFunctionData;
92}
93
94
95
96
97
98G4bool G4VLowEnergyDiscretePhotonProcess::IsApplicable(const G4ParticleDefinition& particleDefinition)
99{
100  return (&particleDefinition)==G4Gamma::Gamma(); 
101}
102
103
104
105void G4VLowEnergyDiscretePhotonProcess::BuildPhysicsTable(const G4ParticleDefinition&  /*photon*/)
106{
107  crossSectionHandler->Clear();
108  crossSectionHandler->LoadData(crossSectionFileName);
109
110  if (meanFreePathTable)
111    delete meanFreePathTable; 
112  meanFreePathTable=crossSectionHandler->BuildMeanFreePathForMaterials();
113}
114
115
116
117
118
119G4double G4VLowEnergyDiscretePhotonProcess::GetMeanFreePath(const G4Track& aTrack, G4double /*previousStepSize*/, G4ForceCondition*  /*condition*/)
120{
121  G4double photonEnergy;
122  photonEnergy = aTrack.GetDynamicParticle()->GetKineticEnergy();
123
124  if (photonEnergy < lowEnergyLimit)
125    return DBL_MAX;
126 
127  if (photonEnergy > highEnergyLimit)
128    photonEnergy=highEnergyLimit;
129 
130  size_t materialIndex;
131  materialIndex = aTrack.GetMaterialCutsCouple()->GetIndex(); 
132
133  return meanFreePathTable->FindValue(photonEnergy, materialIndex);
134}
135
136
137
138
139
140G4ThreeVector G4VLowEnergyDiscretePhotonProcess::GetPhotonPolarization(const G4DynamicParticle&  photon)
141{
142  G4ThreeVector photonMomentumDirection;
143  G4ThreeVector photonPolarization;
144
145  photonPolarization = photon.GetPolarization(); 
146  photonMomentumDirection = photon.GetMomentumDirection();
147
148  if ((!photonPolarization.isOrthogonal(photonMomentumDirection, 1e-6)) || photonPolarization.mag()==0.)
149    {
150      // if |photonPolarization|==0. or |photonPolarization * photonDirection0| > 1e-6 * |photonPolarization ^ photonDirection0|
151      // then polarization is choosen randomly.
152 
153      G4ThreeVector e1(photonMomentumDirection.orthogonal().unit());
154      G4ThreeVector e2(photonMomentumDirection.cross(e1).unit());
155 
156      G4double angle(G4UniformRand() * twopi);
157 
158      e1*=std::cos(angle);
159      e2*=std::sin(angle);
160 
161      photonPolarization=e1+e2;
162    }
163  else if (photonPolarization.howOrthogonal(photonMomentumDirection) != 0.)
164    {
165      // if |photonPolarization * photonDirection0| != 0.
166      // then polarization is made orthonormal;
167 
168      photonPolarization=photonPolarization.perpPart(photonMomentumDirection);
169    }
170 
171  return photonPolarization.unit();
172}
Note: See TracBrowser for help on using the repository browser.