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