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

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. More...

Functions/Subroutines

real(wp) function f_log_cosh (x)
 Stable computation of ln(cosh(x)).
real(wp) function 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 f_thinc_integral_1d (a, b)
 Analytical 1-D integral of the THINC function.
real(wp) function 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 f_mthinc_volume_integral_dd (n1, n2, n3, d, beta, ndim)
 Derivative dV/dd of the volume integral (for Newton iteration).
real(wp) function f_mthinc_solve_d (n1, n2, n3, beta, alpha_cell, ndim)
 Solve for the interface-position parameter d.
real(wp) function 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 s_initialize_thinc_module ()
subroutine, public s_compute_mthinc_normals (v_vf)
 Compute the unit normal and interface-position parameter d at each interface cell.
subroutine, public 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 s_finalize_thinc_module ()

Variables

real(wp), dimension(3) 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) gq3_wts = [5._wp/18._wp, 8._wp/18._wp, 5._wp/18._wp]
 Weights: 5/18, 8/18, 5/18.
real(wpln2 = 0.6931471805599453_wp
 ln(2)
real(wp), dimension(:,:,:,:), allocatable mthinc_nhat
real(wp), dimension(:,:,:,:), allocatable unit
real(wp), dimension(:,:,:,:), allocatable normal
real(wp), dimension(:,:,:,:), allocatable vector
real(wp), dimension(:,:,:), allocatable mthinc_d
real(wp), dimension(:,:,:), allocatable interface
real(wp), dimension(:,:,:), allocatable position
real(wp), dimension(:,:,:), allocatable parameter

Detailed Description

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.

Function/Subroutine Documentation

◆ f_log_cosh()

real(wp) function m_thinc::f_log_cosh ( real(wp), intent(in) x)

Stable computation of ln(cosh(x)).

Definition at line 373 of file m_thinc.fpp.f90.

◆ f_log_cosh_diff()

real(wp) function m_thinc::f_log_cosh_diff ( real(wp), intent(in) a,
real(wp), intent(in) 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.

Definition at line 404 of file m_thinc.fpp.f90.

Here is the caller graph for this function:

◆ f_mthinc_face_average()

real(wp) function m_thinc::f_mthinc_face_average ( real(wp), intent(in) n1,
real(wp), intent(in) n2,
real(wp), intent(in) n3,
real(wp), intent(in) d,
real(wp), intent(in) beta,
integer, intent(in) face_dir,
real(wp), intent(in) face_pos,
integer, intent(in) 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).

Definition at line 589 of file m_thinc.fpp.f90.

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

◆ f_mthinc_solve_d()

real(wp) function m_thinc::f_mthinc_solve_d ( real(wp), intent(in) n1,
real(wp), intent(in) n2,
real(wp), intent(in) n3,
real(wp), intent(in) beta,
real(wp), intent(in) alpha_cell,
integer, intent(in) ndim )

Solve for the interface-position parameter d.

Definition at line 552 of file m_thinc.fpp.f90.

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

◆ f_mthinc_volume_integral()

real(wp) function m_thinc::f_mthinc_volume_integral ( real(wp), intent(in) n1,
real(wp), intent(in) n2,
real(wp), intent(in) n3,
real(wp), intent(in) d,
real(wp), intent(in) beta,
integer, intent(in) ndim )

Volume integral of H(xi) = 0.5*(1 + tanh(beta*(n.xi + d))) over the cell [-1/2, 1/2]^ndim.

Definition at line 459 of file m_thinc.fpp.f90.

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

◆ f_mthinc_volume_integral_dd()

real(wp) function m_thinc::f_mthinc_volume_integral_dd ( real(wp), intent(in) n1,
real(wp), intent(in) n2,
real(wp), intent(in) n3,
real(wp), intent(in) d,
real(wp), intent(in) beta,
integer, intent(in) ndim )

Derivative dV/dd of the volume integral (for Newton iteration).

Definition at line 502 of file m_thinc.fpp.f90.

Here is the caller graph for this function:

◆ f_thinc_integral_1d()

real(wp) function m_thinc::f_thinc_integral_1d ( real(wp), intent(in) a,
real(wp), intent(in) b )

Analytical 1-D integral of the THINC function.

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

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

◆ s_compute_mthinc_normals()

subroutine, public m_thinc::s_compute_mthinc_normals ( type(scalar_field), dimension(:), intent(in) v_vf)

Compute the unit normal and interface-position parameter d at each interface cell.

Definition at line 715 of file m_thinc.fpp.f90.

Here is the call graph for this function:

◆ s_finalize_thinc_module()

subroutine, public m_thinc::s_finalize_thinc_module

Definition at line 1253 of file m_thinc.fpp.f90.

Here is the caller graph for this function:

◆ s_initialize_thinc_module()

subroutine, public m_thinc::s_initialize_thinc_module

Definition at line 640 of file m_thinc.fpp.f90.

Here is the caller graph for this function:

◆ s_thinc_compression()

subroutine, public m_thinc::s_thinc_compression ( real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(in) v_rs_ws,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) vl_rs_vf_x,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) vl_rs_vf_y,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) vl_rs_vf_z,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) vr_rs_vf_x,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) vr_rs_vf_y,
real(wp), dimension(idwbuff(1)%beg:,idwbuff(2)%beg:,idwbuff(3)%beg:,1:), intent(inout) vr_rs_vf_z,
integer, intent(in) recon_dir,
type(int_bounds_info), intent(in) is1_d,
type(int_bounds_info), intent(in) is2_d,
type(int_bounds_info), intent(in) 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.

Definition at line 860 of file m_thinc.fpp.f90.

Here is the call graph for this function:

Variable Documentation

◆ gq3_pts

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

Definition at line 339 of file m_thinc.fpp.f90.

◆ gq3_wts

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.

Definition at line 341 of file m_thinc.fpp.f90.

◆ interface

real(wp), dimension(:,:,:), allocatable m_thinc::interface

Definition at line 357 of file m_thinc.fpp.f90.

◆ ln2

real(wp) m_thinc::ln2 = 0.6931471805599453_wp

ln(2)

Definition at line 343 of file m_thinc.fpp.f90.

◆ mthinc_d

real(wp), dimension(:,:,:), allocatable m_thinc::mthinc_d

Definition at line 357 of file m_thinc.fpp.f90.

◆ mthinc_nhat

real(wp), dimension(:,:,:,:), allocatable m_thinc::mthinc_nhat

Definition at line 356 of file m_thinc.fpp.f90.

◆ normal

real(wp), dimension(:,:,:,:), allocatable m_thinc::normal

Definition at line 356 of file m_thinc.fpp.f90.

◆ parameter

real(wp), dimension(:,:,:), allocatable m_thinc::parameter

Definition at line 357 of file m_thinc.fpp.f90.

◆ position

real(wp), dimension(:,:,:), allocatable m_thinc::position

Definition at line 357 of file m_thinc.fpp.f90.

◆ unit

real(wp), dimension(:,:,:,:), allocatable m_thinc::unit

Definition at line 356 of file m_thinc.fpp.f90.

◆ vector

real(wp), dimension(:,:,:,:), allocatable m_thinc::vector

Definition at line 356 of file m_thinc.fpp.f90.