version 0.6.0
mod_integrate_cell_advection_term_explicit_generic Module Reference

Integration methods for the advection equation. More...

Data Types

interface  divut_term_computer_interface
 The abstract function for computing div(U phi) More...
 
interface  divut_term_computer_interface_fd
 The abstract function for computing div(U phi) for the FD WENO schemes. More...
 

Functions/Subroutines

subroutine integrate_cell_advection_term_explicit_generic_euler (divut, time_step, dt_nm1, dt_n, velocity_nm1, velocity_n, velocity_np1, scalar, divut_term_computer, flux_type, cfl_coefficient, equation_has_immersed_boundaries, equation_isd_target, ibc_variable, boundary_condition)
 Adds the advection term \( c \nabla \cdot ( \mathbf{u} \phi ) \) with a generic spatial scheme euler (in time) version. More...
 
subroutine integrate_cell_advection_term_explicit_generic_euler_split (divut, time_step, dt_nm1, dt_n, velocity_nm1, velocity_n, velocity_np1, scalar, divut_term_computer, flux_type, cfl_coefficient, equation_has_immersed_boundaries, equation_isd_target, ibc_variable, boundary_condition)
 Adds the advection term \( c \nabla \cdot ( \mathbf{u} \phi ) \) with a generic spatial scheme euler (in time) version. More...
 
subroutine integrate_cell_advection_term_explicit_generic_sspso (divut, time_step, dt_nm1, dt_n, velocity_nm1, velocity_n, velocity_np1, scalar, divut_term_computer, flux_type, cfl_coefficient, aki, bi, ci, boundary_condition)
 Adds the advection term \( c \nabla \cdot ( \mathbf{u} \phi ) \) with a generic spatial scheme SSP 2 stage 2nd order without velocity extrapolation. More...
 
subroutine integrate_cell_advection_term_explicit_generic_sspso_split (divut, time_step, dt_nm1, dt_n, velocity_nm1, velocity_n, velocity_np1, scalar, divut_term_computer, flux_type, cfl_coefficient, aki, bi, ci, boundary_condition)
 Adds the advection term \( c \nabla \cdot ( \mathbf{u} \phi ) \) with a generic spatial scheme in a splitted manner SSP 2 stage 2nd order without velocity extrapolation. More...
 
subroutine integrate_cell_advection_term_explicit_generic_nssp32 (divut, time_step, dt_nm1, dt_n, velocity_nm1, velocity_n, velocity_np1, scalar, divut_term_computer, flux_type, cfl_coefficient, boundary_condition)
 Adds the advection term \( c \nabla \cdot ( \mathbf{u} \phi ) \) with a generic spatial scheme. More...
 
subroutine integrate_cell_advection_term_explicit_generic_nssp32_split (divut, time_step, dt_nm1, dt_n, velocity_nm1, velocity_n, velocity_np1, scalar, divut_term_computer, flux_type, cfl_coefficient, boundary_condition)
 Integrate the advection term \( c \nabla \cdot ( \mathbf{u} \phi ) \) with a generic spatial scheme. The temporal scheme is a 3 steps NSSP at 2nd order (NSSP32). The integration is done in a split manner (x then y then z)
 
subroutine integrate_cell_advection_term_explicit_generic_nssp53 (divut, time_step, dt_nm1, dt_n, velocity_nm1, velocity_n, velocity_np1, scalar, divut_term_computer, flux_type, cfl_coefficient, boundary_condition)
 Adds the advection term \( c \nabla \cdot ( \mathbf{u} \phi ) \) with a generic spatial scheme. More...
 
subroutine integrate_cell_advection_term_explicit_generic_nssp53_split (divut, time_step, dt_nm1, dt_n, velocity_nm1, velocity_n, velocity_np1, scalar, divut_term_computer, flux_type, cfl_coefficient, boundary_condition)
 Integrate the advection term \( c \nabla \cdot ( \mathbf{u} \phi ) \) with a generic spatial scheme. The temporal scheme is a 5 steps NSSP at 3nd order (NSSP53). The integration is done in a split manner (x then y then z)
 

Detailed Description

This module furnishes various temporal integration schemes (mainly RK schemes) for the advection equation. As the explicit resolution of the advection equation is subjected to a stability condition (the CFL), sub iterations can occur during the integration process.

Note
A generic abstract function divut_term_computer_interface is defined for simplifying the process. With this manner, any spatial scheme can be used, may it be first order or WENO.

Function/Subroutine Documentation

◆ integrate_cell_advection_term_explicit_generic_euler()

subroutine mod_integrate_cell_advection_term_explicit_generic::integrate_cell_advection_term_explicit_generic_euler ( double precision, dimension(:,:,:), intent(inout)  divut,
double precision, intent(in)  time_step,
double precision, intent(in)  dt_nm1,
double precision, intent(in)  dt_n,
type(t_face_field), intent(in)  velocity_nm1,
type(t_face_field), intent(in)  velocity_n,
type(t_face_field), intent(in)  velocity_np1,
double precision, dimension(:,:,:), intent(in)  scalar,
procedure(divut_term_computer_interface divut_term_computer,
type(t_fv_flux), intent(in)  flux_type,
double precision, intent(in)  cfl_coefficient,
logical, intent(in)  equation_has_immersed_boundaries,
integer, dimension(:), intent(in)  equation_isd_target,
type(t_immersed_boundary_condition), dimension(:), intent(in)  ibc_variable,
type(t_boundary_condition), intent(in)  boundary_condition 
)
Parameters
[in,out]divutthe resulting term
[in]time_stepthe time step on which we integrate (might be different from dt_n for sub integration steps!)
[in]velocity_nm1velocity field at \( t^{n-1} \)
[in]velocity_nvelocity field at \( t^n \)
[in]velocity_np1velocity field at \( t^{n+1} \)
[in]dt_nm1the previous time step \( \delta t^{n-1} = t^{n} - t^{n-1} \)
[in]dt_nthe current time step \( \delta t^{n} = t^{n+1} - t^{n} \)
[in]flux_typethe flux type to use for the explicit advection (type_fv_flux)
[in]boundary_conditionthe boundary conditions of the cell scalar
[in]scalarthe scalar to advect
[in]divut_term_computerthe procedure that will compute the div(u phi) term giving u and phi
[in]divut_term_computer_fdthe procedure that will compute the div(u phi) term giving u and phi for FD WENO schemes
[in]cfl_coefficientthe cfl coefficient to use (depending on the time and space schemes)
[in]equation_isd_targetSubset of immersed boundaries.
[in]ibc_variableboundary conditions on immersed boundaries.
[in]boundary_condition

◆ integrate_cell_advection_term_explicit_generic_euler_split()

subroutine mod_integrate_cell_advection_term_explicit_generic::integrate_cell_advection_term_explicit_generic_euler_split ( double precision, dimension(:,:,:), intent(inout)  divut,
double precision, intent(in)  time_step,
double precision, intent(in)  dt_nm1,
double precision, intent(in)  dt_n,
type(t_face_field), intent(in)  velocity_nm1,
type(t_face_field), intent(in)  velocity_n,
type(t_face_field), intent(in)  velocity_np1,
double precision, dimension(:,:,:), intent(in)  scalar,
procedure(divut_term_computer_interface divut_term_computer,
type(t_fv_flux), intent(in)  flux_type,
double precision, intent(in)  cfl_coefficient,
logical, intent(in)  equation_has_immersed_boundaries,
integer, dimension(:), intent(in)  equation_isd_target,
type(t_immersed_boundary_condition), dimension(:), intent(in)  ibc_variable,
type(t_boundary_condition), intent(in)  boundary_condition 
)
Parameters
[in,out]divutthe resulting term
[in]time_stepthe time step on which we integrate (might be different from dt_n for sub integration steps!)
[in]velocity_nm1velocity field at \( t^{n-1} \)
[in]velocity_nvelocity field at \( t^n \)
[in]velocity_np1velocity field at \( t^{n+1} \)
[in]dt_nm1the previous time step \( \delta t^{n-1} = t^{n} - t^{n-1} \)
[in]dt_nthe current time step \( \delta t^{n} = t^{n+1} - t^{n} \)
[in]flux_typethe flux type to use for the explicit advection (type_fv_flux)
[in]boundary_conditionthe boundary conditions of the cell scalar
[in]scalarthe scalar to advect
[in]divut_term_computerthe procedure that will compute the div(u phi) term giving u and phi
[in]divut_term_computer_fdthe procedure that will compute the div(u phi) term giving u and phi for FD WENO schemes
[in]cfl_coefficientthe cfl coefficient to use (depending on the time and space schemes)

◆ integrate_cell_advection_term_explicit_generic_nssp32()

subroutine mod_integrate_cell_advection_term_explicit_generic::integrate_cell_advection_term_explicit_generic_nssp32 ( double precision, dimension(:,:,:), intent(inout)  divut,
double precision, intent(in)  time_step,
double precision, intent(in)  dt_nm1,
double precision, intent(in)  dt_n,
type(t_face_field), intent(in)  velocity_nm1,
type(t_face_field), intent(in)  velocity_n,
type(t_face_field), intent(in)  velocity_np1,
double precision, dimension(:,:,:), intent(in)  scalar,
procedure(divut_term_computer_interface divut_term_computer,
type(t_fv_flux), intent(in)  flux_type,
double precision, intent(in)  cfl_coefficient,
type(t_boundary_condition), intent(in)  boundary_condition 
)

Optimized version for 3 stages at order 2

Parameters
[in,out]divutthe resulting field
[in]time_stepthe time step on which we integrate (might be different from dt_n for sub integration steps!)
[in]velocity_nm1velocity field at \( t^{n-1} \)
[in]velocity_nvelocity field at \( t^n \)
[in]velocity_np1velocity field at \( t^{n+1} \)
[in]dt_nm1the previous time step \( \delta t^{n-1} = t^{n} - t^{n-1} \)
[in]dt_nthe current time step \( \delta t^{n} = t^{n+1} - t^{n} \)
[in]scalarthe scalar to advect
[in]divut_term_computerthe procedure that will compute the div(u phi) term giving u and phi
[in]flux_typethe flux type to use for the explicit advection (type_fv_flux)
[in]cfl_coefficientthe cfl coefficient to use (depending on the time and space schemes)
[in]boundary_conditionthe boundary conditions of the cell scalar

◆ integrate_cell_advection_term_explicit_generic_nssp53()

subroutine mod_integrate_cell_advection_term_explicit_generic::integrate_cell_advection_term_explicit_generic_nssp53 ( double precision, dimension(:,:,:), intent(inout)  divut,
double precision, intent(in)  time_step,
double precision, intent(in)  dt_nm1,
double precision, intent(in)  dt_n,
type(t_face_field), intent(in)  velocity_nm1,
type(t_face_field), intent(in)  velocity_n,
type(t_face_field), intent(in)  velocity_np1,
double precision, dimension(:,:,:), intent(in)  scalar,
procedure(divut_term_computer_interface divut_term_computer,
type(t_fv_flux), intent(in)  flux_type,
double precision, intent(in)  cfl_coefficient,
type(t_boundary_condition), intent(in)  boundary_condition 
)

Optimized version for 5 stages at order 3

Parameters
[in,out]divutthe resulting field
[in]time_stepthe time step on which we integrate (might be different from dt_n for sub integration steps!)
[in]velocity_nm1velocity field at \( t^{n-1} \)
[in]velocity_nvelocity field at \( t^n \)
[in]velocity_np1velocity field at \( t^{n+1} \)
[in]dt_nm1the previous time step \( \delta t^{n-1} = t^{n} - t^{n-1} \)
[in]dt_nthe current time step \( \delta t^{n} = t^{n+1} - t^{n} \)
[in]scalarthe scalar to advect
[in]divut_term_computerthe procedure that will compute the div(u phi) term giving u and phi
[in]divut_term_computer_fdthe procedure that will compute the div(u phi) term giving u and phi for FD WENO schemes
[in]flux_typethe flux type to use for the explicit advection (type_fv_flux)
[in]cfl_coefficientthe cfl coefficient to use (depending on the time and space schemes)
[in]boundary_conditionthe boundary conditions of the cell scalar

◆ integrate_cell_advection_term_explicit_generic_sspso()

subroutine mod_integrate_cell_advection_term_explicit_generic::integrate_cell_advection_term_explicit_generic_sspso ( double precision, dimension(:,:,:), intent(inout)  divut,
double precision, intent(in)  time_step,
double precision, intent(in)  dt_nm1,
double precision, intent(in)  dt_n,
type(t_face_field), intent(in)  velocity_nm1,
type(t_face_field), intent(in)  velocity_n,
type(t_face_field), intent(in)  velocity_np1,
double precision, dimension(:,:,:), intent(in)  scalar,
procedure(divut_term_computer_interface divut_term_computer,
type(t_fv_flux), intent(in)  flux_type,
double precision, intent(in)  cfl_coefficient,
double precision, dimension(:,:), intent(in)  aki,
double precision, dimension(:), intent(in)  bi,
double precision, dimension(:), intent(in)  ci,
type(t_boundary_condition), intent(in)  boundary_condition 
)
Parameters
[in,out]divutthe resulting field
[in]time_stepthe time step on which we integrate (might be different from dt_n for sub integration steps!)
[in]velocity_nm1velocity field at \( t^{n-1} \)
[in]velocity_nvelocity field at \( t^n \)
[in]velocity_np1velocity field at \( t^{n+1} \)
[in]dt_nm1the previous time step \( \delta t^{n-1} = t^{n} - t^{n-1} \)
[in]dt_nthe current time step \( \delta t^{n} = t^{n+1} - t^{n} \)
[in]scalarthe scalar to advect
[in]divut_term_computerthe procedure that will compute the div(u phi) term giving u and phi
[in]divut_term_computer_fdthe procedure that will compute the div(u phi) term giving u and phi for FD WENO schemes
[in]flux_typethe flux type to use for the explicit advection (type_fv_flux)
[in]cfl_coefficientthe cfl coefficient to use (depending on the time and space schemes)
[in]boundary_conditionthe boundary conditions of the cell scalar
[in]aki,bi,cithe given coefficients of the Butcher tableau

◆ integrate_cell_advection_term_explicit_generic_sspso_split()

subroutine mod_integrate_cell_advection_term_explicit_generic::integrate_cell_advection_term_explicit_generic_sspso_split ( double precision, dimension(:,:,:), intent(inout)  divut,
double precision, intent(in)  time_step,
double precision, intent(in)  dt_nm1,
double precision, intent(in)  dt_n,
type(t_face_field), intent(in)  velocity_nm1,
type(t_face_field), intent(in)  velocity_n,
type(t_face_field), intent(in)  velocity_np1,
double precision, dimension(:,:,:), intent(in)  scalar,
procedure(divut_term_computer_interface divut_term_computer,
type(t_fv_flux), intent(in)  flux_type,
double precision, intent(in)  cfl_coefficient,
double precision, dimension(:,:), intent(in)  aki,
double precision, dimension(:), intent(in)  bi,
double precision, dimension(:), intent(in)  ci,
type(t_boundary_condition), intent(in)  boundary_condition 
)
Parameters
[in,out]divutthe resulting field
[in]time_stepthe time step on which we integrate (might be different from dt_n for sub integration steps!)
[in]velocity_nm1velocity field at \( t^{n-1} \)
[in]velocity_nvelocity field at \( t^n \)
[in]velocity_np1velocity field at \( t^{n+1} \)
[in]dt_nm1the previous time step \( \delta t^{n-1} = t^{n} - t^{n-1} \)
[in]dt_nthe current time step \( \delta t^{n} = t^{n+1} - t^{n} \)
[in]scalarthe scalar to advect
[in]divut_term_computerthe procedure that will compute the div(u phi) term giving u and phi
[in]flux_typethe flux type to use for the explicit advection (type_fv_flux)
[in]cfl_coefficientthe cfl coefficient to use (depending on the time and space schemes)
[in]boundary_conditionthe boundary conditions of the cell scalar
[in]aki,bi,cithe given coefficients of the Butcher tableau