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
32 integer, private :: flg !< Dimensionality flag: 1 = 3D dataset, 0 = otherwise
33
34contains
35
36 !> Computation of parameters, allocation procedures, and/or any other tasks needed to properly setup the module
38
39 ! Allocate density gradient magnitude if Schlieren output requested
40 if (schlieren_wrt) then
41 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))
42 end if
43
44 ! Allocate FD coefficients (up to 4th order; higher orders need extension)
45
46 if (omega_wrt(2) .or. omega_wrt(3) .or. schlieren_wrt .or. liutex_wrt) then
47 allocate (fd_coeff_x(-fd_number:fd_number,-offset_x%beg:m + offset_x%end))
48 end if
49
50 if (omega_wrt(1) .or. omega_wrt(3) .or. liutex_wrt .or. (n > 0 .and. schlieren_wrt)) then
51 allocate (fd_coeff_y(-fd_number:fd_number,-offset_y%beg:n + offset_y%end))
52 end if
53
54 if (omega_wrt(1) .or. omega_wrt(2) .or. liutex_wrt .or. (p > 0 .and. schlieren_wrt)) then
55 allocate (fd_coeff_z(-fd_number:fd_number,-offset_z%beg:p + offset_z%end))
56 end if
57
58 if (p > 0) then
59 flg = 1
60 else
61 flg = 0
62 end if
63
65
66 !> Derive the specific heat ratio from the specific heat ratio function gamma_sf. The latter is stored in the derived flow
67 !! quantity storage variable, q_sf.
69
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
72
73 integer :: i, j, k
74 do k = -offset_z%beg, p + offset_z%end
75 do j = -offset_y%beg, n + offset_y%end
76 do i = -offset_x%beg, m + offset_x%end
77 q_sf(i, j, k) = 1._wp + 1._wp/gamma_sf(i, j, k)
78 end do
79 end do
80 end do
81
82 end subroutine s_derive_specific_heat_ratio
83
84 !> Compute the liquid stiffness from the specific heat ratio function gamma_sf and the liquid stiffness function pi_inf_sf,
85 !! respectively. These are used to calculate the values of the liquid stiffness, which are stored in the derived flow quantity
86 !! storage variable, q_sf.
88
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
91
92 integer :: i, j, k
93 do k = -offset_z%beg, p + offset_z%end
94 do j = -offset_y%beg, n + offset_y%end
95 do i = -offset_x%beg, m + offset_x%end
96 q_sf(i, j, k) = pi_inf_sf(i, j, k)/(gamma_sf(i, j, k) + 1._wp)
97 end do
98 end do
99 end do
100
101 end subroutine s_derive_liquid_stiffness
102
103 !> Compute the speed of sound from the primitive variables, density, specific heat ratio function, and liquid stiffness
104 !! function. It then computes from those variables the values of the speed of sound, which are stored in the derived flow
105 !! quantity storage variable, q_sf.
106 subroutine s_derive_sound_speed(q_prim_vf, q_sf)
107
108 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
109
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
112
113 integer :: i, j, k
114 real(wp) :: blkmod1, blkmod2
115
116 do k = -offset_z%beg, p + offset_z%end
117 do j = -offset_y%beg, n + offset_y%end
118 do i = -offset_x%beg, m + offset_x%end
119 if (alt_soundspeed .neqv. .true.) then
120 q_sf(i, j, k) = (((gamma_sf(i, j, k) + 1._wp)*q_prim_vf(e_idx)%sf(i, j, k) + pi_inf_sf(i, j, &
121 & k))/(gamma_sf(i, j, k)*rho_sf(i, j, k)))
122 else
123 blkmod1 = ((gammas(1) + 1._wp)*q_prim_vf(e_idx)%sf(i, j, k) + pi_infs(1))/gammas(1)
124 blkmod2 = ((gammas(2) + 1._wp)*q_prim_vf(e_idx)%sf(i, j, k) + pi_infs(2))/gammas(2)
125 q_sf(i, j, k) = (1._wp/(rho_sf(i, j, k)*(q_prim_vf(adv_idx%beg)%sf(i, j, &
126 & k)/blkmod1 + (1._wp - q_prim_vf(adv_idx%beg)%sf(i, j, k))/blkmod2)))
127 end if
128
129 if (mixture_err .and. q_sf(i, j, k) < 0._wp) then
130 q_sf(i, j, k) = 1.e-16_wp
131 else
132 q_sf(i, j, k) = sqrt(q_sf(i, j, k))
133 end if
134 end do
135 end do
136 end do
137
138 end subroutine s_derive_sound_speed
139
140 !> Derive the flux limiter at cell boundary i+1/2. This is an approximation because the velocity used to determine the upwind
141 !! direction is the velocity at the cell center i instead of the contact velocity at the cell boundary from the Riemann solver.
142 subroutine s_derive_flux_limiter(i, q_prim_vf, q_sf)
143
144 integer, intent(in) :: i
145 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
146
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
149
150 real(wp) :: top, bottom, slope
151 integer :: j, k, l
152 do l = -offset_z%beg, p + offset_z%end
153 do k = -offset_y%beg, n + offset_y%end
154 do j = -offset_x%beg, m + offset_x%end
155 if (i == 1) then
156 if (q_prim_vf(cont_idx%end + i)%sf(j, k, l) >= 0._wp) then
157 top = q_prim_vf(adv_idx%beg)%sf(j, k, l) - q_prim_vf(adv_idx%beg)%sf(j - 1, k, l)
158 bottom = q_prim_vf(adv_idx%beg)%sf(j + 1, k, l) - q_prim_vf(adv_idx%beg)%sf(j, k, l)
159 else
160 top = q_prim_vf(adv_idx%beg)%sf(j + 2, k, l) - q_prim_vf(adv_idx%beg)%sf(j + 1, k, l)
161 bottom = q_prim_vf(adv_idx%beg)%sf(j + 1, k, l) - q_prim_vf(adv_idx%beg)%sf(j, k, l)
162 end if
163 else if (i == 2) then
164 if (q_prim_vf(cont_idx%end + i)%sf(j, k, l) >= 0._wp) then
165 top = q_prim_vf(adv_idx%beg)%sf(j, k, l) - q_prim_vf(adv_idx%beg)%sf(j, k - 1, l)
166 bottom = q_prim_vf(adv_idx%beg)%sf(j, k + 1, l) - q_prim_vf(adv_idx%beg)%sf(j, k, l)
167 else
168 top = q_prim_vf(adv_idx%beg)%sf(j, k + 2, l) - q_prim_vf(adv_idx%beg)%sf(j, k + 1, l)
169 bottom = q_prim_vf(adv_idx%beg)%sf(j, k + 1, l) - q_prim_vf(adv_idx%beg)%sf(j, k, l)
170 end if
171 else
172 if (q_prim_vf(cont_idx%end + i)%sf(j, k, l) >= 0._wp) then
173 top = q_prim_vf(adv_idx%beg)%sf(j, k, l) - q_prim_vf(adv_idx%beg)%sf(j, k, l - 1)
174 bottom = q_prim_vf(adv_idx%beg)%sf(j, k, l + 1) - q_prim_vf(adv_idx%beg)%sf(j, k, l)
175 else
176 top = q_prim_vf(adv_idx%beg)%sf(j, k, l + 2) - q_prim_vf(adv_idx%beg)%sf(j, k, l + 1)
177 bottom = q_prim_vf(adv_idx%beg)%sf(j, k, l + 1) - q_prim_vf(adv_idx%beg)%sf(j, k, l)
178 end if
179 end if
180
181 if (abs(top) < 1.e-8_wp) top = 0._wp
182 if (abs(bottom) < 1.e-8_wp) bottom = 0._wp
183
184 if (f_approx_equal(top, bottom)) then
185 slope = 1._wp
186 else
187 slope = (top*bottom)/(bottom**2._wp + 1.e-16_wp)
188 end if
189
190 if (flux_lim == 1) then ! MINMOD (MM)
191 q_sf(j, k, l) = max(0._wp, min(1._wp, slope))
192 else if (flux_lim == 2) then ! MUSCL (MC)
193 q_sf(j, k, l) = max(0._wp, min(2._wp*slope, 5.e-1_wp*(1._wp + slope), 2._wp))
194 else if (flux_lim == 3) then ! OSPRE (OP)
195 q_sf(j, k, l) = (15.e-1_wp*(slope**2._wp + slope))/(slope**2._wp + slope + 1._wp)
196 else if (flux_lim == 4) then ! SUPERBEE (SB)
197 q_sf(j, k, l) = max(0._wp, min(1._wp, 2._wp*slope), min(slope, 2._wp))
198 else if (flux_lim == 5) then ! SWEBY (SW) (beta = 1.5)
199 q_sf(j, k, l) = max(0._wp, min(15.e-1_wp*slope, 1._wp), min(slope, 15.e-1_wp))
200 else if (flux_lim == 6) then ! VAN ALBADA (VA)
201 q_sf(j, k, l) = (slope**2._wp + slope)/(slope**2._wp + 1._wp)
202 else if (flux_lim == 7) then ! VAN LEER (VL)
203 q_sf(j, k, l) = (abs(slope) + slope)/(1._wp + abs(slope))
204 end if
205 end do
206 end do
207 end do
208
209 end subroutine s_derive_flux_limiter
210
211 !> Solve Ax=b via Gaussian elimination with partial pivoting
212 subroutine s_solve_linear_system(A, b, sol, ndim)
213
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
218 integer :: i, j, k
219
220 ! Forward elimination with partial pivoting
221
222 do i = 1, ndim
223 j = i - 1 + maxloc(abs(a(i:ndim,i)), 1)
224 sol = a(i,:)
225 a(i,:) = a(j,:)
226 a(j,:) = sol
227 sol(1) = b(i)
228 b(i) = b(j)
229 b(j) = sol(1)
230 b(i) = b(i)/a(i, i)
231 a(i,:) = a(i,:)/a(i, i)
232 do k = i + 1, ndim
233 b(k) = b(k) - a(k, i)*b(i)
234 a(k,:) = a(k,:) - a(k, i)*a(i,:)
235 end do
236 end do
237
238 do i = ndim, 1, -1
239 sol(i) = b(i)
240 do k = i + 1, ndim
241 sol(i) = sol(i) - a(i, k)*sol(k)
242 end do
243 end do
244
245 end subroutine s_solve_linear_system
246
247 !> Compute the specified component of the vorticity from the primitive variables. From those inputs, it proceeds to calculate
248 !! values of the desired vorticity component, which are subsequently stored in derived flow quantity storage variable, q_sf.
249 subroutine s_derive_vorticity_component(i, q_prim_vf, q_sf)
250
251 integer, intent(in) :: i
252 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
253
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
256
257 integer :: j, k, l, r
258 if (i == 1) then
259 do l = -offset_z%beg, p + offset_z%end
260 do k = -offset_y%beg, n + offset_y%end
261 do j = -offset_x%beg, m + offset_x%end
262 q_sf(j, k, l) = 0._wp
263
264 do r = -fd_number, fd_number
265 if (grid_geometry == 3) then
266 q_sf(j, k, l) = q_sf(j, k, l) + 1._wp/y_cc(k)*(fd_coeff_y(r, &
267 & k)*y_cc(r + k)*q_prim_vf(mom_idx%end)%sf(j, r + k, l) - fd_coeff_z(r, &
268 & l)*q_prim_vf(mom_idx%beg + 1)%sf(j, k, r + l))
269 else
270 q_sf(j, k, l) = q_sf(j, k, l) + fd_coeff_y(r, k)*q_prim_vf(mom_idx%end)%sf(j, r + k, &
271 & l) - fd_coeff_z(r, l)*q_prim_vf(mom_idx%beg + 1)%sf(j, k, r + l)
272 end if
273 end do
274 end do
275 end do
276 end do
277 else if (i == 2) then
278 do l = -offset_z%beg, p + offset_z%end
279 do k = -offset_y%beg, n + offset_y%end
280 do j = -offset_x%beg, m + offset_x%end
281 q_sf(j, k, l) = 0._wp
282
283 do r = -fd_number, fd_number
284 if (grid_geometry == 3) then
285 q_sf(j, k, l) = q_sf(j, k, l) + fd_coeff_z(r, l)/y_cc(k)*q_prim_vf(mom_idx%beg)%sf(j, k, &
286 & r + l) - fd_coeff_x(r, j)*q_prim_vf(mom_idx%end)%sf(r + j, k, l)
287 else
288 q_sf(j, k, l) = q_sf(j, k, l) + fd_coeff_z(r, l)*q_prim_vf(mom_idx%beg)%sf(j, k, &
289 & r + l) - fd_coeff_x(r, j)*q_prim_vf(mom_idx%end)%sf(r + j, k, l)
290 end if
291 end do
292 end do
293 end do
294 end do
295 else
296 do l = -offset_z%beg, p + offset_z%end
297 do k = -offset_y%beg, n + offset_y%end
298 do j = -offset_x%beg, m + offset_x%end
299 q_sf(j, k, l) = 0._wp
300
301 do r = -fd_number, fd_number
302 q_sf(j, k, l) = q_sf(j, k, l) + fd_coeff_x(r, j)*q_prim_vf(mom_idx%beg + 1)%sf(r + j, k, &
303 & l) - fd_coeff_y(r, k)*q_prim_vf(mom_idx%beg)%sf(j, r + k, l)
304 end do
305 end do
306 end do
307 end do
308 end if
309
310 end subroutine s_derive_vorticity_component
311
312 !> Compute the Q_M criterion from the primitive variables. The Q_M function, which are subsequently stored in the derived flow
313 !! quantity storage variable, q_sf.
314 subroutine s_derive_qm(q_prim_vf, q_sf)
315
316 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
317
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
320
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
324 do l = -offset_z%beg, p + offset_z%end
325 do k = -offset_y%beg, n + offset_y%end
326 do j = -offset_x%beg, m + offset_x%end
327 ! Get velocity gradient tensor
328 q_jacobian_sf(:,:) = 0._wp
329
330 do r = -fd_number, fd_number
331 do jj = 1, 3
332 ! d()/dx
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)
335 ! d()/dy
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, &
337 & r + k, l)
338 ! d()/dz
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, &
340 & k, r + l)
341 end do
342 end do
343
344 ! Decompose velocity gradient into symmetric strain-rate S and skew-symmetric rotation-rate O
345 do jj = 1, 3
346 do kk = 1, 3
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))
349 end do
350 end do
351
352 do jj = 1, 3
353 do kk = 1, 3
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)
356 end do
357 end do
358
359 ! Q-criterion: Q = (||O||^2 - ||S||^2)/2, Hunt et al. CTR (1988)
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)
362 ! Second invariant of strain-rate tensor
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
365 end do
366 end do
367 end do
368
369 end subroutine s_derive_qm
370
371 !> Compute the Liutex vector and its magnitude based on Xu et al. (2019).
372 impure subroutine s_derive_liutex(q_prim_vf, liutex_mag, liutex_axis)
373
374 ! Liutex vortex identification via real eigenvector of velocity gradient, Xu et al. PoF (2019)
375
376 integer, parameter :: nm = 3
377 type(scalar_field), dimension(sys_size), intent(in) :: q_prim_vf
378
379 !> Liutex magnitude
380
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
383 !> Liutex rigid rotation axis
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' !< compute left eigenvectors
387 character, parameter :: ivr = 'V' !< compute right eigenvectors
388 real(wp), dimension(nm, nm) :: vgt !< velocity gradient tensor
389 real(wp), dimension(nm) :: lr, li !< real and imaginary parts of eigenvalues
390 real(wp), dimension(nm, nm) :: vl, vr !< left and right eigenvectors
391 integer, parameter :: lwork = 4*nm !< size of work array (4*nm recommended)
392 real(wp), dimension(lwork) :: work !< work array
393 integer :: info
394 real(wp), dimension(nm) :: eigvec !< real eigenvector
395 real(wp) :: eigvec_mag !< magnitude of real eigenvector
396 real(wp) :: omega_proj !< projection of vorticity on real eigenvector
397 real(wp) :: lci !< imaginary part of complex eigenvalue
398 real(wp) :: alpha
399 integer :: j, k, l, r, i
400 integer :: idx
401
402 do l = -offset_z%beg, p + offset_z%end
403 do k = -offset_y%beg, n + offset_y%end
404 do j = -offset_x%beg, m + offset_x%end
405 ! Get velocity gradient tensor (VGT)
406 vgt(:,:) = 0._wp
407
408 do r = -fd_number, fd_number
409 do i = 1, 3
410 ! d()/dx
411 vgt(i, 1) = vgt(i, 1) + fd_coeff_x(r, j)*q_prim_vf(mom_idx%beg + i - 1)%sf(r + j, k, l)
412 ! d()/dy
413 vgt(i, 2) = vgt(i, 2) + fd_coeff_y(r, k)*q_prim_vf(mom_idx%beg + i - 1)%sf(j, r + k, l)
414 ! d()/dz
415 vgt(i, 3) = vgt(i, 3) + fd_coeff_z(r, l)*q_prim_vf(mom_idx%beg + i - 1)%sf(j, k, r + l)
416 end do
417 end do
418
419 ! Call appropriate LAPACK routine based on precision
420#ifdef MFC_SINGLE_PRECISION
421 call sgeev(ivl, ivr, nm, vgt, nm, lr, li, vl, nm, vr, nm, work, lwork, info)
422#else
423 call dgeev(ivl, ivr, nm, vgt, nm, lr, li, vl, nm, vr, nm, work, lwork, info)
424#endif
425
426 ! Find eigenvector with smallest imaginary eigenvalue (real eigenvector of VGT)
427 idx = 1
428 do r = 2, 3
429 if (abs(li(r)) < abs(li(idx))) then
430 idx = r
431 end if
432 end do
433 eigvec = vr(:,idx)
434
435 ! Normalize real eigenvector if it is effectively non-zero
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
439 else
440 eigvec = 0._wp
441 end if
442
443 ! Compute vorticity projected on the eigenvector
444 omega_proj = (vgt(3, 2) - vgt(2, 3))*eigvec(1) + (vgt(1, 3) - vgt(3, 1))*eigvec(2) + (vgt(2, 1) - vgt(1, &
445 & 2))*eigvec(3)
446
447 ! As eigenvector can have +/- signs, we can choose the sign so that omega_proj is positive
448 if (omega_proj < 0._wp) then
449 eigvec = -eigvec
450 omega_proj = -omega_proj
451 end if
452
453 ! Imaginary eigenvalue of the complex conjugate pair (cyclic index selection)
454 lci = li(mod(idx, 3) + 1)
455
456 ! Discriminant: determines whether rotation dominates strain
457 alpha = omega_proj**2._wp - 4._wp*lci**2._wp
458 ! Liutex magnitude = omega_proj - sqrt(discriminant) when rotation dominates
459 if (alpha > 0._wp) then
460 liutex_mag(j, k, l) = omega_proj - sqrt(alpha)
461 else
462 liutex_mag(j, k, l) = omega_proj
463 end if
464
465 ! Compute Liutex axis
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)
469 end do
470 end do
471 end do
472
473 end subroutine s_derive_liutex
474
475 !> Compute the values of the numerical Schlieren function, which are subsequently stored in the derived flow quantity storage
476 !! variable, q_sf.
477 impure subroutine s_derive_numerical_schlieren_function(q_cons_vf, q_sf)
478
479 type(scalar_field), dimension(sys_size), intent(in) :: q_cons_vf
480
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
483
484 real(wp) :: drho_dx, drho_dy, drho_dz !< Spatial derivatives of the density in the x-, y- and z-directions
485 real(wp), dimension(2) :: gm_rho_max !< Global (max gradient magnitude, rank) pair for density
486 integer :: i, j, k, l
487
488 do l = -offset_z%beg, p + offset_z%end
489 do k = -offset_y%beg, n + offset_y%end
490 do j = -offset_x%beg, m + offset_x%end
491 drho_dx = 0._wp
492 drho_dy = 0._wp
493
494 do i = -fd_number, fd_number
495 drho_dx = drho_dx + fd_coeff_x(i, j)*rho_sf(i + j, k, l)
496 drho_dy = drho_dy + fd_coeff_y(i, k)*rho_sf(j, i + k, l)
497 end do
498
499 gm_rho_sf(j, k, l) = drho_dx*drho_dx + drho_dy*drho_dy
500 end do
501 end do
502 end do
503
504 if (p > 0) then
505 do l = -offset_z%beg, p + offset_z%end
506 do k = -offset_y%beg, n + offset_y%end
507 do j = -offset_x%beg, m + offset_x%end
508 drho_dz = 0._wp
509
510 do i = -fd_number, fd_number
511 if (grid_geometry == 3) then
512 drho_dz = drho_dz + fd_coeff_z(i, l)/y_cc(k)*rho_sf(j, k, i + l)
513 else
514 drho_dz = drho_dz + fd_coeff_z(i, l)*rho_sf(j, k, i + l)
515 end if
516 end do
517
518 gm_rho_sf(j, k, l) = gm_rho_sf(j, k, l) + drho_dz*drho_dz
519 end do
520 end do
521 end do
522 end if
523
524 gm_rho_sf = sqrt(gm_rho_sf)
525
526 gm_rho_max = (/maxval(gm_rho_sf), real(proc_rank, wp)/)
527
528 if (num_procs > 1) call s_mpi_reduce_maxloc(gm_rho_max)
529
530 ! The form of the numerical Schlieren function depends on the choice of the multicomponent flow model. For the gamma/pi_inf
531 ! model, the exponential of the negative, normalized, gradient magnitude of the density is computed. For the volume fraction
532 ! model, the amplitude of the exponential's inside is also modulated with respect to the identity of the fluid in which the
533 ! function is evaluated. For more information, refer to Marquina and Mulet (2003).
534
535 if (model_eqns == 1) then ! Gamma/pi_inf model
536 q_sf = -gm_rho_sf/gm_rho_max(1)
537 else ! Volume fraction model
538 do l = -offset_z%beg, p + offset_z%end
539 do k = -offset_y%beg, n + offset_y%end
540 do j = -offset_x%beg, m + offset_x%end
541 q_sf(j, k, l) = 0._wp
542
543 do i = 1, adv_idx%end - e_idx
544 q_sf(j, k, l) = q_sf(j, k, l) - schlieren_alpha(i)*q_cons_vf(i + e_idx)%sf(j, k, l)*gm_rho_sf(j, k, &
545 & l)/gm_rho_max(1)
546 end do
547 end do
548 end do
549 end do
550 end if
551
552 ! Up until now, only the inside of the exponential of the numerical Schlieren function has been evaluated and stored. Then,
553 ! to finish the computation, the exponential of the inside quantity is taken.
554 q_sf = exp(q_sf)
555
557
558 !> Deallocation procedures for the module
560
561 ! Deallocating the variable containing the gradient magnitude of the density field provided that the numerical Schlieren
562 ! function was was outputted during the post-process
563 if (schlieren_wrt) deallocate (gm_rho_sf)
564
565 ! Deallocating the variables that might have been used to bookkeep the finite-difference coefficients in the x-, y- and
566 ! z-directions
567 if (allocated(fd_coeff_x)) deallocate (fd_coeff_x)
568 if (allocated(fd_coeff_y)) deallocate (fd_coeff_y)
569 if (allocated(fd_coeff_z)) deallocate (fd_coeff_z)
570
572
573end 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...
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).