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