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
32 integer,
private ::
flg
70 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
71 & intent(inout) :: q_sf
89 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
90 & intent(inout) :: q_sf
108 type(
scalar_field),
dimension(sys_size),
intent(in) :: q_prim_vf
110 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
111 & intent(inout) :: q_sf
114 real(wp) :: blkmod1, blkmod2
126 &
k)/blkmod1 + (1._wp - q_prim_vf(
adv_idx%beg)%sf(i,
j,
k))/blkmod2)))
130 q_sf(i,
j,
k) = 1.e-16_wp
132 q_sf(i,
j,
k) = sqrt(q_sf(i,
j,
k))
144 integer,
intent(in) :: i
145 type(
scalar_field),
dimension(sys_size),
intent(in) :: q_prim_vf
147 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
148 & intent(inout) :: q_sf
150 real(wp) :: top, bottom, slope
156 if (q_prim_vf(
cont_idx%end + i)%sf(
j,
k,
l) >= 0._wp)
then
163 else if (i == 2)
then
164 if (q_prim_vf(
cont_idx%end + i)%sf(
j,
k,
l) >= 0._wp)
then
172 if (q_prim_vf(
cont_idx%end + i)%sf(
j,
k,
l) >= 0._wp)
then
181 if (abs(top) < 1.e-8_wp) top = 0._wp
182 if (abs(bottom) < 1.e-8_wp) bottom = 0._wp
187 slope = (top*bottom)/(bottom**2._wp + 1.e-16_wp)
191 q_sf(
j,
k,
l) = max(0._wp, min(1._wp, slope))
193 q_sf(
j,
k,
l) = max(0._wp, min(2._wp*slope, 5.e-1_wp*(1._wp + slope), 2._wp))
195 q_sf(
j,
k,
l) = (15.e-1_wp*(slope**2._wp + slope))/(slope**2._wp + slope + 1._wp)
197 q_sf(
j,
k,
l) = max(0._wp, min(1._wp, 2._wp*slope), min(slope, 2._wp))
199 q_sf(
j,
k,
l) = max(0._wp, min(15.e-1_wp*slope, 1._wp), min(slope, 15.e-1_wp))
201 q_sf(
j,
k,
l) = (slope**2._wp + slope)/(slope**2._wp + 1._wp)
203 q_sf(
j,
k,
l) = (abs(slope) + slope)/(1._wp + abs(slope))
214 integer,
intent(in) :: ndim
215 real(wp),
dimension(ndim, ndim),
intent(inout) :: A
216 real(wp),
dimension(ndim),
intent(inout) :: b
217 real(wp),
dimension(ndim),
intent(out) :: sol
223 j = i - 1 + maxloc(abs(a(i:ndim,i)), 1)
231 a(i,:) = a(i,:)/a(i, i)
233 b(k) = b(k) - a(k, i)*b(i)
234 a(k,:) = a(k,:) - a(k, i)*a(i,:)
241 sol(i) = sol(i) - a(i, k)*sol(k)
251 integer,
intent(in) :: i
252 type(
scalar_field),
dimension(sys_size),
intent(in) :: q_prim_vf
254 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
255 & intent(inout) :: q_sf
257 integer ::
j,
k,
l, r
262 q_sf(
j,
k,
l) = 0._wp
277 else if (i == 2)
then
281 q_sf(
j,
k,
l) = 0._wp
299 q_sf(
j,
k,
l) = 0._wp
316 type(
scalar_field),
dimension(sys_size),
intent(in) :: q_prim_vf
318 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
319 & intent(inout) :: q_sf
321 real(wp),
dimension(1:3,1:3) :: q_jacobian_sf, s, s2, o, o2
322 real(wp) :: trs, q, iis
323 integer ::
j,
k,
l, r, jj, kk
328 q_jacobian_sf(:,:) = 0._wp
333 q_jacobian_sf(jj, 1) = q_jacobian_sf(jj, 1) +
fd_coeff_x(r, &
334 &
j)*q_prim_vf(
mom_idx%beg + jj - 1)%sf(r +
j,
k,
l)
336 q_jacobian_sf(jj, 2) = q_jacobian_sf(jj, 2) +
fd_coeff_y(r,
k)*q_prim_vf(
mom_idx%beg + jj - 1)%sf(
j, &
339 q_jacobian_sf(jj, 3) = q_jacobian_sf(jj, 3) +
fd_coeff_z(r,
l)*q_prim_vf(
mom_idx%beg + jj - 1)%sf(
j, &
347 s(jj, kk) = 0.5_wp*(q_jacobian_sf(jj, kk) + q_jacobian_sf(kk, jj))
348 o(jj, kk) = 0.5_wp*(q_jacobian_sf(jj, kk) - q_jacobian_sf(kk, jj))
354 o2(jj, kk) = o(jj, 1)*o(kk, 1) + o(jj, 2)*o(kk, 2) + o(jj, 3)*o(kk, 3)
355 s2(jj, kk) = s(jj, 1)*s(kk, 1) + s(jj, 2)*s(kk, 2) + s(jj, 3)*s(kk, 3)
360 q = 0.5_wp*((o2(1, 1) + o2(2, 2) + o2(3, 3)) - (s2(1, 1) + s2(2, 2) + s2(3, 3)))
361 trs = s(1, 1) + s(2, 2) + s(3, 3)
363 iis = 0.5_wp*((s(1, 1) + s(2, 2) + s(3, 3))**2 - (s2(1, 1) + s2(2, 2) + s2(3, 3)))
364 q_sf(
j,
k,
l) = q + iis
376 integer,
parameter :: nm = 3
377 type(
scalar_field),
dimension(sys_size),
intent(in) :: q_prim_vf
381 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
382 & intent(out) :: liutex_mag
384 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), &
385 & intent(out) :: liutex_axis
386 character,
parameter :: ivl =
'N'
387 character,
parameter :: ivr =
'V'
388 real(wp),
dimension(nm, nm) :: vgt
389 real(wp),
dimension(nm) :: lr, li
390 real(wp),
dimension(nm, nm) :: vl, vr
391 integer,
parameter :: lwork = 4*nm
392 real(wp),
dimension(lwork) :: work
394 real(wp),
dimension(nm) :: eigvec
395 real(wp) :: eigvec_mag
396 real(wp) :: omega_proj
399 integer ::
j,
k,
l, r, i
420#ifdef MFC_SINGLE_PRECISION
421 call sgeev(ivl, ivr, nm, vgt, nm, lr, li, vl, nm, vr, nm, work, lwork, info)
423 call dgeev(ivl, ivr, nm, vgt, nm, lr, li, vl, nm, vr, nm, work, lwork, info)
429 if (abs(li(r)) < abs(li(idx)))
then
436 eigvec_mag = sqrt(eigvec(1)**2._wp + eigvec(2)**2._wp + eigvec(3)**2._wp)
437 if (eigvec_mag > sgm_eps)
then
438 eigvec = eigvec/eigvec_mag
444 omega_proj = (vgt(3, 2) - vgt(2, 3))*eigvec(1) + (vgt(1, 3) - vgt(3, 1))*eigvec(2) + (vgt(2, 1) - vgt(1, &
448 if (omega_proj < 0._wp)
then
450 omega_proj = -omega_proj
454 lci = li(mod(idx, 3) + 1)
457 alpha = omega_proj**2._wp - 4._wp*lci**2._wp
459 if (alpha > 0._wp)
then
460 liutex_mag(
j,
k,
l) = omega_proj - sqrt(alpha)
462 liutex_mag(
j,
k,
l) = omega_proj
466 liutex_axis(
j,
k,
l, 1) = eigvec(1)
467 liutex_axis(
j,
k,
l, 2) = eigvec(2)
468 liutex_axis(
j,
k,
l, 3) = eigvec(3)
481 real(wp),
dimension(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end), &
482 & intent(inout) :: q_sf
484 real(wp) :: drho_dx, drho_dy, drho_dz
485 real(wp),
dimension(2) :: gm_rho_max
486 integer :: i,
j,
k,
l
528 if (
num_procs > 1)
call s_mpi_reduce_maxloc(gm_rho_max)
541 q_sf(
j,
k,
l) = 0._wp
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...
subroutine s_solve_linear_system(a, b, sol, ndim)
Solve Ax=b via Gaussian elimination with partial pivoting.
real(wp), dimension(:,:,:), allocatable gm_rho_sf
Density gradient magnitude for numerical Schlieren.
integer, private flg
Dimensionality flag: 1 = 3D dataset, 0 = otherwise.
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
type(int_bounds_info) mom_idx
Indexes of first & last momentum eqns.
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.
type(int_bounds_info) cont_idx
Indexes of first & last continuity eqns.
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) adv_idx
Indexes of first & last advection eqns.
type(int_bounds_info) offset_z
integer e_idx
Index of energy equation.
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).