#include <NonLTELineGasMix.hpp>
Classes | |
struct | ColPartner |
Public Types | |
enum class | Species : int { Test , Hydroxyl , HydroxylHFS , Formyl , CarbonMonoxide , AtomicCarbon , IonizedCarbon , MolecularHydrogen } |
![]() | |
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 ¶ms) 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< SnapshotParameter > | parameterInfo () 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< StateVariable > | specificStateVariableInfo () const override |
bool | storeMeanIntensities () const |
UpdateStatus | updateSpecificState (MaterialState *state, const Array &Jv) const override |
![]() | |
MaterialType | materialType () const override |
double | sourceWeight () const |
double | wavelengthBias () const |
WavelengthDistribution * | wavelengthBiasDistribution () const |
Range | wavelengthRange () const override |
![]() | |
virtual double | asymmpar (double lambda) const |
virtual Array | emissionSpectrum (const MaterialState *state, const Array &Jv) const |
virtual DisjointWavelengthGrid * | emissionWavelengthGrid () 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 ¶ms) 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< SnapshotParameter > | parameterInfo () 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 Array & | sectionsAbs (double lambda) const |
virtual const Array & | sectionsAbspol (double lambda) const |
virtual double | sectionSca (double lambda) const =0 |
virtual vector< StateVariable > | specificStateVariableInfo () const =0 |
virtual const Array & | thetaGrid () const |
virtual UpdateStatus | updateSpecificState (MaterialState *state, const Array &Jv) const |
![]() | |
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 |
![]() | |
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 Item * | getItemProperty (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) |
Item & | operator= (const Item &)=delete |
Item * | parent () 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 |
![]() | |
virtual | ~WavelengthRangeInterface () |
virtual Range | wavelengthRange () const =0 |
Protected Member Functions | |
NonLTELineGasMix () | |
void | setupSelfBefore () override |
![]() | |
EmittingGasMix () | |
![]() | |
MaterialMix () | |
Configuration * | config () const |
Random * | random () const |
void | setupSelfBefore () override |
![]() | |
SimulationItem () | |
virtual bool | offersInterface (const std::type_info &interfaceTypeInfo) const |
virtual void | setupSelfAfter () |
virtual void | setupSelfBefore () |
![]() | |
Item () | |
![]() | |
SourceWavelengthRangeInterface () | |
![]() | |
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 |
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
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:
Two-level
test
molecule (TT): a fictive test molecule (called TT for our purposes) that has just two rotational energy levels with a corresponding transition line at 1666.67 Hydroxyl
radical (OH): includes rotational energy levels up to Hydroxyl
radical (OHhfs): includes rotational energy levels up to Formyl
cation (HCO+): includes rotational energy levels up to Carbon
monoxide
(CO): includes rotational energy levels up to Atomic
carbon
(C): includes three hyperfine split levels at the electronic ground level of 3P. The corresponding transition lines are at wavelengths 230.3, 370.4 and 609.1 Ionized
carbon
(C+): includes four electronic levels (2p, 4p, 2D, and 2S) and the hyperfine split levels in the electronic levels of 2p and 4p. The corresponding fine-structure transition lines for levels 2p and 4p are at wavelengths 157.74, 198.9, 353.57, 454.67 Molecular
hydrogen
(H2): includes rotational energy levels up to 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
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.
Thermal motion and turbulence
The thermal velocity in a medium of particles with mass
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
We can artificially combine the effect of both thermal motion and turbulence into an effective temperature,
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
where
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
where
The Einstein coefficients
and
where
Emission
The integrated line luminosity
where
Absorption
The absorption opacity
where the sum runs over all supported transitions,
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
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.
|
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)" .
|
inlineprotected |
Default constructor for concrete Item subclass NonLTELineGasMix : "A gas mix supporting rotational transitions in specific molecules and atoms" .
|
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.
|
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.
|
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.
|
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.
|
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.
|
overridevirtual |
This function returns true, indicating that this material supports secondary line emission from gas.
Reimplemented from MaterialMix.
|
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.
|
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.
|
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.
|
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.
|
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.
|
overridevirtual |
This function returns a list with the line centers of the supported transitions.
Reimplemented from MaterialMix.
|
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.
|
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.
|
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.
|
overridevirtual |
This function returns the mass of a molecule or atom of the species under consideration.
Implements MaterialMix.
|
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.
|
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.
|
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.
|
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.
|
overridevirtual |
This function returns the absorption opacity
Implements MaterialMix.
|
overridevirtual |
This function returns the extinction opacity
Implements MaterialMix.
|
overridevirtual |
This function returns the scattering opacity
Implements MaterialMix.
|
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.
|
overridevirtual |
This function does nothing because the lines under consideration do not scatter.
Implements MaterialMix.
|
overridevirtual |
This function does nothing because the lines under consideration do not scatter.
Implements MaterialMix.
|
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.
|
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.
|
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.
|
overrideprotectedvirtual |
This function loads the required information on the configured species from resource files:
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.
|
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" .
|
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.
|
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.
|
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.