Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GammaXTRadiator.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
29#include <complex>
30
31#include "G4GammaXTRadiator.hh"
32#include "Randomize.hh"
33
34#include "G4Gamma.hh"
35
36////////////////////////////////////////////////////////////////////////////
37//
38// Constructor, destructor
39
41 G4double alphaPlate,
42 G4double alphaGas,
43 G4Material* foilMat,G4Material* gasMat,
44 G4double a, G4double b, G4int n,
45 const G4String& processName) :
46 G4VXTRenergyLoss(anEnvelope,foilMat,gasMat,a,b,n,processName)
47{
48 G4cout<<"Gamma distributed X-ray TR radiator model is called"<<G4endl ;
49
50 // Build energy and angular integral spectra of X-ray TR photons from
51 // a radiator
52
53 fAlphaPlate = alphaPlate ;
54 fAlphaGas = alphaGas ;
55 G4cout<<"fAlphaPlate = "<<fAlphaPlate<<" ; fAlphaGas = "<<fAlphaGas<<G4endl ;
56
57 // BuildTable() ;
58}
59
60///////////////////////////////////////////////////////////////////////////
61
63{
64 ;
65}
66
67
68
69///////////////////////////////////////////////////////////////////////////
70//
71// Rough approximation for radiator interference factor for the case of
72// fully GamDistr radiator. The plate and gas gap thicknesses are distributed
73// according to exponent. The mean values of the plate and gas gap thicknesses
74// are supposed to be about XTR formation zones but much less than
75// mean absorption length of XTR photons in coresponding material.
76
79 G4double gamma, G4double varAngle )
80{
81 G4double result, Za, Zb, Ma, Mb ;
82
83 Za = GetPlateFormationZone(energy,gamma,varAngle) ;
84 Zb = GetGasFormationZone(energy,gamma,varAngle) ;
85
86 Ma = GetPlateLinearPhotoAbs(energy) ;
87 Mb = GetGasLinearPhotoAbs(energy) ;
88
89
92
93 G4complex Ha = std::pow(Ca,-fAlphaPlate) ;
94 G4complex Hb = std::pow(Cb,-fAlphaGas) ;
95 G4complex H = Ha*Hb ;
96
97 G4complex F1 = (1.0 - Ha)*(1.0 - Hb )/(1.0 - H)
99
100 G4complex F2 = (1.0-Ha)*(1.0-Ha)*Hb/(1.0-H)/(1.0-H)
101 * (1.0 - std::pow(H,fPlateNumber)) ;
102
103 G4complex R = (F1 + F2)*OneInterfaceXTRdEdx(energy,gamma,varAngle) ;
104
105 result = 2.0*std::real(R) ;
106
107 return result ;
108}
109
110
111//
112//
113////////////////////////////////////////////////////////////////////////////
114
115
116
117
118
119
120
121
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
G4GammaXTRadiator(G4LogicalVolume *anEnvelope, G4double, G4double, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="XTRgammaRadiator")
G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle) override
G4double GetPlateLinearPhotoAbs(G4double)
G4double GetGasFormationZone(G4double, G4double, G4double)
G4complex OneInterfaceXTRdEdx(G4double energy, G4double gamma, G4double varAngle)
G4double GetPlateFormationZone(G4double, G4double, G4double)
G4double GetGasLinearPhotoAbs(G4double)