Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4StrawTubeXTRadiator.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//
27//
28
31#include "G4SystemOfUnits.hh"
32#include "Randomize.hh"
33#include "G4Gamma.hh"
34
35////////////////////////////////////////////////////////////////////////////
36//
37// Constructor, destructor
38
40 G4Material* foilMat,G4Material* gasMat,
41 G4double a, G4double b, G4Material* mediumMat,
42 G4bool unishut,
43 const G4String& processName) :
44 G4VXTRenergyLoss(anEnvelope,foilMat,gasMat,a,b,1,processName)
45{
46 if(verboseLevel > 0)
47 G4cout<<"Straw tube X-ray TR radiator EM process is called"<<G4endl;
48
49 if( unishut )
50 {
51 fAlphaPlate = 1./3.;
52 fAlphaGas = 12.4;
53 if(verboseLevel > 0)
54 G4cout<<"straw uniform shooting: "<<"fAlphaPlate = "
55 <<fAlphaPlate<<" ; fAlphaGas = "<<fAlphaGas<<G4endl;
56
57 }
58 else
59 {
60 fAlphaPlate = 0.5;
61 fAlphaGas = 5.;
62 if(verboseLevel > 0)
63 G4cout<<"straw isotropical shooting: "<<"fAlphaPlate = "
64 <<fAlphaPlate<<" ; fAlphaGas = "<<fAlphaGas<<G4endl;
65
66
67 }
68 // index of medium material
69
70 fMatIndex3 = mediumMat->GetIndex();
71 if(verboseLevel > 0)
72 G4cout<<"medium material = "<<mediumMat->GetName()<<G4endl;
73
74 // plasma energy squared for plate material
75
77 if(verboseLevel > 0)
78 G4cout<<"medium plasma energy = "<<std::sqrt(fSigma3)/eV<<" eV"<<G4endl;
79
80 // Compute cofs for preparation of linear photo absorption in external medium
81
83
84 // Build energy and angular integral spectra of X-ray TR photons from
85 // a radiator
86
87 // BuildTable();
88}
89
90///////////////////////////////////////////////////////////////////////////
91
93{
94}
95
96///////////////////////////////////////////////////////////////////////////
97//
98// Approximation for radiator interference factor for the case of
99// straw tube radiator. The plate (window, straw wall) and gas (inside straw)
100// gap thicknesses are gamma distributed.
101// The mean values of the plate and gas gap thicknesses
102// are supposed to be about XTR formation zone.
103
106 G4double gamma, G4double varAngle )
107{
108
109
110 G4double result, L2, L3, M2, M3;
111
112 L2 = GetPlateFormationZone(energy,gamma,varAngle);
113 L3 = GetGasFormationZone(energy,gamma,varAngle);
114
115 M2 = GetPlateLinearPhotoAbs(energy);
116 M3 = GetGasLinearPhotoAbs(energy);
117
120
121 G4complex H2 = std::pow(C2,-fAlphaPlate);
122 G4complex H3 = std::pow(C3,-fAlphaGas);
123 G4complex H = H2*H3;
124
125 G4complex Z1 = GetMediumComplexFZ(energy,gamma,varAngle);
126 G4complex Z2 = GetPlateComplexFZ(energy,gamma,varAngle);
127 G4complex Z3 = GetGasComplexFZ(energy,gamma,varAngle);
128
129
130 G4complex R = ( Z1 - Z2 )*( Z1 - Z2 )*( 1. - H2*H ) +
131 ( Z2 - Z3 )*( Z2 - Z3 )*( 1. - H3 ) +
132 2.*( Z1 - Z2 )*( Z2 - Z3 )*H2*( 1. - H3 ) ;
133
134 result = 2.0*std::real(R)*(varAngle*energy/hbarc/hbarc);
135
136 return result;
137
138}
139
140
141//////////////////////////////////////////////////////////////////////
142//////////////////////////////////////////////////////////////////////
143//////////////////////////////////////////////////////////////////////
144//
145// Calculates formation zone for external medium. Omega is energy !!!
146
148 G4double gamma ,
149 G4double varAngle )
150{
151 G4double cof, lambda;
152 lambda = 1.0/gamma/gamma + varAngle + fSigma3/omega/omega;
153 cof = 2.0*hbarc/omega/lambda ;
154 return cof ;
155}
156
157//////////////////////////////////////////////////////////////////////
158//
159// Calculates complex formation zone for external medium. Omega is energy !!!
160
162 G4double gamma ,
163 G4double varAngle )
164{
165 G4double cof, length,delta, real_v, image_v;
166
167 length = 0.5*GetMediumFormationZone(omega,gamma,varAngle);
168 delta = length*GetMediumLinearPhotoAbs(omega);
169 cof = 1.0/(1.0 + delta*delta);
170
171 real_v = length*cof;
172 image_v = real_v*delta;
173
174 G4complex zone(real_v,image_v);
175 return zone;
176}
177
178////////////////////////////////////////////////////////////////////////
179//
180// Computes matrix of Sandia photo absorption cross section coefficients for
181// medium material
182
184{
185 const G4MaterialTable* theMaterialTable = G4Material::GetMaterialTable();
186 const G4Material* mat = (*theMaterialTable)[fMatIndex3];
188}
189
190//////////////////////////////////////////////////////////////////////
191//
192// Returns the value of linear photo absorption coefficient (in reciprocal
193// length) for medium for given energy of X-ray photon omega
194
196{
197 G4double omega2, omega3, omega4;
198
199 omega2 = omega*omega;
200 omega3 = omega2*omega;
201 omega4 = omega2*omega2;
202
203 const G4double* SandiaCof = fMediumPhotoAbsCof->GetSandiaCofForMaterial(omega);
204
205 G4double cross = SandiaCof[0]/omega + SandiaCof[1]/omega2 +
206 SandiaCof[2]/omega3 + SandiaCof[3]/omega4;
207 return cross;
208}
209
210//
211//
212////////////////////////////////////////////////////////////////////////////
213
214
215
216
217
218
219
220
std::vector< G4Material * > G4MaterialTable
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
std::complex< G4double > G4complex
Definition: G4Types.hh:88
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const double C2
#define C3
G4SandiaTable * GetSandiaTable() const
Definition: G4Material.hh:227
G4double GetElectronDensity() const
Definition: G4Material.hh:215
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:637
const G4String & GetName() const
Definition: G4Material.hh:175
size_t GetIndex() const
Definition: G4Material.hh:258
G4double GetSandiaCofForMaterial(G4int, G4int) const
G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle) override
G4SandiaTable * fMediumPhotoAbsCof
G4StrawTubeXTRadiator(G4LogicalVolume *anEnvelope, G4Material *, G4Material *, G4double, G4double, G4Material *, G4bool unishut=false, const G4String &processName="StrawTubeXTRadiator")
G4double GetMediumFormationZone(G4double, G4double, G4double)
G4double GetMediumLinearPhotoAbs(G4double)
G4complex GetMediumComplexFZ(G4double, G4double, G4double)
G4int verboseLevel
Definition: G4VProcess.hh:356
G4double GetPlateLinearPhotoAbs(G4double)
G4double GetGasFormationZone(G4double, G4double, G4double)
G4complex GetPlateComplexFZ(G4double, G4double, G4double)
G4double GetPlateFormationZone(G4double, G4double, G4double)
G4double GetGasLinearPhotoAbs(G4double)
G4complex GetGasComplexFZ(G4double, G4double, G4double)