G4TransitionRadiation Class Reference

#include <G4TransitionRadiation.hh>

Inheritance diagram for G4TransitionRadiation:

G4VDiscreteProcess G4VProcess G4ForwardXrayTR

Public Member Functions

 G4TransitionRadiation (const G4String &processName="TR", G4ProcessType type=fElectromagnetic)
virtual ~G4TransitionRadiation ()
G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *condition)
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
virtual G4double SpectralAngleTRdensity (G4double energy, G4double varAngle) const =0
G4double IntegralOverEnergy (G4double energy1, G4double energy2, G4double varAngle) const
G4double IntegralOverAngle (G4double energy, G4double varAngle1, G4double varAngle2) const
G4double AngleIntegralDistribution (G4double varAngle1, G4double varAngle2) const
G4double EnergyIntegralDistribution (G4double energy1, G4double energy2) const

Protected Attributes

G4int fMatIndex1
G4int fMatIndex2
G4double fGamma
G4double fEnergy
G4double fVarAngle
G4double fMinEnergy
G4double fMaxEnergy
G4double fMaxTheta
G4double fSigma1
G4double fSigma2

Static Protected Attributes

static const G4int fSympsonNumber = 100
static const G4int fGammaNumber = 15
static const G4int fPointNumber = 100

Detailed Description

Definition at line 57 of file G4TransitionRadiation.hh.


Constructor & Destructor Documentation

G4TransitionRadiation::G4TransitionRadiation ( const G4String processName = "TR",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 60 of file G4TransitionRadiation.cc.

References fEnergy, fGamma, fMatIndex1, fMatIndex2, fMaxEnergy, fMaxTheta, fMinEnergy, fSigma1, fSigma2, fTransitionRadiation, fVarAngle, and G4VProcess::SetProcessSubType().

00061   : G4VDiscreteProcess(processName, type)
00062 {
00063   SetProcessSubType(fTransitionRadiation);
00064   fMatIndex1 = fMatIndex2 = 0;
00065 
00066   fGamma = fEnergy = fVarAngle = fMinEnergy = fMaxEnergy = fMaxTheta = fSigma1 = fSigma2 = 0.0;
00067 }

G4TransitionRadiation::~G4TransitionRadiation (  )  [virtual]

Definition at line 74 of file G4TransitionRadiation.cc.

00075 {}


Member Function Documentation

G4double G4TransitionRadiation::AngleIntegralDistribution ( G4double  varAngle1,
G4double  varAngle2 
) const

Definition at line 156 of file G4TransitionRadiation.cc.

References fMaxEnergy, fMinEnergy, fSympsonNumber, and IntegralOverEnergy().

00158 {
00159   G4int i ;
00160   G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
00161   h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ;
00162   for(i=1;i<fSympsonNumber;i++)
00163   {
00164    sumEven += IntegralOverEnergy(fMinEnergy,
00165                                  fMinEnergy +0.3*(fMaxEnergy-fMinEnergy),
00166                                  varAngle1 + 2*i*h)
00167             + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00168                                  fMaxEnergy,
00169                                  varAngle1 + 2*i*h);
00170    sumOdd  += IntegralOverEnergy(fMinEnergy,
00171                                  fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00172                                  varAngle1 + (2*i - 1)*h)
00173             + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00174                                  fMaxEnergy,
00175                                  varAngle1 + (2*i - 1)*h) ;
00176   }
00177   sumOdd += IntegralOverEnergy(fMinEnergy,
00178                                fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00179                                varAngle1 + (2*fSympsonNumber - 1)*h)
00180           + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00181                                fMaxEnergy,
00182                                varAngle1 + (2*fSympsonNumber - 1)*h) ;
00183 
00184   return h*(IntegralOverEnergy(fMinEnergy,
00185                                fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00186                                varAngle1)
00187           + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00188                                fMaxEnergy,
00189                                varAngle1)
00190           + IntegralOverEnergy(fMinEnergy,
00191                                fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00192                                varAngle2)
00193           + IntegralOverEnergy(fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
00194                                fMaxEnergy,
00195                                varAngle2)
00196             + 4.0*sumOdd + 2.0*sumEven    )/3.0 ;
00197 }

G4double G4TransitionRadiation::EnergyIntegralDistribution ( G4double  energy1,
G4double  energy2 
) const

Definition at line 206 of file G4TransitionRadiation.cc.

References fMaxTheta, fSympsonNumber, and IntegralOverAngle().

00208 {
00209   G4int i ;
00210   G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
00211   h = 0.5*(energy2 - energy1)/fSympsonNumber ;
00212   for(i=1;i<fSympsonNumber;i++)
00213   {
00214    sumEven += IntegralOverAngle(energy1 + 2*i*h,0.0,0.01*fMaxTheta )
00215             + IntegralOverAngle(energy1 + 2*i*h,0.01*fMaxTheta,fMaxTheta);
00216    sumOdd  += IntegralOverAngle(energy1 + (2*i - 1)*h,0.0,0.01*fMaxTheta)
00217             + IntegralOverAngle(energy1 + (2*i - 1)*h,0.01*fMaxTheta,fMaxTheta) ;
00218   }
00219   sumOdd += IntegralOverAngle(energy1 + (2*fSympsonNumber - 1)*h,
00220                               0.0,0.01*fMaxTheta)
00221           + IntegralOverAngle(energy1 + (2*fSympsonNumber - 1)*h,
00222                               0.01*fMaxTheta,fMaxTheta) ;
00223 
00224   return h*(IntegralOverAngle(energy1,0.0,0.01*fMaxTheta)
00225           + IntegralOverAngle(energy1,0.01*fMaxTheta,fMaxTheta)
00226           + IntegralOverAngle(energy2,0.0,0.01*fMaxTheta)
00227           + IntegralOverAngle(energy2,0.01*fMaxTheta,fMaxTheta)
00228             + 4.0*sumOdd + 2.0*sumEven    )/3.0 ;
00229 }

G4double G4TransitionRadiation::GetMeanFreePath ( const G4Track ,
G4double  ,
G4ForceCondition condition 
) [virtual]

Implements G4VDiscreteProcess.

Reimplemented in G4ForwardXrayTR.

Definition at line 83 of file G4TransitionRadiation.cc.

References DBL_MAX, and Forced.

00086 {
00087   *condition = Forced;
00088   return DBL_MAX;      // so TR doesn't limit mean free path
00089 }

G4double G4TransitionRadiation::IntegralOverAngle ( G4double  energy,
G4double  varAngle1,
G4double  varAngle2 
) const

Definition at line 130 of file G4TransitionRadiation.cc.

References fSympsonNumber, and SpectralAngleTRdensity().

Referenced by EnergyIntegralDistribution().

00133 {
00134   G4int i ;
00135   G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
00136   h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ;
00137   for(i=1;i<fSympsonNumber;i++)
00138   {
00139     sumEven += SpectralAngleTRdensity(energy,varAngle1 + 2*i*h)  ;
00140     sumOdd  += SpectralAngleTRdensity(energy,varAngle1 + (2*i - 1)*h) ;
00141   }
00142   sumOdd += SpectralAngleTRdensity(energy,varAngle1 + (2*fSympsonNumber - 1)*h) ;
00143 
00144   return h*(  SpectralAngleTRdensity(energy,varAngle1)
00145             + SpectralAngleTRdensity(energy,varAngle2)
00146             + 4.0*sumOdd + 2.0*sumEven    )/3.0 ;
00147 }

G4double G4TransitionRadiation::IntegralOverEnergy ( G4double  energy1,
G4double  energy2,
G4double  varAngle 
) const

Definition at line 104 of file G4TransitionRadiation.cc.

References fSympsonNumber, and SpectralAngleTRdensity().

Referenced by AngleIntegralDistribution().

00107 {
00108   G4int i ;
00109   G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
00110   h = 0.5*(energy2 - energy1)/fSympsonNumber ;
00111   for(i=1;i<fSympsonNumber;i++)
00112   {
00113     sumEven += SpectralAngleTRdensity(energy1 + 2*i*h,varAngle)  ;
00114     sumOdd  += SpectralAngleTRdensity(energy1 + (2*i - 1)*h,varAngle) ;
00115   }
00116   sumOdd += SpectralAngleTRdensity(energy1 + (2*fSympsonNumber - 1)*h,varAngle) ;
00117   return h*(  SpectralAngleTRdensity(energy1,varAngle)
00118             + SpectralAngleTRdensity(energy2,varAngle)
00119             + 4.0*sumOdd + 2.0*sumEven    )/3.0 ;
00120 }

G4bool G4TransitionRadiation::IsApplicable ( const G4ParticleDefinition aParticleType  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 78 of file G4TransitionRadiation.cc.

References G4ParticleDefinition::GetPDGCharge().

00079 {
00080   return ( aParticleType.GetPDGCharge() != 0.0 );
00081 }

G4VParticleChange * G4TransitionRadiation::PostStepDoIt ( const G4Track ,
const G4Step  
) [virtual]

Reimplemented from G4VDiscreteProcess.

Reimplemented in G4ForwardXrayTR.

Definition at line 91 of file G4TransitionRadiation.cc.

References G4VProcess::aParticleChange, and G4VProcess::ClearNumberOfInteractionLengthLeft().

00093 {
00094   ClearNumberOfInteractionLengthLeft();
00095   return &aParticleChange;
00096 }

virtual G4double G4TransitionRadiation::SpectralAngleTRdensity ( G4double  energy,
G4double  varAngle 
) const [pure virtual]

Implemented in G4ForwardXrayTR.

Referenced by IntegralOverAngle(), and IntegralOverEnergy().


Field Documentation

G4double G4TransitionRadiation::fEnergy [protected]

Definition at line 105 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().

G4double G4TransitionRadiation::fGamma [protected]

Reimplemented in G4ForwardXrayTR.

Definition at line 104 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().

const G4int G4TransitionRadiation::fGammaNumber = 15 [static, protected]

Definition at line 110 of file G4TransitionRadiation.hh.

G4int G4TransitionRadiation::fMatIndex1 [protected]

Definition at line 99 of file G4TransitionRadiation.hh.

Referenced by G4ForwardXrayTR::BuildXrayTRtables(), G4ForwardXrayTR::G4ForwardXrayTR(), G4TransitionRadiation(), and G4ForwardXrayTR::PostStepDoIt().

G4int G4TransitionRadiation::fMatIndex2 [protected]

Definition at line 100 of file G4TransitionRadiation.hh.

Referenced by G4ForwardXrayTR::BuildXrayTRtables(), G4ForwardXrayTR::G4ForwardXrayTR(), G4TransitionRadiation(), and G4ForwardXrayTR::PostStepDoIt().

G4double G4TransitionRadiation::fMaxEnergy [protected]

Definition at line 114 of file G4TransitionRadiation.hh.

Referenced by AngleIntegralDistribution(), and G4TransitionRadiation().

G4double G4TransitionRadiation::fMaxTheta [protected]

Definition at line 115 of file G4TransitionRadiation.hh.

Referenced by EnergyIntegralDistribution(), and G4TransitionRadiation().

G4double G4TransitionRadiation::fMinEnergy [protected]

Definition at line 113 of file G4TransitionRadiation.hh.

Referenced by AngleIntegralDistribution(), and G4TransitionRadiation().

const G4int G4TransitionRadiation::fPointNumber = 100 [static, protected]

Definition at line 111 of file G4TransitionRadiation.hh.

G4double G4TransitionRadiation::fSigma1 [protected]

Reimplemented in G4ForwardXrayTR.

Definition at line 117 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().

G4double G4TransitionRadiation::fSigma2 [protected]

Reimplemented in G4ForwardXrayTR.

Definition at line 118 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().

const G4int G4TransitionRadiation::fSympsonNumber = 100 [static, protected]

Reimplemented in G4ForwardXrayTR.

Definition at line 109 of file G4TransitionRadiation.hh.

Referenced by AngleIntegralDistribution(), EnergyIntegralDistribution(), IntegralOverAngle(), and IntegralOverEnergy().

G4double G4TransitionRadiation::fVarAngle [protected]

Definition at line 106 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:53:33 2013 for Geant4 by  doxygen 1.4.7