Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4TransparentRegXTRadiator.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 // $Id: G4TransparentRegXTRadiator.cc 68037 2013-03-13 14:15:08Z gcosmo $
28 //
29 
30 #include <complex>
31 
33 #include "G4PhysicalConstants.hh"
34 #include "Randomize.hh"
35 #include "G4Integrator.hh"
36 #include "G4Gamma.hh"
37 
38 ////////////////////////////////////////////////////////////////////////////
39 //
40 // Constructor, destructor
41 
43  G4Material* foilMat,G4Material* gasMat,
45  const G4String& processName) :
46  G4VXTRenergyLoss(anEnvelope,foilMat,gasMat,a,b,n,processName)
47 {
48  if(verboseLevel > 0)
49  G4cout<<"Regular transparent X-ray TR radiator EM process is called"<<G4endl;
50 
51  // Build energy and angular integral spectra of X-ray TR photons from
52  // a radiator
53 
54  fAlphaPlate = 10000;
55  fAlphaGas = 1000;
56 
57  // BuildTable();
58 }
59 
60 ///////////////////////////////////////////////////////////////////////////
61 
63 {
64  ;
65 }
66 
67 ///////////////////////////////////////////////////////////////////////////
68 //
69 //
70 
72 {
73  G4double result, sum = 0., tmp, cof1, cof2, cofMin, cofPHC, theta2, theta2k /*, aMa, bMb ,sigma*/;
74  G4int k, kMax, kMin;
75 
76  //aMa = fPlateThick*GetPlateLinearPhotoAbs(energy);
77  //bMb = fGasThick*GetGasLinearPhotoAbs(energy);
78  //sigma = aMa + bMb;
79 
80  cofPHC = 4*pi*hbarc;
81  tmp = (fSigma1 - fSigma2)/cofPHC/energy;
82  cof1 = fPlateThick*tmp;
83  cof2 = fGasThick*tmp;
84 
85  cofMin = energy*(fPlateThick + fGasThick)/fGamma/fGamma;
86  cofMin += (fPlateThick*fSigma1 + fGasThick*fSigma2)/energy;
87  cofMin /= cofPHC;
88 
89  theta2 = cofPHC/(energy*(fPlateThick + fGasThick));
90 
91  // if (fGamma < 1200) kMin = G4int(cofMin); // 1200 ?
92  // else kMin = 1;
93 
94 
95  kMin = G4int(cofMin);
96  if (cofMin > kMin) kMin++;
97 
98  // tmp = (fPlateThick + fGasThick)*energy*fMaxThetaTR;
99  // tmp /= cofPHC;
100  // kMax = G4int(tmp);
101  // if(kMax < 0) kMax = 0;
102  // kMax += kMin;
103 
104 
105  kMax = kMin + 49; // 19; // kMin + G4int(tmp);
106 
107  // tmp /= fGamma;
108  // if( G4int(tmp) < kMin ) kMin = G4int(tmp);
109 
110  if(verboseLevel > 2)
111  {
112  G4cout<<cof1<<" "<<cof2<<" "<<cofMin<<G4endl;
113  G4cout<<"kMin = "<<kMin<<"; kMax = "<<kMax<<G4endl;
114  }
115  for( k = kMin; k <= kMax; k++ )
116  {
117  tmp = pi*fPlateThick*(k + cof2)/(fPlateThick + fGasThick);
118  result = (k - cof1)*(k - cof1)*(k + cof2)*(k + cof2);
119  // tmp = std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
120  if( k == kMin && kMin == G4int(cofMin) )
121  {
122  sum += 0.5*std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
123  }
124  else
125  {
126  sum += std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
127  }
128  theta2k = std::sqrt(theta2*std::abs(k-cofMin));
129 
130  if(verboseLevel > 2)
131  {
132  // G4cout<<"k = "<<k<<"; sqrt(theta2k) = "<<theta2k<<"; tmp = "<<std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result
133  // <<"; sum = "<<sum<<G4endl;
134  G4cout<<k<<" "<<theta2k<<" "<<std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result
135  <<" "<<sum<<G4endl;
136  }
137  }
138  result = 4*( cof1 + cof2 )*( cof1 + cof2 )*sum/energy;
139  // result *= ( 1 - std::exp(-0.5*fPlateNumber*sigma) )/( 1 - std::exp(-0.5*sigma) );
140  // fPlateNumber;
141  result *= fPlateNumber; // *std::exp(-0.5*fPlateNumber*sigma);
142  // +1-std::exp(-0.5*fPlateNumber*sigma);
143  /*
144  fEnergy = energy;
145  // G4Integrator<G4VXTRenergyLoss,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
146  G4Integrator<G4TransparentRegXTRadiator,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
147 
148  tmp = integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
149  0.0,0.3*fMaxThetaTR) +
150  integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
151  0.3*fMaxThetaTR,0.6*fMaxThetaTR) +
152  integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
153  0.6*fMaxThetaTR,fMaxThetaTR) ;
154  result += tmp;
155  */
156  return result;
157 }
158 
159 
160 ///////////////////////////////////////////////////////////////////////////
161 //
162 // Approximation for radiator interference factor for the case of
163 // fully Regular radiator. The plate and gas gap thicknesses are fixed .
164 // The mean values of the plate and gas gap thicknesses
165 // are supposed to be about XTR formation zones but much less than
166 // mean absorption length of XTR photons in coresponding material.
167 
168 G4double
170  G4double gamma, G4double varAngle )
171 {
172  /*
173  G4double result, Za, Zb, Ma, Mb, sigma;
174 
175  Za = GetPlateFormationZone(energy,gamma,varAngle);
176  Zb = GetGasFormationZone(energy,gamma,varAngle);
177  Ma = GetPlateLinearPhotoAbs(energy);
178  Mb = GetGasLinearPhotoAbs(energy);
179  sigma = Ma*fPlateThick + Mb*fGasThick;
180 
181  G4complex Ca(1.0+0.5*fPlateThick*Ma/fAlphaPlate,fPlateThick/Za/fAlphaPlate);
182  G4complex Cb(1.0+0.5*fGasThick*Mb/fAlphaGas,fGasThick/Zb/fAlphaGas);
183 
184  G4complex Ha = std::pow(Ca,-fAlphaPlate);
185  G4complex Hb = std::pow(Cb,-fAlphaGas);
186  G4complex H = Ha*Hb;
187  G4complex F1 = (1.0 - Ha)*(1.0 - Hb )/(1.0 - H)
188  * G4double(fPlateNumber) ;
189  G4complex F2 = (1.0-Ha)*(1.0-Ha)*Hb/(1.0-H)/(1.0-H)
190  * (1.0 - std::exp(-0.5*fPlateNumber*sigma)) ;
191  // *(1.0 - std::pow(H,fPlateNumber)) ;
192  G4complex R = (F1 + F2)*OneInterfaceXTRdEdx(energy,gamma,varAngle);
193  // G4complex R = F2*OneInterfaceXTRdEdx(energy,gamma,varAngle);
194  result = 2.0*std::real(R);
195  return result;
196  */
197  // numerically unstable result
198 
199  G4double result, Qa, Qb, Q, aZa, bZb, aMa, bMb, D, sigma;
200 
201  aZa = fPlateThick/GetPlateFormationZone(energy,gamma,varAngle);
202  bZb = fGasThick/GetGasFormationZone(energy,gamma,varAngle);
203  aMa = fPlateThick*GetPlateLinearPhotoAbs(energy);
204  bMb = fGasThick*GetGasLinearPhotoAbs(energy);
205  sigma = aMa*fPlateThick + bMb*fGasThick;
206  Qa = std::exp(-0.5*aMa);
207  Qb = std::exp(-0.5*bMb);
208  Q = Qa*Qb;
209 
210  G4complex Ha( Qa*std::cos(aZa), -Qa*std::sin(aZa) );
211  G4complex Hb( Qb*std::cos(bZb), -Qb*std::sin(bZb) );
212  G4complex H = Ha*Hb;
213  G4complex Hs = conj(H);
214  D = 1.0 /( (1 - Q)*(1 - Q) +
215  4*Q*std::sin(0.5*(aZa + bZb))*std::sin(0.5*(aZa + bZb)) );
216  G4complex F1 = (1.0 - Ha)*(1.0 - Hb)*(1.0 - Hs)
217  * G4double(fPlateNumber)*D;
218  G4complex F2 = (1.0 - Ha)*(1.0 - Ha)*Hb*(1.0 - Hs)*(1.0 - Hs)
219  // * (1.0 - std::pow(H,fPlateNumber)) * D*D;
220  * (1.0 - std::exp(-0.5*fPlateNumber*sigma)) * D*D;
221  G4complex R = (F1 + F2)*OneInterfaceXTRdEdx(energy,gamma,varAngle);
222  result = 2.0*std::real(R);
223  return result;
224 
225 }
226 
227 
228 //
229 //
230 ////////////////////////////////////////////////////////////////////////////
231 
232 
233 
234 
235 
236 
237 
238 
G4double GetGasFormationZone(G4double, G4double, G4double)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double GetPlateLinearPhotoAbs(G4double)
G4double GetGasLinearPhotoAbs(G4double)
int G4int
Definition: G4Types.hh:78
G4double GetPlateFormationZone(G4double, G4double, G4double)
G4double GetStackFactor(G4double energy, G4double gamma, G4double varAngle)
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
std::complex< G4double > G4complex
Definition: G4Types.hh:81
G4GLOB_DLL std::ostream G4cout
const G4int n
G4double SpectralXTRdEdx(G4double energy)
G4TransparentRegXTRadiator(G4LogicalVolume *anEnvelope, G4Material *, G4Material *, G4double, G4double, G4int, const G4String &processName="TransparentRegXTRadiator")
G4complex OneInterfaceXTRdEdx(G4double energy, G4double gamma, G4double varAngle)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76