Geant4 11.2.2
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 ()
 
 G4TransitionRadiation (const G4TransitionRadiation &right)=delete
 
G4TransitionRadiationoperator= (const G4TransitionRadiation &right)=delete
 
G4bool IsApplicable (const G4ParticleDefinition &aParticleType) override
 
virtual G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *condition) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
virtual void ProcessDescription (std::ostream &) const override
 
virtual void DumpInfo () const 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 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 &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- 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
 
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 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 const G4VProcessGetCreatorProcess () const
 
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
 
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

G4double fGamma
 
G4double fEnergy
 
G4double fVarAngle
 
G4double fMinEnergy
 
G4double fMaxEnergy
 
G4double fMaxTheta
 
G4double fSigma1
 
G4double fSigma2
 
G4int fMatIndex1
 
G4int fMatIndex2
 
- 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 constexpr G4int fSympsonNumber = 100
 
static constexpr G4int fGammaNumber = 15
 
static constexpr G4int fPointNumber = 100
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VDiscreteProcess
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 52 of file G4TransitionRadiation.hh.

Constructor & Destructor Documentation

◆ G4TransitionRadiation() [1/2]

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

◆ ~G4TransitionRadiation()

G4TransitionRadiation::~G4TransitionRadiation ( )
virtualdefault

◆ G4TransitionRadiation() [2/2]

G4TransitionRadiation::G4TransitionRadiation ( const G4TransitionRadiation & right)
delete

Member Function Documentation

◆ AngleIntegralDistribution()

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

Definition at line 137 of file G4TransitionRadiation.cc.

139{
140 G4int i;
141 G4double h, sumEven = 0.0, sumOdd = 0.0;
142 h = 0.5 * (varAngle2 - varAngle1) / fSympsonNumber;
143 for(i = 1; i < fSympsonNumber; ++i)
144 {
146 fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
147 varAngle1 + 2 * i * h) +
149 fMaxEnergy, varAngle1 + 2 * i * h);
151 fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
152 varAngle1 + (2 * i - 1) * h) +
154 fMaxEnergy, varAngle1 + (2 * i - 1) * h);
155 }
156 sumOdd +=
158 varAngle1 + (2 * fSympsonNumber - 1) * h) +
160 varAngle1 + (2 * fSympsonNumber - 1) * h);
161
162 return h *
164 fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
165 varAngle1) +
167 fMaxEnergy, varAngle1) +
169 fMinEnergy + 0.3 * (fMaxEnergy - fMinEnergy),
170 varAngle2) +
172 fMaxEnergy, varAngle2) +
173 4.0 * sumOdd + 2.0 * sumEven) /
174 3.0;
175}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
static constexpr G4int fSympsonNumber
G4double IntegralOverEnergy(G4double energy1, G4double energy2, G4double varAngle) const

◆ DumpInfo()

virtual void G4TransitionRadiation::DumpInfo ( ) const
inlineoverridevirtual

Reimplemented from G4VProcess.

Reimplemented in G4ForwardXrayTR.

Definition at line 74 of file G4TransitionRadiation.hh.

G4GLOB_DLL std::ostream G4cout
virtual void ProcessDescription(std::ostream &) const override

◆ EnergyIntegralDistribution()

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

Definition at line 180 of file G4TransitionRadiation.cc.

182{
183 G4int i;
184 G4double h, sumEven = 0.0, sumOdd = 0.0;
185 h = 0.5 * (energy2 - energy1) / fSympsonNumber;
186 for(i = 1; i < fSympsonNumber; ++i)
187 {
188 sumEven +=
189 IntegralOverAngle(energy1 + 2 * i * h, 0.0, 0.01 * fMaxTheta) +
190 IntegralOverAngle(energy1 + 2 * i * h, 0.01 * fMaxTheta, fMaxTheta);
191 sumOdd +=
192 IntegralOverAngle(energy1 + (2 * i - 1) * h, 0.0, 0.01 * fMaxTheta) +
193 IntegralOverAngle(energy1 + (2 * i - 1) * h, 0.01 * fMaxTheta, fMaxTheta);
194 }
195 sumOdd += IntegralOverAngle(energy1 + (2 * fSympsonNumber - 1) * h, 0.0,
196 0.01 * fMaxTheta) +
197 IntegralOverAngle(energy1 + (2 * fSympsonNumber - 1) * h,
198 0.01 * fMaxTheta, fMaxTheta);
199
200 return h *
201 (IntegralOverAngle(energy1, 0.0, 0.01 * fMaxTheta) +
202 IntegralOverAngle(energy1, 0.01 * fMaxTheta, fMaxTheta) +
203 IntegralOverAngle(energy2, 0.0, 0.01 * fMaxTheta) +
204 IntegralOverAngle(energy2, 0.01 * fMaxTheta, fMaxTheta) +
205 4.0 * sumOdd + 2.0 * sumEven) /
206 3.0;
207}
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 71 of file G4TransitionRadiation.cc.

73{
75 return DBL_MAX; // so TR doesn't limit mean free path
76}
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 112 of file G4TransitionRadiation.cc.

115{
116 G4int i;
117 G4double h, sumEven = 0.0, sumOdd = 0.0;
118 h = 0.5 * (varAngle2 - varAngle1) / fSympsonNumber;
119 for(i = 1; i < fSympsonNumber; ++i)
120 {
121 sumEven += SpectralAngleTRdensity(energy, varAngle1 + 2 * i * h);
122 sumOdd += SpectralAngleTRdensity(energy, varAngle1 + (2 * i - 1) * h);
123 }
124 sumOdd +=
125 SpectralAngleTRdensity(energy, varAngle1 + (2 * fSympsonNumber - 1) * h);
126
127 return h *
128 (SpectralAngleTRdensity(energy, varAngle1) +
129 SpectralAngleTRdensity(energy, varAngle2) + 4.0 * sumOdd +
130 2.0 * sumEven) /
131 3.0;
132}
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 88 of file G4TransitionRadiation.cc.

91{
92 G4int i;
93 G4double h, sumEven = 0.0, sumOdd = 0.0;
94 h = 0.5 * (energy2 - energy1) / fSympsonNumber;
95 for(i = 1; i < fSympsonNumber; i++)
96 {
97 sumEven += SpectralAngleTRdensity(energy1 + 2 * i * h, varAngle);
98 sumOdd += SpectralAngleTRdensity(energy1 + (2 * i - 1) * h, varAngle);
99 }
100 sumOdd +=
101 SpectralAngleTRdensity(energy1 + (2 * fSympsonNumber - 1) * h, varAngle);
102 return h *
103 (SpectralAngleTRdensity(energy1, varAngle) +
104 SpectralAngleTRdensity(energy2, varAngle) + 4.0 * sumOdd +
105 2.0 * sumEven) /
106 3.0;
107}

Referenced by AngleIntegralDistribution().

◆ IsApplicable()

G4bool G4TransitionRadiation::IsApplicable ( const G4ParticleDefinition & aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 65 of file G4TransitionRadiation.cc.

67{
68 return (aParticleType.GetPDGCharge() != 0.0);
69}

◆ operator=()

G4TransitionRadiation & G4TransitionRadiation::operator= ( const G4TransitionRadiation & right)
delete

◆ PostStepDoIt()

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

Reimplemented from G4VDiscreteProcess.

Reimplemented in G4ForwardXrayTR.

Definition at line 78 of file G4TransitionRadiation.cc.

80{
82 return &aParticleChange;
83}
G4ParticleChange aParticleChange
void ClearNumberOfInteractionLengthLeft()

◆ ProcessDescription()

void G4TransitionRadiation::ProcessDescription ( std::ostream & out) const
overridevirtual

Reimplemented from G4VProcess.

Reimplemented in G4ForwardXrayTR.

Definition at line 60 of file G4TransitionRadiation.cc.

61{
62 out << "Base class for simulation of x-ray transition radiation.\n";
63}

Referenced by DumpInfo().

◆ 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 98 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().

◆ fGamma

G4double G4TransitionRadiation::fGamma
protected

Definition at line 97 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().

◆ fGammaNumber

G4int G4TransitionRadiation::fGammaNumber = 15
staticconstexprprotected

Definition at line 94 of file G4TransitionRadiation.hh.

◆ fMatIndex1

◆ fMatIndex2

◆ fMaxEnergy

G4double G4TransitionRadiation::fMaxEnergy
protected

Definition at line 102 of file G4TransitionRadiation.hh.

Referenced by AngleIntegralDistribution(), and G4TransitionRadiation().

◆ fMaxTheta

G4double G4TransitionRadiation::fMaxTheta
protected

Definition at line 103 of file G4TransitionRadiation.hh.

Referenced by EnergyIntegralDistribution(), and G4TransitionRadiation().

◆ fMinEnergy

G4double G4TransitionRadiation::fMinEnergy
protected

Definition at line 101 of file G4TransitionRadiation.hh.

Referenced by AngleIntegralDistribution(), and G4TransitionRadiation().

◆ fPointNumber

G4int G4TransitionRadiation::fPointNumber = 100
staticconstexprprotected

Definition at line 95 of file G4TransitionRadiation.hh.

◆ fSigma1

G4double G4TransitionRadiation::fSigma1
protected

Definition at line 105 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().

◆ fSigma2

G4double G4TransitionRadiation::fSigma2
protected

Definition at line 106 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().

◆ fSympsonNumber

G4int G4TransitionRadiation::fSympsonNumber = 100
staticconstexprprotected

◆ fVarAngle

G4double G4TransitionRadiation::fVarAngle
protected

Definition at line 99 of file G4TransitionRadiation.hh.

Referenced by G4TransitionRadiation().


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