Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TransparentRegXTRadiator.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
30
31////////////////////////////////////////////////////////////////////////////
32// Constructor, destructor
34 G4LogicalVolume* anEnvelope, G4Material* foilMat, G4Material* gasMat,
35 G4double a, G4double b, G4int n, const G4String& processName)
36 : G4VXTRenergyLoss(anEnvelope, foilMat, gasMat, a, b, n, processName)
37{
38 if(verboseLevel > 0)
39 G4cout << "Regular transparent X-ray TR radiator EM process is called"
40 << G4endl;
41
42 // Build energy and angular integral spectra of X-ray TR photons from
43 // a radiator
44
45 fAlphaPlate = 10000;
46 fAlphaGas = 1000;
47}
48
49///////////////////////////////////////////////////////////////////////////
51
52///////////////////////////////////////////////////////////////////////////
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// Approximation for radiator interference factor for the case of
118// fully Regular radiator. The plate and gas gap thicknesses are fixed.
119// The mean values of the plate and gas gap thicknesses
120// are supposed to be about XTR formation zones but much less than
121// mean absorption length of XTR photons in corresponding material.
123 G4double gamma,
124 G4double varAngle)
125{
126 G4double result, Qa, Qb, Q, aZa, bZb, aMa, bMb, D, sigma;
127
128 aZa = fPlateThick / GetPlateFormationZone(energy, gamma, varAngle);
129 bZb = fGasThick / GetGasFormationZone(energy, gamma, varAngle);
130 aMa = fPlateThick * GetPlateLinearPhotoAbs(energy);
131 bMb = fGasThick * GetGasLinearPhotoAbs(energy);
132 sigma = aMa * fPlateThick + bMb * fGasThick;
133 Qa = std::exp(-0.5 * aMa);
134 Qb = std::exp(-0.5 * bMb);
135 Q = Qa * Qb;
136
137 G4complex Ha(Qa * std::cos(aZa), -Qa * std::sin(aZa));
138 G4complex Hb(Qb * std::cos(bZb), -Qb * std::sin(bZb));
139 G4complex H = Ha * Hb;
140 G4complex Hs = conj(H);
141 D = 1.0 / ((1 - Q) * (1 - Q) +
142 4 * Q * std::sin(0.5 * (aZa + bZb)) * std::sin(0.5 * (aZa + bZb)));
143 G4complex F1 =
144 (1.0 - Ha) * (1.0 - Hb) * (1.0 - Hs) * G4double(fPlateNumber) * D;
145 G4complex F2 = (1.0 - Ha) * (1.0 - Ha) * Hb * (1.0 - Hs) * (1.0 - Hs) *
146 (1.0 - std::exp(-0.5 * fPlateNumber * sigma)) * D * D;
147 G4complex R = (F1 + F2) * OneInterfaceXTRdEdx(energy, gamma, varAngle);
148 result = 2.0 * std::real(R);
149 return result;
150}
G4double D(G4double temp)
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
G4double SpectralXTRdEdx(G4double energy) override
G4TransparentRegXTRadiator(G4LogicalVolume *anEnvelope, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="TransparentRegXTRadiator")
void ProcessDescription(std::ostream &) const override
G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle) override
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)