|
MFC: Simulation
High-fidelity multiphase flow simulation
|
This module features a database of approximate and exact Riemann problem solvers for the Navier-Stokes system of equations, which is supplemented by appropriate advection equations that are used to capture the material interfaces. The closure of the system is achieved by the stiffened gas equation of state and any required mixture relations. Surface tension effects are accounted for and are modeled by means of a volume force acting across the diffuse material interface region. The implementation details of viscous and capillary effects, into the Riemann solvers, may be found in Perigaud and Saurel (2005). Note that both effects are available only in the volume fraction model. At this time, the approximate and exact Riemann solvers that are listed below are available: 1) Harten-Lax-van Leer (HLL) 2) Harten-Lax-van Leer-Contact (HLLC) 3) Exact 4) Harten-Lax-van Leer Discontinuities (HLLD) - for MHD only. 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_exact_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) |
| 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) |
| 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) |
| This procedure is the implementation of the Harten, Lax, van Leer, and contact (HLLC) approximate Riemann solver, see Toro (1999) and Johnsen (2007). The viscous and the surface tension effects have been included by modifying the exact Riemann solver of Perigaud and Saurel (2005). | |
| 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 resolves 5 of the 7 waves of MHD equations: 1 entropy wave, 2 Alfvén waves, 2 fast magnetosonic waves. | |
| impure subroutine, public | s_initialize_riemann_solvers_module |
| The computation of parameters, the allocation of memory, the association of pointers and/or the execution of any other procedures that are necessary to setup the 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) |
| The purpose of this subroutine is to populate the buffers of the left and right Riemann states variables, depending on the boundary conditions. | |
| subroutine | s_initialize_riemann_solver (flux_src_vf, norm_dir) |
| The computation of parameters, the allocation of memory, the association of pointers and/or the execution of any other procedures needed to configure the chosen Riemann solver algorithm. | |
| 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) |
Computes cylindrical viscous source flux contributions for momentum and energy. Calculates Cartesian components of the stress tensor using averaged velocity derivatives and cylindrical geometric factors, then updates flux_src_vf. Assumes x-dir is axial (z_cyl), y-dir is radial (r_cyl), z-dir is azimuthal (theta_cyl for derivatives). | |
| 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) |
Computes Cartesian viscous source flux contributions for momentum and energy. Calculates averaged velocity gradients, gets Re and interface velocities, calls helpers for shear/bulk stress, then updates flux_src_vf. | |
| subroutine | s_calculate_shear_stress_tensor (vel_grad_avg, re_shear, divergence_v, tau_shear_out) |
| Calculates shear stress tensor components. tau_ij_shear = ( (dui/dxj + duj/dxi) - (2/3)*(div_v)*delta_ij ) / Re_shear. | |
| subroutine | s_calculate_bulk_stress_tensor (re_bulk, divergence_v, tau_bulk_out) |
| Calculates bulk stress tensor components (diagonal only). tau_ii_bulk = (div_v) / Re_bulk. Off-diagonals are zero. | |
| 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_info) | is1 |
| type(int_bounds_info) | is2 |
| type(int_bounds_info) | is3 |
| type(int_bounds_info) | isx |
| type(int_bounds_info) | isy |
| type(int_bounds_info) | isz |
This module features a database of approximate and exact Riemann problem solvers for the Navier-Stokes system of equations, which is supplemented by appropriate advection equations that are used to capture the material interfaces. The closure of the system is achieved by the stiffened gas equation of state and any required mixture relations. Surface tension effects are accounted for and are modeled by means of a volume force acting across the diffuse material interface region. The implementation details of viscous and capillary effects, into the Riemann solvers, may be found in Perigaud and Saurel (2005). Note that both effects are available only in the volume fraction model. At this time, the approximate and exact Riemann solvers that are listed below are available: 1) Harten-Lax-van Leer (HLL) 2) Harten-Lax-van Leer-Contact (HLLC) 3) Exact 4) Harten-Lax-van Leer Discontinuities (HLLD) - for MHD only.
| 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 ) |
Calculates bulk stress tensor components (diagonal only). tau_ii_bulk = (div_v) / Re_bulk. Off-diagonals are zero.
| [in] | Re_bulk | Bulk Reynolds number. |
| [in] | divergence_v | Velocity divergence (du/dx + dv/dy + dw/dz). |
| [out] | tau_bulk_out | Calculated bulk stress tensor (stress on i-face, i-direction). |
| 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 ) |
Calculates shear stress tensor components. tau_ij_shear = ( (dui/dxj + duj/dxi) - (2/3)*(div_v)*delta_ij ) / Re_shear.
| [in] | vel_grad_avg | Averaged velocity gradient tensor (d(vel_i)/d(coord_j)). |
| [in] | Re_shear | Shear Reynolds number. |
| [in] | divergence_v | Velocity divergence (du/dx + dv/dy + dw/dz). |
| [out] | tau_shear_out | Calculated shear stress tensor (stress on i-face, j-direction). |
| 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 ) |
Computes Cartesian viscous source flux contributions for momentum and energy. Calculates averaged velocity gradients, gets Re and interface velocities, calls helpers for shear/bulk stress, then updates flux_src_vf.
| [in] | velL_vf | Left boundary velocity (num_dims scalar_field). |
| [in] | dvelL_dx_vf | Left boundary d(vel)/dx (num_dims scalar_field). |
| [in] | dvelL_dy_vf | Left boundary d(vel)/dy (num_dims scalar_field). |
| [in] | dvelL_dz_vf | Left boundary d(vel)/dz (num_dims scalar_field). |
| [in] | velR_vf | Right boundary velocity (num_dims scalar_field). |
| [in] | dvelR_dx_vf | Right boundary d(vel)/dx (num_dims scalar_field). |
| [in] | dvelR_dy_vf | Right boundary d(vel)/dy (num_dims scalar_field). |
| [in] | dvelR_dz_vf | Right boundary d(vel)/dz (num_dims scalar_field). |
| [in,out] | flux_src_vf | Intercell source flux array to update (sys_size scalar_field). |
| [in] | norm_dir | Interface normal direction (1=x, 2=y, 3=z). |
| [in] | ix | X-direction loop bounds (int_bounds_info). |
| [in] | iy | Y-direction loop bounds (int_bounds_info). |
| [in] | iz | Z-direction loop bounds (int_bounds_info). |
| 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 ) |
Computes cylindrical viscous source flux contributions for momentum and energy. Calculates Cartesian components of the stress tensor using averaged velocity derivatives and cylindrical geometric factors, then updates flux_src_vf. Assumes x-dir is axial (z_cyl), y-dir is radial (r_cyl), z-dir is azimuthal (theta_cyl for derivatives).
| [in] | velL_vf | Left boundary velocity ($v_x, v_y, v_z$) (num_dims scalar_field). |
| [in] | dvelL_dx_vf | Left boundary $\partial v_i/\partial x$ (num_dims scalar_field). |
| [in] | dvelL_dy_vf | Left boundary $\partial v_i/\partial y$ (num_dims scalar_field). |
| [in] | dvelL_dz_vf | Left boundary $\partial v_i/\partial z$ (num_dims scalar_field). |
| [in] | velR_vf | Right boundary velocity ($v_x, v_y, v_z$) (num_dims scalar_field). |
| [in] | dvelR_dx_vf | Right boundary $\partial v_i/\partial x$ (num_dims scalar_field). |
| [in] | dvelR_dy_vf | Right boundary $\partial v_i/\partial y$ (num_dims scalar_field). |
| [in] | dvelR_dz_vf | Right boundary $\partial v_i/\partial z$ (num_dims scalar_field). |
| [in,out] | flux_src_vf | Intercell source flux array to update (sys_size scalar_field). |
| [in] | norm_dir | Interface normal direction (1=x-face, 2=y-face, 3=z-face). |
| [in] | ix | Global X-direction loop bounds (int_bounds_info). |
| [in] | iy | Global Y-direction loop bounds (int_bounds_info). |
| [in] | iz | Global Z-direction loop bounds (int_bounds_info). |
| 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.
| 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.
| flux_vf | Intercell fluxes |
| flux_src_vf | Intercell source fluxes |
| flux_gsrc_vf | Intercell geometric source fluxes |
| norm_dir | Dimensional splitting coordinate direction |
| impure subroutine, public m_riemann_solvers::s_finalize_riemann_solvers_module |
Module deallocation and/or disassociation procedures.
| 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 ) |
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'
| 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 ) |
This procedure is the implementation of the Harten, Lax, van Leer, and contact (HLLC) approximate Riemann solver, see Toro (1999) and Johnsen (2007). The viscous and the surface tension effects have been included by modifying the exact Riemann solver of Perigaud and Saurel (2005).
| qL_prim_vf | The left WENO-reconstructed cell-boundary values of the cell-average primitive variables |
| qR_prim_vf | The right WENO-reconstructed cell-boundary values of the cell-average primitive variables |
| dqL_prim_dx_vf | The left WENO-reconstructed cell-boundary values of the first-order x-dir spatial derivatives |
| dqL_prim_dy_vf | The left WENO-reconstructed cell-boundary values of the first-order y-dir spatial derivatives |
| dqL_prim_dz_vf | The left WENO-reconstructed cell-boundary values of the first-order z-dir spatial derivatives |
| dqR_prim_dx_vf | The right WENO-reconstructed cell-boundary values of the first-order x-dir spatial derivatives |
| dqR_prim_dy_vf | The right WENO-reconstructed cell-boundary values of the first-order y-dir spatial derivatives |
| dqR_prim_dz_vf | The right WENO-reconstructed cell-boundary values of the first-order z-dir spatial derivatives |
| gm_alphaL_vf | Left averaged gradient magnitude |
| gm_alphaR_vf | Right averaged gradient magnitude |
| flux_vf | Intra-cell fluxes |
| flux_src_vf | Intra-cell fluxes sources |
| flux_gsrc_vf | Intra-cell geometric fluxes sources |
| norm_dir | Dir. splitting direction |
| ix | Index bounds in the x-dir |
| iy | Index bounds in the y-dir |
| iz | Index bounds in the z-dir |
| q_prim_vf | Cell-averaged primitive variables |
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'
| 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 resolves 5 of the 7 waves of MHD equations: 1 entropy wave, 2 Alfvén waves, 2 fast magnetosonic waves.
| subroutine m_riemann_solvers::s_initialize_riemann_solver | ( | type(scalar_field), dimension(sys_size), intent(inout) | flux_src_vf, |
| integer, intent(in) | norm_dir ) |
The computation of parameters, the allocation of memory, the association of pointers and/or the execution of any other procedures needed to configure the chosen Riemann solver algorithm.
| qL_prim_vf | The left WENO-reconstructed cell-boundary values of the cell-average primitive variables |
| qR_prim_vf | The right WENO-reconstructed cell-boundary values of the cell-average primitive variables |
| flux_vf | Intra-cell fluxes |
| flux_src_vf | Intra-cell fluxes sources |
| flux_gsrc_vf | Intra-cell geometric fluxes sources |
| norm_dir | Dir. splitting direction |
| ix | Index bounds in the x-dir |
| iy | Index bounds in the y-dir |
| iz | Index bounds in the z-dir |
| q_prim_vf | Cell-averaged primitive variables |
| impure subroutine, public m_riemann_solvers::s_initialize_riemann_solvers_module |
The computation of parameters, the allocation of memory, the association of pointers and/or the execution of any other procedures that are necessary to setup the module.
| 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 ) |
| 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 ) |
The purpose of this subroutine is to populate the buffers of the left and right Riemann states variables, depending on the boundary conditions.
| qL_prim_vf | The left WENO-reconstructed cell-boundary values of the cell-average primitive variables |
| qR_prim_vf | The right WENO-reconstructed cell-boundary values of the cell-average primitive variables |
| dqL_prim_dx_vf | The left WENO-reconstructed cell-boundary values of the first-order x-dir spatial derivatives |
| dqL_prim_dy_vf | The left WENO-reconstructed cell-boundary values of the first-order y-dir spatial derivatives |
| dqL_prim_dz_vf | The left WENO-reconstructed cell-boundary values of the first-order z-dir spatial derivatives |
| dqR_prim_dx_vf | The right WENO-reconstructed cell-boundary values of the first-order x-dir spatial derivatives |
| dqR_prim_dy_vf | The right WENO-reconstructed cell-boundary values of the first-order y-dir spatial derivatives |
| dqR_prim_dz_vf | The right WENO-reconstructed cell-boundary values of the first-order z-dir spatial derivatives |
| gm_alphaL_vf | Left averaged gradient magnitude |
| gm_alphaR_vf | Right averaged gradient magnitude |
| norm_dir | Dir. splitting direction |
| ix | Index bounds in the x-dir |
| iy | Index bounds in the y-dir |
| iz | Index bounds in the z-dir |
| 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_exact_riemann_solver 4) s_hlld_riemann_solver.
| qL_prim_vf | The left WENO-reconstructed cell-boundary values of the cell-average primitive variables |
| qR_prim_vf | The right WENO-reconstructed cell-boundary values of the cell-average primitive variables |
| dqL_prim_dx_vf | The left WENO-reconstructed cell-boundary values of the first-order x-dir spatial derivatives |
| dqL_prim_dy_vf | The left WENO-reconstructed cell-boundary values of the first-order y-dir spatial derivatives |
| dqL_prim_dz_vf | The left WENO-reconstructed cell-boundary values of the first-order z-dir spatial derivatives |
| dqR_prim_dx_vf | The right WENO-reconstructed cell-boundary values of the first-order x-dir spatial derivatives |
| dqR_prim_dy_vf | The right WENO-reconstructed cell-boundary values of the first-order y-dir spatial derivatives |
| dqR_prim_dz_vf | The right WENO-reconstructed cell-boundary values of the first-order z-dir spatial derivatives |
| gm_alphaL_vf | Left averaged gradient magnitude |
| gm_alphaR_vf | Right averaged gradient magnitude |
| flux_vf | Intra-cell fluxes |
| flux_src_vf | Intra-cell fluxes sources |
| flux_gsrc_vf | Intra-cell geometric fluxes sources |
| norm_dir | Dir. splitting direction |
| ix | Index bounds in the x-dir |
| iy | Index bounds in the y-dir |
| iz | Index bounds in the z-dir |
| q_prim_vf | Cell-averaged primitive variables |
| 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.
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::flux_gsrc_rsy_vf |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::flux_gsrc_rsz_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.
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::flux_rsy_vf |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::flux_rsz_vf |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::flux_src_rsx_vf |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::flux_src_rsy_vf |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::flux_src_rsz_vf |
| real(wp), dimension(:), allocatable m_riemann_solvers::gs_rs |
| type(int_bounds_info) m_riemann_solvers::is1 |
| type(int_bounds_info) m_riemann_solvers::is2 |
| type(int_bounds_info) m_riemann_solvers::is3 |
| type(int_bounds_info) m_riemann_solvers::isx |
| type(int_bounds_info) m_riemann_solvers::isy |
| type(int_bounds_info) m_riemann_solvers::isz |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::mom_sp_rsx_vf |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::mom_sp_rsy_vf |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::mom_sp_rsz_vf |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::re_avg_rsx_vf |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::re_avg_rsy_vf |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::re_avg_rsz_vf |
| real(wp), dimension(:, :), allocatable m_riemann_solvers::res_gs |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::vel_src_rsx_vf |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::vel_src_rsy_vf |
| real(wp), dimension(:, :, :, :), allocatable m_riemann_solvers::vel_src_rsz_vf |