Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4RegularXTRadiator.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
28
29#include "G4Gamma.hh"
31
32////////////////////////////////////////////////////////////////////////////
33// Constructor, destructor
35 G4Material* foilMat,
36 G4Material* gasMat, G4double a,
37 G4double b, G4int n,
38 const G4String& processName)
39 : G4VXTRenergyLoss(anEnvelope, foilMat, gasMat, a, b, n, processName)
40{
41 G4cout << "Regular X-ray TR radiator EM process is called" << G4endl;
42
43 // Build energy and angular integral spectra of X-ray TR photons from
44 // a radiator
45
46 fAlphaPlate = 10000;
47 fAlphaGas = 1000;
48 G4cout << "fAlphaPlate = " << fAlphaPlate << " ; fAlphaGas = " << fAlphaGas
49 << G4endl;
50}
51
52///////////////////////////////////////////////////////////////////////////
54
55void G4RegularXTRadiator::ProcessDescription(std::ostream& out) const
56{
57 out << "Simulation of X-ray transition radiation generated by\n"
58 "relativistic charged particles crossing the interface between\n"
59 "two materials. Thicknesses of plates and gaps are fixed.\n";
60}
61
62///////////////////////////////////////////////////////////////////////////
64{
65 G4double result, sum = 0., tmp, cof1, cof2, cofMin, cofPHC, theta2, theta2k;
66 G4double aMa, bMb, sigma, dump;
67 G4int k, kMax, kMin;
68
69 aMa = fPlateThick * GetPlateLinearPhotoAbs(energy);
70 bMb = fGasThick * GetGasLinearPhotoAbs(energy);
71 sigma = 0.5 * (aMa + bMb);
72 dump = std::exp(-fPlateNumber * sigma);
73 if(verboseLevel > 2)
74 G4cout << " dump = " << dump << G4endl;
75 cofPHC = 4 * pi * hbarc;
76 tmp = (fSigma1 - fSigma2) / cofPHC / energy;
77 cof1 = fPlateThick * tmp;
78 cof2 = fGasThick * tmp;
79
80 cofMin = energy * (fPlateThick + fGasThick) / fGamma / fGamma;
81 cofMin += (fPlateThick * fSigma1 + fGasThick * fSigma2) / energy;
82 cofMin /= cofPHC;
83
84 theta2 = cofPHC / (energy * (fPlateThick + fGasThick));
85
86 kMin = G4int(cofMin);
87 if(cofMin > kMin)
88 kMin++;
89
90 kMax = kMin + 49;
91
92 if(verboseLevel > 2)
93 {
94 G4cout << cof1 << " " << cof2 << " " << cofMin << G4endl;
95 G4cout << "kMin = " << kMin << "; kMax = " << kMax << G4endl;
96 }
97 for(k = kMin; k <= kMax; ++k)
98 {
99 tmp = pi * fPlateThick * (k + cof2) / (fPlateThick + fGasThick);
100 result = (k - cof1) * (k - cof1) * (k + cof2) * (k + cof2);
101 if(k == kMin && kMin == G4int(cofMin))
102 {
103 sum +=
104 0.5 * std::sin(tmp) * std::sin(tmp) * std::abs(k - cofMin) / result;
105 }
106 else
107 {
108 sum += std::sin(tmp) * std::sin(tmp) * std::abs(k - cofMin) / result;
109 }
110 theta2k = std::sqrt(theta2 * std::abs(k - cofMin));
111
112 if(verboseLevel > 2)
113 {
114 G4cout << k << " " << theta2k << " "
115 << std::sin(tmp) * std::sin(tmp) * std::abs(k - cofMin) / result
116 << " " << sum << G4endl;
117 }
118 }
119 result = 2 * (cof1 + cof2) * (cof1 + cof2) * sum / energy;
120 result *= (1 - dump + 2 * dump * fPlateNumber);
121
122 return result;
123}
124
125///////////////////////////////////////////////////////////////////////////
126// Approximation for radiator interference factor for the case of
127// fully Regular radiator. The plate and gas gap thicknesses are fixed.
128// The mean values of the plate and gas gap thicknesses
129// are supposed to be about XTR formation zones but much less than
130// mean absorption length of XTR photons in corresponding material.
131
133 G4double varAngle)
134{
135 // some gamma (10000/1000) like algorithm
136
137 G4double result, Za, Zb, Ma, Mb;
138
139 Za = GetPlateFormationZone(energy, gamma, varAngle);
140 Zb = GetGasFormationZone(energy, gamma, varAngle);
141
142 Ma = GetPlateLinearPhotoAbs(energy);
143 Mb = GetGasLinearPhotoAbs(energy);
144
145 G4complex Ca(1.0 + 0.5 * fPlateThick * Ma / fAlphaPlate,
146 fPlateThick / Za / fAlphaPlate);
147 G4complex Cb(1.0 + 0.5 * fGasThick * Mb / fAlphaGas,
148 fGasThick / Zb / fAlphaGas);
149
150 G4complex Ha = std::pow(Ca, -fAlphaPlate);
151 G4complex Hb = std::pow(Cb, -fAlphaGas);
152 G4complex H = Ha * Hb;
153
154 G4complex F1 = (1.0 - Ha) * (1.0 - Hb) / (1.0 - H) * G4double(fPlateNumber);
155
156 G4complex F2 = (1.0 - Ha) * (1.0 - Ha) * Hb / (1.0 - H) / (1.0 - H) *
157 (1.0 - std::pow(H, fPlateNumber));
158
159 G4complex R = (F1 + F2) * OneInterfaceXTRdEdx(energy, gamma, varAngle);
160
161 result = 2.0 * std::real(R);
162
163 return result;
164}
double G4double
Definition G4Types.hh:83
std::complex< G4double > G4complex
Definition G4Types.hh:88
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4RegularXTRadiator(G4LogicalVolume *anEnvelope, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="XTRegularRadiator")
G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle) override
G4double SpectralXTRdEdx(G4double energy) override
void ProcessDescription(std::ostream &) const override
G4int verboseLevel
G4double GetPlateLinearPhotoAbs(G4double)
G4double GetGasFormationZone(G4double, G4double, G4double)
G4complex OneInterfaceXTRdEdx(G4double energy, G4double gamma, G4double varAngle)
G4double GetPlateFormationZone(G4double, G4double, G4double)
G4double GetGasLinearPhotoAbs(G4double)