MFC
Exascale flow solver
Loading...
Searching...
No Matches
m_derived_variables.fpp.f90
Go to the documentation of this file.
1# 1 "/home/runner/work/MFC/MFC/src/post_process/m_derived_variables.fpp"
2!>
3!! @file
4!! @brief Contains module m_derived_variables
5
6!> @brief Computes derived flow quantities (sound speed, vorticity, Schlieren, etc.) from conservative and primitive variables
7
9
12 use m_mpi_proxy
15
16 implicit none
17
21
22 real(wp), allocatable, dimension(:,:,:) :: gm_rho_sf !< Density gradient magnitude for numerical Schlieren
23 !> @name Finite-difference (fd) coefficients in x-, y- and z-coordinate directions. Note that because sufficient boundary
24 !! information is available for all the active coordinate directions, the centered family of the finite-difference schemes is
25 !! used.
26 !> @{
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
30 !> @}
31
32contains
33
34 !> Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the module
36
37 ! Allocate density gradient magnitude if Schlieren output requested
38 if (schlieren_wrt) then
39 allocate (gm_rho_sf(-offset_x%beg:m + offset_x%end,-offset_y%beg:n + offset_y%end,-offset_z%beg:p + offset_z%end))
40 end if
41
42 ! Allocate FD coefficients (up to 4th order; higher orders need extension)
43
44 if (omega_wrt(2) .or. omega_wrt(3) .or. schlieren_wrt .or. liutex_wrt) then
45 allocate (fd_coeff_x(-fd_number:fd_number,-offset_x%beg:m + offset_x%end))
46 end if
47
48 if (omega_wrt(1) .or. omega_wrt(3) .or. liutex_wrt .or. (n > 0 .and. schlieren_wrt)) then
49 allocate (fd_coeff_y(-fd_number:fd_number,-offset_y%beg:n + offset_y%end))
50 end if
51
52 if (omega_wrt(1) .or. omega_wrt(2) .or. liutex_wrt .or. (p > 0 .and. schlieren_wrt)) then
53 allocate (fd_coeff_z(-fd_number:fd_number,-offset_z%beg:p + offset_z%end))
54 end if
55
57
58 !> Derive the specific heat ratio from the specific heat ratio function gamma_sf. The latter is stored in the derived flow
59 !! quantity storage variable, q_sf.
61
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
64
65 integer :: i, j, k
66 do k = -offset_z%beg, p + offset_z%end
67 do j = -offset_y%beg, n + offset_y%end
68 do i = -offset_x%beg, m + offset_x%end
69 q_sf(i, j, k) = 1._wp + 1._wp/gamma_sf(i, j, k)
70 end do
71 end do
72 end do
73
74 end subroutine s_derive_specific_heat_ratio
75
76 !> Compute the liquid stiffness from the specific heat ratio function gamma_sf and the liquid stiffness function pi_inf_sf,
77 !! respectively. These are used to calculate the values of the liquid stiffness, which are stored in the derived flow quantity
78 !! storage variable, q_sf.
80
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
83
84 integer :: i, j, k
85 do k = -offset_z%beg, p + offset_z%end
86 do j = -offset_y%beg, n + offset_y%end
87 do i = -offset_x%beg, m + offset_x%end
88 q_sf(i, j, k) = pi_inf_sf(i, j, k)/(gamma_sf(i, j, k) + 1._wp)
89 end do
90 end do
91 end do
92
93 end subroutine s_derive_liquid_stiffness
94
95 !> Compute the speed of sound from the primitive variables, density, specific heat ratio function, and liquid stiffness
96 !! function. It then computes from those variables the values of the speed of sound, which are stored in the derived flow
97 !! quantity storage variable, q_sf.
98 subroutine s_derive_sound_speed(q_prim_vf, q_sf)
99
100 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
101
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
104
105 integer :: i, j, k
106 real(wp) :: blkmod1, blkmod2
107
108 do k = -offset_z%beg, p + offset_z%end
109 do j = -offset_y%beg, n + offset_y%end
110 do i = -offset_x%beg, m + offset_x%end
111 if (alt_soundspeed .neqv. .true.) then
112 q_sf(i, j, k) = (((gamma_sf(i, j, k) + 1._wp)*q_prim_vf(eqn_idx%E)%sf(i, j, k) + pi_inf_sf(i, j, &
113 & k))/(gamma_sf(i, j, k)*rho_sf(i, j, k)))
114 else
115 blkmod1 = ((gammas(1) + 1._wp)*q_prim_vf(eqn_idx%E)%sf(i, j, k) + pi_infs(1))/gammas(1)
116 blkmod2 = ((gammas(2) + 1._wp)*q_prim_vf(eqn_idx%E)%sf(i, j, k) + pi_infs(2))/gammas(2)
117 q_sf(i, j, k) = (1._wp/(rho_sf(i, j, k)*(q_prim_vf(eqn_idx%adv%beg)%sf(i, j, &
118 & k)/blkmod1 + (1._wp - q_prim_vf(eqn_idx%adv%beg)%sf(i, j, k))/blkmod2)))
119 end if
120
121 if (mixture_err .and. q_sf(i, j, k) < 0._wp) then
122 q_sf(i, j, k) = 1.e-16_wp
123 else
124 q_sf(i, j, k) = sqrt(q_sf(i, j, k))
125 end if
126 end do
127 end do
128 end do
129
130 end subroutine s_derive_sound_speed
131
132 !> Derive the flux limiter at cell boundary i+1/2. This is an approximation because the velocity used to determine the upwind
133 !! direction is the velocity at the cell center i instead of the contact velocity at the cell boundary from the Riemann solver.
134 subroutine s_derive_flux_limiter(i, q_prim_vf, q_sf)
135
136 integer, intent(in) :: i
137 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
138
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
141
142 real(wp) :: top, bottom, slope
143 integer :: j, k, l
144 do l = -offset_z%beg, p + offset_z%end
145 do k = -offset_y%beg, n + offset_y%end
146 do j = -offset_x%beg, m + offset_x%end
147 if (i == 1) then
148 if (q_prim_vf(eqn_idx%cont%end + i)%sf(j, k, l) >= 0._wp) then
149 top = q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l) - q_prim_vf(eqn_idx%adv%beg)%sf(j - 1, k, l)
150 bottom = q_prim_vf(eqn_idx%adv%beg)%sf(j + 1, k, l) - q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l)
151 else
152 top = q_prim_vf(eqn_idx%adv%beg)%sf(j + 2, k, l) - q_prim_vf(eqn_idx%adv%beg)%sf(j + 1, k, l)
153 bottom = q_prim_vf(eqn_idx%adv%beg)%sf(j + 1, k, l) - q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l)
154 end if
155 else if (i == 2) then
156 if (q_prim_vf(eqn_idx%cont%end + i)%sf(j, k, l) >= 0._wp) then
157 top = q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l) - q_prim_vf(eqn_idx%adv%beg)%sf(j, k - 1, l)
158 bottom = q_prim_vf(eqn_idx%adv%beg)%sf(j, k + 1, l) - q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l)
159 else
160 top = q_prim_vf(eqn_idx%adv%beg)%sf(j, k + 2, l) - q_prim_vf(eqn_idx%adv%beg)%sf(j, k + 1, l)
161 bottom = q_prim_vf(eqn_idx%adv%beg)%sf(j, k + 1, l) - q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l)
162 end if
163 else
164 if (q_prim_vf(eqn_idx%cont%end + i)%sf(j, k, l) >= 0._wp) then
165 top = q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l) - q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l - 1)
166 bottom = q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l + 1) - q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l)
167 else
168 top = q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l + 2) - q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l + 1)
169 bottom = q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l + 1) - q_prim_vf(eqn_idx%adv%beg)%sf(j, k, l)
170 end if
171 end if
172
173 if (abs(top) < 1.e-8_wp) top = 0._wp
174 if (abs(bottom) < 1.e-8_wp) bottom = 0._wp
175
176 if (f_approx_equal(top, bottom)) then
177 slope = 1._wp
178 else
179 slope = (top*bottom)/(bottom**2._wp + 1.e-16_wp)
180 end if
181
182 if (flux_lim == 1) then ! MINMOD (MM)
183 q_sf(j, k, l) = max(0._wp, min(1._wp, slope))
184 else if (flux_lim == 2) then ! MUSCL (MC)
185 q_sf(j, k, l) = max(0._wp, min(2._wp*slope, 5.e-1_wp*(1._wp + slope), 2._wp))
186 else if (flux_lim == 3) then ! OSPRE (OP)
187 q_sf(j, k, l) = (15.e-1_wp*(slope**2._wp + slope))/(slope**2._wp + slope + 1._wp)
188 else if (flux_lim == 4) then ! SUPERBEE (SB)
189 q_sf(j, k, l) = max(0._wp, min(1._wp, 2._wp*slope), min(slope, 2._wp))
190 else if (flux_lim == 5) then ! SWEBY (SW) (beta = 1.5)
191 q_sf(j, k, l) = max(0._wp, min(15.e-1_wp*slope, 1._wp), min(slope, 15.e-1_wp))
192 else if (flux_lim == 6) then ! VAN ALBADA (VA)
193 q_sf(j, k, l) = (slope**2._wp + slope)/(slope**2._wp + 1._wp)
194 else if (flux_lim == 7) then ! VAN LEER (VL)
195 q_sf(j, k, l) = (abs(slope) + slope)/(1._wp + abs(slope))
196 end if
197 end do
198 end do
199 end do
200
201 end subroutine s_derive_flux_limiter
202
203 !> Compute the specified component of the vorticity from the primitive variables. From those inputs, it proceeds to calculate
204 !! values of the desired vorticity component, which are subsequently stored in derived flow quantity storage variable, q_sf.
205 subroutine s_derive_vorticity_component(i, q_prim_vf, q_sf)
206
207 integer, intent(in) :: i
208 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
209
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
212
213 integer :: j, k, l, r
214 if (i == 1) then
215 do l = -offset_z%beg, p + offset_z%end
216 do k = -offset_y%beg, n + offset_y%end
217 do j = -offset_x%beg, m + offset_x%end
218 q_sf(j, k, l) = 0._wp
219
220 do r = -fd_number, fd_number
221 if (grid_geometry == 3) then
222 q_sf(j, k, l) = q_sf(j, k, l) + 1._wp/y_cc(k)*(fd_coeff_y(r, &
223 & k)*y_cc(r + k)*q_prim_vf(eqn_idx%mom%end)%sf(j, r + k, l) - fd_coeff_z(r, &
224 & l)*q_prim_vf(eqn_idx%mom%beg + 1)%sf(j, k, r + l))
225 else
226 q_sf(j, k, l) = q_sf(j, k, l) + fd_coeff_y(r, k)*q_prim_vf(eqn_idx%mom%end)%sf(j, r + k, &
227 & l) - fd_coeff_z(r, l)*q_prim_vf(eqn_idx%mom%beg + 1)%sf(j, k, r + l)
228 end if
229 end do
230 end do
231 end do
232 end do
233 else if (i == 2) then
234 do l = -offset_z%beg, p + offset_z%end
235 do k = -offset_y%beg, n + offset_y%end
236 do j = -offset_x%beg, m + offset_x%end
237 q_sf(j, k, l) = 0._wp
238
239 do r = -fd_number, fd_number
240 if (grid_geometry == 3) then
241 q_sf(j, k, l) = q_sf(j, k, l) + fd_coeff_z(r, l)/y_cc(k)*q_prim_vf(eqn_idx%mom%beg)%sf(j, k, &
242 & r + l) - fd_coeff_x(r, j)*q_prim_vf(eqn_idx%mom%end)%sf(r + j, k, l)
243 else
244 q_sf(j, k, l) = q_sf(j, k, l) + fd_coeff_z(r, l)*q_prim_vf(eqn_idx%mom%beg)%sf(j, k, &
245 & r + l) - fd_coeff_x(r, j)*q_prim_vf(eqn_idx%mom%end)%sf(r + j, k, l)
246 end if
247 end do
248 end do
249 end do
250 end do
251 else
252 do l = -offset_z%beg, p + offset_z%end
253 do k = -offset_y%beg, n + offset_y%end
254 do j = -offset_x%beg, m + offset_x%end
255 q_sf(j, k, l) = 0._wp
256
257 do r = -fd_number, fd_number
258 q_sf(j, k, l) = q_sf(j, k, l) + fd_coeff_x(r, j)*q_prim_vf(eqn_idx%mom%beg + 1)%sf(r + j, k, &
259 & l) - fd_coeff_y(r, k)*q_prim_vf(eqn_idx%mom%beg)%sf(j, r + k, l)
260 end do
261 end do
262 end do
263 end do
264 end if
265
266 end subroutine s_derive_vorticity_component
267
268 !> Compute the Q_M criterion from the primitive variables. The Q_M function, which are subsequently stored in the derived flow
269 !! quantity storage variable, q_sf.
270 subroutine s_derive_qm(q_prim_vf, q_sf)
271
272 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
273
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
276
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
280 do l = -offset_z%beg, p + offset_z%end
281 do k = -offset_y%beg, n + offset_y%end
282 do j = -offset_x%beg, m + offset_x%end
283 ! Get velocity gradient tensor
284 q_jacobian_sf(:,:) = 0._wp
285
286 do r = -fd_number, fd_number
287 do jj = 1, 3
288 ! d()/dx
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)
291 ! d()/dy
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)
294 ! d()/dz
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)
297 end do
298 end do
299
300 ! Decompose velocity gradient into symmetric strain-rate S and skew-symmetric rotation-rate O
301 do jj = 1, 3
302 do kk = 1, 3
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))
305 end do
306 end do
307
308 do jj = 1, 3
309 do kk = 1, 3
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)
312 end do
313 end do
314
315 ! Q-criterion: Q = (||O||^2 - ||S||^2)/2, Hunt et al. CTR (1988)
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)
318 ! Second invariant of strain-rate tensor
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
321 end do
322 end do
323 end do
324
325 end subroutine s_derive_qm
326
327 !> Compute the Liutex vector and its magnitude based on Xu et al. (2019).
328 impure subroutine s_derive_liutex(q_prim_vf, liutex_mag, liutex_axis)
329
330 ! Liutex vortex identification via real eigenvector of velocity gradient, Xu et al. PoF (2019)
331
332 integer, parameter :: nm = 3
333 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
334
335 !> Liutex magnitude
336
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
339 !> Liutex rigid rotation axis
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' !< compute left eigenvectors
343 character, parameter :: ivr = 'V' !< compute right eigenvectors
344 real(wp), dimension(nm, nm) :: vgt !< velocity gradient tensor
345 real(wp), dimension(nm) :: lr, li !< real and imaginary parts of eigenvalues
346 real(wp), dimension(nm, nm) :: vl, vr !< left and right eigenvectors
347 integer, parameter :: lwork = 4*nm !< size of work array (4*nm recommended)
348 real(wp), dimension(lwork) :: work !< work array
349 integer :: info
350 real(wp), dimension(nm) :: eigvec !< real eigenvector
351 real(wp) :: eigvec_mag !< magnitude of real eigenvector
352 real(wp) :: omega_proj !< projection of vorticity on real eigenvector
353 real(wp) :: lci !< imaginary part of complex eigenvalue
354 real(wp) :: alpha
355 integer :: j, k, l, r, i
356 integer :: idx
357
358 do l = -offset_z%beg, p + offset_z%end
359 do k = -offset_y%beg, n + offset_y%end
360 do j = -offset_x%beg, m + offset_x%end
361 ! Get velocity gradient tensor (VGT)
362 vgt(:,:) = 0._wp
363
364 do r = -fd_number, fd_number
365 do i = 1, 3
366 ! d()/dx
367 vgt(i, 1) = vgt(i, 1) + fd_coeff_x(r, j)*q_prim_vf(eqn_idx%mom%beg + i - 1)%sf(r + j, k, l)
368 ! d()/dy
369 vgt(i, 2) = vgt(i, 2) + fd_coeff_y(r, k)*q_prim_vf(eqn_idx%mom%beg + i - 1)%sf(j, r + k, l)
370 ! d()/dz
371 vgt(i, 3) = vgt(i, 3) + fd_coeff_z(r, l)*q_prim_vf(eqn_idx%mom%beg + i - 1)%sf(j, k, r + l)
372 end do
373 end do
374
375 ! Call appropriate LAPACK routine based on precision
376#ifdef MFC_SINGLE_PRECISION
377 call sgeev(ivl, ivr, nm, vgt, nm, lr, li, vl, nm, vr, nm, work, lwork, info)
378#else
379 call dgeev(ivl, ivr, nm, vgt, nm, lr, li, vl, nm, vr, nm, work, lwork, info)
380#endif
381
382 ! Find eigenvector with smallest imaginary eigenvalue (real eigenvector of VGT)
383 idx = 1
384 do r = 2, 3
385 if (abs(li(r)) < abs(li(idx))) then
386 idx = r
387 end if
388 end do
389 eigvec = vr(:,idx)
390
391 ! Normalize real eigenvector if it is effectively non-zero
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
395 else
396 eigvec = 0._wp
397 end if
398
399 ! Compute vorticity projected on the eigenvector
400 omega_proj = (vgt(3, 2) - vgt(2, 3))*eigvec(1) + (vgt(1, 3) - vgt(3, 1))*eigvec(2) + (vgt(2, 1) - vgt(1, &
401 & 2))*eigvec(3)
402
403 ! As eigenvector can have +/- signs, we can choose the sign so that omega_proj is positive
404 if (omega_proj < 0._wp) then
405 eigvec = -eigvec
406 omega_proj = -omega_proj
407 end if
408
409 ! Imaginary eigenvalue of the complex conjugate pair (cyclic index selection)
410 lci = li(mod(idx, 3) + 1)
411
412 ! Discriminant: determines whether rotation dominates strain
413 alpha = omega_proj**2._wp - 4._wp*lci**2._wp
414 ! Liutex magnitude = omega_proj - sqrt(discriminant) when rotation dominates
415 if (alpha > 0._wp) then
416 liutex_mag(j, k, l) = omega_proj - sqrt(alpha)
417 else
418 liutex_mag(j, k, l) = omega_proj
419 end if
420
421 ! Compute Liutex axis
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)
425 end do
426 end do
427 end do
428
429 end subroutine s_derive_liutex
430
431 !> Compute the values of the numerical Schlieren function, which are subsequently stored in the derived flow quantity storage
432 !! variable, q_sf.
433 impure subroutine s_derive_numerical_schlieren_function(q_cons_vf, q_sf)
434
435 type(scalar_field), dimension(sys_size), intent(in) :: q_cons_vf
436
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
439
440 real(wp) :: drho_dx, drho_dy, drho_dz !< Spatial derivatives of the density in the x-, y- and z-directions
441 real(wp), dimension(2) :: gm_rho_max !< Global (max gradient magnitude, rank) pair for density
442 integer :: i, j, k, l
443
444 do l = -offset_z%beg, p + offset_z%end
445 do k = -offset_y%beg, n + offset_y%end
446 do j = -offset_x%beg, m + offset_x%end
447 drho_dx = 0._wp
448 drho_dy = 0._wp
449
450 do i = -fd_number, fd_number
451 drho_dx = drho_dx + fd_coeff_x(i, j)*rho_sf(i + j, k, l)
452 drho_dy = drho_dy + fd_coeff_y(i, k)*rho_sf(j, i + k, l)
453 end do
454
455 gm_rho_sf(j, k, l) = drho_dx*drho_dx + drho_dy*drho_dy
456 end do
457 end do
458 end do
459
460 if (p > 0) then
461 do l = -offset_z%beg, p + offset_z%end
462 do k = -offset_y%beg, n + offset_y%end
463 do j = -offset_x%beg, m + offset_x%end
464 drho_dz = 0._wp
465
466 do i = -fd_number, fd_number
467 if (grid_geometry == 3) then
468 drho_dz = drho_dz + fd_coeff_z(i, l)/y_cc(k)*rho_sf(j, k, i + l)
469 else
470 drho_dz = drho_dz + fd_coeff_z(i, l)*rho_sf(j, k, i + l)
471 end if
472 end do
473
474 gm_rho_sf(j, k, l) = gm_rho_sf(j, k, l) + drho_dz*drho_dz
475 end do
476 end do
477 end do
478 end if
479
480 gm_rho_sf = sqrt(gm_rho_sf)
481
482 gm_rho_max = (/maxval(gm_rho_sf), real(proc_rank, wp)/)
483
484 if (num_procs > 1) call s_mpi_reduce_maxloc(gm_rho_max)
485
486 ! The form of the numerical Schlieren function depends on the choice of the multicomponent flow model. For the gamma/pi_inf
487 ! model, the exponential of the negative, normalized, gradient magnitude of the density is computed. For the volume fraction
488 ! model, the amplitude of the exponential's inside is also modulated with respect to the identity of the fluid in which the
489 ! function is evaluated. For more information, refer to Marquina and Mulet (2003).
490
491 if (model_eqns == 1) then ! Gamma/pi_inf model
492 q_sf = -gm_rho_sf/gm_rho_max(1)
493 else ! Volume fraction model
494 do l = -offset_z%beg, p + offset_z%end
495 do k = -offset_y%beg, n + offset_y%end
496 do j = -offset_x%beg, m + offset_x%end
497 q_sf(j, k, l) = 0._wp
498
499 do i = 1, eqn_idx%adv%end - eqn_idx%E
500 q_sf(j, k, l) = q_sf(j, k, l) - schlieren_alpha(i)*q_cons_vf(i + eqn_idx%E)%sf(j, k, l)*gm_rho_sf(j, &
501 & k, l)/gm_rho_max(1)
502 end do
503 end do
504 end do
505 end do
506 end if
507
508 ! Up until now, only the inside of the exponential of the numerical Schlieren function has been evaluated and stored. Then,
509 ! to finish the computation, the exponential of the inside quantity is taken.
510 q_sf = exp(q_sf)
511
513
514 !> Deallocation procedures for the module
516
517 ! Deallocating the variable containing the gradient magnitude of the density field provided that the numerical Schlieren
518 ! function was was outputted during the post-process
519 if (schlieren_wrt) deallocate (gm_rho_sf)
520
521 ! Deallocating the variables that might have been used to bookkeep the finite-difference coefficients in the x-, y- and
522 ! z-directions
523 if (allocated(fd_coeff_x)) deallocate (fd_coeff_x)
524 if (allocated(fd_coeff_y)) deallocate (fd_coeff_y)
525 if (allocated(fd_coeff_z)) deallocate (fd_coeff_z)
526
528
529end module m_derived_variables
type(scalar_field), dimension(sys_size), intent(inout) q_cons_vf
integer, intent(in) k
integer, intent(in) j
integer, intent(in) l
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).