atl_modg_2d_maxwell_kernel_module Module

Module for routines and datatypes of MOdal Discontinuous Galerkin (MODG) scheme. This scheme is a spectral scheme for linear, purley hyperbolic partial differential equation systems.


Uses

Used by

  • module~~atl_modg_2d_maxwell_kernel_module~~UsedByGraph module~atl_modg_2d_maxwell_kernel_module atl_modg_2d_maxwell_kernel_module proc~compute_rhs_cubes_modg_2d compute_rhs_cubes_modg_2d proc~compute_rhs_cubes_modg_2d->module~atl_modg_2d_maxwell_kernel_module

Contents


Functions

private function atl_physFluxMaxwell_2d(state, material) result(flux)

Function for physical flux of the Maxwell equations in terms of D and B.

Arguments

TypeIntentOptionalAttributesName
real(kind=rk), intent(in) :: state(7)

State to compute the fluxes from (D,B).

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

Material parameters (mu, epsilon) the flux calculation

Return Value real(kind=rk)(7)

The resulting flux in x direction


Subroutines

public subroutine atl_modg_maxwell_2d_numFlux(equation, facedata, scheme, poly_proj, material)

Calculate the numerical flux for Maxwell equation and MODG scheme

Arguments

TypeIntentOptionalAttributesName
type(atl_Equations_type), intent(in) :: equation

The equation you solve.

type(atl_facedata_type), intent(inout) :: facedata

The face representation of the state.

type(atl_modg_2d_scheme_type), intent(inout) :: scheme

Parameters of the modal dg scheme

type(ply_poly_project_type), intent(inout) :: poly_proj

Data for projection method

type(atl_material_type), intent(inout) :: material

Material description for the faces on the current level

public subroutine atl_modg_maxwell_2d_physFlux_const(equation, res, state, iElem, iDir, penalizationData, poly_proj, material, nodal_data, nodal_GradData, nodal_res, ElemLength, scheme_min, scheme_current)

Calculate the physical flux for the MODG scheme and Maxwell equation. used for both space P and Q

Arguments

TypeIntentOptionalAttributesName
type(atl_Equations_type), intent(in) :: equation

The equation you solve.

real(kind=rk), intent(inout) :: res(:,:)

To store the resulting phy flux in modal form

real(kind=rk), intent(in), optional :: state(:,:)

The state of the equation

integer, intent(in) :: iElem

The current Element

integer, intent(in) :: iDir

The current Direction

type(atl_penalizationData_type), intent(inout) :: penalizationData

The penalization data

type(ply_poly_project_type), intent(inout) :: poly_proj

Poly project

type(atl_material_type), intent(inout) :: material

Material description for the faces on the current level

real(kind=rk), intent(in), optional :: nodal_data(:,:)

The state in nodal form

real(kind=rk), intent(in), optional :: nodal_GradData(:,:,:)
real(kind=rk), intent(inout) :: nodal_res(:,:)
real(kind=rk), intent(in) :: ElemLength

Length of the element

type(atl_scheme_type), intent(inout) :: scheme_min

The scheme information

type(atl_scheme_type), intent(inout) :: scheme_current

public subroutine atl_modg_maxwell_2d_physFlux_NonConst(equation, res, state, iElem, iDir, penalizationData, poly_proj, material, nodal_data, nodal_GradData, nodal_res, ElemLength, scheme_min, scheme_current)

Calculate the physical flux for the MODG scheme and Maxwell equation. used for both space P and Q

Arguments

TypeIntentOptionalAttributesName
type(atl_Equations_type), intent(in) :: equation

The equation you solve.

real(kind=rk), intent(inout) :: res(:,:)

To store the resulting phy flux in modal form

real(kind=rk), intent(in), optional :: state(:,:)

The state of the equation

integer, intent(in) :: iElem

The current Element

integer, intent(in) :: iDir

The current Direction

type(atl_penalizationData_type), intent(inout) :: penalizationData

The penalization data

type(ply_poly_project_type), intent(inout) :: poly_proj

Poly project

type(atl_material_type), intent(inout) :: material

Material description for the faces on the current level

real(kind=rk), intent(in), optional :: nodal_data(:,:)

The state in nodal form

real(kind=rk), intent(in), optional :: nodal_GradData(:,:,:)
real(kind=rk), intent(inout) :: nodal_res(:,:)
real(kind=rk), intent(in) :: ElemLength

Length of the element

type(atl_scheme_type), intent(inout) :: scheme_min

The scheme information

type(atl_scheme_type), intent(inout) :: scheme_current

public subroutine atl_modg_maxwell_2d_penalization_NonConst(equation, poly_proj, nodal_data, scheme_min, penalizationData, iElem, material)

Read more…

Arguments

TypeIntentOptionalAttributesName
type(atl_Equations_type), intent(in) :: equation

The equation you solve.

type(ply_poly_project_type), intent(inout) :: poly_proj

Poly project

real(kind=rk), intent(in), optional :: nodal_data(:,:)

The state in nodal form

type(atl_scheme_type), intent(inout) :: scheme_min

The scheme information

type(atl_penalizationData_type), intent(inout) :: penalizationData

The penalization data

integer, intent(in) :: iElem

The current Element

type(atl_material_type), intent(inout) :: material

Material description for the faces on the current level

public subroutine atl_modg_maxwell_2d_penalization_Const(equation, poly_proj, nodal_data, scheme_min, penalizationData, iElem, material)

Arguments

TypeIntentOptionalAttributesName
type(atl_Equations_type), intent(in) :: equation

The equation you solve.

type(ply_poly_project_type), intent(inout) :: poly_proj

Poly project

real(kind=rk), intent(in), optional :: nodal_data(:,:)

The state in nodal form

type(atl_scheme_type), intent(inout) :: scheme_min

The scheme information

type(atl_penalizationData_type), intent(inout) :: penalizationData

The penalization data

integer, intent(in) :: iElem

The current Element

type(atl_material_type), intent(inout) :: material

Material description for the faces on the current level

private subroutine compute_physFlux_2d(nDofs, nScalars, state_der, state, rot, inv_mu, inv_epsi)

Compute the physical flux in x direction.

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: nDofs

dimensions

integer, intent(in) :: nScalars

dimensions

real(kind=rk), intent(inout) :: state_der(:,:)

Array to store the fluxes in.

real(kind=rk), intent(in) :: state(nDofs,nScalars)

State to compute the fluxes from.

integer, intent(in) :: rot(7)

Rotationing to index the variables.

real(kind=rk), intent(in) :: inv_mu
real(kind=rk), intent(in) :: inv_epsi

private subroutine compute_physFlux_nonConst_2d(nDofs, nScalars, iElem, state_der, state, rot, nElems, elems, material, poly_proj, modalCoeffs, nodalPhysFlux)

Compute the physical flux in x direction.

Read more…

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: nDofs

dimensions

integer, intent(in) :: nScalars

dimensions

integer, intent(in) :: iElem
real(kind=rk), intent(inout) :: state_der(:,:)

Array to store the fluxes in.

real(kind=rk), intent(in) :: state(nDofs,nScalars)

State to compute the fluxes from.

integer, intent(in) :: rot(7)

Rotationing to index the variables.

integer, intent(in) :: nElems

Number of elements.

integer, intent(in) :: elems(nElems)

Element positions in the total state vector

real(kind=rk), intent(in) :: material(nElems,nDofs,3)

Material parameters (mu, epsilon) for all elements

type(ply_poly_project_type) :: poly_proj

Data for projection method

real(kind=rk), intent(inout) :: modalCoeffs(poly_proj%body_2D%oversamp_dofs,size(state,2),1)

Working array for modal coefficients of the current element in the loop.

real(kind=rk), intent(inout) :: nodalPhysFlux(poly_proj%body_2D%nquadpoints,size(state,2),2)

Working array for nodal representation of the physical flux along the 3 spatial directions.