atl_acoustic_2d_flux_module Module

module that holds all routines to calculate the flux for hyperbolic linearzied gas dynamic equations.

VK use atl_laxFriedrichFlux_module, only: atl_laxFriedAcoustic



Interfaces

public interface atl_acoustic_2d_numflux

Interface for fluxes of acoustic equations.

  • private subroutine atl_acoustic_2d_numflux_cube_vec(nTotalFaces, nSides, nFaceDofs, faceRep, faceFlux, leftPos, rightPos, var, acoustic, idir)

    VK ! call lax friedrich flux VK call atl_laxFriedAcoustic(left = leftstate, & VK & right = rightstate, & VK & acoustic = acoustic, & VK & flux = flux, & VK & iDir = idir ) VK VK faceFlux(left,iDof,var,2) = flux

    Arguments

    Type IntentOptional Attributes Name
    integer, intent(in) :: nTotalFaces
    integer, intent(in) :: nSides
    integer, intent(in) :: nFaceDofs
    real(kind=rk), intent(in) :: faceRep(nTotalFaces,nFaceDofs,3,2)
    real(kind=rk), intent(inout) :: faceFlux(nTotalFaces,nFaceDofs,3,2)
    integer, intent(in) :: leftPos(nSides)
    integer, intent(in) :: rightPos(nsides)
    integer, intent(in) :: var(3)
    type(atl_acoustic_type), intent(in) :: acoustic

    Datatype for acoustic equation include all background data

    integer, intent(in) :: idir

    Direction of the flow, used for background velocity


Functions

public function atl_acoustic_2d_physFlux(state, acoustic, iDir) result(flux)

Function for physical flux of the acoustic equation F, 1D? Since it is 1d, there need to be passed the correct background velocity (u0 for F - flux in x direction, v0 for G - flux in y direction, w0 for H - flux in z direction)

Arguments

Type IntentOptional Attributes Name
real(kind=rk), intent(in) :: state(acoustic%ndims+1)

State to compute the fluxes (rho, u, v, w)

type(atl_acoustic_type), intent(in) :: acoustic

Datatype for acoustic equation include all background data

integer, intent(in) :: iDir

Direction of flux, used fot background velocity

Return Value real(kind=rk), (acoustic%ndims+1)

The resulting flux in x direction

private function atl_acoustic_2d_numFlux_oneDir(left, right, acoustic, iDir) result(flux)

Arguments

Type IntentOptional Attributes Name
real(kind=rk), intent(in) :: left(3)

State to compute the fluxes (rho, u, v, w)

real(kind=rk), intent(in) :: right(3)

State to compute the fluxes (rho, u, v, w)

type(atl_acoustic_type), intent(in) :: acoustic

Datatype for acoustic equation include all background data

integer, intent(in) :: iDir

Direction of flux, used fot background velocity

Return Value real(kind=rk), (3)

The resulting flux in x direction


Subroutines

private subroutine atl_acoustic_2d_numflux_cube_vec(nTotalFaces, nSides, nFaceDofs, faceRep, faceFlux, leftPos, rightPos, var, acoustic, idir)

calculate flux of pure acoustic equation directly on the face-vector

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: nTotalFaces
integer, intent(in) :: nSides
integer, intent(in) :: nFaceDofs
real(kind=rk), intent(in) :: faceRep(nTotalFaces,nFaceDofs,3,2)
real(kind=rk), intent(inout) :: faceFlux(nTotalFaces,nFaceDofs,3,2)
integer, intent(in) :: leftPos(nSides)
integer, intent(in) :: rightPos(nsides)
integer, intent(in) :: var(3)
type(atl_acoustic_type), intent(in) :: acoustic

Datatype for acoustic equation include all background data

integer, intent(in) :: idir

Direction of the flow, used for background velocity