Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4TransitionRadiation Class Referenceabstract

#include <G4TransitionRadiation.hh>

+ Inheritance diagram for G4TransitionRadiation:

Public Member Functions

 G4TransitionRadiation (const G4String &processName="TR", G4ProcessType type=fElectromagnetic)
 
virtual ~G4TransitionRadiation ()
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
 
virtual G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *condition) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
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
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
virtual void ProcessDescription (std::ostream &outfile) const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Attributes

G4int fMatIndex1
 
G4int fMatIndex2
 
G4double fGamma
 
G4double fEnergy
 
G4double fVarAngle
 
G4double fMinEnergy
 
G4double fMaxEnergy
 
G4double fMaxTheta
 
G4double fSigma1
 
G4double fSigma2
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Static Protected Attributes

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

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 55 of file G4TransitionRadiation.hh.

Constructor & Destructor Documentation

◆ G4TransitionRadiation()

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

◆ ~G4TransitionRadiation()

G4TransitionRadiation::~G4TransitionRadiation ( )
virtual

Definition at line 73 of file G4TransitionRadiation.cc.

74{}

Member Function Documentation

◆ AngleIntegralDistribution()

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

Definition at line 154 of file G4TransitionRadiation.cc.

157{
158 G4int i ;
159 G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
160 h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ;
161 for(i=1;i<fSympsonNumber;i++)
162 {
165 varAngle1 + 2*i*h)
168 varAngle1 + 2*i*h);
171 varAngle1 + (2*i - 1)*h)
174 varAngle1 + (2*i - 1)*h) ;
175 }
178 varAngle1 + (2*fSympsonNumber - 1)*h)
181 varAngle1 + (2*fSympsonNumber - 1)*h) ;
182
185 varAngle1)
188 varAngle1)
191 varAngle2)
194 varAngle2)
195 + 4.0*sumOdd + 2.0*sumEven )/3.0 ;
196}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
G4double IntegralOverEnergy(G4double energy1, G4double energy2, G4double varAngle) const
static const G4int fSympsonNumber

◆ EnergyIntegralDistribution()

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

Definition at line 204 of file G4TransitionRadiation.cc.

207{
208 G4int i ;
209 G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
210 h = 0.5*(energy2 - energy1)/fSympsonNumber ;
211 for(i=1;i<fSympsonNumber;i++)
212 {
213 sumEven += IntegralOverAngle(energy1 + 2*i*h,0.0,0.01*fMaxTheta )
214 + IntegralOverAngle(energy1 + 2*i*h,0.01*fMaxTheta,fMaxTheta);
215 sumOdd += IntegralOverAngle(energy1 + (2*i - 1)*h,0.0,0.01*fMaxTheta)
216 + IntegralOverAngle(energy1 + (2*i - 1)*h,0.01*fMaxTheta,fMaxTheta) ;
217 }
218 sumOdd += IntegralOverAngle(energy1 + (2*fSympsonNumber - 1)*h,
219 0.0,0.01*fMaxTheta)
220 + IntegralOverAngle(energy1 + (2*fSympsonNumber - 1)*h,
221 0.01*fMaxTheta,fMaxTheta) ;
222
223 return h*(IntegralOverAngle(energy1,0.0,0.01*fMaxTheta)
225 + IntegralOverAngle(energy2,0.0,0.01*fMaxTheta)
227 + 4.0*sumOdd + 2.0*sumEven )/3.0 ;
228}
G4double IntegralOverAngle(G4double energy, G4double varAngle1, G4double varAngle2) const

◆ GetMeanFreePath()

G4double G4TransitionRadiation::GetMeanFreePath ( const G4Track ,
G4double  ,
G4ForceCondition condition 
)
overridevirtual

Implements G4VDiscreteProcess.

Reimplemented in G4ForwardXrayTR.

Definition at line 82 of file G4TransitionRadiation.cc.

85{
87 return DBL_MAX; // so TR doesn't limit mean free path
88}
G4double condition(const G4ErrorSymMatrix &m)
@ Forced
#define DBL_MAX
Definition: templates.hh:62

◆ IntegralOverAngle()

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

Definition at line 129 of file G4TransitionRadiation.cc.

132{
133 G4int i ;
134 G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
135 h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ;
136 for(i=1;i<fSympsonNumber;i++)
137 {
138 sumEven += SpectralAngleTRdensity(energy,varAngle1 + 2*i*h) ;
139 sumOdd += SpectralAngleTRdensity(energy,varAngle1 + (2*i - 1)*h) ;
140 }
141 sumOdd += SpectralAngleTRdensity(energy,varAngle1 + (2*fSympsonNumber - 1)*h) ;
142
143 return h*( SpectralAngleTRdensity(energy,varAngle1)
144 + SpectralAngleTRdensity(energy,varAngle2)
145 + 4.0*sumOdd + 2.0*sumEven )/3.0 ;
146}
virtual G4double SpectralAngleTRdensity(G4double energy, G4double varAngle) const =0

Referenced by EnergyIntegralDistribution().

◆ IntegralOverEnergy()

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

Definition at line 103 of file G4TransitionRadiation.cc.

106{
107 G4int i ;
108 G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
109 h = 0.5*(energy2 - energy1)/fSympsonNumber ;
110 for(i=1;i<fSympsonNumber;i++)
111 {
112 sumEven += SpectralAngleTRdensity(energy1 + 2*i*h,varAngle) ;
113 sumOdd += SpectralAngleTRdensity(energy1 + (2*i - 1)*h,varAngle) ;
114 }
115 sumOdd += SpectralAngleTRdensity(energy1 + (2*fSympsonNumber - 1)*h,varAngle) ;
116 return h*( SpectralAngleTRdensity(energy1,varAngle)
117 + SpectralAngleTRdensity(energy2,varAngle)
118 + 4.0*sumOdd + 2.0*sumEven )/3.0 ;
119}

Referenced by AngleIntegralDistribution().

◆ IsApplicable()

G4bool G4TransitionRadiation::IsApplicable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 77 of file G4TransitionRadiation.cc.

78{
79 return ( aParticleType.GetPDGCharge() != 0.0 );
80}
G4double GetPDGCharge() const

◆ PostStepDoIt()

G4VParticleChange * G4TransitionRadiation::PostStepDoIt ( const G4Track ,
const G4Step  
)
overridevirtual

Reimplemented from G4VDiscreteProcess.

Reimplemented in G4ForwardXrayTR.

Definition at line 90 of file G4TransitionRadiation.cc.

92{
94 return &aParticleChange;
95}
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:327
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:424

◆ SpectralAngleTRdensity()

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

Implemented in G4ForwardXrayTR.

Referenced by IntegralOverAngle(), and IntegralOverEnergy().

Member Data Documentation

◆ fEnergy

G4double G4TransitionRadiation::fEnergy
protected

Definition at line 104 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().

◆ fGamma

G4double G4TransitionRadiation::fGamma
protected

Definition at line 103 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().

◆ fGammaNumber

const G4int G4TransitionRadiation::fGammaNumber = 15
staticprotected

Definition at line 109 of file G4TransitionRadiation.hh.

◆ fMatIndex1

◆ fMatIndex2

◆ fMaxEnergy

G4double G4TransitionRadiation::fMaxEnergy
protected

Definition at line 113 of file G4TransitionRadiation.hh.

Referenced by AngleIntegralDistribution(), and G4TransitionRadiation().

◆ fMaxTheta

G4double G4TransitionRadiation::fMaxTheta
protected

Definition at line 114 of file G4TransitionRadiation.hh.

Referenced by EnergyIntegralDistribution(), and G4TransitionRadiation().

◆ fMinEnergy

G4double G4TransitionRadiation::fMinEnergy
protected

Definition at line 112 of file G4TransitionRadiation.hh.

Referenced by AngleIntegralDistribution(), and G4TransitionRadiation().

◆ fPointNumber

const G4int G4TransitionRadiation::fPointNumber = 100
staticprotected

Definition at line 110 of file G4TransitionRadiation.hh.

◆ fSigma1

G4double G4TransitionRadiation::fSigma1
protected

Definition at line 116 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().

◆ fSigma2

G4double G4TransitionRadiation::fSigma2
protected

Definition at line 117 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().

◆ fSympsonNumber

const G4int G4TransitionRadiation::fSympsonNumber = 100
staticprotected

◆ fVarAngle

G4double G4TransitionRadiation::fVarAngle
protected

Definition at line 105 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().


The documentation for this class was generated from the following files: