First level routines to discretize PDE term by term on faces. More...
Topics | |
Stress tensor divergence term | |
Namespaces | |
module | mod_momentum_boundary_conditions |
Boundary conditions related routines for momentum. | |
module | enum_navier_nonlinear_term_scheme |
Discretization schemes for the nonlinear momentum term. | |
module | mod_face_vector_implicit_discretization_collection |
Module that encapsulates all face vector term by term discretization modules. | |
Functions | |
subroutine | mod_discretize_face_temporal_term::discretize_face_temporal_term (matrix, rhs, density_face_n, density_face_star, velocity_n, velocity_nm1, equation_time_step, equation_time_step_n, equation_time_order_discretization, is_momentum_preserving_method, equation_stencil, equation_ls_map) |
Discretize the advection term. | |
subroutine | mod_discretize_nonlinear_momentum_term::discretize_nonlinear_momentum_term (matrix, navier_nonlinear_term_discretization_type, navier_nonlinear_term_scheme, navier_has_div_u_advection_term, is_navier_momentum_preserving_method, navier_has_immersed_boundaries, navier_ib_has_one_sided_inner_discretization, navier_isd_target, navier_ib_inner_discretization_order, density_face_n, density_face_star, divergence_extrapolated, time_step_nm1, time_step, cfl_navier_advection, velocity_extrapolated, velocity_nm1, velocity_n, explicit_time_order_discretization, explicit_specify_temporal_stability_factor, explicit_temporal_stability_factor, flux_type, explicit_nonlinear_term_n, navier_stencil, navier_ls_map, boundary_condition) |
Add nonlinear momentum term to the given linear system. | |
First level routines to discretize PDE term by term on faces.
In this directory you will find the first level of routines to solve PDE term by term on face based variables.
subroutine mod_discretize_face_temporal_term::discretize_face_temporal_term | ( | double precision, dimension(:), intent(inout) | matrix, |
double precision, dimension(:), intent(inout) | rhs, | ||
type(t_face_field), intent(in) | density_face_n, | ||
type(t_face_field), intent(in) | density_face_star, | ||
type(t_face_field), intent(in) | velocity_n, | ||
type(t_face_field), intent(in) | velocity_nm1, | ||
double precision, intent(in) | equation_time_step, | ||
double precision, intent(in) | equation_time_step_n, | ||
integer, intent(in) | equation_time_order_discretization, | ||
logical, intent(in) | is_momentum_preserving_method, | ||
type(t_face_stencil), intent(in) | equation_stencil, | ||
type(t_face_ls_map), intent(in) | equation_ls_map ) |
Discretize the advection term.
This term is treated implicitly thanks to backward first order Euler scheme or a second order Backward Differentiation Formula.
The diagonal of the matrix and the right-hand-side are modified accordingly to the scheme chosen.
[in,out] | matrix,rhs | linear system |
[in] | density_face_n | |
[in] | density_face_star | |
[in] | velocity_n | |
[in] | velocity_nm1 | |
[in] | equation_time_step | |
[in] | equation_time_step_n | |
[in] | equation_time_order_discretization | |
[in] | is_momentum_preserving_method | |
[in] | equation_stencil | |
[in] | equation_ls_map |
subroutine mod_discretize_nonlinear_momentum_term::discretize_nonlinear_momentum_term | ( | double precision, dimension(:), intent(inout) | matrix, |
integer, intent(in) | navier_nonlinear_term_discretization_type, | ||
integer, intent(in) | navier_nonlinear_term_scheme, | ||
logical, intent(in) | navier_has_div_u_advection_term, | ||
logical, intent(in) | is_navier_momentum_preserving_method, | ||
logical, intent(in) | navier_has_immersed_boundaries, | ||
logical, intent(in) | navier_ib_has_one_sided_inner_discretization, | ||
integer, dimension(:), intent(in) | navier_isd_target, | ||
integer, dimension(:), intent(in) | navier_ib_inner_discretization_order, | ||
type(t_face_field), intent(in) | density_face_n, | ||
type(t_face_field), intent(in) | density_face_star, | ||
double precision, dimension(:,:,:), intent(in) | divergence_extrapolated, | ||
double precision, intent(in) | time_step_nm1, | ||
double precision, intent(in) | time_step, | ||
double precision, intent(out) | cfl_navier_advection, | ||
type(t_face_field), intent(inout) | velocity_extrapolated, | ||
type(t_face_field), intent(inout) | velocity_nm1, | ||
type(t_face_field), intent(inout) | velocity_n, | ||
integer, intent(in) | explicit_time_order_discretization, | ||
logical, intent(in) | explicit_specify_temporal_stability_factor, | ||
double precision, intent(in) | explicit_temporal_stability_factor, | ||
type(t_fv_flux), intent(in) | flux_type, | ||
type(t_face_field), intent(inout) | explicit_nonlinear_term_n, | ||
type(t_face_stencil), intent(in) | navier_stencil, | ||
type(t_face_ls_map), intent(in) | navier_ls_map, | ||
type(t_boundary_condition_face), intent(in) | boundary_condition ) |
Add nonlinear momentum term to the given linear system.
The nonlinear momentum term is written as follows:
\[ \rho \mathbf{u} \mathbf{\nabla} \mathbf{u} = \rho \nabla \cdot (\mathbf{u} \otimes \mathbf{u}) - \rho \mathbf{u} \nabla \cdot(\mathbf{u}), \]
where each term is discretized separately. The term \(\nabla \cdot (\mathbf{u} \otimes \mathbf{u})\) is a vector that can be decomposed as:
\[ \nabla \cdot (\mathbf{u} \otimes \mathbf{u}) = \left[ \begin{matrix} \nabla \cdot (u \times \mathbf{u}) \\ \nabla \cdot (v \times \mathbf{u}) \\ \nabla \cdot (w \times \mathbf{u}) \end{matrix} \right] \]