MFC
Exascale flow solver
Loading...
Searching...
No Matches
m_thinc.fpp.f90 File Reference

Contains module m_thinc. More...

Go to the source code of this file.

Modules

module  m_thinc
 THINC and MTHINC interface compression for volume fraction sharpening. THINC (int_comp=1): 1D directional interface compression applied after MUSCL/WENO reconstruction. Uses hyperbolic tangent profile per dimension. MTHINC (int_comp=2): Multi-dimensional THINC that reconstructs a tanh profile oriented along the interface normal (computed from the gradient of alpha), then integrates that profile over each cell face using Gaussian quadrature. A Newton iteration enforces the conservation constraint (cell-averaged alpha is preserved). Reference: B. Xie and F. Xiao, "Toward efficient and accurate interface capturing on arbitrary hybrid unstructured grids: The THINC method with quadratic surface representation and Gaussian quadrature," Journal of Computational Physics, vol. 349, pp. 415-440, 2017.

Functions/Subroutines

real(wp) function m_thinc::f_log_cosh (x)
 Stable computation of ln(cosh(x)).
real(wp) function m_thinc::f_log_cosh_diff (a, h)
 Stable difference: log_cosh(a+h) - log_cosh(a-h) = 2*atanh(tanh(a)*tanh(h)). Avoids catastrophic cancellation when h is small relative to a.
real(wp) function m_thinc::f_thinc_integral_1d (a, b)
 Analytical 1-D integral of the THINC function.
real(wp) function m_thinc::f_mthinc_volume_integral (n1, n2, n3, d, beta, ndim)
 Volume integral of H(xi) = 0.5*(1 + tanh(beta*(n.xi + d))) over the cell [-1/2, 1/2]^ndim.
real(wp) function m_thinc::f_mthinc_volume_integral_dd (n1, n2, n3, d, beta, ndim)
 Derivative dV/dd of the volume integral (for Newton iteration).
real(wp) function m_thinc::f_mthinc_solve_d (n1, n2, n3, beta, alpha_cell, ndim)
 Solve for the interface-position parameter d.
real(wp) function m_thinc::f_mthinc_face_average (n1, n2, n3, d, beta, face_dir, face_pos, ndim)
 Face-averaged THINC function at a cell face. face_dir: 1=x, 2=y, 3=z. face_pos: -0.5 (low) or +0.5 (high). The face coordinate in face_dir is fixed; remaining directions are integrated over [-1/2, 1/2] (analytically along one, Gauss for others).
subroutine, public m_thinc::s_initialize_thinc_module ()
subroutine, public m_thinc::s_compute_mthinc_normals (v_vf)
 Compute the unit normal and interface-position parameter d at each interface cell.
subroutine, public m_thinc::s_thinc_compression (v_rs_ws, vl_rs_vf_x, vl_rs_vf_y, vl_rs_vf_z, vr_rs_vf_x, vr_rs_vf_y, vr_rs_vf_z, recon_dir, is1_d, is2_d, is3_d)
 Applies THINC (int_comp=1) or MTHINC (int_comp=2) interface compression to sharpen volume-fraction and density reconstructions at material interfaces. Called after WENO/MUSCL reconstruction per direction.
subroutine, public m_thinc::s_finalize_thinc_module ()

Variables

real(wp), dimension(3) m_thinc::gq3_pts = [-5e-1_wp*0.7745966692414834_wp, 0._wp, 5e-1_wp*0.7745966692414834_wp]
 3-point Gauss-Legendre quadrature on [-1/2, 1/2] Node locations: +-sqrt(3/5)/2, 0
real(wp), dimension(3) m_thinc::gq3_wts = [5._wp/18._wp, 8._wp/18._wp, 5._wp/18._wp]
 Weights: 5/18, 8/18, 5/18.
real(wpm_thinc::ln2 = 0.6931471805599453_wp
 ln(2)
real(wp), dimension(:,:,:,:), allocatable m_thinc::mthinc_nhat
real(wp), dimension(:,:,:,:), allocatable m_thinc::unit
real(wp), dimension(:,:,:,:), allocatable m_thinc::normal
real(wp), dimension(:,:,:,:), allocatable m_thinc::vector
real(wp), dimension(:,:,:), allocatable m_thinc::mthinc_d
real(wp), dimension(:,:,:), allocatable m_thinc::interface
real(wp), dimension(:,:,:), allocatable m_thinc::position
real(wp), dimension(:,:,:), allocatable m_thinc::parameter

Detailed Description

Contains module m_thinc.

Definition in file m_thinc.fpp.f90.