Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4TransitionRadiation Class Referenceabstract

#include <G4TransitionRadiation.hh>

Inheritance diagram for G4TransitionRadiation:
Collaboration 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 &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
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 &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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 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
 
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
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

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)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

Definition at line 57 of file G4TransitionRadiation.hh.

Constructor & Destructor Documentation

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

Definition at line 74 of file G4TransitionRadiation.cc.

75 {}

Member Function Documentation

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

Definition at line 156 of file G4TransitionRadiation.cc.

158 {
159  G4int i ;
160  G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
161  h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ;
162  for(i=1;i<fSympsonNumber;i++)
163  {
164  sumEven += IntegralOverEnergy(fMinEnergy,
166  varAngle1 + 2*i*h)
168  fMaxEnergy,
169  varAngle1 + 2*i*h);
170  sumOdd += IntegralOverEnergy(fMinEnergy,
171  fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
172  varAngle1 + (2*i - 1)*h)
174  fMaxEnergy,
175  varAngle1 + (2*i - 1)*h) ;
176  }
177  sumOdd += IntegralOverEnergy(fMinEnergy,
178  fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
179  varAngle1 + (2*fSympsonNumber - 1)*h)
181  fMaxEnergy,
182  varAngle1 + (2*fSympsonNumber - 1)*h) ;
183 
184  return h*(IntegralOverEnergy(fMinEnergy,
185  fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
186  varAngle1)
188  fMaxEnergy,
189  varAngle1)
191  fMinEnergy + 0.3*(fMaxEnergy - fMinEnergy),
192  varAngle2)
194  fMaxEnergy,
195  varAngle2)
196  + 4.0*sumOdd + 2.0*sumEven )/3.0 ;
197 }
static const G4int fSympsonNumber
int G4int
Definition: G4Types.hh:78
double G4double
Definition: G4Types.hh:76
G4double IntegralOverEnergy(G4double energy1, G4double energy2, G4double varAngle) const

Here is the call graph for this function:

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

Definition at line 206 of file G4TransitionRadiation.cc.

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

Here is the call graph for this function:

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

Implements G4VDiscreteProcess.

Reimplemented in G4ForwardXrayTR.

Definition at line 83 of file G4TransitionRadiation.cc.

86 {
87  *condition = Forced;
88  return DBL_MAX; // so TR doesn't limit mean free path
89 }
G4double condition(const G4ErrorSymMatrix &m)
#define DBL_MAX
Definition: templates.hh:83
G4double G4TransitionRadiation::IntegralOverAngle ( G4double  energy,
G4double  varAngle1,
G4double  varAngle2 
) const

Definition at line 130 of file G4TransitionRadiation.cc.

133 {
134  G4int i ;
135  G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
136  h = 0.5*(varAngle2 - varAngle1)/fSympsonNumber ;
137  for(i=1;i<fSympsonNumber;i++)
138  {
139  sumEven += SpectralAngleTRdensity(energy,varAngle1 + 2*i*h) ;
140  sumOdd += SpectralAngleTRdensity(energy,varAngle1 + (2*i - 1)*h) ;
141  }
142  sumOdd += SpectralAngleTRdensity(energy,varAngle1 + (2*fSympsonNumber - 1)*h) ;
143 
144  return h*( SpectralAngleTRdensity(energy,varAngle1)
145  + SpectralAngleTRdensity(energy,varAngle2)
146  + 4.0*sumOdd + 2.0*sumEven )/3.0 ;
147 }
static const G4int fSympsonNumber
int G4int
Definition: G4Types.hh:78
virtual G4double SpectralAngleTRdensity(G4double energy, G4double varAngle) const =0
G4double energy(const ThreeVector &p, const G4double m)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 104 of file G4TransitionRadiation.cc.

107 {
108  G4int i ;
109  G4double h , sumEven = 0.0 , sumOdd = 0.0 ;
110  h = 0.5*(energy2 - energy1)/fSympsonNumber ;
111  for(i=1;i<fSympsonNumber;i++)
112  {
113  sumEven += SpectralAngleTRdensity(energy1 + 2*i*h,varAngle) ;
114  sumOdd += SpectralAngleTRdensity(energy1 + (2*i - 1)*h,varAngle) ;
115  }
116  sumOdd += SpectralAngleTRdensity(energy1 + (2*fSympsonNumber - 1)*h,varAngle) ;
117  return h*( SpectralAngleTRdensity(energy1,varAngle)
118  + SpectralAngleTRdensity(energy2,varAngle)
119  + 4.0*sumOdd + 2.0*sumEven )/3.0 ;
120 }
static const G4int fSympsonNumber
int G4int
Definition: G4Types.hh:78
virtual G4double SpectralAngleTRdensity(G4double energy, G4double varAngle) const =0
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

G4bool G4TransitionRadiation::IsApplicable ( const G4ParticleDefinition aParticleType)
overridevirtual

Reimplemented from G4VProcess.

Definition at line 78 of file G4TransitionRadiation.cc.

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

Here is the call graph for this function:

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

Reimplemented from G4VDiscreteProcess.

Reimplemented in G4ForwardXrayTR.

Definition at line 91 of file G4TransitionRadiation.cc.

93 {
95  return &aParticleChange;
96 }
void ClearNumberOfInteractionLengthLeft()
Definition: G4VProcess.hh:447
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289

Here is the call graph for this function:

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

Implemented in G4ForwardXrayTR.

Here is the caller graph for this function:

Member Data Documentation

G4double G4TransitionRadiation::fEnergy
protected

Definition at line 106 of file G4TransitionRadiation.hh.

G4double G4TransitionRadiation::fGamma
protected

Definition at line 105 of file G4TransitionRadiation.hh.

const G4int G4TransitionRadiation::fGammaNumber = 15
staticprotected

Definition at line 111 of file G4TransitionRadiation.hh.

G4int G4TransitionRadiation::fMatIndex1
protected

Definition at line 100 of file G4TransitionRadiation.hh.

G4int G4TransitionRadiation::fMatIndex2
protected

Definition at line 101 of file G4TransitionRadiation.hh.

G4double G4TransitionRadiation::fMaxEnergy
protected

Definition at line 115 of file G4TransitionRadiation.hh.

G4double G4TransitionRadiation::fMaxTheta
protected

Definition at line 116 of file G4TransitionRadiation.hh.

G4double G4TransitionRadiation::fMinEnergy
protected

Definition at line 114 of file G4TransitionRadiation.hh.

const G4int G4TransitionRadiation::fPointNumber = 100
staticprotected

Definition at line 112 of file G4TransitionRadiation.hh.

G4double G4TransitionRadiation::fSigma1
protected

Definition at line 118 of file G4TransitionRadiation.hh.

G4double G4TransitionRadiation::fSigma2
protected

Definition at line 119 of file G4TransitionRadiation.hh.

const G4int G4TransitionRadiation::fSympsonNumber = 100
staticprotected

Definition at line 110 of file G4TransitionRadiation.hh.

G4double G4TransitionRadiation::fVarAngle
protected

Definition at line 107 of file G4TransitionRadiation.hh.


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