Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VLowEnergyDiscretePhotonProcess.cc
Go to the documentation of this file.
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$
27//
28// --------------------------------------------------------------
29//
30// File name: G4VLowEnergyDiscretePhotonProcess.cc
31//
32// Author: Capra Riccardo
33//
34// Creation date: May 2005
35//
36// History:
37// -----------
38// 02 May 2005 R. Capra 1st implementation
39//
40//----------------------------------------------------------------
41
42
43
45
47#include "G4SystemOfUnits.hh"
48#include "G4String.hh"
51#include "G4Gamma.hh"
52#include "G4Track.hh"
53#include "G4DynamicParticle.hh"
54#include "G4ThreeVector.hh"
55#include "Randomize.hh" // G4UniformRand
56
57G4VLowEnergyDiscretePhotonProcess :: G4VLowEnergyDiscretePhotonProcess(const G4String& processName,
58 const G4String& aCrossSectionFileName,
59 const G4String& aScatterFileName,
60 G4VDataSetAlgorithm* aScatterInterpolation,
61 G4double aLowEnergyLimit,
62 G4double aHighEnergyLimit)
63 :
65 lowEnergyLimit(aLowEnergyLimit),
66 highEnergyLimit(aHighEnergyLimit),
67 crossSectionFileName(aCrossSectionFileName),
68 meanFreePathTable(0)
69{
70 crossSectionHandler = new G4CrossSectionHandler();
71 scatterFunctionData = new G4CompositeEMDataSet(aScatterInterpolation, 1., 1.);
72 scatterFunctionData->LoadData(aScatterFileName);
73
74 if (verboseLevel > 0)
75 {
76 G4cout << GetProcessName() << " is created " << G4endl
77 << "Energy range: "
78 << lowEnergyLimit / keV << " keV - "
79 << highEnergyLimit / GeV << " GeV"
80 << G4endl;
81 }
82}
83
84
85
87{
88 if (meanFreePathTable)
89 delete meanFreePathTable;
90
91 delete crossSectionHandler;
92 delete scatterFunctionData;
93}
94
95
96
97
98
100{
101 return (&particleDefinition)==G4Gamma::Gamma();
102}
103
104
105
107{
108 crossSectionHandler->Clear();
109 crossSectionHandler->LoadData(crossSectionFileName);
110
111 if (meanFreePathTable)
112 delete meanFreePathTable;
113 meanFreePathTable=crossSectionHandler->BuildMeanFreePathForMaterials();
114}
115
116
117
118
119
121{
122 G4double photonEnergy;
123 photonEnergy = aTrack.GetDynamicParticle()->GetKineticEnergy();
124
125 if (photonEnergy < lowEnergyLimit)
126 return DBL_MAX;
127
128 if (photonEnergy > highEnergyLimit)
129 photonEnergy=highEnergyLimit;
130
131 size_t materialIndex;
132 materialIndex = aTrack.GetMaterialCutsCouple()->GetIndex();
133
134 return meanFreePathTable->FindValue(photonEnergy, materialIndex);
135}
136
137
138
139
140
142{
143 G4ThreeVector photonMomentumDirection;
144 G4ThreeVector photonPolarization;
145
146 photonPolarization = photon.GetPolarization();
147 photonMomentumDirection = photon.GetMomentumDirection();
148
149 if ((!photonPolarization.isOrthogonal(photonMomentumDirection, 1e-6)) || photonPolarization.mag()==0.)
150 {
151 // if |photonPolarization|==0. or |photonPolarization * photonDirection0| > 1e-6 * |photonPolarization ^ photonDirection0|
152 // then polarization is choosen randomly.
153
154 G4ThreeVector e1(photonMomentumDirection.orthogonal().unit());
155 G4ThreeVector e2(photonMomentumDirection.cross(e1).unit());
156
157 G4double angle(G4UniformRand() * twopi);
158
159 e1*=std::cos(angle);
160 e2*=std::sin(angle);
161
162 photonPolarization=e1+e2;
163 }
164 else if (photonPolarization.howOrthogonal(photonMomentumDirection) != 0.)
165 {
166 // if |photonPolarization * photonDirection0| != 0.
167 // then polarization is made orthonormal;
168
169 photonPolarization=photonPolarization.perpPart(photonMomentumDirection);
170 }
171
172 return photonPolarization.unit();
173}
@ photon
G4ForceCondition
double G4double
Definition: G4Types.hh:64
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
Hep3Vector unit() const
Hep3Vector orthogonal() const
Hep3Vector cross(const Hep3Vector &) const
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
Definition: SpaceVector.cc:237
Hep3Vector perpPart() const
double mag() const
double howOrthogonal(const Hep3Vector &v) const
Definition: SpaceVector.cc:219
G4double GetKineticEnergy() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4DynamicParticle * GetDynamicParticle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4VEMDataSet * BuildMeanFreePathForMaterials(const G4DataVector *energyCuts=0)
void LoadData(const G4String &dataFile)
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
virtual G4bool LoadData(const G4String &fileName)=0
virtual G4bool IsApplicable(const G4ParticleDefinition &particleDefinition)
G4ThreeVector GetPhotonPolarization(const G4DynamicParticle &photon)
virtual G4double GetMeanFreePath(const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
virtual void BuildPhysicsTable(const G4ParticleDefinition &photon)
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4String & GetProcessName() const
Definition: G4VProcess.hh:379
#define DBL_MAX
Definition: templates.hh:83