1# 1 "/home/runner/work/MFC/MFC/src/post_process/m_derived_variables.fpp"
22 real(wp),
allocatable,
dimension(:,:,:) ::
gm_rho_sf
27 real(wp),
allocatable,
dimension(:,:),
public ::
fd_coeff_x
28 real(wp),
allocatable,
dimension(:,:),
public ::
fd_coeff_y
29 real(wp),
allocatable,
dimension(:,:),
public ::
fd_coeff_z
62 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
63 & intent(inout) :: q_sf
81 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
82 & intent(inout) :: q_sf
100 type(
scalar_field),
dimension(sys_size),
intent(in) :: q_prim_vf
102 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
103 & intent(inout) :: q_sf
106 real(wp) :: blkmod1, blkmod2
118 &
k)/blkmod1 + (1._wp - q_prim_vf(
eqn_idx%adv%beg)%sf(i,
j,
k))/blkmod2)))
122 q_sf(i,
j,
k) = 1.e-16_wp
124 q_sf(i,
j,
k) = sqrt(q_sf(i,
j,
k))
136 integer,
intent(in) :: i
137 type(
scalar_field),
dimension(sys_size),
intent(in) :: q_prim_vf
139 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
140 & intent(inout) :: q_sf
142 real(wp) :: top, bottom, slope
148 if (q_prim_vf(
eqn_idx%cont%end + i)%sf(
j,
k,
l) >= 0._wp)
then
155 else if (i == 2)
then
156 if (q_prim_vf(
eqn_idx%cont%end + i)%sf(
j,
k,
l) >= 0._wp)
then
164 if (q_prim_vf(
eqn_idx%cont%end + i)%sf(
j,
k,
l) >= 0._wp)
then
173 if (abs(top) < 1.e-8_wp) top = 0._wp
174 if (abs(bottom) < 1.e-8_wp) bottom = 0._wp
179 slope = (top*bottom)/(bottom**2._wp + 1.e-16_wp)
183 q_sf(
j,
k,
l) = max(0._wp, min(1._wp, slope))
185 q_sf(
j,
k,
l) = max(0._wp, min(2._wp*slope, 5.e-1_wp*(1._wp + slope), 2._wp))
187 q_sf(
j,
k,
l) = (15.e-1_wp*(slope**2._wp + slope))/(slope**2._wp + slope + 1._wp)
189 q_sf(
j,
k,
l) = max(0._wp, min(1._wp, 2._wp*slope), min(slope, 2._wp))
191 q_sf(
j,
k,
l) = max(0._wp, min(15.e-1_wp*slope, 1._wp), min(slope, 15.e-1_wp))
193 q_sf(
j,
k,
l) = (slope**2._wp + slope)/(slope**2._wp + 1._wp)
195 q_sf(
j,
k,
l) = (abs(slope) + slope)/(1._wp + abs(slope))
207 integer,
intent(in) :: i
208 type(
scalar_field),
dimension(sys_size),
intent(in) :: q_prim_vf
210 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
211 & intent(inout) :: q_sf
213 integer ::
j,
k,
l, r
218 q_sf(
j,
k,
l) = 0._wp
224 &
l)*q_prim_vf(
eqn_idx%mom%beg + 1)%sf(
j,
k, r +
l))
233 else if (i == 2)
then
237 q_sf(
j,
k,
l) = 0._wp
255 q_sf(
j,
k,
l) = 0._wp
272 type(
scalar_field),
dimension(sys_size),
intent(in) :: q_prim_vf
274 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
275 & intent(inout) :: q_sf
277 real(wp),
dimension(1:3,1:3) :: q_jacobian_sf, s, s2, o, o2
278 real(wp) :: trs, q, iis
279 integer ::
j,
k,
l, r, jj, kk
284 q_jacobian_sf(:,:) = 0._wp
289 q_jacobian_sf(jj, 1) = q_jacobian_sf(jj, 1) +
fd_coeff_x(r, &
290 &
j)*q_prim_vf(
eqn_idx%mom%beg + jj - 1)%sf(r +
j,
k,
l)
292 q_jacobian_sf(jj, 2) = q_jacobian_sf(jj, 2) +
fd_coeff_y(r, &
293 &
k)*q_prim_vf(
eqn_idx%mom%beg + jj - 1)%sf(
j, r +
k,
l)
295 q_jacobian_sf(jj, 3) = q_jacobian_sf(jj, 3) +
fd_coeff_z(r, &
296 &
l)*q_prim_vf(
eqn_idx%mom%beg + jj - 1)%sf(
j,
k, r +
l)
303 s(jj, kk) = 0.5_wp*(q_jacobian_sf(jj, kk) + q_jacobian_sf(kk, jj))
304 o(jj, kk) = 0.5_wp*(q_jacobian_sf(jj, kk) - q_jacobian_sf(kk, jj))
310 o2(jj, kk) = o(jj, 1)*o(kk, 1) + o(jj, 2)*o(kk, 2) + o(jj, 3)*o(kk, 3)
311 s2(jj, kk) = s(jj, 1)*s(kk, 1) + s(jj, 2)*s(kk, 2) + s(jj, 3)*s(kk, 3)
316 q = 0.5_wp*((o2(1, 1) + o2(2, 2) + o2(3, 3)) - (s2(1, 1) + s2(2, 2) + s2(3, 3)))
317 trs = s(1, 1) + s(2, 2) + s(3, 3)
319 iis = 0.5_wp*((s(1, 1) + s(2, 2) + s(3, 3))**2 - (s2(1, 1) + s2(2, 2) + s2(3, 3)))
320 q_sf(
j,
k,
l) = q + iis
332 integer,
parameter :: nm = 3
333 type(
scalar_field),
dimension(sys_size),
intent(in) :: q_prim_vf
337 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
338 & intent(out) :: liutex_mag
340 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end,nm), &
341 & intent(out) :: liutex_axis
342 character,
parameter :: ivl =
'N'
343 character,
parameter :: ivr =
'V'
344 real(wp),
dimension(nm, nm) :: vgt
345 real(wp),
dimension(nm) :: lr, li
346 real(wp),
dimension(nm, nm) :: vl, vr
347 integer,
parameter :: lwork = 4*nm
348 real(wp),
dimension(lwork) :: work
350 real(wp),
dimension(nm) :: eigvec
351 real(wp) :: eigvec_mag
352 real(wp) :: omega_proj
355 integer ::
j,
k,
l, r, i
376#ifdef MFC_SINGLE_PRECISION
377 call sgeev(ivl, ivr, nm, vgt, nm, lr, li, vl, nm, vr, nm, work, lwork, info)
379 call dgeev(ivl, ivr, nm, vgt, nm, lr, li, vl, nm, vr, nm, work, lwork, info)
385 if (abs(li(r)) < abs(li(idx)))
then
392 eigvec_mag = sqrt(eigvec(1)**2._wp + eigvec(2)**2._wp + eigvec(3)**2._wp)
393 if (eigvec_mag > sgm_eps)
then
394 eigvec = eigvec/eigvec_mag
400 omega_proj = (vgt(3, 2) - vgt(2, 3))*eigvec(1) + (vgt(1, 3) - vgt(3, 1))*eigvec(2) + (vgt(2, 1) - vgt(1, &
404 if (omega_proj < 0._wp)
then
406 omega_proj = -omega_proj
410 lci = li(mod(idx, 3) + 1)
413 alpha = omega_proj**2._wp - 4._wp*lci**2._wp
415 if (alpha > 0._wp)
then
416 liutex_mag(
j,
k,
l) = omega_proj - sqrt(alpha)
418 liutex_mag(
j,
k,
l) = omega_proj
422 liutex_axis(
j,
k,
l, 1) = eigvec(1)
423 liutex_axis(
j,
k,
l, 2) = eigvec(2)
424 liutex_axis(
j,
k,
l, 3) = eigvec(3)
437 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
438 & intent(inout) :: q_sf
440 real(wp) :: drho_dx, drho_dy, drho_dz
441 real(wp),
dimension(2) :: gm_rho_max
442 integer :: i,
j,
k,
l
484 if (
num_procs > 1)
call s_mpi_reduce_maxloc(gm_rho_max)
497 q_sf(
j,
k,
l) = 0._wp
501 &
k,
l)/gm_rho_max(1)
type(scalar_field), dimension(sys_size), intent(inout) q_cons_vf
Shared derived types for field data, patch geometry, bubble dynamics, and MPI I/O structures.
Computes derived flow quantities (sound speed, vorticity, Schlieren, etc.) from conservative and prim...
real(wp), dimension(:,:,:), allocatable gm_rho_sf
Density gradient magnitude for numerical Schlieren.
impure subroutine, public s_derive_liutex(q_prim_vf, liutex_mag, liutex_axis)
Compute the Liutex vector and its magnitude based on Xu et al. (2019).
subroutine, public s_derive_specific_heat_ratio(q_sf)
Derive the specific heat ratio from the specific heat ratio function gamma_sf. The latter is stored i...
subroutine, public s_derive_liquid_stiffness(q_sf)
Compute the liquid stiffness from the specific heat ratio function gamma_sf and the liquid stiffness ...
real(wp), dimension(:,:), allocatable, public fd_coeff_z
real(wp), dimension(:,:), allocatable, public fd_coeff_x
impure subroutine, public s_initialize_derived_variables_module
Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the...
real(wp), dimension(:,:), allocatable, public fd_coeff_y
subroutine, public s_derive_vorticity_component(i, q_prim_vf, q_sf)
Compute the specified component of the vorticity from the primitive variables. From those inputs,...
subroutine, public s_derive_sound_speed(q_prim_vf, q_sf)
Compute the speed of sound from the primitive variables, density, specific heat ratio function,...
impure subroutine, public s_derive_numerical_schlieren_function(q_cons_vf, q_sf)
Compute the values of the numerical Schlieren function, which are subsequently stored in the derived ...
subroutine, public s_derive_flux_limiter(i, q_prim_vf, q_sf)
Derive the flux limiter at cell boundary i+1/2. This is an approximation because the velocity used to...
impure subroutine, public s_finalize_derived_variables_module
Deallocation procedures for the module.
subroutine, public s_derive_qm(q_prim_vf, q_sf)
Compute the Q_M criterion from the primitive variables. The Q_M function, which are subsequently stor...
Global parameters for the post-process: domain geometry, equation of state, and output database setti...
type(int_bounds_info) offset_y
real(wp), dimension(num_fluids_max) schlieren_alpha
Per-fluid Schlieren intensity amplitude coefficients.
real(wp), dimension(:), allocatable y_cc
integer proc_rank
Rank of the local processor.
logical mixture_err
Mixture error limiter.
logical alt_soundspeed
Alternate sound speed.
integer fd_number
Finite-difference half-stencil size: MAX(1, fd_order/2).
integer model_eqns
Multicomponent flow model.
type(int_bounds_info) offset_x
logical, dimension(3) omega_wrt
integer num_procs
Number of processors.
type(int_bounds_info) offset_z
type(eqn_idx_info) eqn_idx
All conserved-variable equation index ranges and scalars.
Basic floating-point utilities: approximate equality, default detection, and coordinate bounds.
logical elemental function, public f_approx_equal(a, b, tol_input)
Check if two floating point numbers of wp are within tolerance.
MPI gather and scatter operations for distributing post-process grid and flow-variable data.
Conservative-to-primitive variable conversion, mixture property evaluation, and pressure computation.
real(wp), dimension(:), allocatable, public gammas
subroutine s_compute_speed_of_sound(pres, rho, gamma, pi_inf, h, adv, vel_sum, c_c, c, qv)
Compute the speed of sound from thermodynamic state variables, supporting multiple equation-of-state ...
real(wp), dimension(:,:,:), allocatable, public pi_inf_sf
Scalar liquid stiffness function.
real(wp), dimension(:,:,:), allocatable, public gamma_sf
Scalar sp. heat ratio function.
real(wp), dimension(:,:,:), allocatable, public rho_sf
Scalar density function.
real(wp), dimension(:), allocatable, public pi_infs
Derived type annexing a scalar field (SF).