Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4XTRTransparentRegRadModel.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 G4cout << "Regular transparent X-ray TR radiator EM process is called"
39 << G4endl;
40
41 fExitFlux = true;
42 fAlphaPlate = 10000;
43 fAlphaGas = 1000;
44}
45
46///////////////////////////////////////////////////////////////////////////
48
49///////////////////////////////////////////////////////////////////////////
51{
52 out << "Process describing radiator of X-ray transition radiation.\n";
53}
54
55///////////////////////////////////////////////////////////////////////////
57{
58 static constexpr G4double cofPHC = 4. * pi * hbarc;
59 G4double result, sum = 0., tmp, cof1, cof2, cofMin, aMa, bMb, sigma;
60 G4int k, kMax, kMin;
61
62 aMa = GetPlateLinearPhotoAbs(energy);
63 bMb = GetGasLinearPhotoAbs(energy);
64
65 if(fCompton)
66 {
67 aMa += GetPlateCompton(energy);
68 bMb += GetGasCompton(energy);
69 }
70 aMa *= fPlateThick;
71 bMb *= fGasThick;
72
73 sigma = aMa + bMb;
74
75 tmp = (fSigma1 - fSigma2) / cofPHC / energy;
76 cof1 = fPlateThick * tmp;
77 cof2 = fGasThick * tmp;
78
79 cofMin = energy * (fPlateThick + fGasThick) / fGamma / fGamma;
80 cofMin += (fPlateThick * fSigma1 + fGasThick * fSigma2) / energy;
81 cofMin /= cofPHC;
82
83 kMin = G4int(cofMin);
84 if(cofMin > kMin)
85 kMin++;
86
87 kMax = kMin + 19;
88
89 for(k = kMin; k <= kMax; k++)
90 {
91 tmp = pi * fPlateThick * (k + cof2) / (fPlateThick + fGasThick);
92 result = (k - cof1) * (k - cof1) * (k + cof2) * (k + cof2);
93
94 if(k == kMin && kMin == G4int(cofMin))
95 {
96 sum +=
97 0.5 * std::sin(tmp) * std::sin(tmp) * std::abs(k - cofMin) / result;
98 }
99 else
100 {
101 sum += std::sin(tmp) * std::sin(tmp) * std::abs(k - cofMin) / result;
102 }
103 }
104 result = 4. * (cof1 + cof2) * (cof1 + cof2) * sum / energy;
105 result *= (1. - std::exp(-fPlateNumber * sigma)) / (1. - std::exp(-sigma));
106 return result;
107}
108
109///////////////////////////////////////////////////////////////////////////
110// Approximation for radiator interference factor for the case of
111// fully Regular radiator. The plate and gas gap thicknesses are fixed.
112// The mean values of the plate and gas gap thicknesses
113// are supposed to be about XTR formation zones but much less than
114// mean absorption length of XTR photons in corresponding material.
116 G4double gamma,
117 G4double varAngle)
118{
119 G4double aZa = fPlateThick / GetPlateFormationZone(energy, gamma, varAngle);
120 G4double bZb = fGasThick / GetGasFormationZone(energy, gamma, varAngle);
123 G4double sigma = aMa * fPlateThick + bMb * fGasThick;
124 G4double Qa = std::exp(-0.5 * aMa);
125 G4double Qb = std::exp(-0.5 * bMb);
126 G4double Q = Qa * Qb;
127
128 G4complex Ha(Qa * std::cos(aZa), -Qa * std::sin(aZa));
129 G4complex Hb(Qb * std::cos(bZb), -Qb * std::sin(bZb));
130 G4complex H = Ha * Hb;
131 G4complex Hs = conj(H);
132 G4double D =
133 1.0 / ((1. - Q) * (1. - Q) +
134 4. * Q * std::sin(0.5 * (aZa + bZb)) * std::sin(0.5 * (aZa + bZb)));
135 G4complex F1 =
136 (1.0 - Ha) * (1.0 - Hb) * (1.0 - Hs) * G4double(fPlateNumber) * D;
137 G4complex F2 = (1.0 - Ha) * (1.0 - Ha) * Hb * (1.0 - Hs) * (1.0 - Hs) *
138 (1.0 - std::exp(-0.5 * fPlateNumber * sigma)) * D * D;
139 G4complex R = (F1 + F2) * OneInterfaceXTRdEdx(energy, gamma, varAngle);
140 return 2.0 * std::real(R);
141}
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 GetPlateLinearPhotoAbs(G4double)
G4double GetGasFormationZone(G4double, G4double, G4double)
G4complex OneInterfaceXTRdEdx(G4double energy, G4double gamma, G4double varAngle)
G4double GetGasCompton(G4double)
G4double GetPlateFormationZone(G4double, G4double, G4double)
G4double GetGasLinearPhotoAbs(G4double)
G4double GetPlateCompton(G4double)
G4XTRTransparentRegRadModel(G4LogicalVolume *anEnvelope, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="XTRTransparentRegRadModel")
G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle) override
void ProcessDescription(std::ostream &) const override
G4double SpectralXTRdEdx(G4double energy) override