geometry_class Module

Houses geometry object responsible for storing cell widths and Jacobians


Used by


Interfaces

interface

  • private pure module function getCellWidths(this, dualGrid, extendedBoundaryCells) result(dx)

    Getter for cellWidths, if dualGrid is true returns dual grid values based off of cellWidths. If extendedBoundaryCells is true will extend dual grid cells at boundaries to the boundaries themselves if not periodic (defaults to true)

    Arguments

    Type IntentOptional Attributes Name
    class(Geometry), intent(in) :: this
    logical, intent(in), optional :: dualGrid
    logical, intent(in), optional :: extendedBoundaryCells

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

interface

  • private pure module function getJacobianLeft(this, dualGrid) result(jLeft)

    Getter for jacobianLeft, if dualGrid is true returns dual grid values

    Arguments

    Type IntentOptional Attributes Name
    class(Geometry), intent(in) :: this
    logical, intent(in), optional :: dualGrid

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

interface

  • private pure module function getJacobianRight(this, dualGrid, extendedBoundaryCells) result(jRight)

    Getter for jacobianRight, if dualGrid is true returns dual grid values

    Arguments

    Type IntentOptional Attributes Name
    class(Geometry), intent(in) :: this
    logical, intent(in), optional :: dualGrid
    logical, intent(in), optional :: extendedBoundaryCells

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

interface

  • private pure module function getJacobianCentre(this, dualGrid, extendedBoundaryCells) result(jCentre)

    Getter for jacobianCentre, if dualGrid is true returns dual grid values

    Arguments

    Type IntentOptional Attributes Name
    class(Geometry), intent(in) :: this
    logical, intent(in), optional :: dualGrid
    logical, intent(in), optional :: extendedBoundaryCells

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

interface

  • private pure module function getLinInterp(this, dualGrid) result(linInterp)

    Getter for linInterp, if dualGrid is true returns dual grid values (0.5 everywhere)

    Arguments

    Type IntentOptional Attributes Name
    class(Geometry), intent(in) :: this
    logical, intent(in), optional :: dualGrid

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

interface

  • private pure module function isPeriodic(this) result(periodic)

    Getter for periodicGrid

    Arguments

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

    Return Value logical

interface

  • private pure module subroutine initGeometry(this, cellWidths, jLeft, jRight, periodicGrid)

    Geometry initialization routine

    Arguments

    Type IntentOptional Attributes Name
    class(Geometry), intent(inout) :: this
    real(kind=rk), intent(in), dimension(:) :: cellWidths

    Spatial cell widths

    real(kind=rk), intent(in), dimension(:) :: jLeft

    Left face jacobian values

    real(kind=rk), intent(in), dimension(:) :: jRight

    Right face jacobian values

    logical, intent(in), optional :: periodicGrid

    Set to true if the grid is periodic


Derived Types

type, public, extends(Object) ::  Geometry

Geometry object storing spatial cell widths, cell face and cell centre jacobians, as well as linear interpolation coefficients

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(:) :: cellWidths

x cell widths

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

Jacobian of left cell faces

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

Jacobian of right cell faces

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

Jacobian at cell centre (avg of faces)

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

Linear interpolation coefficient for inner faces

logical, private :: periodicGrid

Type-Bound Procedures

procedure, public :: isDefined => isDefinedObject
procedure, public :: makeDefined => makeDefinedObject
procedure, public :: makeUndefined => makeUndefinedObject
procedure, public :: getCellWidths
procedure, public :: getJacobianLeft
procedure, public :: getJacobianRight
procedure, public :: getJacobianCentre
procedure, public :: isPeriodic
procedure, public :: getLinInterp
procedure, public :: init => initGeometry