Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PhotoElectricAngularGeneratorPolarized.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 // -------------------------------------------------------------------
28 //
29 // GEANT4 Class file
30 //
31 //
32 // File name: G4PhotoElectricAngularGeneratorPolarized
33 //
34 // Author: A. C. Farinha, L. Peralta, P. Rodrigues and A. Trindade
35 //
36 // Creation date:
37 //
38 // Modifications:
39 // 10 January 2006
40 // 06 May 2011, Replace FILE with std::ifstream
41 //
42 // Class Description:
43 //
44 // Concrete class for PhotoElectric Electron Angular Polarized Distribution Generation
45 //
46 // Class Description:
47 // PhotoElectric Electron Angular Generator based on the general Gavrila photoelectron angular distribution.
48 // Includes polarization effects for K and L1 atomic shells, according to Gavrila (1959, 1961).
49 // For higher shells the L1 cross-section is used.
50 //
51 // The Gavrila photoelectron angular distribution is a complex function which can not be sampled using
52 // the inverse-transform method (James 1980). Instead a more general approach based on the one already
53 // used to sample bremsstrahlung 2BN cross section (G4Generator2BN, Peralta, 2005) was used.
54 //
55 // M. Gavrila, "Relativistic K-Shell Photoeffect", Phys. Rev. 113, 514-526 (1959)
56 // M. Gavrila, "Relativistic L-Shell Photoeffect", Phys. Rev. 124, 1132-1141 (1961)
57 // F. James, Rept. on Prog. in Phys. 43, 1145 (1980)
58 // L. Peralta et al., "A new low-energy bremsstrahlung generator for GEANT4", Radiat. Prot. Dosimetry. 116, 59-64 (2005)
59 //
60 //
61 // -------------------------------------------------------------------
62 //
63 //
64 
66 #include "G4PhysicalConstants.hh"
67 #include "G4RotationMatrix.hh"
68 #include "Randomize.hh"
69 
71  :G4VEmAngularDistribution("AngularGenSauterGavrilaPolarized")
72 {
73  const G4int arrayDim = 980;
74 
75  //minimum electron beta parameter allowed
76  betaArray[0] = 0.02;
77  //beta step
78  betaArray[1] = 0.001;
79  //maximum index array for a and c tables
80  betaArray[2] = arrayDim - 1;
81 
82  // read Majorant Surface Parameters. This are required in order to generate Gavrila angular photoelectron distribution
83  for(G4int level = 0; level < 2; level++){
84 
85  char nameChar0[100] = "ftab0.dat"; // K-shell Majorant Surface Parameters
86  char nameChar1[100] = "ftab1.dat"; // L-shell Majorant Surface Parameters
87 
88  G4String filename;
89  if(level == 0) filename = nameChar0;
90  if(level == 1) filename = nameChar1;
91 
92  char* path = getenv("G4LEDATA");
93  if (!path)
94  {
95  G4String excep = "G4EMDataSet - G4LEDATA environment variable not set";
96  G4Exception("G4PhotoElectricAngularGeneratorPolarized::G4PhotoElectricAngularGeneratorPolarized",
97  "em0006",FatalException,"G4LEDATA environment variable not set");
98  return;
99  }
100 
101  G4String pathString(path);
102  G4String dirFile = pathString + "/photoelectric_angular/" + filename;
103  std::ifstream infile(dirFile);
104  if (!infile.is_open())
105  {
106  G4String excep = "data file: " + dirFile + " not found";
107  G4Exception("G4PhotoElectricAngularGeneratorPolarized::G4PhotoElectricAngularGeneratorPolarized",
108  "em0003",FatalException,excep);
109  return;
110  }
111 
112  // Read parameters into tables. The parameters are function of incident electron energy and shell level
113  G4float aRead=0,cRead=0, beta=0;
114  for(G4int i=0 ; i<arrayDim ;i++){
115  //fscanf(infile,"%f\t %e\t %e",&beta,&aRead,&cRead);
116  infile >> beta >> aRead >> cRead;
117  aMajorantSurfaceParameterTable[i][level] = aRead;
118  cMajorantSurfaceParameterTable[i][level] = cRead;
119  }
120  infile.close();
121  }
122 }
123 
125 {}
126 
129  const G4DynamicParticle* dp,
130  G4double eKinEnergy,
131  G4int shellId,
132  const G4Material*)
133 {
134  // (shellId == 0) - K-shell - Polarized model for K-shell
135  // (shellId > 0) - L1-shell - Polarized model for L1 and higher shells
136 
137  // Calculate Lorentz term (gamma) and beta parameters
138  G4double gamma = 1 + eKinEnergy/electron_mass_c2;
139  G4double beta = std::sqrt((gamma - 1)*(gamma + 1))/gamma;
140 
141  const G4ThreeVector& direction = dp->GetMomentumDirection();
142  const G4ThreeVector& polarization = dp->GetPolarization();
143 
144  G4double theta, phi = 0;
145  // Majorant surface parameters
146  // function of the outgoing electron kinetic energy
147  G4double aBeta = 0;
148  G4double cBeta = 0;
149 
150  // For the outgoing kinetic energy
151  // find the current majorant surface parameters
152  PhotoElectronGetMajorantSurfaceAandCParameters(shellId, beta, &aBeta, &cBeta);
153 
154  // Generate pho and theta according to the shell level
155  // and beta parameter of the electron
156  PhotoElectronGeneratePhiAndTheta(shellId, beta, aBeta, cBeta, &phi, &theta);
157 
158  // Determine the rotation matrix
159  const G4RotationMatrix rotation =
160  PhotoElectronRotationMatrix(direction, polarization);
161 
162  // Compute final direction of the outgoing electron
163  fLocalDirection = PhotoElectronComputeFinalDirection(rotation, theta, phi);
164 
165  return fLocalDirection;
166 }
167 
168 void
169 G4PhotoElectricAngularGeneratorPolarized::PhotoElectronGeneratePhiAndTheta(
170  G4int shellLevel, G4double beta, G4double aBeta, G4double cBeta,
171  G4double *pphi, G4double *ptheta) const
172 {
173  G4double rand1, rand2, rand3 = 0;
174  G4double phi = 0;
175  G4double theta = 0;
176  G4double crossSectionValue = 0;
177  G4double crossSectionMajorantFunctionValue = 0;
178  G4double maxBeta = 0;
179 
180  //G4cout << "shell= " << shellLevel << " beta= " << beta
181  // << " aBeta= " << aBeta << " cBeta= " << cBeta << G4endl;
182 
183  do {
184 
185  rand1 = G4UniformRand();
186  rand2 = G4UniformRand();
187  rand3 = G4UniformRand();
188 
189  phi=2*pi*rand1;
190 
191  if(shellLevel == 0){
192 
193  // Polarized Gavrila Cross-Section for K-shell (1959)
194  theta=std::sqrt(((std::exp(rand2*std::log(1+cBeta*pi*pi)))-1)/cBeta);
195  crossSectionMajorantFunctionValue =
196  CrossSectionMajorantFunction(theta, cBeta);
197  crossSectionValue = DSigmaKshellGavrila1959(beta, theta, phi);
198 
199  } else {
200 
201  // Polarized Gavrila Cross-Section for other shells (L1-shell) (1961)
202  theta = std::sqrt(((std::exp(rand2*std::log(1+cBeta*pi*pi)))-1)/cBeta);
203  crossSectionMajorantFunctionValue =
204  CrossSectionMajorantFunction(theta, cBeta);
205  crossSectionValue = DSigmaL1shellGavrila(beta, theta, phi);
206 
207  }
208 
209  maxBeta=rand3*aBeta*crossSectionMajorantFunctionValue;
210  //G4cout << " crossSectionValue= " << crossSectionValue
211  // << " max= " << maxBeta << G4endl;
212  if(crossSectionValue < 0.0) { crossSectionValue = maxBeta; }
213 
214  } while(maxBeta > crossSectionValue || theta > CLHEP::pi);
215 
216  *pphi = phi;
217  *ptheta = theta;
218 }
219 
220 G4double
221 G4PhotoElectricAngularGeneratorPolarized::CrossSectionMajorantFunction(
222  G4double theta, G4double cBeta) const
223 {
224  // Compute Majorant Function
225  G4double crossSectionMajorantFunctionValue = 0;
226  crossSectionMajorantFunctionValue = theta/(1+cBeta*theta*theta);
227  return crossSectionMajorantFunctionValue;
228 }
229 
230 G4double
231 G4PhotoElectricAngularGeneratorPolarized::DSigmaKshellGavrila1959(
232  G4double beta, G4double theta, G4double phi) const
233 {
234  //Double differential K shell cross-section (Gavrila 1959)
235 
236  G4double beta2 = beta*beta;
237  G4double oneBeta2 = 1 - beta2;
238  G4double sqrtOneBeta2 = std::sqrt(oneBeta2);
239  G4double oneBeta2_to_3_2 = std::pow(oneBeta2,1.5);
240  G4double cosTheta = std::cos(theta);
241  G4double sinTheta2 = std::sin(theta)*std::sin(theta);
242  G4double cosPhi2 = std::cos(phi)*std::cos(phi);
243  G4double oneBetaCosTheta = 1-beta*cosTheta;
244  G4double dsigma = 0;
245  G4double firstTerm = 0;
246  G4double secondTerm = 0;
247 
248  firstTerm = sinTheta2*cosPhi2/std::pow(oneBetaCosTheta,4)-(1 - sqrtOneBeta2)/(2*oneBeta2) *
249  (sinTheta2 * cosPhi2)/std::pow(oneBetaCosTheta,3) + (1-sqrtOneBeta2)*
250  (1-sqrtOneBeta2)/(4*oneBeta2_to_3_2) * sinTheta2/std::pow(oneBetaCosTheta,3);
251 
252  secondTerm = std::sqrt(1 - sqrtOneBeta2)/(std::pow(2.,3.5)*beta2*std::pow(oneBetaCosTheta,2.5)) *
253  (4*beta2/sqrtOneBeta2 * sinTheta2*cosPhi2/oneBetaCosTheta + 4*beta/oneBeta2 * cosTheta * cosPhi2
254  - 4*(1-sqrtOneBeta2)/oneBeta2 *(1+cosPhi2) - beta2 * (1-sqrtOneBeta2)/oneBeta2 * sinTheta2/oneBetaCosTheta
255  + 4*beta2*(1-sqrtOneBeta2)/oneBeta2_to_3_2 - 4*beta*(1-sqrtOneBeta2)*(1-sqrtOneBeta2)/oneBeta2_to_3_2 * cosTheta)
256  + (1-sqrtOneBeta2)/(4*beta2*oneBetaCosTheta*oneBetaCosTheta) * (beta/oneBeta2 - 2/oneBeta2 * cosTheta * cosPhi2 +
257  (1-sqrtOneBeta2)/oneBeta2_to_3_2 * cosTheta - beta * (1-sqrtOneBeta2)/oneBeta2_to_3_2);
258 
259  dsigma = ( firstTerm*(1-pi*fine_structure_const/beta) + secondTerm*(pi*fine_structure_const) );
260 
261  return dsigma;
262 }
263 
264 //
265 
266 G4double
267 G4PhotoElectricAngularGeneratorPolarized::DSigmaL1shellGavrila(
268  G4double beta, G4double theta, G4double phi) const
269 {
270  //Double differential L1 shell cross-section (Gavrila 1961)
271 
272  G4double beta2 = beta*beta;
273  G4double oneBeta2 = 1-beta2;
274  G4double sqrtOneBeta2 = std::sqrt(oneBeta2);
275  G4double oneBeta2_to_3_2=std::pow(oneBeta2,1.5);
276  G4double cosTheta = std::cos(theta);
277  G4double sinTheta2 =std::sin(theta)*std::sin(theta);
278  G4double cosPhi2 = std::cos(phi)*std::cos(phi);
279  G4double oneBetaCosTheta = 1-beta*cosTheta;
280 
281  G4double dsigma = 0;
282  G4double firstTerm = 0;
283  G4double secondTerm = 0;
284 
285  firstTerm = sinTheta2*cosPhi2/std::pow(oneBetaCosTheta,4)-(1 - sqrtOneBeta2)/(2*oneBeta2)
286  * (sinTheta2 * cosPhi2)/std::pow(oneBetaCosTheta,3) + (1-sqrtOneBeta2)*
287  (1-sqrtOneBeta2)/(4*oneBeta2_to_3_2) * sinTheta2/std::pow(oneBetaCosTheta,3);
288 
289  secondTerm = std::sqrt(1 - sqrtOneBeta2)/(std::pow(2.,3.5)*beta2*std::pow(oneBetaCosTheta,2.5)) *
290  (4*beta2/sqrtOneBeta2 * sinTheta2*cosPhi2/oneBetaCosTheta + 4*beta/oneBeta2 * cosTheta * cosPhi2
291  - 4*(1-sqrtOneBeta2)/oneBeta2 *(1+cosPhi2) - beta2 * (1-sqrtOneBeta2)/oneBeta2 * sinTheta2/oneBetaCosTheta
292  + 4*beta2*(1-sqrtOneBeta2)/oneBeta2_to_3_2 - 4*beta*(1-sqrtOneBeta2)*(1-sqrtOneBeta2)/oneBeta2_to_3_2 * cosTheta)
293  + (1-sqrtOneBeta2)/(4*beta2*oneBetaCosTheta*oneBetaCosTheta) * (beta/oneBeta2 - 2/oneBeta2 * cosTheta * cosPhi2 +
294  (1-sqrtOneBeta2)/oneBeta2_to_3_2*cosTheta - beta*(1-sqrtOneBeta2)/oneBeta2_to_3_2);
295 
296  dsigma = ( firstTerm*(1-pi*fine_structure_const/beta) + secondTerm*(pi*fine_structure_const) );
297 
298  return dsigma;
299 }
300 
302 G4PhotoElectricAngularGeneratorPolarized::PhotoElectronRotationMatrix(
303  const G4ThreeVector& direction,
304  const G4ThreeVector& polarization)
305 {
306  G4double mK = direction.mag();
307  G4double mS = polarization.mag();
308  G4ThreeVector polarization2 = polarization;
309  const G4double kTolerance = 1e-6;
310 
311  if(!(polarization.isOrthogonal(direction,kTolerance)) || mS == 0){
312  G4ThreeVector d0 = direction.unit();
314  G4ThreeVector a0 = a1.unit();
315  G4double rand1 = G4UniformRand();
316  G4double angle = twopi*rand1;
317  G4ThreeVector b0 = d0.cross(a0);
319  c.setX(std::cos(angle)*(a0.x())+std::sin(angle)*b0.x());
320  c.setY(std::cos(angle)*(a0.y())+std::sin(angle)*b0.y());
321  c.setZ(std::cos(angle)*(a0.z())+std::sin(angle)*b0.z());
322  polarization2 = c.unit();
323  mS = polarization2.mag();
324  }else
325  {
326  if ( polarization.howOrthogonal(direction) != 0)
327  {
328  polarization2 = polarization
329  - polarization.dot(direction)/direction.dot(direction) * direction;
330  }
331  }
332 
333  G4ThreeVector direction2 = direction/mK;
334  polarization2 = polarization2/mS;
335 
336  G4ThreeVector y = direction2.cross(polarization2);
337 
338  G4RotationMatrix R(polarization2,y,direction2);
339  return R;
340 }
341 
342 void
343 G4PhotoElectricAngularGeneratorPolarized::PhotoElectronGetMajorantSurfaceAandCParameters(G4int shellId, G4double beta, G4double *majorantSurfaceParameterA, G4double *majorantSurfaceParameterC) const
344 {
345  // This member function finds for a given shell and beta value
346  // of the outgoing electron the correct Majorant Surface parameters
347 
348  G4double aBeta,cBeta;
349  G4double bMin,bStep;
350  G4int indexMax;
351  G4int level = 0;
352  if(shellId > 0) { level = 1; }
353 
354  bMin = betaArray[0];
355  bStep = betaArray[1];
356  indexMax = (G4int)betaArray[2];
357  const G4double kBias = 1e-9;
358 
359  G4int k = (G4int)((beta-bMin+kBias)/bStep);
360 
361  if(k < 0)
362  k = 0;
363  if(k > indexMax)
364  k = indexMax;
365 
366  if(k == 0)
367  aBeta = std::max(aMajorantSurfaceParameterTable[k][level],aMajorantSurfaceParameterTable[k+1][level]);
368  else if(k==indexMax)
369  aBeta = std::max(aMajorantSurfaceParameterTable[k-1][level],aMajorantSurfaceParameterTable[k][level]);
370  else{
371  aBeta = std::max(aMajorantSurfaceParameterTable[k-1][level],aMajorantSurfaceParameterTable[k][level]);
372  aBeta = std::max(aBeta,aMajorantSurfaceParameterTable[k+1][level]);
373  }
374 
375  if(k == 0)
376  cBeta = std::max(cMajorantSurfaceParameterTable[k][level],cMajorantSurfaceParameterTable[k+1][level]);
377  else if(k == indexMax)
378  cBeta = std::max(cMajorantSurfaceParameterTable[k-1][level],cMajorantSurfaceParameterTable[k][level]);
379  else{
380  cBeta = std::max(cMajorantSurfaceParameterTable[k-1][level],cMajorantSurfaceParameterTable[k][level]);
381  cBeta = std::max(cBeta,cMajorantSurfaceParameterTable[k+1][level]);
382  }
383 
384  *majorantSurfaceParameterA = aBeta;
385  *majorantSurfaceParameterC = cBeta;
386 }
387 
388 G4ThreeVector G4PhotoElectricAngularGeneratorPolarized::PhotoElectronComputeFinalDirection(const G4RotationMatrix& rotation, G4double theta, G4double phi) const
389 {
390  //computes the photoelectron momentum unitary vector
391  G4double sint = std::sin(theta);
392  G4double px = std::cos(phi)*sint;
393  G4double py = std::sin(phi)*sint;
394  G4double pz = std::cos(theta);
395 
396  G4ThreeVector samplingDirection(px,py,pz);
397 
398  G4ThreeVector outgoingDirection = rotation*samplingDirection;
399  return outgoingDirection;
400 }
401 
403 {
404  G4cout << "\n" << G4endl;
405  G4cout << "Polarized Photoelectric Angular Generator" << G4endl;
406  G4cout << "PhotoElectric Electron Angular Generator based on the general Gavrila photoelectron angular distribution" << G4endl;
407  G4cout << "Includes polarization effects for K and L1 atomic shells, according to Gavrilla (1959, 1961)." << G4endl;
408  G4cout << "For higher shells the L1 cross-section is used." << G4endl;
409  G4cout << "(see Physics Reference Manual) \n" << G4endl;
410 }
411 
414  const G4ThreeVector& a) const
415 {
416  G4double dx = a.x();
417  G4double dy = a.y();
418  G4double dz = a.z();
419  G4double x = dx < 0.0 ? -dx : dx;
420  G4double y = dy < 0.0 ? -dy : dy;
421  G4double z = dz < 0.0 ? -dz : dz;
422  if (x < y) {
423  return x < z ? G4ThreeVector(-dy,dx,0) : G4ThreeVector(0,-dz,dy);
424  }else{
425  return y < z ? G4ThreeVector(dz,0,-dx) : G4ThreeVector(-dy,dx,0);
426  }
427 }
CLHEP::Hep3Vector G4ThreeVector
double x() const
double dot(const Hep3Vector &) const
G4double z
Definition: TRTMaterials.hh:39
float G4float
Definition: G4Types.hh:77
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
Definition: SpaceVector.cc:237
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double eKinEnergy, G4int shellId, const G4Material *mat=0)
int G4int
Definition: G4Types.hh:78
void setY(double)
double howOrthogonal(const Hep3Vector &v) const
Definition: SpaceVector.cc:219
double z() const
void setZ(double)
void setX(double)
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
function dsigma(XP)
Definition: leptonew.f:6713
float electron_mass_c2
Definition: hepunit.py:274
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
G4ThreeVector PerpendicularVector(const G4ThreeVector &a) const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
Hep3Vector unit() const
double y() const
const G4ThreeVector & GetPolarization() const
#define G4endl
Definition: G4ios.hh:61
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
double mag() const