source: trunk/examples/extended/electromagnetic/TestEm10/src/Em10XTRTransparentRegRadModel.cc@ 1036

Last change on this file since 1036 was 807, checked in by garnier, 17 years ago

update

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