modelbound_CRM_data_procedures Submodule

Contains module procedures associated with the modelbound CRM data


Uses


Module Functions

pure module function getFixedW(this, ind) result(wMat)

Return fixed mapping matrix with given index if inelastic data component is allocated

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
integer(kind=ik), intent(in) :: ind

Return Value type(SparseRowData)

pure module function getFixedEmissionVector(this, ind) result(emit)

Return fixed emission vector for mapping matrix with given index if inelastic data component is allocated

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
integer(kind=ik), intent(in) :: ind

Return Value real(kind=rk), allocatable, dimension(:)

pure module function getInterpolatedEmissionVector(this, E) result(emit)

Return interpolated emission vector for given input energy E if inelastic data component is allocated

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
real(kind=rk), intent(in) :: E

Return Value real(kind=rk), allocatable, dimension(:)

pure module function getTransitionIngoingStates(this, ind) result(inStates)

Get ingoing states of transition with given index

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
integer(kind=ik), intent(in) :: ind

Return Value integer(kind=ik), allocatable, dimension(:)

pure module function getTransitionOutgoingStates(this, ind) result(outStates)

Get outgoing states of transition with given index

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
integer(kind=ik), intent(in) :: ind

Return Value integer(kind=ik), allocatable, dimension(:)

pure module function getTransitionRate(this, ind) result(rate)

Get transition rate from transition with given index

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
integer(kind=ik), intent(in) :: ind

Return Value real(kind=rk), allocatable, dimension(:)

pure module function getTransitionRateMomentum(this, ind) result(rate)

Get momentum transfer rate from transition with given index

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
integer(kind=ik), intent(in) :: ind

Return Value real(kind=rk), allocatable, dimension(:)

pure module function getTransitionRateEnergy(this, ind) result(rate)

Get energy transfer rate from transition with given index

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
integer(kind=ik), intent(in) :: ind

Return Value real(kind=rk), allocatable, dimension(:)

pure module function getTransitionCrossSection(this, col, ind) result(crossSection)

Get cross-section values from column col of the cross-section data of transition with given index

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
integer(kind=ik), intent(in) :: col
integer(kind=ik), intent(in) :: ind

Return Value real(kind=rk), allocatable, dimension(:)

pure module function getTransitionEnergy(this, ind) result(energy)

Get transition energy of transition with given index

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
integer(kind=ik), intent(in) :: ind

Return Value real(kind=rk), allocatable, dimension(:)

pure module function getNumTransitions(this) result(numTrans)

Get number of transitions registered in this object

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this

Return Value integer(kind=ik)

pure module function ratesIncludeElDensity(this, ind) result(includesElDens)

Check whether given transition includes electron density in the rate by default

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
integer(kind=ik), intent(in) :: ind

Return Value logical

pure module function getPopulationChangeMatrix(this, ids, transitionIndices) result(popChangeMat)

For a given set of species indices and transition indices returns a matrix whose entries are the change in particle number of a given species (by index) in transition processes determined by passed indices. Will not provide warnings if any particular ID is not found in any transition.

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
integer(kind=ik), intent(in), dimension(:) :: ids
integer(kind=ik), intent(in), dimension(:) :: transitionIndices

Return Value integer(kind=ik), dimension(:,:), allocatable

pure module function getRequiredDensityData(this, transitionIndex, removeLastState) result(densDataMat)

For a given transition index returns ingoingState data as a (:,2) matrix, where the first column is the list of participating states and the second column the number of times that state appears. If an ID is zero (electrons) and the transition includes electron density in the rate, the second column value is reduced by one. If removeLastState is true, the corresponding number in the second column of the result is reduced by one (useful in implicit terms). If any value in the second column drops to 0, the corresponding row will be removed. The result can then be used to determine density variables and their powers required to convert the transitions rate values into a source.

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
integer(kind=ik), intent(in) :: transitionIndex
logical(kind=ik), intent(in), optional :: removeLastState

Return Value integer(kind=ik), dimension(:,:), allocatable

module function getDataDimCRM(this, name) result(dim)

Get data dimensionality - currently always returns 1, assuming that the name is associated with a rate

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
character(len=*), intent(in) :: name

Name of data

Return Value integer(kind=ik)

pure module function getElState(this) result(ID)

Getter for electronStateID

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this

Return Value integer(kind=ik)


Module Subroutines

pure module subroutine initCRMData(this, numTransitions)

Collisional-radiative modelbound data initialization

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(inout) :: this
integer(kind=ik), intent(in) :: numTransitions

Expected number of transitions

pure module subroutine addTransition(this, tr)

Add transition to CRM modelbound data

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(inout) :: this
class(Transition), intent(in) :: tr

pure module subroutine setInelData(this, inelData)

Setter for inelData

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(inout) :: this
type(InelasticGridData), intent(in) :: inelData

pure module subroutine interpolateW(this, E, wRes)

Interpolate mapping matrices for given energy and store in passed triangular matrix if inelastic data component is allocated. Assumes upper triangular structure for wRes if E is positive and lower if it's negative.

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
real(kind=rk), intent(in) :: E

Transition energy to interpolate for

type(SparseRowData), intent(inout) :: wRes

Lower/upper triangular matrix to store the interpolated weights

module subroutine updateCRMData(this, hostModel, inputVars, updatePriority)

Update modelbound data based on input variable container

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(inout) :: this
class(ModelSurrogate), intent(in) :: hostModel

Host model passed to transitions

class(VariableContainer), intent(in) :: inputVars

Variable container used to calculate modelbound data

integer(kind=ik), intent(in), optional :: updatePriority

Priority for this update call

module subroutine crmCopy(this, name, container)

Retrieves data based on name format - assumes name of format "rate"//dataSpec//"index"//transIndex, where dataSpec is 0,1,2 and corresponds to rate type (0-particles,1-momentum,2-energy), and transIndex is the transition index in the data or "cssl//dataSpec//"index"//transIndex for the cross-section single harmonics where dataSpec is now the cross-section harmoinic index/column

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(in) :: this
character(len=*), intent(in) :: name
real(kind=rk), intent(inout), allocatable, dimension(..) :: container

module subroutine setElState(this, ID)

Setter for electronStateID

Arguments

Type IntentOptional Attributes Name
class(ModelboundCRMData), intent(inout) :: this
integer(kind=ik), intent(in) :: ID