hrweno_fluxes Module

This module contains basic flux schemes for scalar problems. They are mostly intented to help test the other modules. The WENO schemes themselves are applicable to scalar and multiple component problems. Source: ICASE 97-65 by Shu, 1997.


Uses

  • module~~hrweno_fluxes~~UsesGraph module~hrweno_fluxes hrweno_fluxes module~hrweno_kinds hrweno_kinds module~hrweno_fluxes->module~hrweno_kinds iso_fortran_env iso_fortran_env module~hrweno_kinds->iso_fortran_env

Used by

  • module~~hrweno_fluxes~~UsedByGraph module~hrweno_fluxes hrweno_fluxes program~example1_burgers_1d_fv example1_burgers_1d_fv program~example1_burgers_1d_fv->module~hrweno_fluxes program~example_pbe_2d_fv example_pbe_2d_fv program~example_pbe_2d_fv->module~hrweno_fluxes

Functions

public pure function lax_friedrichs(f, vm, vp, x, t, alpha) result(res)

Monotone Lax-Friedrichs flux. It is more dissipative than the Godunov method, but computationally less demanding. Source: Equation 2.72, page 21.

Read more…

Arguments

Type IntentOptional Attributes Name
procedure(flux) :: f

flux function,

real(kind=rk), intent(in) :: vm

left (minus) reconstruction,

real(kind=rk), intent(in) :: vp

right (plus) reconstruction,

real(kind=rk), intent(in) :: x(:)

at flux interface,

real(kind=rk), intent(in) :: t

time,

real(kind=rk), intent(in) :: alpha

in the domain on the problem

Return Value real(kind=rk)

public pure function godunov(f, vm, vp, x, t) result(res)

Monotone Godunov flux. It is less dissipative than the Lax-Friedrichs method, but computationally more demanding because of the if constructs. Source: Equation 2.70, page 21.

Read more…

Arguments

Type IntentOptional Attributes Name
procedure(flux) :: f

flux function,

real(kind=rk), intent(in) :: vm

left (minus) reconstruction,

real(kind=rk), intent(in) :: vp

right (plus) reconstruction,

real(kind=rk), intent(in) :: x(:)

at flux interface,

real(kind=rk), intent(in) :: t

time,

Return Value real(kind=rk)