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