Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RToEConvForGamma.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// G4RToEConvForGamma class implementation
27//
28// Author: H.Kurashige, 05 October 2002 - First implementation
29// --------------------------------------------------------------------
30
31#include "G4RToEConvForGamma.hh"
33#include "G4ParticleTable.hh"
34#include "G4Material.hh"
35#include "G4PhysicsLogVector.hh"
36
37#include "G4ios.hh"
38#include "G4SystemOfUnits.hh"
39
40// --------------------------------------------------------------------
43{
45 if (theParticle == nullptr)
46 {
47#ifdef G4VERBOSE
48 if (GetVerboseLevel()>0)
49 {
50 G4cout << " G4RToEConvForGamma::G4RToEConvForGamma() - ";
51 G4cout << "Gamma is not defined !!" << G4endl;
52 }
53#endif
54 }
55}
56
58{
59}
60
61// ***********************************************************************
62// ******************* BuildAbsorptionLengthVector ***********************
63// ***********************************************************************
65 const G4Material* aMaterial,
66 G4RangeVector* absorptionLengthVector )
67{
68 // fill the absorption length vector for this material
69 // absorption length is defined here as:
70 // absorption length = 5./ macroscopic absorption cross section
71 //
72 const G4CrossSectionTable* aCrossSectionTable
74 const G4ElementVector* elementVector
75 = aMaterial->GetElementVector();
76 const G4double* atomicNumDensityVector
77 = aMaterial->GetAtomicNumDensityVector();
78
79 // fill absorption length vector
80 G4int NumEl = aMaterial->GetNumberOfElements();
81 G4double absorptionLengthMax = 0.0;
82 for (std::size_t ibin=0; ibin<size_t(TotBin); ++ibin)
83 {
84 G4double SIGMA = 0.;
85 for (std::size_t iel=0; iel<size_t(NumEl); ++iel)
86 {
87 G4int IndEl = (*elementVector)[iel]->GetIndex();
88 SIGMA += atomicNumDensityVector[iel]
89 * (*((*aCrossSectionTable)[IndEl]))[ibin];
90 }
91 // absorption length=5./SIGMA
92 absorptionLengthVector->PutValue(ibin, 5./SIGMA);
93 if (absorptionLengthMax < 5./SIGMA )
94 absorptionLengthMax = 5./SIGMA;
95 }
96}
97
98// ***********************************************************************
99// ********************** ComputeCrossSection ****************************
100// ***********************************************************************
102 G4double KineticEnergy)
103{
104 // Compute the "absorption" cross-section of the photon "absorption".
105 // Cross-section means here the sum of the cross-sections of the
106 // pair production, Compton scattering and photoelectric processes
107
108 const G4double t1keV = 1.*keV;
109 const G4double t200keV = 200.*keV;
110 const G4double t100MeV = 100.*MeV;
111
112 // Compute Z dependent quantities in the case of a new AtomicNumber
113 if(std::abs(AtomicNumber-Z)>0.1)
114 {
115 Z = AtomicNumber;
116 G4double Zsquare = Z*Z;
117 G4double Zlog = std::log(Z);
118 G4double Zlogsquare = Zlog*Zlog;
119
120 s200keV = (0.2651-0.1501*Zlog+0.02283*Zlogsquare)*Zsquare;
121 tmin = (0.552+218.5/Z+557.17/Zsquare)*MeV;
122 smin = (0.01239+0.005585*Zlog-0.000923*Zlogsquare)*std::exp(1.5*Zlog);
123 cmin = std::log(s200keV/smin)
124 /(std::log(tmin/t200keV)*std::log(tmin/t200keV));
125 tlow = 0.2*std::exp(-7.355/std::sqrt(Z))*MeV;
126 slow = s200keV
127 * std::exp(0.042*Z*std::log(t200keV/tlow)*std::log(t200keV/tlow));
128 s1keV = 300.*Zsquare;
129 clow = std::log(s1keV/slow)/std::log(tlow/t1keV);
130 chigh = (7.55e-5-0.0542e-5*Z)*Zsquare*Z/std::log(t100MeV/tmin);
131 }
132
133 // Calculate the cross-section (using an approximate empirical formula)
134 G4double xs;
135 if ( KineticEnergy<tlow )
136 {
137 if(KineticEnergy<t1keV) xs = slow*std::exp(clow*std::log(tlow/t1keV));
138 else xs = slow*std::exp(clow*std::log(tlow/KineticEnergy));
139 }
140 else if ( KineticEnergy<t200keV )
141 {
142 xs = s200keV
143 * std::exp(0.042*Z*std::log(t200keV/KineticEnergy)
144 *std::log(t200keV/KineticEnergy));
145 }
146 else if( KineticEnergy<tmin )
147 {
148 xs = smin
149 * std::exp(cmin*std::log(tmin/KineticEnergy)
150 *std::log(tmin/KineticEnergy));
151 }
152 else
153 {
154 xs = smin + chigh*std::log(KineticEnergy/tmin);
155 }
156 return xs * barn;
157}
std::vector< G4Element * > G4ElementVector
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
void PutValue(std::size_t index, G4double theValue)
G4double ComputeCrossSection(G4double AtomicNumber, G4double KineticEnergy)
void BuildAbsorptionLengthVector(const G4Material *aMaterial, G4RangeVector *rangeVector)
const G4ParticleDefinition * theParticle