source: trunk/source/processes/electromagnetic/xrays/src/G4TransparentRegXTRadiator.cc

Last change on this file was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 8.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// $Id: G4TransparentRegXTRadiator.cc,v 1.12 2010/06/16 15:34:15 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-04-beta-01 $
29//
30
31#include <complex>
32
33#include "G4TransparentRegXTRadiator.hh"
34#include "Randomize.hh"
35#include "G4Integrator.hh"
36#include "G4Gamma.hh"
37
38////////////////////////////////////////////////////////////////////////////
39//
40// Constructor, destructor
41
42G4TransparentRegXTRadiator::G4TransparentRegXTRadiator(G4LogicalVolume *anEnvelope,
43                                         G4Material* foilMat,G4Material* gasMat, 
44                                         G4double a, G4double b, G4int n,
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
62G4TransparentRegXTRadiator::~G4TransparentRegXTRadiator()
63{
64  ;
65}
66
67///////////////////////////////////////////////////////////////////////////
68//
69//
70
71G4double G4TransparentRegXTRadiator::SpectralXTRdEdx(G4double energy)
72{
73  G4double result, sum = 0., tmp, cof1, cof2, cofMin, cofPHC,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  //  if (fGamma < 1200) kMin = G4int(cofMin);  // 1200 ?
90  // else               kMin = 1;
91
92
93  kMin = G4int(cofMin);
94  if (cofMin > kMin) kMin++;
95
96  // tmp  = (fPlateThick + fGasThick)*energy*fMaxThetaTR;
97  // tmp /= cofPHC;
98  // kMax = G4int(tmp);
99  // if(kMax < 0) kMax = 0;
100  // kMax += kMin;
101 
102
103  kMax = kMin + 19; //  9; // kMin + G4int(tmp);
104
105  // tmp /= fGamma;
106  // if( G4int(tmp) < kMin ) kMin = G4int(tmp);
107  // G4cout<<"kMin = "<<kMin<<";    kMax = "<<kMax<<G4endl;
108
109  for( k = kMin; k <= kMax; k++ )
110  {
111    tmp    = pi*fPlateThick*(k + cof2)/(fPlateThick + fGasThick);
112    result = (k - cof1)*(k - cof1)*(k + cof2)*(k + cof2);
113    // tmp = std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
114    if( k == kMin && kMin == G4int(cofMin) )
115    {
116      sum   += 0.5*std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
117    }
118    else
119    {
120      sum   += std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result;
121    }
122    if(verboseLevel > 2)
123    {   
124      G4cout<<"k = "<<k<<"; tmp = "<<std::sin(tmp)*std::sin(tmp)*std::abs(k-cofMin)/result
125              <<";    sum = "<<sum<<G4endl; 
126    } 
127  }
128  result = 4*( cof1 + cof2 )*( cof1 + cof2 )*sum/energy;
129  // result *= ( 1 - std::exp(-0.5*fPlateNumber*sigma) )/( 1 - std::exp(-0.5*sigma) ); 
130  // fPlateNumber;
131  result *= fPlateNumber; // *std::exp(-0.5*fPlateNumber*sigma);
132                             // +1-std::exp(-0.5*fPlateNumber*sigma);
133  /* 
134  fEnergy = energy;
135  //  G4Integrator<G4VXTRenergyLoss,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
136  G4Integrator<G4TransparentRegXTRadiator,G4double(G4VXTRenergyLoss::*)(G4double)> integral;
137 
138  tmp = integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
139                             0.0,0.3*fMaxThetaTR) +
140      integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
141                             0.3*fMaxThetaTR,0.6*fMaxThetaTR) +         
142      integral.Legendre96(this,&G4VXTRenergyLoss::SpectralAngleXTRdEdx,
143                             0.6*fMaxThetaTR,fMaxThetaTR) ;
144  result += tmp;
145  */
146  return result;
147}
148
149
150///////////////////////////////////////////////////////////////////////////
151//
152// Approximation for radiator interference factor for the case of
153// fully Regular radiator. The plate and gas gap thicknesses are fixed .
154// The mean values of the plate and gas gap thicknesses
155// are supposed to be about XTR formation zones but much less than
156// mean absorption length of XTR photons in coresponding material.
157
158G4double
159G4TransparentRegXTRadiator::GetStackFactor( G4double energy, 
160                                         G4double gamma, G4double varAngle )
161{
162  /*
163  G4double result, Za, Zb, Ma, Mb, sigma;
164 
165  Za = GetPlateFormationZone(energy,gamma,varAngle);
166  Zb = GetGasFormationZone(energy,gamma,varAngle);
167  Ma = GetPlateLinearPhotoAbs(energy);
168  Mb = GetGasLinearPhotoAbs(energy);
169  sigma = Ma*fPlateThick + Mb*fGasThick;
170
171  G4complex Ca(1.0+0.5*fPlateThick*Ma/fAlphaPlate,fPlateThick/Za/fAlphaPlate);
172  G4complex Cb(1.0+0.5*fGasThick*Mb/fAlphaGas,fGasThick/Zb/fAlphaGas);
173
174  G4complex Ha = std::pow(Ca,-fAlphaPlate); 
175  G4complex Hb = std::pow(Cb,-fAlphaGas);
176  G4complex H  = Ha*Hb;
177  G4complex F1 =   (1.0 - Ha)*(1.0 - Hb )/(1.0 - H)
178                 * G4double(fPlateNumber) ;
179  G4complex F2 =   (1.0-Ha)*(1.0-Ha)*Hb/(1.0-H)/(1.0-H)
180                 * (1.0 - std::exp(-0.5*fPlateNumber*sigma)) ;
181  //    *(1.0 - std::pow(H,fPlateNumber)) ;
182    G4complex R  = (F1 + F2)*OneInterfaceXTRdEdx(energy,gamma,varAngle);
183  // G4complex R  = F2*OneInterfaceXTRdEdx(energy,gamma,varAngle);
184  result       = 2.0*std::real(R); 
185  return      result;
186  */
187   // numerically unstable result
188
189  G4double result, Qa, Qb, Q, aZa, bZb, aMa, bMb, D, sigma; 
190 
191  aZa   = fPlateThick/GetPlateFormationZone(energy,gamma,varAngle);
192  bZb   = fGasThick/GetGasFormationZone(energy,gamma,varAngle);
193  aMa   = fPlateThick*GetPlateLinearPhotoAbs(energy);
194  bMb   = fGasThick*GetGasLinearPhotoAbs(energy);
195  sigma = aMa*fPlateThick + bMb*fGasThick;
196  Qa    = std::exp(-0.5*aMa);
197  Qb    = std::exp(-0.5*bMb);
198  Q     = Qa*Qb;
199
200  G4complex Ha( Qa*std::cos(aZa), -Qa*std::sin(aZa)   ); 
201  G4complex Hb( Qb*std::cos(bZb), -Qb*std::sin(bZb)    );
202  G4complex H  = Ha*Hb;
203  G4complex Hs = conj(H);
204  D            = 1.0 /( (1 - Q)*(1 - Q) + 
205                  4*Q*std::sin(0.5*(aZa + bZb))*std::sin(0.5*(aZa + bZb)) );
206  G4complex F1 = (1.0 - Ha)*(1.0 - Hb)*(1.0 - Hs)
207                 * G4double(fPlateNumber)*D;
208  G4complex F2 = (1.0 - Ha)*(1.0 - Ha)*Hb*(1.0 - Hs)*(1.0 - Hs)
209                   // * (1.0 - std::pow(H,fPlateNumber)) * D*D;
210                 * (1.0 - std::exp(-0.5*fPlateNumber*sigma)) * D*D;
211  G4complex R  = (F1 + F2)*OneInterfaceXTRdEdx(energy,gamma,varAngle);
212  result       = 2.0*std::real(R); 
213  return      result;
214 
215}
216
217
218//
219//
220////////////////////////////////////////////////////////////////////////////
221
222
223
224
225
226
227
228
Note: See TracBrowser for help on using the repository browser.