MFC
Exascale flow solver
Loading...
Searching...
No Matches
m_riemann_solvers Module Reference

Approximate and exact Riemann solvers (HLL, HLLC, HLLD, exact) for the multicomponent Navier–Stokes equations. More...

Functions/Subroutines

subroutine, public s_riemann_solver (ql_prim_rsx_vf, ql_prim_rsy_vf, ql_prim_rsz_vf, dql_prim_dx_vf, dql_prim_dy_vf, dql_prim_dz_vf, ql_prim_vf, qr_prim_rsx_vf, qr_prim_rsy_vf, qr_prim_rsz_vf, dqr_prim_dx_vf, dqr_prim_dy_vf, dqr_prim_dz_vf, qr_prim_vf, q_prim_vf, flux_vf, flux_src_vf, flux_gsrc_vf, norm_dir, ix, iy, iz)
 Dispatch to the subroutines that are utilized to compute the Riemann problem solution. For additional information please reference: 1) s_hll_riemann_solver 2) s_hllc_riemann_solver 3) s_lf_riemann_solver 4) s_hlld_riemann_solver.
subroutine s_compute_viscous_source_flux (vell_vf, dvell_dx_vf, dvell_dy_vf, dvell_dz_vf, velr_vf, dvelr_dx_vf, dvelr_dy_vf, dvelr_dz_vf, flux_src_vf, norm_dir, ix, iy, iz)
 Dispatch to the subroutines that are utilized to compute the viscous source fluxes for either Cartesian or cylindrical geometries. For more information please refer to: 1) s_compute_cartesian_viscous_source_flux 2) s_compute_cylindrical_viscous_source_flux.
subroutine, public s_hll_riemann_solver (ql_prim_rsx_vf, ql_prim_rsy_vf, ql_prim_rsz_vf, dql_prim_dx_vf, dql_prim_dy_vf, dql_prim_dz_vf, ql_prim_vf, qr_prim_rsx_vf, qr_prim_rsy_vf, qr_prim_rsz_vf, dqr_prim_dx_vf, dqr_prim_dy_vf, dqr_prim_dz_vf, qr_prim_vf, q_prim_vf, flux_vf, flux_src_vf, flux_gsrc_vf, norm_dir, ix, iy, iz)
 HLL approximate Riemann solver, Harten et al. SIAM Review (1983).
subroutine, public s_lf_riemann_solver (ql_prim_rsx_vf, ql_prim_rsy_vf, ql_prim_rsz_vf, dql_prim_dx_vf, dql_prim_dy_vf, dql_prim_dz_vf, ql_prim_vf, qr_prim_rsx_vf, qr_prim_rsy_vf, qr_prim_rsz_vf, dqr_prim_dx_vf, dqr_prim_dy_vf, dqr_prim_dz_vf, qr_prim_vf, q_prim_vf, flux_vf, flux_src_vf, flux_gsrc_vf, norm_dir, ix, iy, iz)
 Lax-Friedrichs (Rusanov) approximate Riemann solver.
subroutine, public s_hllc_riemann_solver (ql_prim_rsx_vf, ql_prim_rsy_vf, ql_prim_rsz_vf, dql_prim_dx_vf, dql_prim_dy_vf, dql_prim_dz_vf, ql_prim_vf, qr_prim_rsx_vf, qr_prim_rsy_vf, qr_prim_rsz_vf, dqr_prim_dx_vf, dqr_prim_dy_vf, dqr_prim_dz_vf, qr_prim_vf, q_prim_vf, flux_vf, flux_src_vf, flux_gsrc_vf, norm_dir, ix, iy, iz)
 HLLC Riemann solver with contact restoration, Toro et al. Shock Waves (1994).
subroutine, public s_hlld_riemann_solver (ql_prim_rsx_vf, ql_prim_rsy_vf, ql_prim_rsz_vf, dql_prim_dx_vf, dql_prim_dy_vf, dql_prim_dz_vf, ql_prim_vf, qr_prim_rsx_vf, qr_prim_rsy_vf, qr_prim_rsz_vf, dqr_prim_dx_vf, dqr_prim_dy_vf, dqr_prim_dz_vf, qr_prim_vf, q_prim_vf, flux_vf, flux_src_vf, flux_gsrc_vf, norm_dir, ix, iy, iz)
 HLLD Riemann solver for MHD, Miyoshi & Kusano JCP (2005).
impure subroutine, public s_initialize_riemann_solvers_module
 Initialize the Riemann solvers module.
subroutine s_populate_riemann_states_variables_buffers (ql_prim_rsx_vf, ql_prim_rsy_vf, ql_prim_rsz_vf, dql_prim_dx_vf, dql_prim_dy_vf, dql_prim_dz_vf, qr_prim_rsx_vf, qr_prim_rsy_vf, qr_prim_rsz_vf, dqr_prim_dx_vf, dqr_prim_dy_vf, dqr_prim_dz_vf, norm_dir, ix, iy, iz)
 Populate the left and right Riemann state variable buffers based on boundary conditions.
subroutine s_initialize_riemann_solver (flux_src_vf, norm_dir)
 Set up the chosen Riemann solver algorithm for the current direction.
subroutine s_compute_cylindrical_viscous_source_flux (vell_vf, dvell_dx_vf, dvell_dy_vf, dvell_dz_vf, velr_vf, dvelr_dx_vf, dvelr_dy_vf, dvelr_dz_vf, flux_src_vf, norm_dir, ix, iy, iz)
 Compute cylindrical viscous source flux contributions for momentum and energy.
subroutine s_compute_cartesian_viscous_source_flux (dvell_dx_vf, dvell_dy_vf, dvell_dz_vf, dvelr_dx_vf, dvelr_dy_vf, dvelr_dz_vf, flux_src_vf, norm_dir)
 Compute Cartesian viscous source flux contributions for momentum and energy.
subroutine s_calculate_shear_stress_tensor (vel_grad_avg, re_shear, divergence_v, tau_shear_out)
 Compute shear stress tensor components.
subroutine s_calculate_bulk_stress_tensor (re_bulk, divergence_v, tau_bulk_out)
 Compute bulk stress tensor components (diagonal only).
subroutine s_finalize_riemann_solver (flux_vf, flux_src_vf, flux_gsrc_vf, norm_dir)
 Deallocation and/or disassociation procedures that are needed to finalize the selected Riemann problem solver.
impure subroutine, public s_finalize_riemann_solvers_module
 Module deallocation and/or disassociation procedures.

Variables

real(wp), dimension(:,:,:,:), allocatable vel_src_rsx_vf
real(wp), dimension(:,:,:,:), allocatable vel_src_rsy_vf
real(wp), dimension(:,:,:,:), allocatable vel_src_rsz_vf
real(wp), dimension(:,:,:,:), allocatable mom_sp_rsx_vf
real(wp), dimension(:,:,:,:), allocatable mom_sp_rsy_vf
real(wp), dimension(:,:,:,:), allocatable mom_sp_rsz_vf
real(wp), dimension(:,:,:,:), allocatable re_avg_rsx_vf
real(wp), dimension(:,:,:,:), allocatable re_avg_rsy_vf
real(wp), dimension(:,:,:,:), allocatable re_avg_rsz_vf
real(wp), dimension(:), allocatable gs_rs
real(wp), dimension(:,:), allocatable res_gs
real(wp), dimension(:,:,:,:), allocatable flux_rsx_vf
 The cell-boundary values of the fluxes (src - source) that are computed through the chosen Riemann problem solver, and the direct evaluation of source terms, by using the left and right states given in qK_prim_rs_vf, dqK_prim_ds_vf where ds = dx, dy or dz.
real(wp), dimension(:,:,:,:), allocatable flux_src_rsx_vf
real(wp), dimension(:,:,:,:), allocatable flux_rsy_vf
real(wp), dimension(:,:,:,:), allocatable flux_src_rsy_vf
real(wp), dimension(:,:,:,:), allocatable flux_rsz_vf
real(wp), dimension(:,:,:,:), allocatable flux_src_rsz_vf
real(wp), dimension(:,:,:,:), allocatable flux_gsrc_rsx_vf
 The cell-boundary values of the geometrical source flux that are computed through the chosen Riemann problem solver by using the left and right states given in qK_prim_rs_vf. Currently 2D axisymmetric for inviscid only.
real(wp), dimension(:,:,:,:), allocatable flux_gsrc_rsy_vf
real(wp), dimension(:,:,:,:), allocatable flux_gsrc_rsz_vf
Indical bounds in the s1-, s2- and s3-directions
type(int_bounds_infois1
type(int_bounds_infois2
type(int_bounds_infois3
type(int_bounds_infoisx
type(int_bounds_infoisy
type(int_bounds_infoisz

Detailed Description

Approximate and exact Riemann solvers (HLL, HLLC, HLLD, exact) for the multicomponent Navier–Stokes equations.

Function/Subroutine Documentation

◆ s_calculate_bulk_stress_tensor()

subroutine m_riemann_solvers::s_calculate_bulk_stress_tensor ( real(wp), intent(in) re_bulk,
real(wp), intent(in) divergence_v,
real(wp), dimension(num_dims, num_dims), intent(out) tau_bulk_out )

Compute bulk stress tensor components (diagonal only).

Definition at line 18684 of file m_riemann_solvers.fpp.f90.

Here is the caller graph for this function:

◆ s_calculate_shear_stress_tensor()

subroutine m_riemann_solvers::s_calculate_shear_stress_tensor ( real(wp), dimension(num_dims, num_dims), intent(in) vel_grad_avg,
real(wp), intent(in) re_shear,
real(wp), intent(in) divergence_v,
real(wp), dimension(num_dims, num_dims), intent(out) tau_shear_out )

Compute shear stress tensor components.

Definition at line 18641 of file m_riemann_solvers.fpp.f90.

Here is the caller graph for this function:

◆ s_compute_cartesian_viscous_source_flux()

subroutine m_riemann_solvers::s_compute_cartesian_viscous_source_flux ( type(scalar_field), dimension(num_dims), intent(in) dvell_dx_vf,
type(scalar_field), dimension(num_dims), intent(in) dvell_dy_vf,
type(scalar_field), dimension(num_dims), intent(in) dvell_dz_vf,
type(scalar_field), dimension(num_dims), intent(in) dvelr_dx_vf,
type(scalar_field), dimension(num_dims), intent(in) dvelr_dy_vf,
type(scalar_field), dimension(num_dims), intent(in) dvelr_dz_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_src_vf,
integer, intent(in) norm_dir )

Compute Cartesian viscous source flux contributions for momentum and energy.

Definition at line 18492 of file m_riemann_solvers.fpp.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ s_compute_cylindrical_viscous_source_flux()

subroutine m_riemann_solvers::s_compute_cylindrical_viscous_source_flux ( type(scalar_field), dimension(num_dims), intent(in) vell_vf,
type(scalar_field), dimension(num_dims), intent(in) dvell_dx_vf,
type(scalar_field), dimension(num_dims), intent(in) dvell_dy_vf,
type(scalar_field), dimension(num_dims), intent(in) dvell_dz_vf,
type(scalar_field), dimension(num_dims), intent(in) velr_vf,
type(scalar_field), dimension(num_dims), intent(in) dvelr_dx_vf,
type(scalar_field), dimension(num_dims), intent(in) dvelr_dy_vf,
type(scalar_field), dimension(num_dims), intent(in) dvelr_dz_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_src_vf,
integer, intent(in) norm_dir,
type(int_bounds_info), intent(in) ix,
type(int_bounds_info), intent(in) iy,
type(int_bounds_info), intent(in) iz )

Compute cylindrical viscous source flux contributions for momentum and energy.

Interface velocity ( \(v_1,v_2,v_3\)) (grid directions) for viscous work.

Shear stress vector ( \(\sigma_{N1}, \sigma_{N2}, \sigma_{N3}\)) on N-face (grid directions).

Definition at line 18280 of file m_riemann_solvers.fpp.f90.

Here is the caller graph for this function:

◆ s_compute_viscous_source_flux()

subroutine m_riemann_solvers::s_compute_viscous_source_flux ( type(scalar_field), dimension(num_vels), intent(in) vell_vf,
type(scalar_field), dimension(num_vels), intent(in) dvell_dx_vf,
type(scalar_field), dimension(num_vels), intent(in) dvell_dy_vf,
type(scalar_field), dimension(num_vels), intent(in) dvell_dz_vf,
type(scalar_field), dimension(num_vels), intent(in) velr_vf,
type(scalar_field), dimension(num_vels), intent(in) dvelr_dx_vf,
type(scalar_field), dimension(num_vels), intent(in) dvelr_dy_vf,
type(scalar_field), dimension(num_vels), intent(in) dvelr_dz_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_src_vf,
integer, intent(in) norm_dir,
type(int_bounds_info), intent(in) ix,
type(int_bounds_info), intent(in) iy,
type(int_bounds_info), intent(in) iz )

Dispatch to the subroutines that are utilized to compute the viscous source fluxes for either Cartesian or cylindrical geometries. For more information please refer to: 1) s_compute_cartesian_viscous_source_flux 2) s_compute_cylindrical_viscous_source_flux.

Definition at line 544 of file m_riemann_solvers.fpp.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ s_finalize_riemann_solver()

subroutine m_riemann_solvers::s_finalize_riemann_solver ( type(scalar_field), dimension(sys_size), intent(inout) flux_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_src_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_gsrc_vf,
integer, intent(in) norm_dir )

Deallocation and/or disassociation procedures that are needed to finalize the selected Riemann problem solver.

Definition at line 18720 of file m_riemann_solvers.fpp.f90.

Here is the caller graph for this function:

◆ s_finalize_riemann_solvers_module()

impure subroutine, public m_riemann_solvers::s_finalize_riemann_solvers_module

Module deallocation and/or disassociation procedures.

Definition at line 19198 of file m_riemann_solvers.fpp.f90.

Here is the caller graph for this function:

◆ s_hll_riemann_solver()

subroutine, public m_riemann_solvers::s_hll_riemann_solver ( real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsx_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsy_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dx_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dy_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable ql_prim_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsx_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsy_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dx_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dy_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable qr_prim_vf,
type(scalar_field), dimension(sys_size), intent(in) q_prim_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_src_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_gsrc_vf,
integer, intent(in) norm_dir,
type(int_bounds_info), intent(in) ix,
type(int_bounds_info), intent(in) iy,
type(int_bounds_info), intent(in) iz )

HLL approximate Riemann solver, Harten et al. SIAM Review (1983).

The computation of c_avg does not require all the variables, and therefore the non '_avg'

The computation of c_avg does not require all the variables, and therefore the non '_avg'

The computation of c_avg does not require all the variables, and therefore the non '_avg'

Definition at line 566 of file m_riemann_solvers.fpp.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ s_hllc_riemann_solver()

subroutine, public m_riemann_solvers::s_hllc_riemann_solver ( real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsx_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsy_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dx_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dy_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable ql_prim_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsx_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsy_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dx_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dy_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable qr_prim_vf,
type(scalar_field), dimension(sys_size), intent(in) q_prim_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_src_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_gsrc_vf,
integer, intent(in) norm_dir,
type(int_bounds_info), intent(in) ix,
type(int_bounds_info), intent(in) iy,
type(int_bounds_info), intent(in) iz )

HLLC Riemann solver with contact restoration, Toro et al. Shock Waves (1994).

The computation of c_avg does not require all the variables, and therefore the non '_avg'

The computation of c_avg does not require all the variables, and therefore the non '_avg'

The computation of c_avg does not require all the variables, and therefore the non '_avg'

gamma_method = 1: Ref. Section 2.3.1 Formulation of doi:10.7907/ZKW8-ES97.

gamma_method = 2: c_p / c_v where c_p, c_v are specific heats.

The computation of c_avg does not require all the variables, and therefore the non '_avg'

The computation of c_avg does not require all the variables, and therefore the non '_avg'

The computation of c_avg does not require all the variables, and therefore the non '_avg'

The computation of c_avg does not require all the variables, and therefore the non '_avg'

gamma_method = 1: Ref. Section 2.3.1 Formulation of doi:10.7907/ZKW8-ES97.

gamma_method = 2: c_p / c_v where c_p, c_v are specific heats.

The computation of c_avg does not require all the variables, and therefore the non '_avg'

The computation of c_avg does not require all the variables, and therefore the non '_avg'

The computation of c_avg does not require all the variables, and therefore the non '_avg'

The computation of c_avg does not require all the variables, and therefore the non '_avg'

gamma_method = 1: Ref. Section 2.3.1 Formulation of doi:10.7907/ZKW8-ES97.

gamma_method = 2: c_p / c_v where c_p, c_v are specific heats.

The computation of c_avg does not require all the variables, and therefore the non '_avg'

Definition at line 6428 of file m_riemann_solvers.fpp.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ s_hlld_riemann_solver()

subroutine, public m_riemann_solvers::s_hlld_riemann_solver ( real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsx_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsy_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dx_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dy_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable ql_prim_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsx_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsy_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dx_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dy_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable qr_prim_vf,
type(scalar_field), dimension(sys_size), intent(in) q_prim_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_src_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_gsrc_vf,
integer, intent(in) norm_dir,
type(int_bounds_info), intent(in) ix,
type(int_bounds_info), intent(in) iy,
type(int_bounds_info), intent(in) iz )

HLLD Riemann solver for MHD, Miyoshi & Kusano JCP (2005).

Definition at line 15251 of file m_riemann_solvers.fpp.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ s_initialize_riemann_solver()

subroutine m_riemann_solvers::s_initialize_riemann_solver ( type(scalar_field), dimension(sys_size), intent(inout) flux_src_vf,
integer, intent(in) norm_dir )

Set up the chosen Riemann solver algorithm for the current direction.

Definition at line 17861 of file m_riemann_solvers.fpp.f90.

Here is the caller graph for this function:

◆ s_initialize_riemann_solvers_module()

impure subroutine, public m_riemann_solvers::s_initialize_riemann_solvers_module

Initialize the Riemann solvers module.

Definition at line 16025 of file m_riemann_solvers.fpp.f90.

Here is the caller graph for this function:

◆ s_lf_riemann_solver()

subroutine, public m_riemann_solvers::s_lf_riemann_solver ( real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsx_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsy_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dx_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dy_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable ql_prim_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsx_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsy_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dx_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dy_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable qr_prim_vf,
type(scalar_field), dimension(sys_size), intent(in) q_prim_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_src_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_gsrc_vf,
integer, intent(in) norm_dir,
type(int_bounds_info), intent(in) ix,
type(int_bounds_info), intent(in) iy,
type(int_bounds_info), intent(in) iz )

Lax-Friedrichs (Rusanov) approximate Riemann solver.

Averaged velocity gradient tensor d(vel_i)/d(coord_j).

Definition at line 3636 of file m_riemann_solvers.fpp.f90.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ s_populate_riemann_states_variables_buffers()

subroutine m_riemann_solvers::s_populate_riemann_states_variables_buffers ( real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsx_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsy_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dx_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dy_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dz_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsx_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsy_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dx_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dy_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dz_vf,
integer, intent(in) norm_dir,
type(int_bounds_info), intent(in) ix,
type(int_bounds_info), intent(in) iy,
type(int_bounds_info), intent(in) iz )

Populate the left and right Riemann state variable buffers based on boundary conditions.

Definition at line 16777 of file m_riemann_solvers.fpp.f90.

Here is the caller graph for this function:

◆ s_riemann_solver()

subroutine, public m_riemann_solvers::s_riemann_solver ( real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsx_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsy_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) ql_prim_rsz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dx_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dy_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dql_prim_dz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable ql_prim_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsx_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsy_vf,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) qr_prim_rsz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dx_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dy_vf,
type(scalar_field), dimension(:), intent(inout), allocatable dqr_prim_dz_vf,
type(scalar_field), dimension(:), intent(inout), allocatable qr_prim_vf,
type(scalar_field), dimension(sys_size), intent(in) q_prim_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_src_vf,
type(scalar_field), dimension(sys_size), intent(inout) flux_gsrc_vf,
integer, intent(in) norm_dir,
type(int_bounds_info), intent(in) ix,
type(int_bounds_info), intent(in) iy,
type(int_bounds_info), intent(in) iz )

Dispatch to the subroutines that are utilized to compute the Riemann problem solution. For additional information please reference: 1) s_hll_riemann_solver 2) s_hllc_riemann_solver 3) s_lf_riemann_solver 4) s_hlld_riemann_solver.

Definition at line 493 of file m_riemann_solvers.fpp.f90.

Here is the call graph for this function:

Variable Documentation

◆ flux_gsrc_rsx_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::flux_gsrc_rsx_vf

The cell-boundary values of the geometrical source flux that are computed through the chosen Riemann problem solver by using the left and right states given in qK_prim_rs_vf. Currently 2D axisymmetric for inviscid only.

Definition at line 382 of file m_riemann_solvers.fpp.f90.

◆ flux_gsrc_rsy_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::flux_gsrc_rsy_vf

Definition at line 383 of file m_riemann_solvers.fpp.f90.

◆ flux_gsrc_rsz_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::flux_gsrc_rsz_vf

Definition at line 384 of file m_riemann_solvers.fpp.f90.

◆ flux_rsx_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::flux_rsx_vf

The cell-boundary values of the fluxes (src - source) that are computed through the chosen Riemann problem solver, and the direct evaluation of source terms, by using the left and right states given in qK_prim_rs_vf, dqK_prim_ds_vf where ds = dx, dy or dz.

Definition at line 363 of file m_riemann_solvers.fpp.f90.

◆ flux_rsy_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::flux_rsy_vf

Definition at line 364 of file m_riemann_solvers.fpp.f90.

◆ flux_rsz_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::flux_rsz_vf

Definition at line 365 of file m_riemann_solvers.fpp.f90.

◆ flux_src_rsx_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::flux_src_rsx_vf

Definition at line 363 of file m_riemann_solvers.fpp.f90.

◆ flux_src_rsy_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::flux_src_rsy_vf

Definition at line 364 of file m_riemann_solvers.fpp.f90.

◆ flux_src_rsz_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::flux_src_rsz_vf

Definition at line 365 of file m_riemann_solvers.fpp.f90.

◆ gs_rs

real(wp), dimension(:), allocatable m_riemann_solvers::gs_rs

Definition at line 463 of file m_riemann_solvers.fpp.f90.

◆ is1

type(int_bounds_info) m_riemann_solvers::is1

Definition at line 447 of file m_riemann_solvers.fpp.f90.

◆ is2

type(int_bounds_info) m_riemann_solvers::is2

Definition at line 447 of file m_riemann_solvers.fpp.f90.

◆ is3

type(int_bounds_info) m_riemann_solvers::is3

Definition at line 447 of file m_riemann_solvers.fpp.f90.

◆ isx

type(int_bounds_info) m_riemann_solvers::isx

Definition at line 448 of file m_riemann_solvers.fpp.f90.

◆ isy

type(int_bounds_info) m_riemann_solvers::isy

Definition at line 448 of file m_riemann_solvers.fpp.f90.

◆ isz

type(int_bounds_info) m_riemann_solvers::isz

Definition at line 448 of file m_riemann_solvers.fpp.f90.

◆ mom_sp_rsx_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::mom_sp_rsx_vf

Definition at line 415 of file m_riemann_solvers.fpp.f90.

◆ mom_sp_rsy_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::mom_sp_rsy_vf

Definition at line 416 of file m_riemann_solvers.fpp.f90.

◆ mom_sp_rsz_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::mom_sp_rsz_vf

Definition at line 417 of file m_riemann_solvers.fpp.f90.

◆ re_avg_rsx_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::re_avg_rsx_vf

Definition at line 430 of file m_riemann_solvers.fpp.f90.

◆ re_avg_rsy_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::re_avg_rsy_vf

Definition at line 431 of file m_riemann_solvers.fpp.f90.

◆ re_avg_rsz_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::re_avg_rsz_vf

Definition at line 432 of file m_riemann_solvers.fpp.f90.

◆ res_gs

real(wp), dimension(:,:), allocatable m_riemann_solvers::res_gs

Definition at line 476 of file m_riemann_solvers.fpp.f90.

◆ vel_src_rsx_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::vel_src_rsx_vf

Definition at line 400 of file m_riemann_solvers.fpp.f90.

◆ vel_src_rsy_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::vel_src_rsy_vf

Definition at line 401 of file m_riemann_solvers.fpp.f90.

◆ vel_src_rsz_vf

real(wp), dimension(:,:,:,:), allocatable m_riemann_solvers::vel_src_rsz_vf

Definition at line 402 of file m_riemann_solvers.fpp.f90.