fixed_boltzmann_stencil_gen_class Module

Stencil generator for v stencil of Boltzmann collision operator emission/absorption terms with fixed cross-sections



Interfaces

interface

  • private module subroutine initBoltzGen(this, vspaceObj, transitionIndex, fixedWIndex, lNum, absorptionTerm, dbTerm)

    Boltzmann emission/absorption stencil value generator initialization routine

    Arguments

    Type IntentOptional Attributes Name
    class(FixedBoltzmannStencilGen), intent(inout) :: this
    type(VSpace), intent(in) :: vspaceObj

    VSpace object used to get grid

    integer(kind=ik), intent(in) :: transitionIndex

    Index of transition whose cross-section is required

    integer(kind=ik), intent(in) :: fixedWIndex

    Fixed inelastic mapping index associated with the transition

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

    Harmonic index of cross-section required. Used only if absorption term

    logical, intent(in), optional :: absorptionTerm

    True if this is an absorption term stencil generator. Defaults to false.

    logical, intent(in), optional :: dbTerm

    True if this is a detailed balance term stencil generator. Defaults to false.

interface

  • private module subroutine calcBoltzVals(this, varCont, res, mbData, hostModel)

    Calculate Boltzmann emission/absorption stencil values in place (does not depend on varCont)

    Arguments

    Type IntentOptional Attributes Name
    class(FixedBoltzmannStencilGen), intent(inout) :: this
    type(VariableContainer), intent(in) :: varCont
    type(RealArray), intent(inout), allocatable, dimension(:) :: res
    class(ModelboundData), intent(in), optional :: mbData
    class(ModelSurrogate), intent(in), optional :: hostModel

Derived Types

type, public, extends(JaggedArrayGenerator) ::  FixedBoltzmannStencilGen

JaggedArrayGenerator for calculating velocity component of Boltzmann collision operator emissionabsorption terms with fixed inelastic mappings. Expects fixed vStencil derived from the weights matrix if absorption term, otherwise vStencil = [0]. Generates Wvsigma(v) for absorption term or - v*sigma(v) for emission. If cross-section depends on x, it assumes no halo.

Components

Type Visibility Attributes Name Initial
logical, public :: userDefined = .false.

True only if user explicitly sets it to true

real(kind=rk), private, allocatable, dimension(:) :: vGridCopy

Copy of velocity grid

integer(kind=ik), private :: transitionIndex

Transition index associated with this stencil generator

integer(kind=ik), private :: fixedWIndex

Inelastic mapping index associated with this stencil generator

integer(kind=ik), private :: lNum

Harmonic index of cross-section required, valid only if emission term

type(SparseRowData), private, allocatable :: bufferW
logical, private :: absorptionTerm

True if this is for an absorption term

logical, private :: dbTerm

True if this is for a detailed balance term

Type-Bound Procedures

procedure, public :: isDefined => isDefinedObject
procedure, public :: makeDefined => makeDefinedObject
procedure, public :: makeUndefined => makeUndefinedObject
procedure, public :: calculate
procedure, public :: init => initBoltzGen
procedure, public :: calculateInPlace => calcBoltzVals