atl_ssprk2_module Module

Routines,functions, datatypes related to second order, strong stability preserving explicit Runge Kutta method.

This module implements all routines and datatypes that are necessary for a second order, strong stability preserving Runge-Kutta method. The implementation is based on section 4.1 of: Strong stability-preserving High-Order Time-Discretization Methods, Gottlieb S. , Shu Chi-Wang, Tadmor Etian , SIAM Review, Vol. 43, No. 1, pp. 89 - 112


Uses

Used by

  • module~~atl_ssprk2_module~~UsedByGraph module~atl_ssprk2_module atl_ssprk2_module module~atl_global_time_integration_module atl_global_time_integration_module module~atl_global_time_integration_module->module~atl_ssprk2_module module~atl_container_module atl_container_module module~atl_container_module->module~atl_global_time_integration_module module~atl_program_module atl_program_module module~atl_program_module->module~atl_global_time_integration_module module~atl_program_module->module~atl_container_module module~atl_initialize_module atl_initialize_module module~atl_program_module->module~atl_initialize_module program~ateles ateles program~ateles->module~atl_container_module program~ateles->module~atl_program_module program~atl_harvesting atl_harvesting program~atl_harvesting->module~atl_container_module program~atl_harvesting->module~atl_program_module program~atl_harvesting->module~atl_initialize_module module~atl_initialize_module->module~atl_container_module

Contents


Subroutines

public subroutine atl_init_explicit_ssprk(me, minLevel, maxLevel, steps, statedata_list)

Routine to initialize explicit runge kutta scheme for timestepping.

Arguments

TypeIntentOptionalAttributesName
type(atl_global_timestep_type), intent(inout) :: me

The datatype to initialize.

integer, intent(in) :: minLevel

The minimum of level of the mesh.

integer, intent(in) :: maxLevel

The maximum of level of the mesh.

integer, intent(in) :: steps

The number of steps in the runge kutta procedure

type(atl_statedata_type), intent(in) :: statedata_list(minLevel:maxLevel)

The state list used in your solver, for each level one entry.

private subroutine update_timestep_ssprk2(me, timestepInfo)

Levelwise updating of runge kutta of order 2

Arguments

TypeIntentOptionalAttributesName
class(atl_timestep_type), intent(inout) :: me

The type of your timestepping.

type(atl_local_timestep_type), intent(in) :: timestepInfo

Local timestepping information for that part of the mesh

private subroutine elemental_timestep_ssprk2(me, state, cell, dof, sideFlux)

Elemental operation for timestepping of order 2.

Arguments

TypeIntentOptionalAttributesName
class(atl_timestep_type), intent(inout) :: me

Description of the timestep integration method.

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

The state of all cells on this level. This field will be updated ad the cell position. See kerneldata type for more explanations.

integer, intent(in) :: cell

Position of the cell to update in the state vector.

integer, intent(in) :: dof

The degree of freedom to update

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

The flux for one of the sides of this cell. The length of this array is the number of conservative variables of your equation.

private subroutine elemental_timestep_vec_ssprk2(me, state, kerneldata)

Elemental operation for timestepping of order 2.

Arguments

TypeIntentOptionalAttributesName
class(atl_timestep_type), intent(inout) :: me

Description of the timestep integration method.

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

The state of all cells on this level. This field will be updated ad the cell position. See kerneldata type for more explanations.

type(atl_kerneldata_type), intent(in) :: kerneldata

Complete kerneldata to get the flux from with additional information.

private subroutine compute_vec(nTotal, nDofs, nScalars, state, state_der)

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: nTotal
integer, intent(in) :: nDofs
integer, intent(in) :: nScalars
real(kind=rk), intent(inout) :: state(nTotal,nDofs,nScalars)
real(kind=rk), intent(in) :: state_der(:,:,:)

private subroutine mesh_timestep_ssprk2(minLevel, maxLevel, currentLevel, cubes, tree, timestep_list, nSteps, equation, general, commStateTimer, poly_proj_list)

Subroutine for timestepping with explicit runge kutta of order 2.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: minLevel

The minimum refinement level of the mesh.

integer, intent(in) :: maxLevel

The maximum refinement level of the mesh.

integer, intent(in) :: currentLevel

The level the timestep has to be performed for.

type(atl_cube_container_type), intent(inout) :: cubes

Container for the cubical elements.

type(treelmesh_type), intent(in) :: tree

treelm mesh

type(atl_timestep_type), intent(inout) :: timestep_list(minLevel:)

List of levelwise timestepping algorihtms

integer, intent(in) :: nSteps

The number of steps of the time stepping scheme (assumed to be 2)

type(atl_Equations_type), intent(inout) :: equation

The equation you are operating with.

type(tem_general_type), intent(inout) :: general

General treelm settings

integer, intent(inout) :: commStateTimer

Timer for measuring the communication time inside this routine.

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

unique list for projection methods

private subroutine compute_intermediate(nTotal, nDofs, nScalars, state_tmp, state1, state2, factor)

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: nTotal
integer, intent(in) :: nDofs
integer, intent(in) :: nScalars
real(kind=rk), intent(out) :: state_tmp(nTotal,nDofs,nScalars)
real(kind=rk), intent(in) :: state1(nTotal,nDofs,nScalars)
real(kind=rk), intent(in) :: state2(nTotal,nDofs,nScalars)
real(kind=rk), intent(in) :: factor

private subroutine rk2_update(statedata_list, scheme_list, timestep_list)

Subroutine calculates the final update step of the Runge-Kutta method. It is performing levelwise.

Arguments

TypeIntentOptionalAttributesName
type(atl_statedata_type), intent(inout) :: statedata_list

List of states you want to calc the rhs for. For each level we have one.

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

List of schemes, for each level.

type(atl_timestep_type), intent(inout) :: timestep_list

List of levelwise timestepping algorihtms

private subroutine compute_up(nTotal, nDofs, nScalars, state, state1, state2, dt)

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: nTotal
integer, intent(in) :: nDofs
integer, intent(in) :: nScalars
real(kind=rk), intent(inout) :: state(nTotal,nDofs,nScalars)
real(kind=rk), intent(in) :: state1(nTotal,nDofs,nScalars)
real(kind=rk), intent(in) :: state2(nTotal,nDofs,nScalars)
real(kind=rk), intent(in) :: dt

private recursive subroutine rk2_substep(minLevel, maxLevel, currentLevel, mesh_list, tree, levelPointer, kerneldata_list, statedata_list, facedata_list, source, penalizationdata_list, boundary_list, bc, scheme_list, poly_proj_pos, poly_proj_list, timestep_list, equation, material_list, general, commStateTimer)

Subroutine calculates a substep of the Runge-Kutta timestepping scheme. Calls itself recursively for the finer levels until the finest level is reached.

Arguments

TypeIntentOptionalAttributesName
integer, intent(in) :: minLevel

The minimum refinement level of the mesh.

integer, intent(in) :: maxLevel

The maximum refinement level of the mesh.

integer, intent(in) :: currentLevel

The level the timestep has to be performed for.

type(atl_cube_elem_type), intent(inout) :: mesh_list(minLevel:maxLevel)

List of mesh parts. For each level we have one.

type(treelmesh_type), intent(in) :: tree

treelm mesh

integer, intent(in) :: levelPointer(:)

Pointer for elements from global treeID list index to index in levelwise fluid lists

type(atl_kerneldata_type), intent(inout) :: kerneldata_list(minLevel:maxLevel)

List of kerneldatas. For each level we have one

type(atl_statedata_type), intent(inout) :: statedata_list(minLevel:maxLevel)

List of states you want to calc the rhs for. For each level we have one.

type(atl_facedata_type), intent(inout) :: facedata_list(minLevel:maxLevel)

List of faces states you want to calc the rhs for. For each level we have one.

type(atl_source_type), intent(inout) :: source

List of sources, for each level

type(atl_penalizationData_type), intent(inout) :: penalizationdata_list(minLevel:maxLevel)

List of penalization data, for each level.

type(atl_level_boundary_type), intent(inout) :: boundary_list(minLevel:maxLevel)

List of boundaries, for each level.

type(atl_boundary_type), intent(in) :: bc(:)

Global description of the boundaries

type(atl_scheme_type), intent(inout) :: scheme_list(minLevel:maxLevel)

List of schemes, for each level.

integer, intent(in) :: poly_proj_pos(minLevel:maxLevel)

List of position pointer of projection, for each level.

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

List of projection methods.

type(atl_timestep_type), intent(inout) :: timestep_list(minLevel:maxLevel)

List of levelwise timestepping algorihtms

type(atl_Equations_type), intent(inout) :: equation

The equation you are operating with.

type(atl_material_type), intent(inout) :: material_list(minlevel:maxlevel)

Material parameter description.

type(tem_general_type), intent(inout) :: general

General treelm settings

integer, intent(inout) :: commStateTimer

Timer for measuring the communication time inside this routine.