NullDGModel3D_t Derived Type

type, public, extends(DGModel3D) :: NullDGModel3D_t


Contents

Source Code


Components

TypeVisibilityAttributesNameInitial
type(MappedScalar3D), public :: dSdt
real(kind=prec), public :: dt
real(kind=prec), public :: entropy
type(MappedVector3D), public :: flux
type(MappedScalar3D), public :: fluxDivergence
type(SEMHex), public, pointer:: geometry
logical, public :: gradient_enabled =.false.
integer, public :: ioIterate =0
type(Mesh3D), public, pointer:: mesh
integer, public :: nvar
type(MappedScalar3D), public :: solution
type(MappedVector3D), public :: solutionGradient
type(MappedScalar3D), public :: source
real(kind=prec), public :: t
procedure(SELF_timeIntegrator), public, pointer:: timeIntegrator=> Euler_timeIntegrator
type(MappedScalar3D), public :: workSol

Type-Bound Procedures

procedure, public :: BoundaryFlux => BoundaryFlux_DGModel3D_t

procedure, public :: CalculateEntropy => CalculateEntropy_DGModel3D_t

procedure, public :: CalculateSolutionGradient => CalculateSolutionGradient_DGModel3D_t

procedure, public :: CalculateTendency => CalculateTendency_DGModel3D_t

procedure, public :: Euler_timeIntegrator

  • public subroutine Euler_timeIntegrator(this, tn)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(inout) :: this
    real(kind=prec), intent(in) :: tn

procedure, public :: FluxMethod => fluxmethod_DGModel3D_t

procedure, public :: ForwardStep => ForwardStep_Model

  • public subroutine ForwardStep_Model(this, tn, dt, ioInterval)

    Forward steps the model using the associated tendency procedure and time integrator

    Read more…

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(inout) :: this
    real(kind=prec), intent(in) :: tn
    real(kind=prec), intent(in) :: dt
    real(kind=prec), intent(in) :: ioInterval

procedure, public :: Free => Free_DGModel3D_t

procedure, public :: GetSimulationTime

  • public subroutine GetSimulationTime(this, t)

    Returns the current simulation time stored in the model % t attribute

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(out) :: t

procedure, public :: IncrementIOCounter

  • public subroutine IncrementIOCounter(this)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(inout) :: this

procedure, public :: Init => Init_DGModel3D_t

  • public subroutine Init_DGModel3D_t(this, nvar, mesh, geometry)

    Arguments

    TypeIntentOptionalAttributesName
    class(DGModel3D_t), intent(out) :: this
    integer, intent(in) :: nvar
    type(Mesh3D), intent(in), target:: mesh
    type(SEMHex), intent(in), target:: geometry

procedure, public :: LowStorageRK2_timeIntegrator

  • public subroutine LowStorageRK2_timeIntegrator(this, tn)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(inout) :: this
    real(kind=prec), intent(in) :: tn

procedure, public :: LowStorageRK3_timeIntegrator

  • public subroutine LowStorageRK3_timeIntegrator(this, tn)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(inout) :: this
    real(kind=prec), intent(in) :: tn

procedure, public :: LowStorageRK4_timeIntegrator

  • public subroutine LowStorageRK4_timeIntegrator(this, tn)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(inout) :: this
    real(kind=prec), intent(in) :: tn

procedure, public :: PreTendency => PreTendency_Model

  • public subroutine PreTendency_Model(this)

    PreTendency is a template routine that is used to house any additional calculations that you want to execute at the beginning of the tendency calculation routine. This default PreTendency simply returns back to the caller without executing any instructions

    Read more…

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(inout) :: this

procedure, public :: PrintType => PrintType_Model

  • public subroutine PrintType_Model(this)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this

procedure, public :: ReadModel => Read_DGModel3D_t

  • public subroutine Read_DGModel3D_t(this, fileName)

    Arguments

    TypeIntentOptionalAttributesName
    class(DGModel3D_t), intent(inout) :: this
    character, intent(in) :: fileName

procedure, public :: ReportEntropy => ReportEntropy_Model

  • public subroutine ReportEntropy_Model(this)

    Base method for reporting the entropy of a model to stdout. Only override this procedure if additional reporting is needed. Alternatively, if you think additional reporting would be valuable for all models, open a pull request with modifications to this base method.

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this

procedure, public :: SetBoundaryCondition => setboundarycondition_DGModel3D_t

  • public subroutine setboundarycondition_DGModel3D_t(this)

    Boundary conditions for the solution are set to 0 for the external state to provide radiation type boundary conditions.

    Arguments

    TypeIntentOptionalAttributesName
    class(DGModel3D_t), intent(inout) :: this

procedure, public :: SetGradientBoundaryCondition => setgradientboundarycondition_DGModel3D_t

  • public subroutine setgradientboundarycondition_DGModel3D_t(this)

    Boundary conditions for the solution are set to 0 for the external state to provide radiation type boundary conditions.

    Arguments

    TypeIntentOptionalAttributesName
    class(DGModel3D_t), intent(inout) :: this

procedure, public :: SetMetadata => SetMetadata_DGModel3D_t

procedure, public :: SetSimulationTime

  • public subroutine SetSimulationTime(this, t)

    Sets the model % t attribute with the provided simulation time

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(inout) :: this
    real(kind=prec), intent(in) :: t

generic, public :: SetTimeIntegrator => SetTimeIntegrator_withChar

  • public subroutine SetTimeIntegrator_withChar(this, integrator)

    Sets the time integrator method, using a character input

    Read more…

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(inout) :: this
    character, intent(in) :: integrator

procedure, public :: SourceMethod => sourcemethod_DGModel3D_t

procedure, public :: UpdateGRK2 => UpdateGRK2_DGModel3D_t

procedure, public :: UpdateGRK3 => UpdateGRK3_DGModel3D_t

procedure, public :: UpdateGRK4 => UpdateGRK4_DGModel3D_t

procedure, public :: UpdateSolution => UpdateSolution_DGModel3D_t

  • public subroutine UpdateSolution_DGModel3D_t(this, dt)

    Computes a solution update as , where dt is either provided through the interface or taken as the Model's stored time step size (model % dt)

    Arguments

    TypeIntentOptionalAttributesName
    class(DGModel3D_t), intent(inout) :: this
    real(kind=prec), intent(in), optional :: dt

procedure, public :: WriteModel => Write_DGModel3D_t

  • public subroutine Write_DGModel3D_t(this, fileName)

    Arguments

    TypeIntentOptionalAttributesName
    class(DGModel3D_t), intent(inout) :: this
    character, intent(in), optional :: fileName

procedure, public :: WriteTecplot => WriteTecplot_DGModel3D_t

  • public subroutine WriteTecplot_DGModel3D_t(this, filename)

    Arguments

    TypeIntentOptionalAttributesName
    class(DGModel3D_t), intent(inout) :: this
    character, intent(in), optional :: filename

procedure, public :: entropy_func => entropy_func_Model

  • public pure function entropy_func_Model(this, s) result(e)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)

    Return Value real(kind=prec)

procedure, public :: flux1D => flux1d_Model

  • public pure function flux1d_Model(this, s, dsdx) result(flux)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)
    real(kind=prec), intent(in) :: dsdx(1:this%nvar)

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: flux2D => flux2d_Model

  • public pure function flux2d_Model(this, s, dsdx) result(flux)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)
    real(kind=prec), intent(in) :: dsdx(1:this%nvar,1:2)

    Return Value real(kind=prec)(1:this%nvar,1:2)

procedure, public :: flux3D => flux3d_Model

  • public pure function flux3d_Model(this, s, dsdx) result(flux)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)
    real(kind=prec), intent(in) :: dsdx(1:this%nvar,1:3)

    Return Value real(kind=prec)(1:this%nvar,1:3)

procedure, public :: hbc1d_NoNormalFlow => hbc1d_Generic_Model

  • public pure function hbc1d_Generic_Model(this, s, nhat) result(exts)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)
    real(kind=prec), intent(in) :: nhat

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: hbc1d_Prescribed => hbc1d_Prescribed_Model

  • public pure function hbc1d_Prescribed_Model(this, x, t) result(exts)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: x
    real(kind=prec), intent(in) :: t

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: hbc1d_Radiation => hbc1d_Generic_Model

  • public pure function hbc1d_Generic_Model(this, s, nhat) result(exts)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)
    real(kind=prec), intent(in) :: nhat

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: hbc2d_NoNormalFlow => hbc2d_Generic_Model

  • public pure function hbc2d_Generic_Model(this, s, nhat) result(exts)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)
    real(kind=prec), intent(in) :: nhat(1:2)

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: hbc2d_Prescribed => hbc2d_Prescribed_Model

  • public pure function hbc2d_Prescribed_Model(this, x, t) result(exts)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: x(1:2)
    real(kind=prec), intent(in) :: t

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: hbc2d_Radiation => hbc2d_Generic_Model

  • public pure function hbc2d_Generic_Model(this, s, nhat) result(exts)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)
    real(kind=prec), intent(in) :: nhat(1:2)

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: hbc3d_NoNormalFlow => hbc3d_Generic_Model

  • public pure function hbc3d_Generic_Model(this, s, nhat) result(exts)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)
    real(kind=prec), intent(in) :: nhat(1:3)

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: hbc3d_Prescribed => hbc3d_Prescribed_Model

  • public pure function hbc3d_Prescribed_Model(this, x, t) result(exts)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: x(1:3)
    real(kind=prec), intent(in) :: t

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: hbc3d_Radiation => hbc3d_Generic_Model

  • public pure function hbc3d_Generic_Model(this, s, nhat) result(exts)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)
    real(kind=prec), intent(in) :: nhat(1:3)

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: pbc1d_NoNormalFlow => pbc1d_Generic_Model

  • public pure function pbc1d_Generic_Model(this, dsdx, nhat) result(extDsdx)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: dsdx(1:this%nvar)
    real(kind=prec), intent(in) :: nhat

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: pbc1d_Prescribed => pbc1d_Prescribed_Model

  • public pure function pbc1d_Prescribed_Model(this, x, t) result(extDsdx)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: x
    real(kind=prec), intent(in) :: t

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: pbc1d_Radiation => pbc1d_Generic_Model

  • public pure function pbc1d_Generic_Model(this, dsdx, nhat) result(extDsdx)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: dsdx(1:this%nvar)
    real(kind=prec), intent(in) :: nhat

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: pbc2d_NoNormalFlow => pbc2d_Generic_Model

  • public pure function pbc2d_Generic_Model(this, dsdx, nhat) result(extDsdx)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: dsdx(1:this%nvar,1:2)
    real(kind=prec), intent(in) :: nhat(1:2)

    Return Value real(kind=prec)(1:this%nvar,1:2)

procedure, public :: pbc2d_Prescribed => pbc2d_Prescribed_Model

  • public pure function pbc2d_Prescribed_Model(this, x, t) result(extDsdx)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: x(1:2)
    real(kind=prec), intent(in) :: t

    Return Value real(kind=prec)(1:this%nvar,1:2)

procedure, public :: pbc2d_Radiation => pbc2d_Generic_Model

  • public pure function pbc2d_Generic_Model(this, dsdx, nhat) result(extDsdx)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: dsdx(1:this%nvar,1:2)
    real(kind=prec), intent(in) :: nhat(1:2)

    Return Value real(kind=prec)(1:this%nvar,1:2)

procedure, public :: pbc3d_NoNormalFlow => pbc3d_Generic_Model

  • public pure function pbc3d_Generic_Model(this, dsdx, nhat) result(extDsdx)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: dsdx(1:this%nvar,1:3)
    real(kind=prec), intent(in) :: nhat(1:3)

    Return Value real(kind=prec)(1:this%nvar,1:3)

procedure, public :: pbc3d_Prescribed => pbc3d_Prescribed_Model

  • public pure function pbc3d_Prescribed_Model(this, x, t) result(extDsdx)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: x(1:3)
    real(kind=prec), intent(in) :: t

    Return Value real(kind=prec)(1:this%nvar,1:3)

procedure, public :: pbc3d_Radiation => pbc3d_Generic_Model

  • public pure function pbc3d_Generic_Model(this, dsdx, nhat) result(extDsdx)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: dsdx(1:this%nvar,1:3)
    real(kind=prec), intent(in) :: nhat(1:3)

    Return Value real(kind=prec)(1:this%nvar,1:3)

procedure, public :: riemannflux1d => riemannflux1d_Model

  • public pure function riemannflux1d_Model(this, sL, sR, dsdx, nhat) result(flux)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: sL(1:this%nvar)
    real(kind=prec), intent(in) :: sR(1:this%nvar)
    real(kind=prec), intent(in) :: dsdx(1:this%nvar)
    real(kind=prec), intent(in) :: nhat

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: riemannflux2d => riemannflux2d_Model

  • public pure function riemannflux2d_Model(this, sL, sR, dsdx, nhat) result(flux)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: sL(1:this%nvar)
    real(kind=prec), intent(in) :: sR(1:this%nvar)
    real(kind=prec), intent(in) :: dsdx(1:this%nvar,1:2)
    real(kind=prec), intent(in) :: nhat(1:2)

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: riemannflux3d => riemannflux3d_Model

  • public pure function riemannflux3d_Model(this, sL, sR, dsdx, nhat) result(flux)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: sL(1:this%nvar)
    real(kind=prec), intent(in) :: sR(1:this%nvar)
    real(kind=prec), intent(in) :: dsdx(1:this%nvar,1:3)
    real(kind=prec), intent(in) :: nhat(1:3)

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: source1d => source1d_Model

  • public pure function source1d_Model(this, s, dsdx) result(source)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)
    real(kind=prec), intent(in) :: dsdx(1:this%nvar)

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: source2d => source2d_Model

  • public pure function source2d_Model(this, s, dsdx) result(source)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)
    real(kind=prec), intent(in) :: dsdx(1:this%nvar,1:2)

    Return Value real(kind=prec)(1:this%nvar)

procedure, public :: source3d => source3d_Model

  • public pure function source3d_Model(this, s, dsdx) result(source)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(in) :: this
    real(kind=prec), intent(in) :: s(1:this%nvar)
    real(kind=prec), intent(in) :: dsdx(1:this%nvar,1:3)

    Return Value real(kind=prec)(1:this%nvar)

Source Code

  type,extends(dgmodel3d) :: NullDGModel3D_t
    ! Add any additional attributes here that are specific to your model

  contains
    !   procedure :: hbc3d_Prescribed => hbc3d_Generic_Model
    !   procedure :: hbc3d_Radiation => hbc3d_Generic_Model
    !   procedure :: hbc3d_NoNormalFlow => hbc3d_Generic_Model
    !   procedure :: pbc3d_Prescribed => pbc3d_Generic_Model
    !   procedure :: pbc3d_Radiation => pbc3d_Generic_Model
    !   procedure :: pbc3d_NoNormalFlow => pbc3d_Generic_Model
    !   procedure :: SetMetadata => SetMetadata_NullDGModel3D_t
    !   procedure :: pretendency => pretendency_NullDGModel3D_t
    !   procedure :: entropy_func => entropy_func_NullDGModel3D_t
    !   procedure :: flux3d => flux3d_NullDGModel3D_t
    !   procedure :: riemannflux3d => riemannflux3d_NullDGModel3D_t
    !   procedure :: source3d => source3d_NullDGModel3D_t

  endtype NullDGModel3D_t