DGModel2D_t Derived Type

type, public, extends(Model) :: DGModel2D_t


Inherits

type~~dgmodel2d_t~~InheritsGraph type~dgmodel2d_t DGModel2D_t type~mappedscalar2d MappedScalar2D type~dgmodel2d_t->type~mappedscalar2d solution, source, fluxDivergence, dSdt, workSol type~semquad SEMQuad type~dgmodel2d_t->type~semquad geometry type~mappedvector2d MappedVector2D type~dgmodel2d_t->type~mappedvector2d solutionGradient, flux type~model Model type~dgmodel2d_t->type~model type~mesh2d Mesh2D type~dgmodel2d_t->type~mesh2d mesh type~mappedscalar2d_t MappedScalar2D_t type~mappedscalar2d->type~mappedscalar2d_t type~tensor2d Tensor2D type~semquad->type~tensor2d dxds, dsdx type~scalar2d Scalar2D type~semquad->type~scalar2d nScale, J type~vector2d Vector2D type~semquad->type~vector2d x, nHat type~mappedvector2d_t MappedVector2D_t type~mappedvector2d->type~mappedvector2d_t type~mesh2d_t Mesh2D_t type~mesh2d->type~mesh2d_t type~mappedscalar2d_t->type~semquad geometry type~mappedscalar2d_t->type~scalar2d type~semmesh SEMMesh type~mesh2d_t->type~semmesh type~mappedvector2d_t->type~semquad geometry type~mappedvector2d_t->type~vector2d type~tensor2d_t Tensor2D_t type~tensor2d->type~tensor2d_t type~scalar2d_t Scalar2D_t type~scalar2d->type~scalar2d_t type~vector2d_t Vector2D_t type~vector2d->type~vector2d_t type~domaindecomposition DomainDecomposition type~semmesh->type~domaindecomposition decomp type~self_dataobj SELF_DataObj type~vector2d_t->type~self_dataobj type~tensor2d_t->type~self_dataobj type~scalar2d_t->type~self_dataobj type~domaindecomposition_t DomainDecomposition_t type~domaindecomposition->type~domaindecomposition_t EquationParser EquationParser type~self_dataobj->EquationParser eqn type~metadata Metadata type~self_dataobj->type~metadata meta type~lagrange Lagrange type~self_dataobj->type~lagrange interp type~lagrange_t Lagrange_t type~lagrange->type~lagrange_t c_ptr c_ptr type~lagrange_t->c_ptr blas_handle

Inherited by

type~~dgmodel2d_t~~InheritedByGraph type~dgmodel2d_t DGModel2D_t type~dgmodel2d DGModel2D type~dgmodel2d->type~dgmodel2d_t type~dgmodel2d~2 DGModel2D type~dgmodel2d~2->type~dgmodel2d_t type~advection_diffusion_2d_t advection_diffusion_2d_t type~advection_diffusion_2d_t->type~dgmodel2d type~linearshallowwater2d_t LinearShallowWater2D_t type~linearshallowwater2d_t->type~dgmodel2d type~nulldgmodel2d_t NullDGModel2D_t type~nulldgmodel2d_t->type~dgmodel2d type~lineareuler2d_t LinearEuler2D_t type~lineareuler2d_t->type~dgmodel2d type~linearshallowwater2d LinearShallowWater2D type~linearshallowwater2d->type~linearshallowwater2d_t type~advection_diffusion_2d~2 advection_diffusion_2d type~advection_diffusion_2d~2->type~advection_diffusion_2d_t type~lineareuler2d LinearEuler2D type~lineareuler2d->type~lineareuler2d_t type~advection_diffusion_2d advection_diffusion_2d type~advection_diffusion_2d->type~advection_diffusion_2d_t type~nulldgmodel2d NullDGModel2D type~nulldgmodel2d->type~nulldgmodel2d_t type~linearshallowwater2d~2 LinearShallowWater2D type~linearshallowwater2d~2->type~linearshallowwater2d_t type~nulldgmodel2d~2 NullDGModel2D type~nulldgmodel2d~2->type~nulldgmodel2d_t type~lineareuler2d~2 LinearEuler2D type~lineareuler2d~2->type~lineareuler2d_t

Contents

Source Code


Components

TypeVisibilityAttributesNameInitial
type(MappedScalar2D), public :: dSdt
real(kind=prec), public :: dt
real(kind=prec), public :: entropy
type(MappedVector2D), public :: flux
type(MappedScalar2D), public :: fluxDivergence
type(SEMQuad), public, pointer:: geometry
logical, public :: gradient_enabled =.false.
integer, public :: ioIterate =0
type(Mesh2D), public, pointer:: mesh
integer, public :: nvar
logical, public :: prescribed_bcs_enabled =.true.
type(MappedScalar2D), public :: solution
type(MappedVector2D), public :: solutionGradient
type(MappedScalar2D), public :: source
real(kind=prec), public :: t
logical, public :: tecplot_enabled =.true.
procedure(SELF_timeIntegrator), public, pointer:: timeIntegrator=> Euler_timeIntegrator
type(MappedScalar2D), public :: workSol

Type-Bound Procedures

procedure, public :: AdditionalFree => AdditionalFree_Model

  • public subroutine AdditionalFree_Model(this)

    Arguments

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

procedure, public :: AdditionalInit => AdditionalInit_Model

  • public subroutine AdditionalInit_Model(this)

    Arguments

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

procedure, public :: AdditionalOutput => AdditionalOutput_Model

  • public subroutine AdditionalOutput_Model(this, fileid)

    Arguments

    TypeIntentOptionalAttributesName
    class(Model), intent(inout) :: this
    integer(kind=HID_T), intent(in) :: fileid

procedure, public :: BoundaryFlux => BoundaryFlux_DGModel2D_t

procedure, public :: CalculateEntropy => CalculateEntropy_DGModel2D_t

procedure, public :: CalculateSolutionGradient => CalculateSolutionGradient_DGModel2D_t

procedure, public :: CalculateTendency => CalculateTendency_DGModel2D_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_DGModel2D_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_DGModel2D_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_DGModel2D_t

  • public subroutine Init_DGModel2D_t(this, mesh, geometry)

    Arguments

    TypeIntentOptionalAttributesName
    class(DGModel2D_t), intent(out) :: this
    type(Mesh2D), intent(in), target:: mesh
    type(SEMQuad), 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_DGModel2D_t

  • public subroutine Read_DGModel2D_t(this, fileName)

    Arguments

    TypeIntentOptionalAttributesName
    class(DGModel2D_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 :: ReportMetrics => ReportMetrics_DGModel2D_t

  • public subroutine ReportMetrics_DGModel2D_t(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(DGModel2D_t), intent(inout) :: this

procedure, public :: ReportUserMetrics => ReportUserMetrics_Model

  • public subroutine ReportUserMetrics_Model(this)

    Method that can be overridden by users to report their own custom metrics after file io

    Arguments

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

procedure, public :: SetBoundaryCondition => setboundarycondition_DGModel2D_t

  • public subroutine setboundarycondition_DGModel2D_t(this)

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

    Arguments

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

procedure, public :: SetGradientBoundaryCondition => setgradientboundarycondition_DGModel2D_t

  • public subroutine setgradientboundarycondition_DGModel2D_t(this)

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

    Arguments

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

procedure, public :: SetMetadata => SetMetadata_DGModel2D_t

procedure, public :: SetNumberOfVariables => SetNumberOfVariables_Model

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

procedure, private :: SetSolutionFromChar_DGModel2D_t

procedure, private :: SetSolutionFromEqn_DGModel2D_t

  • public subroutine SetSolutionFromEqn_DGModel2D_t(this, eqn)

    Arguments

    TypeIntentOptionalAttributesName
    class(DGModel2D_t), intent(inout) :: this
    type(EquationParser), intent(in) :: eqn(1:this%solution%nVar)

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_DGModel2D_t

procedure, public :: UpdateGRK2 => UpdateGRK2_DGModel2D_t

procedure, public :: UpdateGRK3 => UpdateGRK3_DGModel2D_t

procedure, public :: UpdateGRK4 => UpdateGRK4_DGModel2D_t

procedure, public :: UpdateSolution => UpdateSolution_DGModel2D_t

  • public subroutine UpdateSolution_DGModel2D_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(DGModel2D_t), intent(inout) :: this
    real(kind=prec), intent(in), optional :: dt

procedure, public :: WriteModel => Write_DGModel2D_t

  • public subroutine Write_DGModel2D_t(this, fileName)

    Arguments

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

procedure, public :: WriteTecplot => WriteTecplot_DGModel2D_t

  • public subroutine WriteTecplot_DGModel2D_t(this, filename)

    Arguments

    TypeIntentOptionalAttributesName
    class(DGModel2D_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(Model) :: DGModel2D_t
    type(MappedScalar2D)   :: solution
    type(MappedVector2D)   :: solutionGradient
    type(MappedVector2D)   :: flux
    type(MappedScalar2D)   :: source
    type(MappedScalar2D)   :: fluxDivergence
    type(MappedScalar2D)   :: dSdt
    type(MappedScalar2D)   :: workSol
    type(Mesh2D),pointer   :: mesh
    type(SEMQuad),pointer  :: geometry

  contains

    procedure :: Init => Init_DGModel2D_t
    procedure :: SetMetadata => SetMetadata_DGModel2D_t
    procedure :: Free => Free_DGModel2D_t

    procedure :: CalculateEntropy => CalculateEntropy_DGModel2D_t
    procedure :: BoundaryFlux => BoundaryFlux_DGModel2D_t
    procedure :: FluxMethod => fluxmethod_DGModel2D_t
    procedure :: SourceMethod => sourcemethod_DGModel2D_t
    procedure :: SetBoundaryCondition => setboundarycondition_DGModel2D_t
    procedure :: SetGradientBoundaryCondition => setgradientboundarycondition_DGModel2D_t
    procedure :: ReportMetrics => ReportMetrics_DGModel2D_t

    procedure :: UpdateSolution => UpdateSolution_DGModel2D_t

    procedure :: UpdateGRK2 => UpdateGRK2_DGModel2D_t
    procedure :: UpdateGRK3 => UpdateGRK3_DGModel2D_t
    procedure :: UpdateGRK4 => UpdateGRK4_DGModel2D_t

    procedure :: CalculateSolutionGradient => CalculateSolutionGradient_DGModel2D_t
    procedure :: CalculateTendency => CalculateTendency_DGModel2D_t

    generic :: SetSolution => SetSolutionFromChar_DGModel2D_t, &
      SetSolutionFromEqn_DGModel2D_t
    procedure,private :: SetSolutionFromChar_DGModel2D_t
    procedure,private :: SetSolutionFromEqn_DGModel2D_t

    procedure :: ReadModel => Read_DGModel2D_t
    procedure :: WriteModel => Write_DGModel2D_t
    procedure :: WriteTecplot => WriteTecplot_DGModel2D_t

  endtype DGModel2D_t