|
MFC: Simulation
High-fidelity multiphase flow simulation
|
Modules | |
| module | m_weno |
| Weighted essentially non-oscillatory (WENO) reconstruction scheme that is supplemented with monotonicity preserving bounds (MPWENO) and a mapping function that boosts the accuracy of the non-linear weights (WENOM). MPWENO, see Balsara and Shu (2000), prevents the reconstructed values to lay outside the range set by the stencil, while WENOM, see Henrick et al. (2005), recovers the formal order of accuracy of the reconstruction at critical points. Please note that the basic WENO approach is implemented according to the work of Jiang and Shu (1996). WENO-Z, which is less dissipative than WENO-JS and WENO-M, is implemented according to the work of Borges, et al. (2008). TENO, which is even less dissipative than WENO-Z but is less robust, is implemented according to the work of Fu et al. (2016). | |
Functions/Subroutines | |
| impure subroutine, public | m_weno::s_initialize_weno_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 | m_weno::s_compute_weno_coefficients (weno_dir, is) |
| The purpose of this subroutine is to compute the grid dependent coefficients of the WENO polynomials, ideal weights and smoothness indicators, provided the order, the coordinate direction and the location of the WENO reconstruction. | |
| subroutine, public | m_weno::s_weno (v_vf, vl_rs_vf_x, vl_rs_vf_y, vl_rs_vf_z, vr_rs_vf_x, vr_rs_vf_y, vr_rs_vf_z, weno_dir, is1_weno_d, is2_weno_d, is3_weno_d) |
| subroutine, public | m_weno::s_initialize_weno (v_vf, weno_dir) |
| The computation of parameters, the allocation of memory, the association of pointers and/or the execution of any other procedures that are required for the setup of the WENO reconstruction. | |
| subroutine | m_weno::s_preserve_monotonicity (v_rs_ws, vl_rs_vf, vr_rs_vf) |
| The goal of this subroutine is to ensure that the WENO reconstruction is monotonic. The latter is achieved by enforcing monotonicity preserving bounds of Suresh and Huynh (1997). The resulting MPWENO reconstruction, see Balsara and Shu (2000), ensures that the reconstructed values do not reside outside the range spanned by WENO stencil. | |
| impure subroutine, public | m_weno::s_finalize_weno_module () |
| Module deallocation and/or disassociation procedures. | |
Variables | |
| integer | m_weno::v_size |
| Number of WENO-reconstructed cell-average variables. | |
The cell-average variables that will be WENO-reconstructed. Formerly, they | |
are stored in v_vf. However, they are transferred to v_rs_wsL and v_rs_wsR as to be reshaped (RS) and/or characteristically decomposed. The reshaping allows the WENO procedure to be independent of the coordinate direction of the reconstruction. Lastly, notice that the left (L) and right (R) results of the characteristic decomposition are stored in custom-constructed WENO- stencils (WS) that are annexed to each position of a given scalar field. | |
| real(wp), dimension(:, :, :, :), allocatable | m_weno::v_rs_ws_x |
| real(wp), dimension(:, :, :, :), allocatable | m_weno::v_rs_ws_y |
| real(wp), dimension(:, :, :, :), allocatable | m_weno::v_rs_ws_z |
Polynomial coefficients at the left and right cell-boundaries (CB) and at | |
the left and right quadrature points (QP), in the x-, y- and z-directions. Note that the first dimension of the array identifies the polynomial, the second dimension identifies the position of its coefficients and the last dimension denotes the cell-location in the relevant coordinate direction. | |
| real(wp), dimension(:, :, :), allocatable, target | m_weno::poly_coef_cbl_x |
| real(wp), dimension(:, :, :), allocatable, target | m_weno::poly_coef_cbl_y |
| real(wp), dimension(:, :, :), allocatable, target | m_weno::poly_coef_cbl_z |
| real(wp), dimension(:, :, :), allocatable, target | m_weno::poly_coef_cbr_x |
| real(wp), dimension(:, :, :), allocatable, target | m_weno::poly_coef_cbr_y |
| real(wp), dimension(:, :, :), allocatable, target | m_weno::poly_coef_cbr_z |
The ideal weights at the left and the right cell-boundaries and at the | |
left and the right quadrature points, in x-, y- and z-directions. Note that the first dimension of the array identifies the weight, while the last denotes the cell-location in the relevant coordinate direction. | |
| real(wp), dimension(:, :), allocatable, target | m_weno::d_cbl_x |
| real(wp), dimension(:, :), allocatable, target | m_weno::d_cbl_y |
| real(wp), dimension(:, :), allocatable, target | m_weno::d_cbl_z |
| real(wp), dimension(:, :), allocatable, target | m_weno::d_cbr_x |
| real(wp), dimension(:, :), allocatable, target | m_weno::d_cbr_y |
| real(wp), dimension(:, :), allocatable, target | m_weno::d_cbr_z |
Smoothness indicator coefficients in the x-, y-, and z-directions. Note | |
that the first array dimension identifies the smoothness indicator, the second identifies the position of its coefficients and the last denotes the cell-location in the relevant coordinate direction. | |
| real(wp), dimension(:, :, :), allocatable, target | m_weno::beta_coef_x |
| real(wp), dimension(:, :, :), allocatable, target | m_weno::beta_coef_y |
| real(wp), dimension(:, :, :), allocatable, target | m_weno::beta_coef_z |
Indical bounds in the s1-, s2- and s3-directions | |
| type(int_bounds_info) | m_weno::is1_weno |
| type(int_bounds_info) | m_weno::is2_weno |
| type(int_bounds_info) | m_weno::is3_weno |