Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GaussXTRadiator.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// 19.09.21 V. Grichine, first version
27//
28
29#include "G4GaussXTRadiator.hh"
30
32
33////////////////////////////////////////////////////////////////////////////
34// Constructor, destructor
35
37 G4LogicalVolume* anEnvelope, G4double alphaPlate, G4double alphaGas, G4Material* foilMat, G4Material* gasMat,
38 G4double a, G4double b, G4int n, const G4String& processName)
39 : G4VXTRenergyLoss(anEnvelope, foilMat, gasMat, a, b, n, processName)
40{
41 if(verboseLevel > 0)
42 G4cout << "Gauss X-ray TR radiator EM process is called"
43 << G4endl;
44
45 fAlphaPlate = alphaPlate;
46 fAlphaGas = alphaGas; // 1000; //
47}
48
49///////////////////////////////////////////////////////////////////////////
51
52///////////////////////////////////////////////////////////////////////////
53void G4GaussXTRadiator::ProcessDescription(std::ostream& out) const
54{
55 out << "Simulation of forward X-ray transition radiation generated by\n"
56 "relativistic charged particles crossing the interface between\n"
57 "two materials.\n";
58}
59
60///////////////////////////////////////////////////////////////////////////
62{
63 G4double result, sum = 0., tmp, cof1, cof2, cofMin, cofPHC, theta2, theta2k;
64 G4int k, kMax, kMin;
65
66 cofPHC = 4. * pi * hbarc;
67 tmp = (fSigma1 - fSigma2) / cofPHC / energy;
68 cof1 = fPlateThick * tmp;
69 cof2 = fGasThick * tmp;
70
71 cofMin = energy * (fPlateThick + fGasThick) / fGamma / fGamma;
72 cofMin += (fPlateThick * fSigma1 + fGasThick * fSigma2) / energy;
73 cofMin /= cofPHC;
74
75 theta2 = cofPHC / (energy * (fPlateThick + fGasThick));
76
77 kMin = G4int(cofMin);
78 if(cofMin > kMin)
79 kMin++;
80
81 kMax = kMin + 49;
82
83 if(verboseLevel > 2)
84 {
85 G4cout << cof1 << " " << cof2 << " " << cofMin << G4endl;
86 G4cout << "kMin = " << kMin << "; kMax = " << kMax << G4endl;
87 }
88 for(k = kMin; k <= kMax; ++k)
89 {
90 tmp = pi * fPlateThick * (k + cof2) / (fPlateThick + fGasThick);
91 result = (k - cof1) * (k - cof1) * (k + cof2) * (k + cof2);
92 if(k == kMin && kMin == G4int(cofMin))
93 {
94 sum +=
95 0.5 * std::sin(tmp) * std::sin(tmp) * std::abs(k - cofMin) / result;
96 }
97 else
98 {
99 sum += std::sin(tmp) * std::sin(tmp) * std::abs(k - cofMin) / result;
100 }
101 theta2k = std::sqrt(theta2 * std::abs(k - cofMin));
102
103 if(verboseLevel > 2)
104 {
105 G4cout << k << " " << theta2k << " "
106 << std::sin(tmp) * std::sin(tmp) * std::abs(k - cofMin) / result
107 << " " << sum << G4endl;
108 }
109 }
110 result = 4. * (cof1 + cof2) * (cof1 + cof2) * sum / energy;
111 result *= fPlateNumber;
112
113 return result;
114}
115
116///////////////////////////////////////////////////////////////////////////
117//
118// Approximation for radiator interference factor for the case of
119// Gauss-distributed regular radiator. The plate and gas gap thicknesses are Gauss distributed with RMS
120// sa and sb for plate and gas, respectively.
121// The mean values of the plate and gas gap thicknesses
122// are supposed to be about XTR formation zones.
123
124
126 G4double gamma,
127 G4double varAngle)
128{
129 G4double result, Qa, Qb, Q, Qn, aZa, bZb, aMa, bMb;
130
131 G4double Ma, Mb, Za, Zb;
132
135
136 Za = GetPlateFormationZone(energy, gamma, varAngle);
137
138 aZa = fPlateThick / Za ;
139
140 Zb = GetGasFormationZone(energy, gamma, varAngle);
141
142 bZb = fGasThick / Zb ;
143
144 Ma = GetPlateLinearPhotoAbs(energy);
145
146 aMa = fPlateThick * Ma;
147
148 Mb = GetGasLinearPhotoAbs(energy);
149
150 bMb = fGasThick * Mb;
151
152 // Gauss fluctuation of gas gaps according to RMS = sb = b/fAlphaGas
153
154 G4double gre, gim, pre, pim;
155
156 pre = -0.5 * aMa - sa * sa * ( 4./ Za / Za - Ma*Ma )/8.;
157
158 gre = -0.5 * bMb - sb * sb * ( 4./ Zb / Zb - Mb*Mb )/8.;
159
160 pim = sa * sa * Ma/2./Za - aZa;
161
162 gim = sb * sb * Mb/2./Zb - bZb;
163
164 Qa = std::exp(pre);
165
166 Qb = std::exp(gre);
167
168 // Q = Qa * Qb;
169
170 G4complex Ha( Qa * std::cos(pim), Qa * std::sin(pim) );
171
172 G4complex Hb( Qb * std::cos(gim), Qb * std::sin(gim) );
173
174 G4double hre, him, hnre, hnim;
175
176 hre = pre + gre;
177
178 him = pim + gim;
179
181
182 hnre = nn*hre;
183
184 hnim = nn*him;
185
186 Q = std::exp(hre);
187
188 Qn = std::exp(hnre);
189
190
191 // G4complex H = Ha * Hb;
192
193
194 G4complex H( Q * std::cos(him), Q * std::sin(him) );
195 G4complex Hn( Qn * std::cos(hnim), Qn * std::sin(hnim) );
196
197 // G4complex Hs = conj(H);
198
199 // G4double sigma, D;
200
201 // sigma = aMa * fPlateThick + bMb * fGasThick;
202
203 // D = 1.0 / ((1 - Q) * (1 - Q) + 4 * Q * std::sin(0.5 * (aZa + bZb)) * std::sin(0.5 * (aZa + bZb)));
204
205 // G4complex F1 = ( 1.0 - Ha ) * ( 1.0 - Hb ) * ( 1.0 - Hs ) * G4double(fPlateNumber) * D;
206
207 G4complex F1 = ( 1.0 - Ha ) * ( 1.0 - Hb ) * nn / ( 1. - H );
208
209 // G4complex F2 = ( 1.0 - Ha ) * ( 1.0 - Ha ) * Hb * ( 1.0 - Hs ) * ( 1.0 - Hs ) * (1.0 - std::exp( -0.5 * fPlateNumber * sigma) ) * D * D;
210
211 G4complex F2 = ( 1.0 - Ha ) * ( 1.0 - Ha ) * Hb * ( 1. - Hn ) / ( 1. - H ) / ( 1. - H );
212
213 G4complex R = (F1 + F2) * OneInterfaceXTRdEdx(energy, gamma, varAngle);
214
215 result = 2.0 * std::real(R);
216
217 return result;
218}
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:57
G4GLOB_DLL std::ostream G4cout
void ProcessDescription(std::ostream &) const override
G4double SpectralXTRdEdx(G4double energy) override
G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle) override
G4GaussXTRadiator(G4LogicalVolume *anEnvelope, G4double, G4double, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="GaussXTRadiator")
G4int verboseLevel
Definition: G4VProcess.hh:360
G4double GetPlateLinearPhotoAbs(G4double)
G4double GetGasFormationZone(G4double, G4double, G4double)
G4complex OneInterfaceXTRdEdx(G4double energy, G4double gamma, G4double varAngle)
G4double GetPlateFormationZone(G4double, G4double, G4double)
G4double GetGasLinearPhotoAbs(G4double)