The SKIRT project
advanced radiative transfer for astrophysics
Public Types | Public Member Functions | Protected Member Functions | Private Types | Private Member Functions | Private Attributes | Friends | List of all members
XRayAtomicGasMix Class Reference

#include <XRayAtomicGasMix.hpp>

Inheritance diagram for XRayAtomicGasMix:
Inheritance graph
[legend]

Public Types

enum class  BoundElectrons : int {
  None , Free , FreeWithPolarization , Good ,
  Exact
}
 
- Public Types inherited from MaterialMix
enum class  DynamicStateType { None , Primary , Secondary , PrimaryIfMergedIterations }
 
enum class  MaterialType { Dust , Electrons , Gas }
 

Public Member Functions

const vector< double > & abundancies () const
 
bool hasPolarizedScattering () const override
 
bool hasScatteringDispersion () const override
 
double indicativeTemperature (const MaterialState *state, const Array &Jv) const override
 
double mass () const override
 
MaterialType materialType () const override
 
double opacityAbs (double lambda, const MaterialState *state, const PhotonPacket *pp) const override
 
double opacityExt (double lambda, const MaterialState *state, const PhotonPacket *pp) const override
 
double opacitySca (double lambda, const MaterialState *state, const PhotonPacket *pp) const override
 
void peeloffScattering (double &I, double &Q, double &U, double &V, double &lambda, Direction bfkobs, Direction bfky, const MaterialState *state, const PhotonPacket *pp) const override
 
void performScattering (double lambda, const MaterialState *state, PhotonPacket *pp) const override
 
BoundElectrons scatterBoundElectrons () const
 
double sectionAbs (double lambda) const override
 
double sectionExt (double lambda) const override
 
double sectionSca (double lambda) const override
 
vector< StateVariablespecificStateVariableInfo () const override
 
double temperature () const
 
- Public Member Functions inherited from MaterialMix
virtual double asymmpar (double lambda) const
 
virtual Array emissionSpectrum (const MaterialState *state, const Array &Jv) const
 
virtual DisjointWavelengthGridemissionWavelengthGrid () const
 
virtual Array emissivity (const Array &Jv) const
 
virtual bool hasContinuumEmission () const
 
virtual DynamicStateType hasDynamicMediumState () const
 
virtual bool hasExtraSpecificState () const
 
virtual bool hasLineEmission () const
 
virtual bool hasNegativeExtinction () const
 
virtual bool hasPolarizedAbsorption () const
 
virtual bool hasPolarizedEmission () const
 
virtual bool hasPolarizedScattering () const
 
virtual bool hasResonantScattering () const
 
virtual bool hasScatteringDispersion () const
 
virtual bool hasStochasticDustEmission () const
 
virtual double indicativeTemperature (const MaterialState *state, const Array &Jv) const
 
virtual void initializeSpecificState (MaterialState *state, double metallicity, double temperature, const Array &params) const
 
bool isDust () const
 
bool isElectrons () const
 
bool isGas () const
 
virtual bool isSpecificStateConverged (int numCells, int numUpdated, int numNotConverged, MaterialState *currentAggregate, MaterialState *previousAggregate) const
 
virtual Array lineEmissionCenters () const
 
virtual Array lineEmissionMasses () const
 
virtual Array lineEmissionSpectrum (const MaterialState *state, const Array &Jv) const
 
virtual double mass () const =0
 
virtual MaterialType materialType () const =0
 
virtual double opacityAbs (double lambda, const MaterialState *state, const PhotonPacket *pp) const =0
 
virtual double opacityExt (double lambda, const MaterialState *state, const PhotonPacket *pp) const =0
 
virtual double opacitySca (double lambda, const MaterialState *state, const PhotonPacket *pp) const =0
 
virtual vector< SnapshotParameterparameterInfo () const
 
virtual void peeloffScattering (double &I, double &Q, double &U, double &V, double &lambda, Direction bfkobs, Direction bfky, const MaterialState *state, const PhotonPacket *pp) const =0
 
virtual void performScattering (double lambda, const MaterialState *state, PhotonPacket *pp) const =0
 
virtual double sectionAbs (double lambda) const =0
 
virtual double sectionExt (double lambda) const =0
 
virtual const ArraysectionsAbs (double lambda) const
 
virtual const ArraysectionsAbspol (double lambda) const
 
virtual double sectionSca (double lambda) const =0
 
virtual vector< StateVariablespecificStateVariableInfo () const =0
 
virtual const ArraythetaGrid () const
 
virtual UpdateStatus updateSpecificState (MaterialState *state, const Array &Jv) const
 
- Public Member Functions inherited from SimulationItem
template<class T >
T * find (bool setup=true) const
 
template<class T >
T * interface (int levels=-999999, bool setup=true) const
 
virtual string itemName () const
 
void setup ()
 
string typeAndName () const
 
- Public Member Functions inherited from Item
 Item (const Item &)=delete
 
virtual ~Item ()
 
void addChild (Item *child)
 
const vector< Item * > & children () const
 
virtual void clearItemListProperty (const PropertyDef *property)
 
void destroyChild (Item *child)
 
virtual bool getBoolProperty (const PropertyDef *property) const
 
virtual vector< double > getDoubleListProperty (const PropertyDef *property) const
 
virtual double getDoubleProperty (const PropertyDef *property) const
 
virtual string getEnumProperty (const PropertyDef *property) const
 
virtual int getIntProperty (const PropertyDef *property) const
 
virtual vector< Item * > getItemListProperty (const PropertyDef *property) const
 
virtual ItemgetItemProperty (const PropertyDef *property) const
 
virtual string getStringProperty (const PropertyDef *property) const
 
int getUtilityProperty (string name) const
 
virtual void insertIntoItemListProperty (const PropertyDef *property, int index, Item *item)
 
Itemoperator= (const Item &)=delete
 
Itemparent () const
 
virtual void removeFromItemListProperty (const PropertyDef *property, int index)
 
virtual void setBoolProperty (const PropertyDef *property, bool value)
 
virtual void setDoubleListProperty (const PropertyDef *property, vector< double > value)
 
virtual void setDoubleProperty (const PropertyDef *property, double value)
 
virtual void setEnumProperty (const PropertyDef *property, string value)
 
virtual void setIntProperty (const PropertyDef *property, int value)
 
virtual void setItemProperty (const PropertyDef *property, Item *item)
 
virtual void setStringProperty (const PropertyDef *property, string value)
 
void setUtilityProperty (string name, int value)
 
virtual string type () const
 

Protected Member Functions

 XRayAtomicGasMix ()
 
 ~XRayAtomicGasMix ()
 
void setupSelfBefore () override
 
- Protected Member Functions inherited from MaterialMix
 MaterialMix ()
 
Configurationconfig () const
 
Randomrandom () const
 
void setupSelfBefore () override
 
- Protected Member Functions inherited from SimulationItem
 SimulationItem ()
 
virtual bool offersInterface (const std::type_info &interfaceTypeInfo) const
 
virtual void setupSelfAfter ()
 
virtual void setupSelfBefore ()
 
- Protected Member Functions inherited from Item
 Item ()
 

Private Types

using BaseType = MaterialMix
 
using ItemType = XRayAtomicGasMix
 

Private Member Functions

int indexForLambda (double lambda) const
 
void setScatteringInfoIfNeeded (PhotonPacket::ScatteringInfo *scatinfo, double lambda) const
 

Private Attributes

vector< double > _abundancies
 
ScatteringHelper * _com
 
ArrayTable< 2 > _cumprobscavv
 
vector< double > _lambdafluov
 
Array _lambdav
 
ScatteringHelper * _ray
 
BoundElectrons _scatterBoundElectrons
 
Array _sigmaextv
 
Array _sigmascav
 
double _temperature
 
vector< double > _vthermscav
 

Friends

class ItemRegistry
 

Detailed Description

The XRayAtomicGasMix class describes the material properties of neutral atomic gas in the X-ray wavelength range, taking into account the effects of photo-absorption, fluorescence, and scattering by bound electrons. To avoid the use of this material mix outside of the regime for which it has been designed, all cross sections are forced to zero below 4.3 eV and above 500 keV (corresponding approximately to a wavelength range from 2.5 pm to 290 nm).

The class assumes a gas containing a mixture of non-ionized elements with atomic numbers from 1 (hydrogen) up to 30 (zinc). The spatial density distribution of the gas is established by setting the hydrogen density. The relative abundances of the 30 elements and the temperature of the gas can be configured by the user as constant properties. In other words, the abundances and the temperature are considered to be spatially constant (for a given medium component), while the overall density can obviously vary across space as usual.

Photo-absorption by an atom is the process where the energy of a photon is used to liberate a bound electron from one of the electron shells of the atom. This class supports photo-absorption by any of the 30 elements in the gas for any of the electron shells, i.e. up to the K, L, M, or N shell depending on the atomic number of the element.

Fluorescence (in this context) is the process where an electron from a higher energy level "falls" into an empty space created by photo-absorption, emitting a new photon with a different energy. For each electron shell and for each possible fluorescence transition towards that shell, the yield defines the probability that such fluorescence event occurs after an electron has been liberated in that shell. This class supports all fluorescent lines that have line energies \(E > 0.1 \, \text{keV}\), for all elements in the gas. These are 20 lines: K \(_{\alpha2}\), K \(_{\alpha1}\), K \(_{\beta3}\), K \(_{\beta1}\), K \(^{II}_{\beta5}\), K \(^{I}_{\beta5}\), L \(_{\beta4}\), L \(_{\beta3}\), L \(_{\beta10}\), L \(_{\beta9}\), L \(_{\eta}\), L \(_{\beta17}\), L \(_{\beta1}\), L \(_{\gamma5}\), L \(_{\ell}\), L \(_{t}\), L \(_{s}\), L \(_{\alpha2}\), L \(_{\alpha1}\), L \(_{\beta6}\) (i.e. all transitions from higher shells towards the K and L shells, except for \(L \to L\) lines).

Because fluorescence only occurs as the result of a photo-absorption event, this class implements fluorescence as a form of scattering (where the wavelength of the photon being scattered changes) as opposed to emission. This allows both photo-absorption and fluorescence to be treated during primary emission. A possible drawback is that the weaker fluorescence lines will be represented by a fairly small number of photon packets.

Electrons bound to the atoms in the gas scatter incoming X-ray photons. This process can be elastic (Rayleigh scattering) or inelastic (bound-Compton scattering). The user can select one of five implementations of bound-electron scattering. In increasing order of accuracy and computational effort, these are:

Configuring the simulation

In addition to a medium component configured with the material mix represented by this class, simulations will usually include primary sources and possibly a dust medium. As indicated above, there is no need to include secondary emission, so the simulation mode can be set to "ExtinctionOnly" and there is no need to store the radiation field. The resulting continuum spectrum and absorption and emission features can be recorded by a single instrument configured with a high-resolution wavelength grid, or separate instruments can be configured with wavelength grids to resolve specific features of interest.

To study fluorescence lines seperated from the background continuum, configure the instrument to record flux components seperately and consider the scattered flux component, which includes the fluorescence lines in addition to any flux scattered from other media components. To ensure that low-intensity lines are properly included in this flux, set the advanced property minWeightReduction in the PhotonPacketOptions section to a value of 1e10 or so. With the default value of 1e4, low-intensity fluorescence photon packets are killed before having a chance to register in the instruments.

The input model must define the spatial distribution of the hydrogen number density \(n_\mathrm{H} = n_\mathrm{HI} + n_\mathrm{HII} + 2\,n_\mathrm{H2}\), i.e. including atomic, ionized, and molecular hydrogen. Note the factor 2 in this equation; \(n_\mathrm{H}\) in fact specifies the number of hydrogen protons per volume rather than the number of hydrogen-like particles per volume.

If this material mix is associated with a geometric medium component, the geometry defines the spatial density distribution \(n_\mathrm{H}\). Alternatively, if the material mix is associated with a subclass of ImportedMedium, the spatial density distribution is read from an input file. In that case, the ski file attributes importMetallicity, importTemperature, and importVariableMixParams must be left at 'false'. For example, if bulk velocities are also imported for this medium component (i.e. importVelocity is 'true'), the column order would be

\[ ..., n_\mathrm{H}, v_\mathrm{x}, v_\mathrm{y}, v_\mathrm{z} \]

The relative abundances of the 30 elements in the gas and the temperature of the gas are configured in the ski file as constant properties. In other words, the abundances and the temperature are considered to be spatially constant (for a given medium component). If the list of abundances is left empty in the ski file, default abundancies are used, taken from Table 2 of Anders & Grevesse (1989), the default abundance table in Xspec. In the default list, the abundance of hydrogen is set to unity. However, it is acceptable to specify a hydrogen abundance lower than one, for example to model an ionized hydrogen fraction.

Photo-absorption cross section

The total photo-absorption cross section per hydrogen atom for this material mix is obtained by accumulating the photo-absorption cross sections for all shells and for all individual elements, weighted by element abundancy, and convolved with a Gaussian profile reflecting each element's thermal velocity. Because the abundancies and the temperature are fixed, this calculation can be performed during setup and the result stored, discretized on a high-resolution wavelength grid for later retrieval.

Verner and Yakovlev (1995, 1996) provide analytic fits to the photo-absorption cross sections \(\sigma_{ph}(E)\) as a function of photon energy \(E\) for the ground-state shells of the first 30 atomic elements:

\[\begin{aligned} \sigma_{ph}(E) &= \begin{cases} 0 & E < E_\mathrm{th} \\ \sigma_0 \, F(y) & E_\mathrm{th} \le E < E_\mathrm{max}\\ 0 & E_\mathrm{max} \le E \end{cases}, \\ x &= \frac{E}{E_0} - y_0, \; \; \; \; y = \sqrt{x^2 + y_1^2},\\ F(y) &= \left[(x-1)^2+y_{\rm w}^2 \right]y^{-Q} \left(1+ \sqrt{(y/y_{\rm a})} \right )^{-P}, \\ Q&=5.5+l-0.5P, \end{aligned} \]

with \(E_\mathrm{th}\) the tabulated ionization threshold energy, \(E_\mathrm{max}\) the tabulated maximum energy for the formula to be valid, \(\sigma_0\), \(E_0\), \(y_{\rm w}\), \(y_{\rm a}\), \(P\), \(y_0\) and \(y_1\) seven tabulated fitting parameters, and \(l\) the subshell orbital quantum number ( \(l=0, 1, 2, 3\) for s, p, d, f orbitals respectively).

Fluorescence cross section

The total fluorescence cross section per hydrogen atom for this material mix is obtained similarly, but now including only the K and L shell photo-absorption cross sections for each element, multiplied by the appropriate fluorescence yields in addition to the element abundancy.

Electron scattering

As described above, this class provides several implementations for the scattering of X-rays photons by the electrons bound to the atoms. These implementations involve four types of scattering: free-electron Compton scattering, bound-electron Compton scattering, smooth Rayleigh scattering, and anomalous Rayleigh scattering. Note that, in all cases, the listed cross sections must be multiplied by the abundance of the corresponding element relative to hydrogen.

For free-electron Compton scattering, the cross section for an element with atomic number \(Z\) is given by \(Z\,\sigma_\mathrm{C}\), where \(\sigma_\mathrm{C}\) is the (wavelength-dependent) Compton cross section for a single free electron. The implementation of the scattering events is delegated to the ComptonPhaseFunction class; see there for more information on the cross section and phase function for free-electron Compton scattering.

The other scattering implementations depend on quantities that are available for each element \(Z\) in tabulated form as a function of the incoming photon energy \(E\) or as a function of the dimensionless momentum transfer parameter \(q\),

\[q= \frac{E}{12.4 \, \mathrm{keV}} \cdot \sin(\theta/2), \]

with \(\theta\) the scattering angle.

For bound-electron Compton scattering, the cross sections \(\sigma_{CS, Z}(E)\) are available as a table. The normalised scattering phase function for element Z is given by

\[ \Phi_{CS, Z}(\theta, E)= \frac{3}{4}\, \frac{\sigma_T}{\sigma_{CS, Z}(E)}\Big[C^3(\theta, E) + C(\theta, E) -C^2(\theta, E)\sin^2\theta\Big] \cdot S_Z(q), \]

with tabulated incoherent scattering functions \(S_Z(q)\) and the Compton factor \(C(\theta, E)\) defined as

\[ C(\theta, E) = {\Big[{1+\frac{E}{m_ec^2}(1-\cos \theta)\Big]}}^{-1}. \]

Also, inelastic bound-Compton scattering will change the photon energy by the Compton factor, just as for free-electron scattering.

For smooth Rayleigh scattering, the cross sections \(\sigma_{RSS, Z}(E)\) are available as a table. The normalised scattering phase function for element Z is given by

\[ \Phi_{RSS, Z}(\theta, E)= \frac{3}{4}\, \frac{\sigma_T}{\sigma_{RSS, Z}(E)}\Big[ 1 + \cos^2\theta \Big] \cdot F_Z^2(q), \]

with tabulated atomic form factors \(F_Z(q)\), which converge to \(Z\) at small \(q\) and decrease to zero for large \(q\).

Similarly, for anomalous Rayleigh scattering, the cross sections \(\sigma_{RSA, Z}(E)\) are available as a table. The normalised scattering phase function for element Z is now given by

\[ \Phi_{RSA, Z}(\theta, E)= \frac{3}{4}\, \frac{\sigma_T}{\sigma_{RSA, Z}(E)}\Big[ 1 + \cos^2\theta \Big] \cdot \Big[\big(F_Z(q) + F'_Z(E)\big)^2 + {F''_Z}^2(E)\Big], \]

with the same atomic form factors \(F_Z(q)\) as before, and tabulated real and imaginary anomalous scattering functions \(F'_Z(E)\) and \(F''_Z(E)\).

Performing scattering

The function performing an actual scattering event randomly selects one of the supported scattering channels (i.e. scattering by an electron bound to one of the supported elements or one of the 20 fluorescent line transition for one of the supported elements). The relative probabilities for these transitions as a function of incoming photon packet wavelength are also calculated during setup. The selected transition determines the scattering mechanism. For bound electrons, Rayleigh or Compton scattering is used. For fluorescence, the emission direction is isotropic, and the outgoing wavelength is the fluorescence wavelength. In both cases, a random Gaussian dispersion reflecting the interacting element's thermal velocity is applied to the outgoing wavelength.

Thermal dispersion

The thermal dispersion appropriate for a given interaction depends on the (constant) temperature configured for this material mix and on the mass of the interacting atom. For electron scattering and fluorescence events, the implementation is straightforward. Once a channel has been randomly selected, the magnitude of the interacting atom's thermal velocity is taken from a precalculated table and a velocity vector is sampled from the Maxwell distribution.

For absorption, the situation is much more involved. In principle, the full cross section curve for each ionization transition must be convolved with a Gaussian kernel of appropriate width. In practice, the cross section curves are very smooth except for the step at the threshold energy. The effect of the convolution is therefore limited to energies near the threshold energy, replacing the infinitely steep step by a sigmoid function. Considering that the convolution of a step function with a Gaussion is given by the error function, this class uses the following approximation. With \(E_\mathrm{th}\) the threshold energy and \(E_\mathrm{s}\) the energy dispersion corresponding to the thermal velocity of the atom being ionized, the cross section near \(E_\mathrm{th}\) is replaced by

\[ \sigma_{ph}'(E) = \frac{\sigma_{ph}(E_\mathrm{th} +2 E_\mathrm{s})} {2} \, \left[ 1+ \mathrm{erf} ( \frac{E - E_\mathrm{th}} {E_\mathrm{s}} ) \right] \qquad \mathrm{for} \; E_\mathrm{th} -2 E_\mathrm{s} < E < E_\mathrm{th} +2 E_\mathrm{s}. \]

The sigmoid is scaled by the value of the actual cross section at the end of the interval, achieving good continuity at that point in most cases, including all important transitions.

When an item of this type is used, the names provided by the conditional value expression "GasMix" are inserted into the name sets used for evaluating Boolean expressions.

Member Enumeration Documentation

◆ BoundElectrons

enum class XRayAtomicGasMix::BoundElectrons : int
strong

The enumeration type indicating the implementation used for scattering by bound electrons.

None : "ignore bound electrons" .

Free : "use free-electron Compton scattering" .

FreeWithPolarization : "use free-electron Compton scattering with support for polarization" .

Good : "use smooth Rayleigh scattering and exact bound-Compton scattering" .

Exact : "use anomalous Rayleigh scattering and exact bound-Compton scattering" .

Constructor & Destructor Documentation

◆ XRayAtomicGasMix()

XRayAtomicGasMix::XRayAtomicGasMix ( )
inlineprotected

Default constructor for concrete Item subclass XRayAtomicGasMix : "A gas mix supporting photo-absorption and fluorescence for X-ray wavelengths" .

◆ ~XRayAtomicGasMix()

XRayAtomicGasMix::~XRayAtomicGasMix ( )
protected

The destructor destructs the phase function helpers that were created during setup.

Member Function Documentation

◆ abundancies()

XRayAtomicGasMix::abundancies ( ) const
inline

This function returns the value of the discoverable double list property abundancies : "the abundancies for the elements with atomic number Z = 1,...,30" .

The minimum value for this property is "[0" .

The maximum value for this property is "1]" .

This property is required only if the Boolean expression "false" evaluates to true after replacing the names by true or false depending on their presence.

This property is displayed only if the Boolean expression "Level2" evaluates to true after replacing the names by true or false depending on their presence.

◆ hasPolarizedScattering()

bool XRayAtomicGasMix::hasPolarizedScattering ( ) const
overridevirtual

This function returns true if this material mix supports polarization during scattering events. For this class, the function returns true if the scatterBoundElectrons property has been set to FreeWithPolarization and false otherwise.

Reimplemented from MaterialMix.

◆ hasScatteringDispersion()

bool XRayAtomicGasMix::hasScatteringDispersion ( ) const
overridevirtual

This function returns true, indicating that a scattering interaction for this material mix may (and usually does) adjust the wavelength of the interacting photon packet.

Reimplemented from MaterialMix.

◆ indexForLambda()

int XRayAtomicGasMix::indexForLambda ( double  lambda) const
private

This function returns the index in the private wavelength grid corresponding to the specified wavelength. The parameters for converting a wavelength to the appropriate index are stored in data members during setup.

◆ indicativeTemperature()

double XRayAtomicGasMix::indicativeTemperature ( const MaterialState state,
const Array Jv 
) const
overridevirtual

This function returns an indicative temperature of the material mix when it would be embedded in a given radiation field. The implementation in this class ignores the radiation field and returns the (spatially constant) temperature configured for this material mix.

Reimplemented from MaterialMix.

◆ mass()

double XRayAtomicGasMix::mass ( ) const
overridevirtual

This function returns the mass of a hydrogen atom.

Implements MaterialMix.

◆ materialType()

MaterialType XRayAtomicGasMix::materialType ( ) const
overridevirtual

This function returns the fundamental material type represented by this material mix, which is MaterialType::Gas.

Implements MaterialMix.

◆ opacityAbs()

double XRayAtomicGasMix::opacityAbs ( double  lambda,
const MaterialState state,
const PhotonPacket pp 
) const
overridevirtual

This function returns the absorption (i.e. extinction minus scattering) opacity at the given wavelength and material state, using the abundances and temperature configured for this material mix. The photon packet properties are not used.

Implements MaterialMix.

◆ opacityExt()

double XRayAtomicGasMix::opacityExt ( double  lambda,
const MaterialState state,
const PhotonPacket pp 
) const
overridevirtual

This function returns the extinction opacity at the given wavelength and material state, using the abundances and temperature configured for this material mix. The photon packet properties are not used.

Implements MaterialMix.

◆ opacitySca()

double XRayAtomicGasMix::opacitySca ( double  lambda,
const MaterialState state,
const PhotonPacket pp 
) const
overridevirtual

This function returns the scattering (i.e. bound electrons plus fluorescence) opacity at the given wavelength and material state, using the abundances and temperature configured for this material mix. The photon packet properties are not used.

Implements MaterialMix.

◆ peeloffScattering()

void XRayAtomicGasMix::peeloffScattering ( double &  I,
double &  Q,
double &  U,
double &  V,
double &  lambda,
Direction  bfkobs,
Direction  bfky,
const MaterialState state,
const PhotonPacket pp 
) const
overridevirtual

This function calculates the contribution of the medium component associated with this material mix to the peel-off photon luminosity, polarization state, and wavelength shift for the given wavelength, geometry, material state, and photon properties. The contributions to the Stokes vector components are stored in the I, Q, U, V arguments, which are guaranteed to be initialized to zero by the caller. If there is wavelength shift, the new wavelength value replaces the incoming value of the lambda argument.

The function first calls the private setScatteringInfoIfNeeded() function to establish a random scattering channel and atom velocity for this event. In case the selected channel is bound electron scattering, the peel-off bias weight and wavelength shift are determined by the Compton phase function and the selected atom velocity. For fluorescence, the peel-off bias weight is trivially one because emission is isotropic and unpolarized, and the outgoing wavelength is determined by Doppler-shifting the rest wavelength of the selected fluorescence transition for the selected atom velocity.

Implements MaterialMix.

◆ performScattering()

void XRayAtomicGasMix::performScattering ( double  lambda,
const MaterialState state,
PhotonPacket pp 
) const
overridevirtual

This function performs a scattering event on the specified photon packet in the spatial cell and medium component represented by the specified material state and the receiving material mix. It first calls the private setScatteringInfoIfNeeded() function to establish a random scattering channel and atom velocity for this event.

In case the selected channel is bound electron scattering, the outgoing direction and adjusted wavelength are determined by the Compton phase function and the selected atom velocity. For fluorescence, emission is isotropic, so the outgoing direction is randomly chosen from the isotropic distribution. The outgoing wavelength is determined by Doppler-shifting the rest wavelength of the selected fluorescence transition for the selected atom velocity.

Implements MaterialMix.

◆ scatterBoundElectrons()

XRayAtomicGasMix::scatterBoundElectrons ( ) const
inline

This function returns the value of the discoverable BoundElectrons enumeration property scatterBoundElectrons : "implementation of scattering by bound electrons" .

The default value for this property is given by the conditional value expression "Good" .

This property is displayed only if the Boolean expression "Level3" evaluates to true after replacing the names by true or false depending on their presence.

◆ sectionAbs()

double XRayAtomicGasMix::sectionAbs ( double  lambda) const
overridevirtual

This function returns the absorption (i.e. extinction minus scattering) cross section per hydrogen atom at the given wavelength and using the abundances and temperature configured for this material mix.

Implements MaterialMix.

◆ sectionExt()

double XRayAtomicGasMix::sectionExt ( double  lambda) const
overridevirtual

This function returns the extinction cross section per hydrogen atom at the given wavelength and using the abundances and temperature configured for this material mix.

Implements MaterialMix.

◆ sectionSca()

double XRayAtomicGasMix::sectionSca ( double  lambda) const
overridevirtual

This function returns the scattering (i.e. bound electrons plus fluorescence) cross section per hydrogen atom at the given wavelength and using the abundances and temperature configured for this material mix.

Implements MaterialMix.

◆ setScatteringInfoIfNeeded()

void XRayAtomicGasMix::setScatteringInfoIfNeeded ( PhotonPacket::ScatteringInfo scatinfo,
double  lambda 
) const
private

This private function draws a random scattering channel and atom velocity and stores this information in the photon packet's scattering information record, unless a previous peel-off stored this already.

◆ setupSelfBefore()

void XRayAtomicGasMix::setupSelfBefore ( )
overrideprotectedvirtual

This function precalculates relevant cross sections and relative contributions over a high-resolution wavelength grid.

Reimplemented from MaterialMix.

◆ specificStateVariableInfo()

vector< StateVariable > XRayAtomicGasMix::specificStateVariableInfo ( ) const
overridevirtual

This function returns a list of StateVariable objects describing the specific state variables used by the receiving material mix. For this class, the function returns just the descriptor for the number density.

Implements MaterialMix.

◆ temperature()

XRayAtomicGasMix::temperature ( ) const
inline

This function returns the value of the discoverable double property temperature : "the temperature of the gas or zero to disable thermal dispersion" .

This property represents a physical quantity of type "temperature" .

The minimum value for this property is "[0" .

The maximum value for this property is "1e9]" .

The default value for this property is given by the conditional value expression "1e4" .

This property is displayed only if the Boolean expression "Level2" evaluates to true after replacing the names by true or false depending on their presence.


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