source: trunk/source/processes/electromagnetic/xrays/src/G4XTRTransparentRegRadModel.cc

Last change on this file was 1337, checked in by garnier, 15 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 7.2 KB
Line 
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#include <complex>
29
30#include "G4XTRTransparentRegRadModel.hh"
31#include "Randomize.hh"
32#include "G4Integrator.hh"
33#include "G4Gamma.hh"
34
35////////////////////////////////////////////////////////////////////////////
36//
37// Constructor, destructor
38
39G4XTRTransparentRegRadModel::G4XTRTransparentRegRadModel(G4LogicalVolume *anEnvelope,
40 G4Material* foilMat,G4Material* gasMat,
41 G4double a, G4double b, G4int n,
42 const G4String& processName) :
43 G4VXTRenergyLoss(anEnvelope,foilMat,gasMat,a,b,n,processName)
44{
45 G4cout<<"Regular transparent X-ray TR radiator EM process is called"<<G4endl;
46
47 // Build energy and angular integral spectra of X-ray TR photons from
48 // a radiator
49 fExitFlux = true;
50 fAlphaPlate = 10000;
51 fAlphaGas = 1000;
52
53 // BuildTable();
54}
55
56///////////////////////////////////////////////////////////////////////////
57
58G4XTRTransparentRegRadModel::~G4XTRTransparentRegRadModel()
59{
60 ;
61}
62
63///////////////////////////////////////////////////////////////////////////
64//
65//
66
67G4double G4XTRTransparentRegRadModel::SpectralXTRdEdx(G4double energy)
68{
69 G4double result, sum = 0., tmp, cof1, cof2, cofMin, cofPHC,aMa, bMb, sigma;
70 G4int k, kMax, kMin;
71
72 aMa = GetPlateLinearPhotoAbs(energy);
73 bMb = GetGasLinearPhotoAbs(energy);
74
75 if(fCompton)
76 {
77 aMa += GetPlateCompton(energy);
78 bMb += GetGasCompton(energy);
79 }
80 aMa *= fPlateThick;
81 bMb *= fGasThick;
82
83 sigma = aMa + bMb;
84
85 cofPHC = 4*pi*hbarc;
86 tmp = (fSigma1 - fSigma2)/cofPHC/energy;
87 cof1 = fPlateThick*tmp;
88 cof2 = fGasThick*tmp;
89
90 cofMin = energy*(fPlateThick + fGasThick)/fGamma/fGamma;
91 cofMin += (fPlateThick*fSigma1 + fGasThick*fSigma2)/energy;
92 cofMin /= cofPHC;
93
94 // if (fGamma < 1200) kMin = G4int(cofMin); // 1200 ?
95 // else kMin = 1;
96
97
98 kMin = G4int(cofMin);
99 if (cofMin > kMin) kMin++;
100
101 // tmp = (fPlateThick + fGasThick)*energy*fMaxThetaTR;
102 // tmp /= cofPHC;
103 // kMax = G4int(tmp);
104 // if(kMax < 0) kMax = 0;
105 // kMax += kMin;
106
107
108 kMax = kMin + 19; // 5; // 9; // kMin + G4int(tmp);
109
110 // tmp /= fGamma;
111 // if( G4int(tmp) < kMin ) kMin = G4int(tmp);
112 // G4cout<<"kMin = "<<kMin<<"; kMax = "<<kMax<<G4endl;
113
114 for( k = kMin; k <= kMax; k++ )
115 {
116 tmp = pi*fPlateThick*(k + cof2)/(fPlateThick + fGasThick);
117 result = (k - cof1)*(k - cof1)*(k + cof2)*(k + cof2);
118
119 if( k == kMin && kMin == G4int(cofMin) )
120 {
121 sum += 0.5*std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
122 }
123 else
124 {
125 sum += std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
126 }
127 // G4cout<<"k = "<<k<<"; sum = "<<sum<<G4endl;
128 }
129 result = 4.*( cof1 + cof2 )*( cof1 + cof2 )*sum/energy;
130 result *= ( 1. - std::exp(-fPlateNumber*sigma) )/( 1. - std::exp(-sigma) );
131 return result;
132}
133
134
135///////////////////////////////////////////////////////////////////////////
136//
137// Approximation for radiator interference factor for the case of
138// fully Regular radiator. The plate and gas gap thicknesses are fixed .
139// The mean values of the plate and gas gap thicknesses
140// are supposed to be about XTR formation zones but much less than
141// mean absorption length of XTR photons in coresponding material.
142
143G4double
144G4XTRTransparentRegRadModel::GetStackFactor( G4double energy,
145 G4double gamma, G4double varAngle )
146{
147 /*
148 G4double result, Za, Zb, Ma, Mb, sigma;
149
150 Za = GetPlateFormationZone(energy,gamma,varAngle);
151 Zb = GetGasFormationZone(energy,gamma,varAngle);
152 Ma = GetPlateLinearPhotoAbs(energy);
153 Mb = GetGasLinearPhotoAbs(energy);
154 sigma = Ma*fPlateThick + Mb*fGasThick;
155
156 G4complex Ca(1.0+0.5*fPlateThick*Ma/fAlphaPlate,fPlateThick/Za/fAlphaPlate);
157 G4complex Cb(1.0+0.5*fGasThick*Mb/fAlphaGas,fGasThick/Zb/fAlphaGas);
158
159 G4complex Ha = std::pow(Ca,-fAlphaPlate);
160 G4complex Hb = std::pow(Cb,-fAlphaGas);
161 G4complex H = Ha*Hb;
162 G4complex F1 = (1.0 - Ha)*(1.0 - Hb )/(1.0 - H)
163 * G4double(fPlateNumber) ;
164 G4complex F2 = (1.0-Ha)*(1.0-Ha)*Hb/(1.0-H)/(1.0-H)
165 * (1.0 - std::exp(-0.5*fPlateNumber*sigma)) ;
166 // *(1.0 - std::pow(H,fPlateNumber)) ;
167 G4complex R = (F1 + F2)*OneInterfaceXTRdEdx(energy,gamma,varAngle);
168 // G4complex R = F2*OneInterfaceXTRdEdx(energy,gamma,varAngle);
169 result = 2.0*std::real(R);
170 return result;
171 */
172 // numerically unstable result
173
174 G4double result, Qa, Qb, Q, aZa, bZb, aMa, bMb, D, sigma;
175
176 aZa = fPlateThick/GetPlateFormationZone(energy,gamma,varAngle);
177 bZb = fGasThick/GetGasFormationZone(energy,gamma,varAngle);
178 aMa = fPlateThick*GetPlateLinearPhotoAbs(energy);
179 bMb = fGasThick*GetGasLinearPhotoAbs(energy);
180 sigma = aMa*fPlateThick + bMb*fGasThick;
181 Qa = std::exp(-0.5*aMa);
182 Qb = std::exp(-0.5*bMb);
183 Q = Qa*Qb;
184
185 G4complex Ha( Qa*std::cos(aZa), -Qa*std::sin(aZa) );
186 G4complex Hb( Qb*std::cos(bZb), -Qb*std::sin(bZb) );
187 G4complex H = Ha*Hb;
188 G4complex Hs = conj(H);
189 D = 1.0 /( (1 - Q)*(1 - Q) +
190 4*Q*std::sin(0.5*(aZa + bZb))*std::sin(0.5*(aZa + bZb)) );
191 G4complex F1 = (1.0 - Ha)*(1.0 - Hb)*(1.0 - Hs)
192 * G4double(fPlateNumber)*D;
193 G4complex F2 = (1.0 - Ha)*(1.0 - Ha)*Hb*(1.0 - Hs)*(1.0 - Hs)
194 // * (1.0 - std::pow(H,fPlateNumber)) * D*D;
195 * (1.0 - std::exp(-0.5*fPlateNumber*sigma)) * D*D;
196 G4complex R = (F1 + F2)*OneInterfaceXTRdEdx(energy,gamma,varAngle);
197 result = 2.0*std::real(R);
198 return result;
199
200}
201
202
203//
204//
205////////////////////////////////////////////////////////////////////////////
206
207
208
209
210
211
212
213
Note: See TracBrowser for help on using the repository browser.