Geant4  10.02
G4PAIySection.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: G4PAIySection.cc 91726 2015-08-03 15:41:36Z gcosmo $
27 //
28 //
29 // G4PAIySection.cc -- class implementation file
30 //
31 // GEANT 4 class implementation file
32 //
33 // For information related to this code, please, contact
34 // the Geant4 Collaboration.
35 //
36 // R&D: Vladimir.Grichine@cern.ch
37 //
38 // History:
39 //
40 // 01.10.07, V.Ivanchenko create using V.Grichine G4PAIxSection class
41 // 26.07.09, V.Ivanchenko added protection for mumerical exceptions for
42 // low-density materials
43 // 21.11.10 V. Grichine bug fixed in Initialise for reading sandia table from
44 // material. Warning: the table is tuned for photo-effect not PAI model.
45 // 23.06.13 V.Grichine arrays->G4DataVectors
46 //
47 
48 #include "G4PAIySection.hh"
49 
50 #include "globals.hh"
51 #include "G4PhysicalConstants.hh"
52 #include "G4SystemOfUnits.hh"
53 #include "G4ios.hh"
54 #include "G4Poisson.hh"
55 #include "G4Material.hh"
56 #include "G4MaterialCutsCouple.hh"
57 #include "G4SandiaTable.hh"
58 
59 using namespace std;
60 
61 // Local class constants
62 
63 const G4double G4PAIySection::fDelta = 0.005; // energy shift from interval border
64 const G4double G4PAIySection::fError = 0.005; // error in lin-log approximation
65 
66 const G4int G4PAIySection::fMaxSplineSize = 500; // Max size of output spline
67  // arrays
68 
69 static const G4double betaBohr = fine_structure_const;
70 static const G4double cofBetaBohr = 4.0;
71 static const G4double betaBohr2 = fine_structure_const*fine_structure_const;
73 
75 //
76 // Constructor
77 //
78 
80 {
81  fSandia = 0;
82  fDensity = fElectronDensity = fNormalizationCof = fLowEnergyCof = 0.0;
83  fIntervalNumber = fSplineNumber = 0;
84  fVerbose = 0;
85 
86  fSplineEnergy = G4DataVector(fMaxSplineSize,0.0);
87  fRePartDielectricConst = G4DataVector(fMaxSplineSize,0.0);
88  fImPartDielectricConst = G4DataVector(fMaxSplineSize,0.0);
89  fIntegralTerm = G4DataVector(fMaxSplineSize,0.0);
90  fDifPAIySection = G4DataVector(fMaxSplineSize,0.0);
91  fdNdxCerenkov = G4DataVector(fMaxSplineSize,0.0);
92  fdNdxPlasmon = G4DataVector(fMaxSplineSize,0.0);
93  fIntegralPAIySection = G4DataVector(fMaxSplineSize,0.0);
94  fIntegralPAIdEdx = G4DataVector(fMaxSplineSize,0.0);
95  fIntegralCerenkov = G4DataVector(fMaxSplineSize,0.0);
96  fIntegralPlasmon = G4DataVector(fMaxSplineSize,0.0);
97 
98  for( G4int i = 0; i < 500; ++i )
99  {
100  for( G4int j = 0; j < 112; ++j ) fPAItable[i][j] = 0.0;
101  }
102 }
103 
105 //
106 // Destructor
107 
109 {}
110 
112 //
113 // Constructor with beta*gamma square value called from G4PAIModel
114 
115 void G4PAIySection::Initialize( const G4Material* material,
116  G4double maxEnergyTransfer,
117  G4double betaGammaSq,
118  G4SandiaTable* sandia)
119 {
120  if(fVerbose > 0)
121  {
122  G4cout<<G4endl;
123  G4cout<<"G4PAIySection::Initialize(...,G4SandiaTable* sandia)"<<G4endl;
124  G4cout<<G4endl;
125  }
126  G4int i, j;
127 
128  fSandia = sandia;
129  fIntervalNumber = sandia->GetMaxInterval();
130  fDensity = material->GetDensity();
131  fElectronDensity = material->GetElectronDensity();
132 
133  // fIntervalNumber--;
134 
135  if( fVerbose > 0 )
136  {
137  G4cout<<"fDensity = "<<fDensity<<"\t"<<fElectronDensity<<"\t fIntervalNumber = "
138  <<fIntervalNumber<< " (beta*gamma)^2= " << betaGammaSq << G4endl;
139  }
140  fEnergyInterval = G4DataVector(fIntervalNumber+2,0.0);
141  fA1 = G4DataVector(fIntervalNumber+2,0.0);
142  fA2 = G4DataVector(fIntervalNumber+2,0.0);
143  fA3 = G4DataVector(fIntervalNumber+2,0.0);
144  fA4 = G4DataVector(fIntervalNumber+2,0.0);
145 
146  for( i = 1; i <= fIntervalNumber; i++ )
147  {
148  if ( sandia->GetSandiaMatTablePAI(i-1,0) < 1.*eV )
149  {
150  fIntervalNumber--;
151  continue;
152  }
153  if( ( sandia->GetSandiaMatTablePAI(i-1,0) >= maxEnergyTransfer )
154  || i >= fIntervalNumber )
155  {
156  fEnergyInterval[i] = maxEnergyTransfer;
157  fIntervalNumber = i;
158  break;
159  }
160  fEnergyInterval[i] = sandia->GetSandiaMatTablePAI(i-1,0);
161  fA1[i] = sandia->GetSandiaMatTablePAI(i-1,1);
162  fA2[i] = sandia->GetSandiaMatTablePAI(i-1,2);
163  fA3[i] = sandia->GetSandiaMatTablePAI(i-1,3);
164  fA4[i] = sandia->GetSandiaMatTablePAI(i-1,4);
165 
166  if( fVerbose > 0 ) {
167  G4cout<<i<<"\t"<<fEnergyInterval[i]/keV<<"\t"<<fA1[i]<<"\t"<<fA2[i]<<"\t"
168  <<fA3[i]<<"\t"<<fA4[i]<<"\t"<<G4endl;
169  }
170  }
171  if( fVerbose > 0 ) {
172  G4cout<<"last i = "<<i<<"; "<<"fIntervalNumber = "
173  <<fIntervalNumber<<G4endl;
174  }
175  if( fEnergyInterval[fIntervalNumber] != maxEnergyTransfer )
176  {
177  fIntervalNumber++;
178  fEnergyInterval[fIntervalNumber] = maxEnergyTransfer;
179  }
180  if( fVerbose > 0 )
181  {
182  for( i = 1; i <= fIntervalNumber; i++ )
183  {
184  G4cout<<i<<"\t"<<fEnergyInterval[i]/keV<<"\t"<<fA1[i]<<"\t"<<fA2[i]<<"\t"
185  <<fA3[i]<<"\t"<<fA4[i]<<"\t"<<G4endl;
186  }
187  }
188  if( fVerbose > 0 ) {
189  G4cout<<"Now checking, if two borders are too close together"<<G4endl;
190  }
191  for( i = 1; i < fIntervalNumber; i++ )
192  {
193  if( fEnergyInterval[i+1]-fEnergyInterval[i] >
194  1.5*fDelta*(fEnergyInterval[i+1]+fEnergyInterval[i]) ) continue;
195  else
196  {
197  for( j = i; j < fIntervalNumber; j++ )
198  {
199  fEnergyInterval[j] = fEnergyInterval[j+1];
200  fA1[j] = fA1[j+1];
201  fA2[j] = fA2[j+1];
202  fA3[j] = fA3[j+1];
203  fA4[j] = fA4[j+1];
204  }
205  fIntervalNumber--;
206  }
207  }
208  if( fVerbose > 0 )
209  {
210  for( i = 1; i <= fIntervalNumber; i++ )
211  {
212  G4cout<<i<<"\t"<<fEnergyInterval[i]/keV<<"\t"<<fA1[i]<<"\t"<<fA2[i]<<"\t"
213  <<fA3[i]<<"\t"<<fA4[i]<<"\t"<<G4endl;
214  }
215  }
216  // Preparation of fSplineEnergy array corresponding to min ionisation, G~4
217 
218  ComputeLowEnergyCof(material);
219 
220  G4double betaGammaSqRef =
221  fLorentzFactor[fRefGammaNumber]*fLorentzFactor[fRefGammaNumber] - 1;
222 
223  NormShift(betaGammaSqRef);
224  SplainPAI(betaGammaSqRef);
225 
226  // Preparation of integral PAI cross section for input betaGammaSq
227 
228  for( i = 1; i <= fSplineNumber; i++ )
229  {
230  fDifPAIySection[i] = DifPAIySection(i,betaGammaSq);
231 
232  if( fVerbose > 0 ) G4cout<<i<<"; dNdxPAI = "<<fDifPAIySection[i]<<G4endl;
233  }
234  IntegralPAIySection();
235 }
236 
238 //
239 // Compute low energy cof. It reduces PAI xsc for Lorentz factors less than 4.
240 //
241 
243 {
244  G4int i, numberOfElements = material->GetNumberOfElements();
245  G4double sumZ = 0., sumCof = 0.;
246 
247  static const G4double p0 = 1.20923e+00;
248  static const G4double p1 = 3.53256e-01;
249  static const G4double p2 = -1.45052e-03;
250 
251  G4double* thisMaterialZ = new G4double[numberOfElements];
252  G4double* thisMaterialCof = new G4double[numberOfElements];
253 
254  for( i = 0; i < numberOfElements; i++ )
255  {
256  thisMaterialZ[i] = material->GetElement(i)->GetZ();
257  sumZ += thisMaterialZ[i];
258  thisMaterialCof[i] = p0+p1*thisMaterialZ[i]+p2*thisMaterialZ[i]*thisMaterialZ[i];
259  }
260  for( i = 0; i < numberOfElements; i++ )
261  {
262  sumCof += thisMaterialCof[i]*thisMaterialZ[i]/sumZ;
263  }
264  fLowEnergyCof = sumCof;
265  delete [] thisMaterialZ;
266  delete [] thisMaterialCof;
267  // G4cout<<"fLowEnergyCof = "<<fLowEnergyCof<<G4endl;
268 }
269 
271 //
272 // General control function for class G4PAIySection
273 //
274 
276 {
277  G4int i;
278  G4double betaGammaSq = fLorentzFactor[fRefGammaNumber]*
279  fLorentzFactor[fRefGammaNumber] - 1;
280 
281  // Preparation of integral PAI cross section for reference gamma
282 
283  NormShift(betaGammaSq);
284  SplainPAI(betaGammaSq);
285 
286  IntegralPAIySection();
287  IntegralCerenkov();
288  IntegralPlasmon();
289 
290  for( i = 0; i<= fSplineNumber; i++)
291  {
292  fPAItable[i][fRefGammaNumber] = fIntegralPAIySection[i];
293 
294  if(i != 0) fPAItable[i][0] = fSplineEnergy[i];
295  }
296  fPAItable[0][0] = fSplineNumber;
297 
298  for( G4int j = 1; j < 112; j++) // for other gammas
299  {
300  if( j == fRefGammaNumber ) continue;
301 
302  betaGammaSq = fLorentzFactor[j]*fLorentzFactor[j] - 1;
303 
304  for(i = 1; i <= fSplineNumber; i++)
305  {
306  fDifPAIySection[i] = DifPAIySection(i,betaGammaSq);
307  fdNdxCerenkov[i] = PAIdNdxCerenkov(i,betaGammaSq);
308  fdNdxPlasmon[i] = PAIdNdxPlasmon(i,betaGammaSq);
309  }
310  IntegralPAIySection();
311  IntegralCerenkov();
312  IntegralPlasmon();
313 
314  for(i = 0; i <= fSplineNumber; i++)
315  {
316  fPAItable[i][j] = fIntegralPAIySection[i];
317  }
318  }
319 }
320 
322 //
323 // Shifting from borders to intervals Creation of first energy points
324 //
325 
327 {
328  G4int i, j;
329 
330  for( i = 1; i <= fIntervalNumber-1; i++ )
331  {
332  for( j = 1; j <= 2; j++ )
333  {
334  fSplineNumber = (i-1)*2 + j;
335 
336  if( j == 1 ) fSplineEnergy[fSplineNumber] = fEnergyInterval[i ]*(1+fDelta);
337  else fSplineEnergy[fSplineNumber] = fEnergyInterval[i+1]*(1-fDelta);
338  // G4cout<<"cn = "<<fSplineNumber<<"; "<<"energy = "
339  // <<fSplineEnergy[fSplineNumber]<<G4endl;
340  }
341  }
342  fIntegralTerm[1]=RutherfordIntegral(1,fEnergyInterval[1],fSplineEnergy[1]);
343 
344  j = 1;
345 
346  for(i=2;i<=fSplineNumber;i++)
347  {
348  if(fSplineEnergy[i]<fEnergyInterval[j+1])
349  {
350  fIntegralTerm[i] = fIntegralTerm[i-1] +
351  RutherfordIntegral(j,fSplineEnergy[i-1],
352  fSplineEnergy[i] );
353  }
354  else
355  {
356  G4double x = RutherfordIntegral(j,fSplineEnergy[i-1],
357  fEnergyInterval[j+1] );
358  j++;
359  fIntegralTerm[i] = fIntegralTerm[i-1] + x +
360  RutherfordIntegral(j,fEnergyInterval[j],
361  fSplineEnergy[i] );
362  }
363  // G4cout<<i<<"\t"<<fSplineEnergy[i]<<"\t"<<fIntegralTerm[i]<<"\n"<<G4endl;
364  }
365  fNormalizationCof = 2*pi*pi*hbarc*hbarc*fine_structure_const/electron_mass_c2;
366  fNormalizationCof *= fElectronDensity/fIntegralTerm[fSplineNumber];
367 
368  // G4cout<<"fNormalizationCof = "<<fNormalizationCof<<G4endl;
369 
370  // Calculation of PAI differrential cross-section (1/(keV*cm))
371  // in the energy points near borders of energy intervals
372 
373  for(G4int k=1;k<=fIntervalNumber-1;k++)
374  {
375  for(j=1;j<=2;j++)
376  {
377  i = (k-1)*2 + j;
378  fImPartDielectricConst[i] = fNormalizationCof*
379  ImPartDielectricConst(k,fSplineEnergy[i]);
380  fRePartDielectricConst[i] = fNormalizationCof*
381  RePartDielectricConst(fSplineEnergy[i]);
382  fIntegralTerm[i] *= fNormalizationCof;
383 
384  fDifPAIySection[i] = DifPAIySection(i,betaGammaSq);
385  fdNdxCerenkov[i] = PAIdNdxCerenkov(i,betaGammaSq);
386  fdNdxPlasmon[i] = PAIdNdxPlasmon(i,betaGammaSq);
387  }
388  }
389 
390 } // end of NormShift
391 
393 //
394 // Creation of new energy points as geometrical mean of existing
395 // one, calculation PAI_cs for them, while the error of logarithmic
396 // linear approximation would be smaller than 'fError'
397 
399 {
400  G4int k = 1;
401  G4int i = 1;
402 
403  while ( (i < fSplineNumber) && (fSplineNumber < fMaxSplineSize-1) )
404  {
405  if(fSplineEnergy[i+1] > fEnergyInterval[k+1])
406  {
407  k++; // Here next energy point is in next energy interval
408  i++;
409  continue;
410  }
411  // Shifting of arrayes for inserting the geometrical
412  // average of 'i' and 'i+1' energy points to 'i+1' place
413  fSplineNumber++;
414 
415  for(G4int j = fSplineNumber; j >= i+2; j-- )
416  {
417  fSplineEnergy[j] = fSplineEnergy[j-1];
418  fImPartDielectricConst[j] = fImPartDielectricConst[j-1];
419  fRePartDielectricConst[j] = fRePartDielectricConst[j-1];
420  fIntegralTerm[j] = fIntegralTerm[j-1];
421 
422  fDifPAIySection[j] = fDifPAIySection[j-1];
423  fdNdxCerenkov[j] = fdNdxCerenkov[j-1];
424  fdNdxPlasmon[j] = fdNdxPlasmon[j-1];
425  }
426  G4double x1 = fSplineEnergy[i];
427  G4double x2 = fSplineEnergy[i+1];
428  G4double yy1 = fDifPAIySection[i];
429  G4double y2 = fDifPAIySection[i+1];
430 
431  G4double en1 = sqrt(x1*x2);
432  fSplineEnergy[i+1] = en1;
433 
434  // Calculation of logarithmic linear approximation
435  // in this (enr) energy point, which number is 'i+1' now
436 
437  G4double a = log10(y2/yy1)/log10(x2/x1);
438  G4double b = log10(yy1) - a*log10(x1);
439  G4double y = a*log10(en1) + b;
440  y = pow(10.,y);
441 
442  // Calculation of the PAI dif. cross-section at this point
443 
444  fImPartDielectricConst[i+1] = fNormalizationCof*
445  ImPartDielectricConst(k,fSplineEnergy[i+1]);
446  fRePartDielectricConst[i+1] = fNormalizationCof*
447  RePartDielectricConst(fSplineEnergy[i+1]);
448  fIntegralTerm[i+1] = fIntegralTerm[i] + fNormalizationCof*
449  RutherfordIntegral(k,fSplineEnergy[i],
450  fSplineEnergy[i+1]);
451 
452  fDifPAIySection[i+1] = DifPAIySection(i+1,betaGammaSq);
453  fdNdxCerenkov[i+1] = PAIdNdxCerenkov(i+1,betaGammaSq);
454  fdNdxPlasmon[i+1] = PAIdNdxPlasmon(i+1,betaGammaSq);
455 
456  // Condition for next division of this segment or to pass
457  // to higher energies
458 
459  G4double x = 2*(fDifPAIySection[i+1] - y)/(fDifPAIySection[i+1] + y);
460 
461  G4double delta = 2.*(fSplineEnergy[i+1]-fSplineEnergy[i])/(fSplineEnergy[i+1]+fSplineEnergy[i]);
462 
463  if( x < 0 )
464  {
465  x = -x;
466  }
467  if( x > fError && fSplineNumber < fMaxSplineSize-1 && delta > 2.*fDelta )
468  {
469  continue; // next division
470  }
471  i += 2; // pass to next segment
472 
473  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
474  } // close 'while'
475 
476 } // end of SplainPAI
477 
478 
480 //
481 // Integration over electrons that could be considered
482 // quasi-free at energy transfer of interest
483 
485  G4double x1,
486  G4double x2 )
487 {
488  G4double c1, c2, c3;
489  // G4cout<<"RI: x1 = "<<x1<<"; "<<"x2 = "<<x2<<G4endl;
490  c1 = (x2 - x1)/x1/x2;
491  c2 = (x2 - x1)*(x2 + x1)/x1/x1/x2/x2;
492  c3 = (x2 - x1)*(x1*x1 + x1*x2 + x2*x2)/x1/x1/x1/x2/x2/x2;
493  // G4cout<<" RI: c1 = "<<c1<<"; "<<"c2 = "<<c2<<"; "<<"c3 = "<<c3<<G4endl;
494 
495  return fA1[k]*log(x2/x1) + fA2[k]*c1 + fA3[k]*c2/2 + fA4[k]*c3/3;
496 
497 } // end of RutherfordIntegral
498 
499 
501 //
502 // Imaginary part of dielectric constant
503 // (G4int k - interval number, G4double en1 - energy point)
504 
506  G4double energy1 )
507 {
508  G4double energy2,energy3,energy4,result;
509 
510  energy2 = energy1*energy1;
511  energy3 = energy2*energy1;
512  energy4 = energy3*energy1;
513 
514  result = fA1[k]/energy1+fA2[k]/energy2+fA3[k]/energy3+fA4[k]/energy4;
515  result *=hbarc/energy1;
516 
517  return result;
518 
519 } // end of ImPartDielectricConst
520 
521 
523 //
524 // Real part of dielectric constant minus unit: epsilon_1 - 1
525 // (G4double enb - energy point)
526 //
527 
529 {
530  G4double x0, x02, x03, x04, x05, x1, x2, xx1 ,xx2 , xx12,
531  c1, c2, c3, cof1, cof2, xln1, xln2, xln3, result;
532 
533  x0 = enb;
534  result = 0;
535 
536  for(G4int i=1;i<=fIntervalNumber-1;i++)
537  {
538  x1 = fEnergyInterval[i];
539  x2 = fEnergyInterval[i+1];
540  xx1 = x1 - x0;
541  xx2 = x2 - x0;
542  xx12 = xx2/xx1;
543 
544  if(xx12<0)
545  {
546  xx12 = -xx12;
547  }
548  xln1 = log(x2/x1);
549  xln2 = log(xx12);
550  xln3 = log((x2 + x0)/(x1 + x0));
551  x02 = x0*x0;
552  x03 = x02*x0;
553  x04 = x03*x0;
554  x05 = x04*x0;
555  c1 = (x2 - x1)/x1/x2;
556  c2 = (x2 - x1)*(x2 +x1)/x1/x1/x2/x2;
557  c3 = (x2 -x1)*(x1*x1 + x1*x2 + x2*x2)/x1/x1/x1/x2/x2/x2;
558 
559  result -= (fA1[i]/x02 + fA3[i]/x04)*xln1;
560  result -= (fA2[i]/x02 + fA4[i]/x04)*c1;
561  result -= fA3[i]*c2/2/x02;
562  result -= fA4[i]*c3/3/x02;
563 
564  cof1 = fA1[i]/x02 + fA3[i]/x04;
565  cof2 = fA2[i]/x03 + fA4[i]/x05;
566 
567  result += 0.5*(cof1 +cof2)*xln2;
568  result += 0.5*(cof1 - cof2)*xln3;
569  }
570  result *= 2*hbarc/pi;
571 
572  return result;
573 
574 } // end of RePartDielectricConst
575 
577 //
578 // PAI differential cross-section in terms of
579 // simplified Allison's equation
580 //
581 
583  G4double betaGammaSq )
584 {
585  G4double beta, be2,cof,x1,x2,x3,x4,x5,x6,x7,x8,result;
586  //G4double beta, be4;
587  //G4double be4;
588  // G4double betaBohr2 = fine_structure_const*fine_structure_const;
589  // G4double betaBohr4 = betaBohr2*betaBohr2*4.0;
590  be2 = betaGammaSq/(1 + betaGammaSq);
591  //be4 = be2*be2;
592  beta = sqrt(be2);
593  cof = 1;
594  x1 = log(2*electron_mass_c2/fSplineEnergy[i]);
595 
596  if( betaGammaSq < 0.01 ) x2 = log(be2);
597  else
598  {
599  x2 = -log( (1/betaGammaSq - fRePartDielectricConst[i])*
600  (1/betaGammaSq - fRePartDielectricConst[i]) +
601  fImPartDielectricConst[i]*fImPartDielectricConst[i] )/2;
602  }
603  if( fImPartDielectricConst[i] == 0.0 ||betaGammaSq < 0.01 )
604  {
605  x6=0;
606  }
607  else
608  {
609  x3 = -fRePartDielectricConst[i] + 1/betaGammaSq;
610  x5 = -1 - fRePartDielectricConst[i] +
611  be2*((1 +fRePartDielectricConst[i])*(1 + fRePartDielectricConst[i]) +
612  fImPartDielectricConst[i]*fImPartDielectricConst[i]);
613 
614  x7 = atan2(fImPartDielectricConst[i],x3);
615  x6 = x5 * x7;
616  }
617  // if(fImPartDielectricConst[i] == 0) x6 = 0;
618 
619  x4 = ((x1 + x2)*fImPartDielectricConst[i] + x6)/hbarc;
620  // if( x4 < 0.0 ) x4 = 0.0;
621  x8 = (1 + fRePartDielectricConst[i])*(1 + fRePartDielectricConst[i]) +
622  fImPartDielectricConst[i]*fImPartDielectricConst[i];
623 
624  result = (x4 + cof*fIntegralTerm[i]/fSplineEnergy[i]/fSplineEnergy[i]);
625  if(result < 1.0e-8) result = 1.0e-8;
626  result *= fine_structure_const/be2/pi;
627  // low energy correction
628 
629  G4double lowCof = fLowEnergyCof; // 6.0 ; // Ar ~ 4.; -> fLowCof as f(Z1,Z2)?
630 
631  result *= (1 - exp(-beta/betaBohr/lowCof));
632 
633  // result *= (1-exp(-beta/betaBohr))*(1-exp(-beta/betaBohr));
634  // result *= (1-exp(-be2/betaBohr2));
635  // result *= (1-exp(-be4/betaBohr4));
636  // if(fDensity >= 0.1)
637  if(x8 > 0.)
638  {
639  result /= x8;
640  }
641  return result;
642 
643 } // end of DifPAIySection
644 
646 //
647 // Calculation od dN/dx of collisions with creation of Cerenkov pseudo-photons
648 
650  G4double betaGammaSq )
651 {
652  G4double logarithm, x3, x5, argument, modul2, dNdxC;
653  G4double be2, be4;
654 
655  //G4double cof = 1.0;
656 
657  be2 = betaGammaSq/(1 + betaGammaSq);
658  be4 = be2*be2;
659 
660  if( betaGammaSq < 0.01 ) logarithm = log(1.0+betaGammaSq); // 0.0;
661  else
662  {
663  logarithm = -log( (1/betaGammaSq - fRePartDielectricConst[i])*
664  (1/betaGammaSq - fRePartDielectricConst[i]) +
665  fImPartDielectricConst[i]*fImPartDielectricConst[i] )*0.5;
666  logarithm += log(1+1.0/betaGammaSq);
667  }
668 
669  if( fImPartDielectricConst[i] == 0.0 || betaGammaSq < 0.01 )
670  {
671  argument = 0.0;
672  }
673  else
674  {
675  x3 = -fRePartDielectricConst[i] + 1.0/betaGammaSq;
676  x5 = -1.0 - fRePartDielectricConst[i] +
677  be2*((1.0 +fRePartDielectricConst[i])*(1.0 + fRePartDielectricConst[i]) +
678  fImPartDielectricConst[i]*fImPartDielectricConst[i]);
679  if( x3 == 0.0 ) argument = 0.5*pi;
680  else argument = atan2(fImPartDielectricConst[i],x3);
681  argument *= x5 ;
682  }
683  dNdxC = ( logarithm*fImPartDielectricConst[i] + argument )/hbarc;
684 
685  if(dNdxC < 1.0e-8) dNdxC = 1.0e-8;
686 
687  dNdxC *= fine_structure_const/be2/pi;
688 
689  dNdxC *= (1-exp(-be4/betaBohr4));
690 
691  // if(fDensity >= 0.1)
692  // {
693  modul2 = (1.0 + fRePartDielectricConst[i])*(1.0 + fRePartDielectricConst[i]) +
694  fImPartDielectricConst[i]*fImPartDielectricConst[i];
695  if(modul2 > 0.)
696  {
697  dNdxC /= modul2;
698  }
699  return dNdxC;
700 
701 } // end of PAIdNdxCerenkov
702 
704 //
705 // Calculation od dN/dx of collisions with creation of longitudinal EM
706 // excitations (plasmons, delta-electrons)
707 
709  G4double betaGammaSq )
710 {
711  G4double cof, resonance, modul2, dNdxP;
712  G4double be2, be4;
713 
714  cof = 1;
715 
716  be2 = betaGammaSq/(1 + betaGammaSq);
717  be4 = be2*be2;
718 
719  resonance = log(2*electron_mass_c2*be2/fSplineEnergy[i]);
720  resonance *= fImPartDielectricConst[i]/hbarc;
721 
722  dNdxP = ( resonance + cof*fIntegralTerm[i]/fSplineEnergy[i]/fSplineEnergy[i] );
723 
724  if( dNdxP < 1.0e-8 ) dNdxP = 1.0e-8;
725 
726  dNdxP *= fine_structure_const/be2/pi;
727  dNdxP *= (1-exp(-be4/betaBohr4));
728 
729 // if( fDensity >= 0.1 )
730 // {
731  modul2 = (1 + fRePartDielectricConst[i])*(1 + fRePartDielectricConst[i]) +
732  fImPartDielectricConst[i]*fImPartDielectricConst[i];
733  if(modul2 > 0.)
734  {
735  dNdxP /= modul2;
736  }
737  return dNdxP;
738 
739 } // end of PAIdNdxPlasmon
740 
742 //
743 // Calculation of the PAI integral cross-section
744 // fIntegralPAIySection[1] = specific primary ionisation, 1/cm
745 // and fIntegralPAIySection[0] = mean energy loss per cm in keV/cm
746 
748 {
749  fIntegralPAIySection[fSplineNumber] = 0;
750  fIntegralPAIdEdx[fSplineNumber] = 0;
751  fIntegralPAIySection[0] = 0;
752  G4int k = fIntervalNumber -1;
753 
754  for(G4int i = fSplineNumber-1; i >= 1; i--)
755  {
756  if(fSplineEnergy[i] >= fEnergyInterval[k])
757  {
758  fIntegralPAIySection[i] = fIntegralPAIySection[i+1] + SumOverInterval(i);
759  fIntegralPAIdEdx[i] = fIntegralPAIdEdx[i+1] + SumOverIntervaldEdx(i);
760  }
761  else
762  {
763  fIntegralPAIySection[i] = fIntegralPAIySection[i+1] +
764  SumOverBorder(i+1,fEnergyInterval[k]);
765  fIntegralPAIdEdx[i] = fIntegralPAIdEdx[i+1] +
766  SumOverBorderdEdx(i+1,fEnergyInterval[k]);
767  k--;
768  }
769  }
770 } // end of IntegralPAIySection
771 
773 //
774 // Calculation of the PAI Cerenkov integral cross-section
775 // fIntegralCrenkov[1] = specific Crenkov ionisation, 1/cm
776 // and fIntegralCerenkov[0] = mean Cerenkov loss per cm in keV/cm
777 
779 {
780  G4int i, k;
781  fIntegralCerenkov[fSplineNumber] = 0;
782  fIntegralCerenkov[0] = 0;
783  k = fIntervalNumber -1;
784 
785  for( i = fSplineNumber-1; i >= 1; i-- )
786  {
787  if(fSplineEnergy[i] >= fEnergyInterval[k])
788  {
789  fIntegralCerenkov[i] = fIntegralCerenkov[i+1] + SumOverInterCerenkov(i);
790  // G4cout<<"int: i = "<<i<<"; sumC = "<<fIntegralCerenkov[i]<<G4endl;
791  }
792  else
793  {
794  fIntegralCerenkov[i] = fIntegralCerenkov[i+1] +
795  SumOverBordCerenkov(i+1,fEnergyInterval[k]);
796  k--;
797  // G4cout<<"bord: i = "<<i<<"; sumC = "<<fIntegralCerenkov[i]<<G4endl;
798  }
799  }
800 
801 } // end of IntegralCerenkov
802 
804 //
805 // Calculation of the PAI Plasmon integral cross-section
806 // fIntegralPlasmon[1] = splasmon primary ionisation, 1/cm
807 // and fIntegralPlasmon[0] = mean plasmon loss per cm in keV/cm
808 
810 {
811  fIntegralPlasmon[fSplineNumber] = 0;
812  fIntegralPlasmon[0] = 0;
813  G4int k = fIntervalNumber -1;
814  for(G4int i=fSplineNumber-1;i>=1;i--)
815  {
816  if(fSplineEnergy[i] >= fEnergyInterval[k])
817  {
818  fIntegralPlasmon[i] = fIntegralPlasmon[i+1] + SumOverInterPlasmon(i);
819  }
820  else
821  {
822  fIntegralPlasmon[i] = fIntegralPlasmon[i+1] +
823  SumOverBordPlasmon(i+1,fEnergyInterval[k]);
824  k--;
825  }
826  }
827 
828 } // end of IntegralPlasmon
829 
831 //
832 // Calculation the PAI integral cross-section inside
833 // of interval of continuous values of photo-ionisation
834 // cross-section. Parameter 'i' is the number of interval.
835 
837 {
838  G4double x0,x1,y0,yy1,a,b,c,result;
839 
840  x0 = fSplineEnergy[i];
841  x1 = fSplineEnergy[i+1];
842 
843  if( std::abs( 2.*(x1-x0)/(x1+x0) ) < 1.e-6) return 0.;
844 
845  y0 = fDifPAIySection[i];
846  yy1 = fDifPAIySection[i+1];
847  //G4cout << "## x0= " << x0 << " x1= " << x1 << G4endl;
848  c = x1/x0;
849  //G4cout << "c= " << c << " y0= " << y0 << " yy1= " << yy1 << G4endl;
850  a = log10(yy1/y0)/log10(c);
851  //G4cout << "a= " << a << G4endl;
852  // b = log10(y0) - a*log10(x0);
853  b = y0/pow(x0,a);
854  a += 1;
855  if(a == 0)
856  {
857  result = b*log(x1/x0);
858  }
859  else
860  {
861  result = y0*(x1*pow(c,a-1) - x0)/a;
862  }
863  a++;
864  if(a == 0)
865  {
866  fIntegralPAIySection[0] += b*log(x1/x0);
867  }
868  else
869  {
870  fIntegralPAIySection[0] += y0*(x1*x1*pow(c,a-2) - x0*x0)/a;
871  }
872  return result;
873 
874 } // end of SumOverInterval
875 
877 
879 {
880  G4double x0,x1,y0,yy1,a,b,c,result;
881 
882  x0 = fSplineEnergy[i];
883  x1 = fSplineEnergy[i+1];
884 
885  if( std::abs( 2.*(x1-x0)/(x1+x0) ) < 1.e-6) return 0.;
886 
887  y0 = fDifPAIySection[i];
888  yy1 = fDifPAIySection[i+1];
889  c = x1/x0;
890  a = log10(yy1/y0)/log10(c);
891  // b = log10(y0) - a*log10(x0);
892  b = y0/pow(x0,a);
893  a += 2;
894  if(a == 0)
895  {
896  result = b*log(x1/x0);
897  }
898  else
899  {
900  result = y0*(x1*x1*pow(c,a-2) - x0*x0)/a;
901  }
902  return result;
903 
904 } // end of SumOverInterval
905 
907 //
908 // Calculation the PAI Cerenkov integral cross-section inside
909 // of interval of continuous values of photo-ionisation Cerenkov
910 // cross-section. Parameter 'i' is the number of interval.
911 
913 {
914  G4double x0,x1,y0,yy1,a,c,result;
915 
916  x0 = fSplineEnergy[i];
917  x1 = fSplineEnergy[i+1];
918 
919  if( std::abs( 2.*(x1-x0)/(x1+x0) ) < 1.e-6) return 0.;
920 
921  y0 = fdNdxCerenkov[i];
922  yy1 = fdNdxCerenkov[i+1];
923  // G4cout<<"SumC, i = "<<i<<"; x0 ="<<x0<<"; x1 = "<<x1
924  // <<"; y0 = "<<y0<<"; yy1 = "<<yy1<<G4endl;
925 
926  c = x1/x0;
927  a = log10(yy1/y0)/log10(c);
928  G4double b = 0.0;
929  if(a < 20.) b = y0/pow(x0,a);
930 
931  a += 1.0;
932  if(a == 0) result = b*log(c);
933  else result = y0*(x1*pow(c,a-1) - x0)/a;
934  a += 1.0;
935 
936  if( a == 0 ) fIntegralCerenkov[0] += b*log(x1/x0);
937  else fIntegralCerenkov[0] += y0*(x1*x1*pow(c,a-2) - x0*x0)/a;
938  // G4cout<<"a = "<<a<<"; b = "<<b<<"; result = "<<result<<G4endl;
939  return result;
940 
941 } // end of SumOverInterCerenkov
942 
944 //
945 // Calculation the PAI Plasmon integral cross-section inside
946 // of interval of continuous values of photo-ionisation Plasmon
947 // cross-section. Parameter 'i' is the number of interval.
948 
950 {
951  G4double x0,x1,y0,yy1,a,c,result;
952 
953  x0 = fSplineEnergy[i];
954  x1 = fSplineEnergy[i+1];
955 
956  if( std::abs( 2.*(x1-x0)/(x1+x0) ) < 1.e-6) return 0.;
957 
958  y0 = fdNdxPlasmon[i];
959  yy1 = fdNdxPlasmon[i+1];
960  c = x1/x0;
961  a = log10(yy1/y0)/log10(c);
962 
963  G4double b = 0.0;
964  if(a < 20.) b = y0/pow(x0,a);
965 
966  a += 1.0;
967  if(a == 0) result = b*log(x1/x0);
968  else result = y0*(x1*pow(c,a-1) - x0)/a;
969  a += 1.0;
970 
971  if( a == 0 ) fIntegralPlasmon[0] += b*log(x1/x0);
972  else fIntegralPlasmon[0] += y0*(x1*x1*pow(c,a-2) - x0*x0)/a;
973 
974  return result;
975 
976 } // end of SumOverInterPlasmon
977 
979 //
980 // Integration of PAI cross-section for the case of
981 // passing across border between intervals
982 
984  G4double en0 )
985 {
986  G4double x0,x1,y0,yy1,a,/*c,*/d,e0,result;
987 
988  e0 = en0;
989  x0 = fSplineEnergy[i];
990  x1 = fSplineEnergy[i+1];
991  y0 = fDifPAIySection[i];
992  yy1 = fDifPAIySection[i+1];
993 
994  //c = x1/x0;
995  d = e0/x0;
996  a = log10(yy1/y0)/log10(x1/x0);
997 
998  G4double b = 0.0;
999  if(a < 20.) b = y0/pow(x0,a);
1000 
1001  a += 1;
1002  if(a == 0)
1003  {
1004  result = b*log(x0/e0);
1005  }
1006  else
1007  {
1008  result = y0*(x0 - e0*pow(d,a-1))/a;
1009  }
1010  a++;
1011  if(a == 0)
1012  {
1013  fIntegralPAIySection[0] += b*log(x0/e0);
1014  }
1015  else
1016  {
1017  fIntegralPAIySection[0] += y0*(x0*x0 - e0*e0*pow(d,a-2))/a;
1018  }
1019  x0 = fSplineEnergy[i - 1];
1020  x1 = fSplineEnergy[i - 2];
1021  y0 = fDifPAIySection[i - 1];
1022  yy1 = fDifPAIySection[i - 2];
1023 
1024  //c = x1/x0;
1025  d = e0/x0;
1026  a = log10(yy1/y0)/log10(x1/x0);
1027  // b0 = log10(y0) - a*log10(x0);
1028  b = y0/pow(x0,a);
1029  a += 1;
1030  if(a == 0)
1031  {
1032  result += b*log(e0/x0);
1033  }
1034  else
1035  {
1036  result += y0*(e0*pow(d,a-1) - x0)/a;
1037  }
1038  a++;
1039  if(a == 0)
1040  {
1041  fIntegralPAIySection[0] += b*log(e0/x0);
1042  }
1043  else
1044  {
1045  fIntegralPAIySection[0] += y0*(e0*e0*pow(d,a-2) - x0*x0)/a;
1046  }
1047  return result;
1048 
1049 }
1050 
1052 
1054  G4double en0 )
1055 {
1056  G4double x0,x1,y0,yy1,a,/*c,*/d,e0,result;
1057 
1058  e0 = en0;
1059  x0 = fSplineEnergy[i];
1060  x1 = fSplineEnergy[i+1];
1061  y0 = fDifPAIySection[i];
1062  yy1 = fDifPAIySection[i+1];
1063 
1064  //c = x1/x0;
1065  d = e0/x0;
1066  a = log10(yy1/y0)/log10(x1/x0);
1067 
1068  G4double b = 0.0;
1069  if(a < 20.) b = y0/pow(x0,a);
1070 
1071  a += 2;
1072  if(a == 0)
1073  {
1074  result = b*log(x0/e0);
1075  }
1076  else
1077  {
1078  result = y0*(x0*x0 - e0*e0*pow(d,a-2))/a;
1079  }
1080  x0 = fSplineEnergy[i - 1];
1081  x1 = fSplineEnergy[i - 2];
1082  y0 = fDifPAIySection[i - 1];
1083  yy1 = fDifPAIySection[i - 2];
1084 
1085  //c = x1/x0;
1086  d = e0/x0;
1087  a = log10(yy1/y0)/log10(x1/x0);
1088 
1089  if(a < 20.) b = y0/pow(x0,a);
1090 
1091  a += 2;
1092  if(a == 0)
1093  {
1094  result += b*log(e0/x0);
1095  }
1096  else
1097  {
1098  result += y0*(e0*e0*pow(d,a-2) - x0*x0)/a;
1099  }
1100  return result;
1101 
1102 }
1103 
1105 //
1106 // Integration of Cerenkov cross-section for the case of
1107 // passing across border between intervals
1108 
1110  G4double en0 )
1111 {
1112  G4double x0,x1,y0,yy1,a,e0,c,d,result;
1113 
1114  e0 = en0;
1115  x0 = fSplineEnergy[i];
1116  x1 = fSplineEnergy[i+1];
1117  y0 = fdNdxCerenkov[i];
1118  yy1 = fdNdxCerenkov[i+1];
1119 
1120  // G4cout<<G4endl;
1121  //G4cout<<"SumBordC, i = "<<i<<"; en0 = "<<en0<<"; x0 ="<<x0<<"; x1 = "<<x1
1122  // <<"; y0 = "<<y0<<"; yy1 = "<<yy1<<G4endl;
1123  c = x1/x0;
1124  d = e0/x0;
1125  a = log10(yy1/y0)/log10(c);
1126 
1127  G4double b = 0.0;
1128  if(a < 20.) b = y0/pow(x0,a);
1129 
1130  a += 1.0;
1131  if( a == 0 ) result = b*log(x0/e0);
1132  else result = y0*(x0 - e0*pow(d,a-1))/a;
1133  a += 1.0;
1134 
1135  if( a == 0 ) fIntegralCerenkov[0] += b*log(x0/e0);
1136  else fIntegralCerenkov[0] += y0*(x0*x0 - e0*e0*pow(d,a-2))/a;
1137 
1138  //G4cout<<"a = "<<a<<"; b = "<<b<<"; result = "<<result<<G4endl;
1139 
1140  x0 = fSplineEnergy[i - 1];
1141  x1 = fSplineEnergy[i - 2];
1142  y0 = fdNdxCerenkov[i - 1];
1143  yy1 = fdNdxCerenkov[i - 2];
1144 
1145  //G4cout<<"x0 ="<<x0<<"; x1 = "<<x1
1146  // <<"; y0 = "<<y0<<"; yy1 = "<<yy1<<G4endl;
1147 
1148  c = x1/x0;
1149  d = e0/x0;
1150  a = log10(yy1/y0)/log10(x1/x0);
1151 
1152  // G4cout << "a= " << a << G4endl;
1153  if(a < 20.) b = y0/pow(x0,a);
1154 
1155  if(a > 20.0) b = 0.0;
1156  else b = y0/pow(x0,a); // pow(10.,b0);
1157 
1158  //G4cout << "b= " << b << G4endl;
1159 
1160  a += 1.0;
1161  if( a == 0 ) result += b*log(e0/x0);
1162  else result += y0*(e0*pow(d,a-1) - x0 )/a;
1163  a += 1.0;
1164  //G4cout << "result= " << result << G4endl;
1165 
1166  if( a == 0 ) fIntegralCerenkov[0] += b*log(e0/x0);
1167  else fIntegralCerenkov[0] += y0*(e0*e0*pow(d,a-2) - x0*x0)/a;
1168 
1169  //G4cout<<"a = "<<a<<"; b = "<<b<<"; result = "<<result<<G4endl;
1170 
1171  return result;
1172 
1173 }
1174 
1176 //
1177 // Integration of Plasmon cross-section for the case of
1178 // passing across border between intervals
1179 
1181  G4double en0 )
1182 {
1183  G4double x0,x1,y0,yy1,a,c,d,e0,result;
1184 
1185  e0 = en0;
1186  x0 = fSplineEnergy[i];
1187  x1 = fSplineEnergy[i+1];
1188  y0 = fdNdxPlasmon[i];
1189  yy1 = fdNdxPlasmon[i+1];
1190 
1191  c = x1/x0;
1192  d = e0/x0;
1193  a = log10(yy1/y0)/log10(c);
1194 
1195  G4double b = 0.0;
1196  if(a < 20.) b = y0/pow(x0,a);
1197 
1198  a += 1.0;
1199  if( a == 0 ) result = b*log(x0/e0);
1200  else result = y0*(x0 - e0*pow(d,a-1))/a;
1201  a += 1.0;
1202 
1203  if( a == 0 ) fIntegralPlasmon[0] += b*log(x0/e0);
1204  else fIntegralPlasmon[0] += y0*(x0*x0 - e0*e0*pow(d,a-2))/a;
1205 
1206  x0 = fSplineEnergy[i - 1];
1207  x1 = fSplineEnergy[i - 2];
1208  y0 = fdNdxPlasmon[i - 1];
1209  yy1 = fdNdxPlasmon[i - 2];
1210 
1211  c = x1/x0;
1212  d = e0/x0;
1213  a = log10(yy1/y0)/log10(c);
1214 
1215  if(a < 20.) b = y0/pow(x0,a);
1216 
1217  a += 1.0;
1218  if( a == 0 ) result += b*log(e0/x0);
1219  else result += y0*(e0*pow(d,a-1) - x0)/a;
1220  a += 1.0;
1221 
1222  if( a == 0 ) fIntegralPlasmon[0] += b*log(e0/x0);
1223  else fIntegralPlasmon[0] += y0*(e0*e0*pow(d,a-2) - x0*x0)/a;
1224 
1225  return result;
1226 
1227 }
1228 
1230 //
1231 //
1232 
1234 {
1235  G4int iTransfer ;
1236  G4long numOfCollisions;
1237  G4double loss = 0.0;
1238  G4double meanNumber, position;
1239 
1240  // G4cout<<" G4PAIySection::GetStepEnergyLoss "<<G4endl;
1241 
1242 
1243 
1244  meanNumber = fIntegralPAIySection[1]*step;
1245  numOfCollisions = G4Poisson(meanNumber);
1246 
1247  // G4cout<<"numOfCollisions = "<<numOfCollisions<<G4endl;
1248 
1249  while(numOfCollisions)
1250  {
1251  position = fIntegralPAIySection[1]*G4UniformRand();
1252 
1253  for( iTransfer=1; iTransfer<=fSplineNumber; iTransfer++ )
1254  {
1255  if( position >= fIntegralPAIySection[iTransfer] ) break;
1256  }
1257  loss += fSplineEnergy[iTransfer] ;
1258  numOfCollisions--;
1259  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
1260  }
1261  // G4cout<<"PAI energy loss = "<<loss/keV<<" keV"<<G4endl;
1262 
1263  return loss;
1264 }
1265 
1267 //
1268 //
1269 
1271 {
1272  G4int iTransfer ;
1273  G4long numOfCollisions;
1274  G4double loss = 0.0;
1275  G4double meanNumber, position;
1276 
1277  // G4cout<<" G4PAIySection::GetStepCreLosnkovs "<<G4endl;
1278 
1279 
1280 
1281  meanNumber = fIntegralCerenkov[1]*step;
1282  numOfCollisions = G4Poisson(meanNumber);
1283 
1284  // G4cout<<"numOfCollisions = "<<numOfCollisions<<G4endl;
1285 
1286  while(numOfCollisions)
1287  {
1288  position = fIntegralCerenkov[1]*G4UniformRand();
1289 
1290  for( iTransfer=1; iTransfer<=fSplineNumber; iTransfer++ )
1291  {
1292  if( position >= fIntegralCerenkov[iTransfer] ) break;
1293  }
1294  loss += fSplineEnergy[iTransfer] ;
1295  numOfCollisions--;
1296  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
1297  }
1298  // G4cout<<"PAI Cerenkov loss = "<<loss/keV<<" keV"<<G4endl;
1299 
1300  return loss;
1301 }
1302 
1304 //
1305 //
1306 
1308 {
1309  G4int iTransfer ;
1310  G4long numOfCollisions;
1311  G4double loss = 0.0;
1312  G4double meanNumber, position;
1313 
1314  // G4cout<<" G4PAIySection::GetStepCreLosnkovs "<<G4endl;
1315 
1316 
1317 
1318  meanNumber = fIntegralPlasmon[1]*step;
1319  numOfCollisions = G4Poisson(meanNumber);
1320 
1321  // G4cout<<"numOfCollisions = "<<numOfCollisions<<G4endl;
1322 
1323  while(numOfCollisions)
1324  {
1325  position = fIntegralPlasmon[1]*G4UniformRand();
1326 
1327  for( iTransfer=1; iTransfer<=fSplineNumber; iTransfer++ )
1328  {
1329  if( position >= fIntegralPlasmon[iTransfer] ) break;
1330  }
1331  loss += fSplineEnergy[iTransfer] ;
1332  numOfCollisions--;
1333  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
1334  }
1335  // G4cout<<"PAI Plasmon loss = "<<loss/keV<<" keV"<<G4endl;
1336 
1337  return loss;
1338 }
1339 
1341 //
1342 
1343 void G4PAIySection::CallError(G4int i, const G4String& methodName) const
1344 {
1345  G4String head = "G4PAIySection::" + methodName + "()";
1347  ed << "Wrong index " << i << " fSplineNumber= " << fSplineNumber;
1348  G4Exception(head,"pai001",FatalException,ed);
1349 }
1350 
1352 //
1353 // Init array of Lorentz factors
1354 //
1355 
1357 
1358 const G4double G4PAIySection::fLorentzFactor[112] = // fNumberOfGammas+1
1359 {
1360 0.0,
1361 1.094989e+00, 1.107813e+00, 1.122369e+00, 1.138890e+00, 1.157642e+00,
1362 1.178925e+00, 1.203082e+00, 1.230500e+00, 1.261620e+00, 1.296942e+00, // 10
1363 1.337032e+00, 1.382535e+00, 1.434181e+00, 1.492800e+00, 1.559334e+00,
1364 1.634850e+00, 1.720562e+00, 1.817845e+00, 1.928263e+00, 2.053589e+00, // 20
1365 2.195835e+00, 2.357285e+00, 2.540533e+00, 2.748522e+00, 2.984591e+00,
1366 3.252533e+00, 3.556649e+00, 3.901824e+00, 4.293602e+00, 4.738274e+00, // 30
1367 5.242981e+00, 5.815829e+00, 6.466019e+00, 7.203990e+00, 8.041596e+00,
1368 8.992288e+00, 1.007133e+01, 1.129606e+01, 1.268614e+01, 1.426390e+01, // 40
1369 1.605467e+01, 1.808721e+01, 2.039417e+01, 2.301259e+01, 2.598453e+01,
1370 2.935771e+01, 3.318630e+01, 3.753180e+01, 4.246399e+01, 4.806208e+01, // 50
1371 5.441597e+01, 6.162770e+01, 6.981310e+01, 7.910361e+01, 8.964844e+01,
1372 1.016169e+02, 1.152013e+02, 1.306197e+02, 1.481198e+02, 1.679826e+02, // 60
1373 1.905270e+02, 2.161152e+02, 2.451581e+02, 2.781221e+02, 3.155365e+02,
1374 3.580024e+02, 4.062016e+02, 4.609081e+02, 5.230007e+02, 5.934765e+02, // 70
1375 6.734672e+02, 7.642575e+02, 8.673056e+02, 9.842662e+02, 1.117018e+03,
1376 1.267692e+03, 1.438709e+03, 1.632816e+03, 1.853128e+03, 2.103186e+03, // 80
1377 2.387004e+03, 2.709140e+03, 3.074768e+03, 3.489760e+03, 3.960780e+03,
1378 4.495394e+03, 5.102185e+03, 5.790900e+03, 6.572600e+03, 7.459837e+03, // 90
1379 8.466860e+03, 9.609843e+03, 1.090714e+04, 1.237959e+04, 1.405083e+04,
1380 1.594771e+04, 1.810069e+04, 2.054434e+04, 2.331792e+04, 2.646595e+04, // 100
1381 3.003901e+04, 3.409446e+04, 3.869745e+04, 4.392189e+04, 4.985168e+04,
1382 5.658206e+04, 6.422112e+04, 7.289153e+04, 8.273254e+04, 9.390219e+04, // 110
1383 1.065799e+05
1384 };
1385 
1387 //
1388 // The number of gamma for creation of spline (near ion-min , G ~ 4 )
1389 //
1390 
1391 const
1393 
1394 
1395 //
1396 // end of G4PAIySection implementation file
1397 //
1399 
static c2_factory< G4double > c2
G4long G4Poisson(G4double mean)
Definition: G4Poisson.hh:51
G4double SumOverIntervaldEdx(G4int intervalNumber)
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
G4double RutherfordIntegral(G4int intervalNumber, G4double limitLow, G4double limitHigh)
static const G4double fDelta
static const G4int fMaxSplineSize
void NormShift(G4double betaGammaSq)
static const G4double betaBohr
static const G4double fLorentzFactor[112]
G4double GetZ() const
Definition: G4Element.hh:131
G4double ImPartDielectricConst(G4int intervalNumber, G4double energy)
G4double SumOverInterCerenkov(G4int intervalNumber)
static const G4double cofBetaBohr
G4double GetDensity() const
Definition: G4Material.hh:180
long G4long
Definition: G4Types.hh:80
void CallError(G4int i, const G4String &methodName) const
G4int GetMaxInterval() const
G4double PAIdNdxCerenkov(G4int intervalNumber, G4double betaGammaSq)
G4double GetSandiaMatTablePAI(G4int, G4int) const
void ComputeLowEnergyCof(const G4Material *material)
G4double a
Definition: TRTMaterials.hh:39
void IntegralCerenkov()
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
G4double SumOverBorderdEdx(G4int intervalNumber, G4double energy)
int G4int
Definition: G4Types.hh:78
void IntegralPlasmon()
static const G4double betaBohr4
static G4int fNumberOfGammas
#define position
Definition: xmlparse.cc:622
G4double SumOverBordPlasmon(G4int intervalNumber, G4double energy)
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4double GetStepEnergyLoss(G4double step)
G4double GetElectronDensity() const
Definition: G4Material.hh:217
void Initialize(const G4Material *material, G4double maxEnergyTransfer, G4double betaGammaSq, G4SandiaTable *)
static const G4double betaBohr2
static const G4double c3
static const G4double c1
G4double GetStepCerenkovLoss(G4double step)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static const double eV
Definition: G4SIunits.hh:212
static const G4int fRefGammaNumber
static const double pi
Definition: G4SIunits.hh:74
G4double SumOverBorder(G4int intervalNumber, G4double energy)
void SplainPAI(G4double betaGammaSq)
static const G4double fError
const G4double x[NPOINTSGL]
#define G4endl
Definition: G4ios.hh:61
G4double SumOverBordCerenkov(G4int intervalNumber, G4double energy)
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
static const double keV
Definition: G4SIunits.hh:213
double G4double
Definition: G4Types.hh:76
G4double RePartDielectricConst(G4double energy)
G4double GetStepPlasmonLoss(G4double step)
G4double DifPAIySection(G4int intervalNumber, G4double betaGammaSq)
G4double SumOverInterval(G4int intervalNumber)
G4double PAIdNdxPlasmon(G4int intervalNumber, G4double betaGammaSq)
void IntegralPAIySection()
G4double SumOverInterPlasmon(G4int intervalNumber)