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

#include <NonLTELineGasMix.hpp>

Inheritance diagram for NonLTELineGasMix:
Inheritance graph
[legend]

Classes

struct  ColPartner
 

Public Types

enum class  Species : int {
  Test , Hydroxyl , HydroxylHFS , Formyl ,
  CarbonMonoxide , AtomicCarbon , IonizedCarbon , MolecularHydrogen
}
 
- Public Types inherited from MaterialMix
enum class  DynamicStateType { None , Primary , Secondary , PrimaryIfMergedIterations }
 
enum class  MaterialType { Dust , Electrons , Gas }
 

Public Member Functions

const vector< double > & defaultCollisionPartnerRatios () const
 
double defaultTemperature () const
 
double defaultTurbulenceVelocity () const
 
DynamicStateType hasDynamicMediumState () const override
 
bool hasExtraSpecificState () const override
 
bool hasLineEmission () const override
 
bool hasNegativeExtinction () const override
 
double indicativeTemperature (const MaterialState *state, const Array &Jv) const override
 
void initializeSpecificState (MaterialState *state, double metallicity, double temperature, const Array &params) const override
 
string initialLevelPopsFilename () const
 
bool isSpecificStateConverged (int numCells, int numUpdated, int numNotConverged, MaterialState *currentAggregate, MaterialState *previousAggregate) const override
 
Array lineEmissionCenters () const override
 
Array lineEmissionMasses () const override
 
Array lineEmissionSpectrum (const MaterialState *state, const Array &Jv) const override
 
double lowestOpticalDepth () const
 
double mass () const override
 
double maxChangeInGlobalLevelPopulations () const
 
double maxChangeInLevelPopulations () const
 
double maxFractionNotConvergedCells () const
 
int numEnergyLevels () const
 
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
 
vector< SnapshotParameterparameterInfo () 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
 
double sectionAbs (double lambda) const override
 
double sectionExt (double lambda) const override
 
double sectionSca (double lambda) const override
 
Species species () const
 
vector< StateVariablespecificStateVariableInfo () const override
 
bool storeMeanIntensities () const
 
UpdateStatus updateSpecificState (MaterialState *state, const Array &Jv) const override
 
- Public Member Functions inherited from EmittingGasMix
MaterialType materialType () const override
 
double sourceWeight () const
 
double wavelengthBias () const
 
WavelengthDistributionwavelengthBiasDistribution () const
 
Range wavelengthRange () const override
 
- 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
 
- Public Member Functions inherited from WavelengthRangeInterface
virtual ~WavelengthRangeInterface ()
 
virtual Range wavelengthRange () const =0
 

Protected Member Functions

 NonLTELineGasMix ()
 
void setupSelfBefore () override
 
- Protected Member Functions inherited from EmittingGasMix
 EmittingGasMix ()
 
- 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 ()
 
- Protected Member Functions inherited from SourceWavelengthRangeInterface
 SourceWavelengthRangeInterface ()
 
- Protected Member Functions inherited from WavelengthRangeInterface
 WavelengthRangeInterface ()
 

Private Types

using BaseType = EmittingGasMix
 
using ItemType = NonLTELineGasMix
 

Private Attributes

Array _center
 
vector< ColPartner_colPartner
 
vector< double > _defaultCollisionPartnerRatios
 
double _defaultTemperature
 
double _defaultTurbulenceVelocity
 
Array _dlambdav
 
vector< double > _einsteinA
 
vector< double > _einsteinBlu
 
vector< double > _einsteinBul
 
vector< double > _energy
 
int _indexFirstColPartnerDensity
 
int _indexFirstLevelPopulation
 
int _indexFirstMeanIntensity
 
int _indexKineticTemperature
 
vector< int > _indexLowRad
 
vector< int > _indexUpRad
 
string _initialLevelPopsFilename
 
vector< Array_initLevelPops
 
Array _lambdav
 
double _lowestOpticalDepth
 
double _mass
 
double _maxChangeInGlobalLevelPopulations
 
double _maxChangeInLevelPopulations
 
double _maxFractionNotConvergedCells
 
int _numColPartners
 
int _numEnergyLevels
 
int _numLevels
 
int _numLines
 
int _numWavelengths
 
Species _species
 
bool _storeMeanIntensities
 
vector< double > _weight
 

Friends

class ItemRegistry
 

Detailed Description

The NonLTELineGasMix class describes the material properties related to selected transitions in selected molecules and atoms. For each supported species, the current implementation includes a number of rotational energy levels (quantum number \(J\)) at the base vibrational level (quantum number \(v=0\)) for molecules and electronic energy levels and hyperfine split levels for atoms, and supports the allowed transitions between these levels. Vibrational energy levels and the corresponding rovibrational transitions may be added later. The class properties allow configuring the species and the number of transitions to be considered.

For each supported transition, the emission luminosity and absorption opacity in a given cell are determined from the gas properties defined in the input model and the local radiation field calculated by the simulation. The class performs an iterative, non-LTE calculation including the effects of collisional transitions (excitation and de-excitation with one or more types of interaction partners) and photonic transitions (spontaneous emission, absorption, and induced emission). This allows establishing the energy level distribution (population) for a wide range of material densities and radiation fields.

Supported species and transitions

The current implementation supports the following molecular or atomic species:

Configuring the simulation

Simulations that include gas represented by the NonLTELineGasMix often also include dust, although this is not a requirement. In any case, the simulation must have one or more primary sources that trigger the molecular lines directly (e.g. the cosmic microwave background) or indirectly (e.g. by heating the dust and thus causing thermal dust emission), or both. The simulationMode should be set to "GasEmission" or "DustAndGasEmission" correspondingly.

Because the secondary emission and the radiation field are calculated self-consistently, iterateSecondaryEmission should always be enabled. If the wavelength range of the primary source(s) significantly overlaps the emission lines of the species under consideration (in other words, if the opacity at the line wavelengths significantly affects the primary radiation field), then includePrimaryEmission in IterationOptions must be enabled as well. If not, this property can be left disabled. In both cases, iteratePrimaryEmission can be left disabled unless iteration is required for other media (for example, to self-consistently calculate radiative dust destruction).

The radiation field wavelength grid should properly resolve the UV, optical, and infrared wavelength range (in case the simulation includes dust) and the wavelength ranges of the supported transition lines. Separate instruments can be configured for the relevant wavelength ranges, e.g. using a logarithmic grid for the continuum spectrum and linear grids for the line profiles.

During the initial primary emission segment, the simulation determines the radiation field resulting from the primary sources, dust attenuation (if present) and molecular line absorption based on default equilibrium level populations. The resulting radiation field allows a first estimation of the level populations for each spatial grid cell; this calculation happens in the updateSpecificState() function.

During secondary emission, the simulation takes into account emission from all configured media in addition to absorption by these same media, including molecular lines in both cases. The previously stored level populations allow calculating the line emission spectrum and the line absorption cross sections. This results in an updated radiation field, which will in turn influence the level populations (and for high optical depths, possibly the dust temperature), which in turn influences the secondary emission spectra. In order to obtain a self-consistent result, the simulation must therefore iterate over secondary emission.

The input model must provide values for the spatial distribution of several medium properties, including the number density of the species under consideration, the number density of any relevant collisional partner species, the kinetic gas temperature, and the turbulence velocity. These values remain constant during the simulation. Most often, this information will be read from an input file by associating the NonLTELineGasMix with a subclass of ImportedMedium. For that medium component, the ski file attribute importTemperature must be set to 'true', and importMetallicity and importVariableMixParams must be left at 'false'. The additional columns required by the material mix are automatically imported and are expected after all other columns. For example, if bulk velocities are also imported for this medium component (i.e. importVelocity is 'true'), the column order would be

\[ ..., n_\mathrm{mol}, T_\mathrm{kin}, v_\mathrm{x}, v_\mathrm{y}, v_\mathrm{z}, n_\mathrm{col1} [, n_\mathrm{col1}, ...], v_\mathrm{turb}\]

For basic testing purposes, the NonLTELineGasMix can also be associated with a geometric medium component. The geometry then defines the spatial density distribution of the species under consideration (i.e. \(n_\mathrm{mol}\)), and the NonLTELineGasMix configuration properties specify a fixed default value for the other properties that will be used across the spatial domain. In this case, the number densities of the collisional partners are defined by a constant multiplier relative to \(n_\mathrm{mol}\).

Thermal motion and turbulence

The thermal velocity in a medium of particles with mass \(m\) at temperature \(T_\mathrm{kin}\) is defined as

\[ v_\mathrm{th} = \sqrt{\frac{2kT_\mathrm{kin}}{m}}. \]

This value corresponds to the most probable particle speed, i.e. the point where the probability distribution of the velocity vector norm reaches its maximum value. One often considers an additional source of line broadening caused by subgrid processes other than those corresponding to the kinetic temperature. This motion is characterized by the turbulent velocity \(v_\mathrm{turb}\). Assuming a Gaussian line profile, the total velocity dispersion \(v_\mathrm{s}\) (the standard deviation of the Gaussian in velocity space) is then defined through

\[ \sqrt{2}\,v_\mathrm{s} = \sqrt{v_\mathrm{th}^2 + v_\mathrm{turb}^2}. \]

We can artificially combine the effect of both thermal motion and turbulence into an effective temperature,

\[ T_\mathrm{eff} = T_\mathrm{kin} + \frac{m v_\mathrm{turb}^2}{2k}. \]

Numerical convergence

As mentioned above, a simulation including this material mix performs iterations to self-consistently calculate the radiation field and the medium state (i.e. the energy level populations that drive the line emission). This class offers three user-configurable properties that specify the convergence criteria for this iterative process.

The first two properties configure a criterion based on statistics per spatial cell. The maxChangeInLevelPopulations property specifies the maximum relative change between consecutive iterations in the level populations for a given spatial cell for that cell to be considered converged. The maxFractionNotConvergedCells property then specifies the maximum fraction of cells (relative to the total number of cells in the simulation) that may be left not converged for the whole spatial domain to be considered converged.

The third property configures a global criterion. The maxChangeInGlobalLevelPopulations property specifies the maximum relative change between consecutive iterations in the global level populations, accumulated over the complete spatial domain, for the spatial domain to be considered converged. In other words, convergence has been reached if

\[ \epsilon = \underset{i}{\mathrm{max}} \left| \frac{N_i^{(k)} - N_i^{(k-1)}} {N_i^{(k-1)}} \right| \le \epsilon_\mathrm{max}\]

where \(N_i^{(k)}\) denotes the total population across the spatial domain for energy level \(i\) at iteration \(k\) and \(\epsilon_\mathrm{max}\) is the user-configured value of the maxChangeInGlobalLevelPopulations property.

Both criteria must be satisfied simultaneously. Howver, the user can effectively disable a criterion by specifying very liberal maximum values.

Level populations

We denote the populations for the \(N\) supported energy levels as \(n_i\), with indices \(i=0,N-1\) increasing from lowest to highest energy state, and with \(\sum_i n_i = n_\mathrm{mol}\). The level populations form the basis to calculate the emission luminosities and the absorption opacities at later stage. Their values are obtained by solving the set of statistical equilibrium equations given by

\[ \sum_{j>i} \Big[n_jA_{ji} + (n_jB_{ji} - n_iB_{ij})J_{\lambda,ji})\Big] - \sum_{j<i} \Big[n_iA_{ij} + (n_iB_{ij} - n_jB_{ji})J_{\lambda,ij})\Big] + \sum_{j\neq i} \Big[n_jC_{ji} - n_iC_{ij}\Big]=0, \quad i=0,N-1\]

where \(J_{\lambda,ul}\) is the mean radiation field intensity integrated over the line profile corresponding to the transition from upper energy level \(u\) to lower energy level \(l\). Furthermore, \(A_{ul}\), \(B_{ul}\), \(B_{lu}\) are the Einstein coefficients for spontaneous emission, induced emission, and absorption, and \(C_{ul}\), \(C_{lu}\) the coefficients for collisional excitation and de-excitation. These coefficients are taken from the literature.

The Einstein coefficients \(A_{ul}\), \(B_{ul}\), and \(B_{lu}\) are constant for each transition, and the \(B_{ul}\) and \(B_{lu}\) coefficients can be obtained from the \(A_{ul}\) coefficients through

\[ B_{ul} = \frac{\lambda_{ul}^5}{2 hc^2}A_{ul} \]

and

\[ B_{lu} = \frac{g_u}{g_l} B_{ul} = \frac{g_u}{g_l} \frac{\lambda_{ul}^5}{2 hc^2}A_{ul}, \]

where \(g_u\) and \(g_l\) represent the degeneracy of the upper and lower energy levels, respectively. The collisional coefficients \(C_{ul}\) and \(C_{lu}\) depend on the number density of the collisional partner and on the kinetic temperature \(T_\mathrm{kin}\) of the gas. These coefficients are related by

\[ \frac{C_{lu}}{C_{ul}} = \frac{g_u}{g_l} \exp(-E_{ul}/kT_\mathrm{kin}).\]

Emission

The integrated line luminosity \(L_{ul}\) corresponding to the transition from upper energy level \(u\) to lower energy level \(l\) for a given spatial cell is given by

\[ L_{ul} = \frac{hc}{\lambda_{ul}} A_{ul} n_u V_\mathrm{cell}, \]

where \(\lambda_{ul}\) is the transition wavelength and \(V_\mathrm{cell}\) is the volume of the cell. The SKIRT framework automatically adds a Gaussian line profile assuming a thermal velocity corresponding to the effective temperature in the cell and the mass of a molecule, in addition to any Doppler shifts caused by the bulk velocity in the cell.

Absorption

The absorption opacity \(k^\text{abs}(\lambda)\) as a function of wavelength \(\lambda\) is given by

\[ k^\text{abs}(\lambda) = \sum_{ul} \, \frac{h c}{4\pi \lambda_{ul}} (n_l\,B_{lu}-n_u\,B_{ul}) \,\phi_{ul}(\lambda) \]

where the sum runs over all supported transitions, \(u\) and \(l\) are the indices of the energy levels before and after the transition, and \(\phi_{ul}(\lambda)\) is the Gaussian profile of the line corresponding to the transition.

The calculation explicitly includes the Gaussian line profile caused by thermal motion of the molecules in the local bulk velocity frame of the cell. Moreover, the absorption profiles of all supported lines are superposed on top of each other. In practice, the implementation includes just the terms that have a significant contribution at any given wavelength.

Limiting negative optical depth

The formula for the absorption opacity given in the previous section yields a negative value in case the medium exhibits stimulated emission. The explicit absorption technique (see the MonteCarloSimulation class) supports the corresponding negative optical depths along a photon path. However, numerical problems arise if the negative optical depth magnitude becomes too high. This may happen, for example, as a result of Monte Carlo noise, or in the early stages when a simulation has not yet converged.

Therefore, this class imposes a user-configurable lower limit on the (negative) optical depth along a cell diagonal before returning an absorption opacity. The default limit is -2.

Storing mean intensities

As discussed above, the mean radiation field intensity \(J_{\lambda,ul}\) is determined for each transition line as a prerequisite for calculating the level populations in a given cell. These values are obtained by integrating the radiation field, as discretized on the simulation's radiation field wavelength grid, over the appropriate Gaussian line profile. There is no need to keep these values around for the regular operation of the code. However, they can form a relevant diagnostic to assess the level of Monte Carlo noise in the simulation or to evaluate iterative convergence. Therefore, if the storeMeanIntensities flag is turned on, the mean radiation field intensity at each transition line is stored in the medium state for each cell, so that the information can be probed using the CustomStateProbe class.

Providing initial level populations

By default, the level populations are initialized in each cell using a Boltzmann distribution as a function of temperature (i.e., assuming LTE conditions). This may be far from the correct non-LTE levels, possibly requiring many iterations to obtain a properly converged solution. It might therefore be meaningful to provide customized initial conditions. If the initialLevelPopsFilename string is nonempty, it specifies the name of a text column file with a column for each energy level and a row for each spatial cell in the simulation. Specifically, the first column lists a cell index that is ignored, and remaining columns list the relative population for each energy level in units of number density (the values are scaled to the total number density in the cell, so the specific units don't really matter). The rows must exactly match the number and ordering of the cells in the simulation's spatial grid (see below).

This input format is designed such that the text file can easily be produced by a CustomStateProbe instance in a previous simulation. For example, assuming 9 energy levels, one could configure a probe as follows:

<CustomStateProbe probeName="populations" indices="2-10" probeAfter="Run">
    <form type="Form">
        <PerCellForm/>
    </form>
</CustomStateProbe>

The requirement to have an identical spatial grid in both simulations is easily met for 1D, 2D and 3D grids with a regular mesh in each spatial direction. For octree and binary tree grids the precise cell hierarchy is usually determined based on random samples of the medium density distribution, resulting in a slightly different structure for each run. To work around this problem, one needs to output the topology of the hierarchical grid in the first simulation using a TreeSpatialGridTopologyProbe instance, and load this topology in subsequent simulations using a FileTreeSpatialGrid instance.

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

Member Enumeration Documentation

◆ Species

enum class NonLTELineGasMix::Species : int
strong

The enumeration type indicating the molecular or atomic species represented by a given NonLTELineGasMix instance. See the class header for more information.

Test : "Fictive two-level test molecule (TT)" .

Hydroxyl : "Hydroxyl radical (OH)" .

HydroxylHFS : "Hydroxyl radical (OH) with hyperfine structure" .

Formyl : "Formyl cation (HCO+)" .

CarbonMonoxide : "Carbon monoxide (CO)" .

AtomicCarbon : "Atomic carbon (C)" .

IonizedCarbon : "Ionized carbon (C+)" .

MolecularHydrogen : "Molecular hydrogen (H2)" .

Constructor & Destructor Documentation

◆ NonLTELineGasMix()

NonLTELineGasMix::NonLTELineGasMix ( )
inlineprotected

Default constructor for concrete Item subclass NonLTELineGasMix : "A gas mix supporting rotational transitions in specific molecules and atoms" .

Member Function Documentation

◆ defaultCollisionPartnerRatios()

NonLTELineGasMix::defaultCollisionPartnerRatios ( ) const
inline

This function returns the value of the discoverable double list property defaultCollisionPartnerRatios : "the default relative abundances of the collisional partners" .

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

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

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.

◆ defaultTemperature()

NonLTELineGasMix::defaultTemperature ( ) const
inline

This function returns the value of the discoverable double property defaultTemperature : "the default temperature of the gas" .

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 "1000" .

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.

◆ defaultTurbulenceVelocity()

NonLTELineGasMix::defaultTurbulenceVelocity ( ) const
inline

This function returns the value of the discoverable double property defaultTurbulenceVelocity : "the default (non-thermal) turbulence velocity" .

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

The minimum value for this property is "[0 km/s" .

The maximum value for this property is "100000 km/s]" .

The default value for this property is given by the conditional value expression "0 km/s" .

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.

◆ hasDynamicMediumState()

DynamicStateType NonLTELineGasMix::hasDynamicMediumState ( ) const
overridevirtual

This function returns DynamicStateType::PrimaryIfMergedIterations, indicating that this material mix has a dynamic medium state with updates that are considered to affect primary emission when the simulation has merged iterations, and only secondary emission if not.

Reimplemented from MaterialMix.

◆ hasExtraSpecificState()

bool NonLTELineGasMix::hasExtraSpecificState ( ) const
overridevirtual

This function returns true, indicating that the cross sections returned by this material mix depend on the values of specific state variables other than the number density.

Reimplemented from MaterialMix.

◆ hasLineEmission()

bool NonLTELineGasMix::hasLineEmission ( ) const
overridevirtual

This function returns true, indicating that this material supports secondary line emission from gas.

Reimplemented from MaterialMix.

◆ hasNegativeExtinction()

bool NonLTELineGasMix::hasNegativeExtinction ( ) const
overridevirtual

This function returns true, indicating that this material may have a negative absorption cross section. This happens when inverted level populations cause net stimulated emission.

Reimplemented from MaterialMix.

◆ indicativeTemperature()

double NonLTELineGasMix::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 temperature stored in the specific state for the relevant spatial cell and medium component. This value corresponds to the temperature defined by the input model at the start of the simulation.

Reimplemented from MaterialMix.

◆ initializeSpecificState()

void NonLTELineGasMix::initializeSpecificState ( MaterialState state,
double  metallicity,
double  temperature,
const Array params 
) const
overridevirtual

This function initializes the specific state variables requested by this material mix through the specificStateVariableInfo() function except for the number density of the species under consideration. For this class, the function uses the imported values, or if not available, the user-configured default values. The level populations are left at zero.

Reimplemented from MaterialMix.

◆ initialLevelPopsFilename()

NonLTELineGasMix::initialLevelPopsFilename ( ) const
inline

This function returns the value of the discoverable string property initialLevelPopsFilename : "the name of the file with initial level populations" .

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 "Level3" evaluates to true after replacing the names by true or false depending on their presence.

◆ isSpecificStateConverged()

bool NonLTELineGasMix::isSpecificStateConverged ( int  numCells,
int  numUpdated,
int  numNotConverged,
MaterialState currentAggregate,
MaterialState previousAggregate 
) const
overridevirtual

This function returns true if the state of the medium component corresponding to this material mix can be considered to be converged based on the given spatial cell statistics and aggregate material states, and false otherwise. For more information on the convergence criteria, see the corresponding section in the class header.

Reimplemented from MaterialMix.

◆ lineEmissionCenters()

Array NonLTELineGasMix::lineEmissionCenters ( ) const
overridevirtual

This function returns a list with the line centers of the supported transitions.

Reimplemented from MaterialMix.

◆ lineEmissionMasses()

Array NonLTELineGasMix::lineEmissionMasses ( ) const
overridevirtual

This function returns a list with the masses of the particles emitting each of the lines, i.e. the mass of a molecule or atom of the species under consideration in each case.

Reimplemented from MaterialMix.

◆ lineEmissionSpectrum()

Array NonLTELineGasMix::lineEmissionSpectrum ( const MaterialState state,
const Array Jv 
) const
overridevirtual

This function returns a list with the line luminosities for the supported transitions in the spatial cell and medium component represented by the specified material state and the receiving material mix when it would be embedded in the specified radiation field.

Reimplemented from MaterialMix.

◆ lowestOpticalDepth()

NonLTELineGasMix::lowestOpticalDepth ( ) const
inline

This function returns the value of the discoverable double property lowestOpticalDepth : "Lower limit of (negative) optical depth along a cell diagonal" .

The minimum value for this property is "[-10" .

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

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

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.

◆ mass()

double NonLTELineGasMix::mass ( ) const
overridevirtual

This function returns the mass of a molecule or atom of the species under consideration.

Implements MaterialMix.

◆ maxChangeInGlobalLevelPopulations()

NonLTELineGasMix::maxChangeInGlobalLevelPopulations ( ) const
inline

This function returns the value of the discoverable double property maxChangeInGlobalLevelPopulations : "the maximum relative change for the global level populations to be considered converged" .

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

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

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

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.

◆ maxChangeInLevelPopulations()

NonLTELineGasMix::maxChangeInLevelPopulations ( ) const
inline

This function returns the value of the discoverable double property maxChangeInLevelPopulations : "the maximum relative change for the level populations in a cell to be considered converged" .

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

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

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

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.

◆ maxFractionNotConvergedCells()

NonLTELineGasMix::maxFractionNotConvergedCells ( ) const
inline

This function returns the value of the discoverable double property maxFractionNotConvergedCells : "the maximum fraction of not-converged cells for all cells to be considered converged" .

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

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

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

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.

◆ numEnergyLevels()

NonLTELineGasMix::numEnergyLevels ( ) const
inline

This function returns the value of the discoverable integer property numEnergyLevels : "the number of energy levels used (or 999 for all supported)" .

The minimum value for this property is "2" .

The maximum value for this property is "999" .

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

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.

◆ opacityAbs()

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

This function returns the absorption opacity \(k^\text{abs}= n_\mathrm{mol} \varsigma^\text{abs}\) for the given wavelength and material state. The photon packet properties are not used.

Implements MaterialMix.

◆ opacityExt()

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

This function returns the extinction opacity \(k^\text{ext}=k^\text{abs}\) for the given wavelength and material state, which equals the absorption opacity because the scattering opacity is zero. The photon packet properties are not used.

Implements MaterialMix.

◆ opacitySca()

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

This function returns the scattering opacity \(k^\text{sca}\) which is trivially zero at all wavelengths.

Implements MaterialMix.

◆ parameterInfo()

vector< SnapshotParameter > NonLTELineGasMix::parameterInfo ( ) const
overridevirtual

This function returns the number and type of import parameters required by this particular material mix as a list of SnapshotParameter objects. For this class, the function returns a descriptor for the number densities of the collisional partners and for the turbulence velocity. Importing the kinetic gas temperature should be enabled through the corresponding standard configuration flag.

Reimplemented from MaterialMix.

◆ peeloffScattering()

void NonLTELineGasMix::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 does nothing because the lines under consideration do not scatter.

Implements MaterialMix.

◆ performScattering()

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

This function does nothing because the lines under consideration do not scatter.

Implements MaterialMix.

◆ sectionAbs()

double NonLTELineGasMix::sectionAbs ( double  lambda) const
overridevirtual

This function should return the absorption cross section using default properties. Because this value is hard to calculate for this material mix, this function returns zero.

Implements MaterialMix.

◆ sectionExt()

double NonLTELineGasMix::sectionExt ( double  lambda) const
overridevirtual

This function should return the extinction cross section using default properties. Because this value is hard to calculate for this material mix, this function returns zero.

Implements MaterialMix.

◆ sectionSca()

double NonLTELineGasMix::sectionSca ( double  lambda) const
overridevirtual

This function should return the scattering cross section using default properties. Because this value is hard to calculate for this material mix, this function returns zero.

Implements MaterialMix.

◆ setupSelfBefore()

void NonLTELineGasMix::setupSelfBefore ( )
overrideprotectedvirtual

This function loads the required information on the configured species from resource files:

  • the molecular weight of the species;
  • the supported energy levels with their corresponding degeneracies;
  • the Einstein A coefficients for the radiative transitions between the energy levels;
  • the gas temperature grid points used for discretizing the collisional coefficients K;
  • the excitation coefficients K for the collisional transitions between the energy levels, given for each of the temperature grid points.

The last two items are repeated for each collisional interaction partner of the species.

In a second step, the function calculates and caches some further values, including for example the Einstein B coefficients.

Reimplemented from MaterialMix.

◆ species()

NonLTELineGasMix::species ( ) const
inline

This function returns the value of the discoverable Species enumeration property species : "the molecular or atomic species being represented" .

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

◆ specificStateVariableInfo()

vector< StateVariable > NonLTELineGasMix::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 a list containing descriptors for the properties defined in the input model and for a number of variables to hold the level populations derived from the radiation field and related information.

Implements MaterialMix.

◆ storeMeanIntensities()

NonLTELineGasMix::storeMeanIntensities ( ) const
inline

This function returns the value of the discoverable Boolean property storeMeanIntensities : "store the mean radiation field intensity at each transition line" .

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

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.

◆ updateSpecificState()

UpdateStatus NonLTELineGasMix::updateSpecificState ( MaterialState state,
const Array Jv 
) const
overridevirtual

Based on the specified radiation field and the input model properties found in the given material state, this function determines the level populations for the supported transitions and stores these results back in the given material state. The function returns the update status as described for the UpdateStatus class.

Reimplemented from MaterialMix.


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