FScalingDerivation Derived Type

type, public, extends(MatDerivation) :: FScalingDerivation

Extrapolates a distribution function to the boundary or last cell centre as a (numH,numV) matrix. Assumes that there are no m>0 harmonics for staggered grids. Expects 1-4 variables. The first is the extrapolated distribution, the second the density vector. The third is the staggered density if the variables are staggered. The value of the density at the boundary (scalar) is either the third or the fourth variable, and should be present only if extrapolating to the boundary.


Components

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

True only if user explicitly sets it to true

integer(kind=ik), public, dimension(2) :: exterpCoords

Coordinates used for extrapolation. Should be the indices of the cells closest to the boundary on the regular and dual grids, respectively.

logical, private :: hasBoundary

True if the processor this derivation is calculated on has the corresponding external boundary. If it doesn't, this derivation returns 0.

logical, private :: leftBoundary
integer(kind=ik), private :: numV

Number of velocity grid points

integer(kind=ik), private :: numH

Number of harmonics

logical, private :: staggeredVars

True if the distribution has staggered harmonics. Defaults to false.

logical, private :: extrapolateToBoundary

True if the extrapolation should be performed to the cell boundary and not the last cell centre before the boundary.


Type-Bound Procedures

procedure, public :: isDefined => isDefinedObject

  • interface

    private pure module function isDefinedObject(this) result(defined)

    Getter for userDefined

    Arguments

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

    Return Value logical

procedure, public :: makeDefined => makeDefinedObject

  • interface

    private pure module subroutine makeDefinedObject(this)

    Set userDefined to .true.

    Arguments

    Type IntentOptional Attributes Name
    class(Object), intent(inout) :: this

procedure, public :: makeUndefined => makeUndefinedObject

  • interface

    private pure module subroutine makeUndefinedObject(this)

    Set userDefined to .false.

    Arguments

    Type IntentOptional Attributes Name
    class(Object), intent(inout) :: this

procedure, public :: init => initFScaling

  • interface

    private module subroutine initFScaling(this, partitionObj, procRank, numV, leftBoundary, staggeredVars, extrapolateToBoundary)

    Initialize distribution scaling extrapolation derivation object

    Arguments

    Type IntentOptional Attributes Name
    class(FScalingDerivation), intent(inout) :: this
    type(Partition), intent(in) :: partitionObj

    Partition object used to determine local processor grid chunk

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

    Current processor rank

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

    Number of expected velocity grid points

    logical, intent(in), optional :: leftBoundary

    True if extrapolating to left boundary. Defaults to false.

    logical, intent(in), optional :: staggeredVars

    True if the distribution has staggered harmonics. Defaults to false.

    logical, intent(in), optional :: extrapolateToBoundary

    True if the extrapolation should be performed to the cell boundary and not the last cell centre before the boundary.

procedure, public :: calculate => calculateFScaling

  • interface

    private module function calculateFScaling(this, inputArray, indices) result(output)

    Arguments

    Type IntentOptional Attributes Name
    class(FScalingDerivation), intent(inout) :: this
    type(RealArray), intent(in), dimension(:) :: inputArray
    integer(kind=ik), intent(in), dimension(:) :: indices

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

procedure, public :: getScalingFactors

  • interface

    private pure module function getScalingFactors(this, inputArray, indices) result(output)

    Arguments

    Type IntentOptional Attributes Name
    class(FScalingDerivation), intent(in) :: this
    type(RealArray), intent(in), dimension(:) :: inputArray
    integer(kind=ik), intent(in), dimension(:) :: indices

    Return Value real(kind=rk), dimension(2)