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