Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fluid_volflow.f90
Go to the documentation of this file.
1! Copyright (c) 2008-2020, UCHICAGO ARGONNE, LLC.
2! Copyright (c) 2025, The Neko Authors
3!
4! The UChicago Argonne, LLC as Operator of Argonne National
5! Laboratory holds copyright in the Software. The copyright holder
6! reserves all rights except those expressly granted to licensees,
7! and U.S. Government license rights.
8!
9! Redistribution and use in source and binary forms, with or without
10! modification, are permitted provided that the following conditions
11! are met:
12!
13! 1. Redistributions of source code must retain the above copyright
14! notice, this list of conditions and the disclaimer below.
15!
16! 2. Redistributions in binary form must reproduce the above copyright
17! notice, this list of conditions and the disclaimer (as noted below)
18! in the documentation and/or other materials provided with the
19! distribution.
20!
21! 3. Neither the name of ANL nor the names of its contributors
22! may be used to endorse or promote products derived from this software
23! without specific prior written permission.
24!
25! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
26! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
27! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
28! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL
29! UCHICAGO ARGONNE, LLC, THE U.S. DEPARTMENT OF
30! ENERGY OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
31! SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED
32! TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
33! DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
34! THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
35! (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
36! OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37!
38! Additional BSD Notice
39! ---------------------
40! 1. This notice is required to be provided under our contract with
41! the U.S. Department of Energy (DOE). This work was produced at
42! Argonne National Laboratory under Contract
43! No. DE-AC02-06CH11357 with the DOE.
44!
45! 2. Neither the United States Government nor UCHICAGO ARGONNE,
46! LLC nor any of their employees, makes any warranty,
47! express or implied, or assumes any liability or responsibility for the
48! accuracy, completeness, or usefulness of any information, apparatus,
49! product, or process disclosed, or represents that its use would not
50! infringe privately-owned rights.
51!
52! 3. Also, reference herein to any specific commercial products, process,
53! or services by trade name, trademark, manufacturer or otherwise does
54! not necessarily constitute or imply its endorsement, recommendation,
55! or favoring by the United States Government or UCHICAGO ARGONNE LLC.
56! The views and opinions of authors expressed
57! herein do not necessarily state or reflect those of the United States
58! Government or UCHICAGO ARGONNE, LLC, and shall
59! not be used for advertising or product endorsement purposes.
60!
62 use operators, only : opgrad, cdtp, rotate_cyc
63 use num_types, only : rp
64 use mathops, only : opchsign
65 use krylov, only : ksp_t, ksp_monitor_t
66 use precon, only : pc_t
67 use dofmap, only : dofmap_t
68 use field, only : field_t
69 use coefs, only : coef_t
71 use math, only : copy, glsc2, glmin, glmax, add2, abscmp
76 use gather_scatter, only : gs_t, gs_op_add
77 use json_module, only : json_file
80 use bc_list, only : bc_list_t
81 use ax_product, only : ax_t
83 use mpi_f08, only : mpi_allreduce, mpi_in_place, mpi_sum
84 implicit none
85 private
86
88 type, public :: fluid_volflow_t
89 integer :: flow_dir
90 logical :: avflow
91 real(kind=rp) :: flow_rate
92 real(kind=rp) :: dtlag = 0d0
93 real(kind=rp) :: bdlag = 0d0
94 type(field_t) :: u_vol, v_vol, w_vol, p_vol
95 real(kind=rp) :: domain_length, base_flow
96 contains
97 procedure, pass(this) :: init => fluid_vol_flow_init
98 procedure, pass(this) :: free => fluid_vol_flow_free
99 procedure, pass(this) :: adjust => fluid_vol_flow
100 procedure, private, pass(this) :: compute => fluid_vol_flow_compute
101 end type fluid_volflow_t
102
103contains
104
105 subroutine fluid_vol_flow_init(this, dm_Xh, params)
106 class(fluid_volflow_t), intent(inout) :: this
107 type(dofmap_t), target, intent(in) :: dm_Xh
108 type(json_file), intent(inout) :: params
109 logical average
110 integer :: direction
111 real(kind=rp) :: rate
112
113 call this%free()
114
115 !Initialize vol_flow (if there is a forced volume flow)
116 call json_get_or_lookup(params, 'case.fluid.flow_rate_force.direction', &
117 direction)
118 call json_get_or_lookup(params, 'case.fluid.flow_rate_force.value', rate)
119 call json_get(params, 'case.fluid.flow_rate_force.use_averaged_flow',&
120 average)
121
122 this%flow_dir = direction
123 this%avflow = average
124 this%flow_rate = rate
125
126 if (this%flow_dir .ne. 0) then
127 call this%u_vol%init(dm_xh, 'u_vol')
128 call this%v_vol%init(dm_xh, 'v_vol')
129 call this%w_vol%init(dm_xh, 'w_vol')
130 call this%p_vol%init(dm_xh, 'p_vol')
131 end if
132
133 end subroutine fluid_vol_flow_init
134
135 subroutine fluid_vol_flow_free(this)
136 class(fluid_volflow_t), intent(inout) :: this
137
138 call this%u_vol%free()
139 call this%v_vol%free()
140 call this%w_vol%free()
141 call this%p_vol%free()
142
143 end subroutine fluid_vol_flow_free
144
148 subroutine fluid_vol_flow_compute(this, u_res, v_res, w_res, p_res, &
149 ext_bdf, gs_Xh, c_Xh, rho, mu, bd, dt, &
150 bclst_dp, bclst_du, bclst_dv, bclst_dw, bclst_vel_res, &
151 Ax_vel, Ax_prs, ksp_prs, ksp_vel, pc_prs, pc_vel, prs_max_iter, &
152 vel_max_iter)
153 class(fluid_volflow_t), intent(inout) :: this
154 type(field_t), intent(inout) :: u_res, v_res, w_res, p_res
155 type(coef_t), intent(inout) :: c_Xh
156 type(gs_t), intent(inout) :: gs_Xh
157 type(time_scheme_controller_t), intent(in) :: ext_bdf
158 type(bc_list_t), intent(inout) :: bclst_dp, bclst_du, bclst_dv, bclst_dw
159 type(bc_list_t), intent(inout) :: bclst_vel_res
160 class(ax_t), intent(in) :: Ax_vel
161 class(ax_t), intent(in) :: Ax_prs
162 class(ksp_t), intent(inout) :: ksp_prs, ksp_vel
163 class(pc_t), intent(inout) :: pc_prs, pc_vel
164 real(kind=rp), intent(in) :: bd
165 real(kind=rp), intent(in) :: rho, dt
166 type(field_t) :: mu
167 integer, intent(in) :: vel_max_iter, prs_max_iter
168 integer :: n, i
169 real(kind=rp) :: xlmin, xlmax
170 real(kind=rp) :: ylmin, ylmax
171 real(kind=rp) :: zlmin, zlmax
172 type(ksp_monitor_t) :: ksp_results(4)
173 type(field_t), pointer :: ta1, ta2, ta3
174 integer :: temp_indices(3)
175
176 call neko_scratch_registry%request_field(ta1, temp_indices(1), .false.)
177 call neko_scratch_registry%request_field(ta2, temp_indices(2), .false.)
178 call neko_scratch_registry%request_field(ta3, temp_indices(3), .false.)
179
180
181 associate(msh => c_xh%msh, p_vol => this%p_vol, &
182 u_vol => this%u_vol, v_vol => this%v_vol, w_vol => this%w_vol)
183
184 n = c_xh%dof%size()
185 xlmin = glmin(c_xh%dof%x, n)
186 xlmax = glmax(c_xh%dof%x, n)
187 ylmin = glmin(c_xh%dof%y, n) ! for Y!
188 ylmax = glmax(c_xh%dof%y, n)
189 zlmin = glmin(c_xh%dof%z, n) ! for Z!
190 zlmax = glmax(c_xh%dof%z, n)
191 if (this%flow_dir .eq. 1) then
192 this%domain_length = xlmax - xlmin
193 end if
194 if (this%flow_dir .eq. 2) then
195 this%domain_length = ylmax - ylmin
196 end if
197 if (this%flow_dir .eq. 3) then
198 this%domain_length = zlmax - zlmin
199 end if
200
201 if (neko_bcknd_device .eq. 1) then
202 call device_cfill(c_xh%h1_d, 1.0_rp/rho, n)
203 call device_rzero(c_xh%h2_d, n)
204 else
205 do i = 1, n
206 c_xh%h1(i,1,1,1) = 1.0_rp / rho
207 c_xh%h2(i,1,1,1) = 0.0_rp
208 end do
209 end if
210 c_xh%ifh2 = .false.
211
212 ! Compute pressure
213
214 if (this%flow_dir .eq. 1) then
215 call cdtp(p_res%x, c_xh%h1, c_xh%drdx, c_xh%dsdx, c_xh%dtdx, c_xh)
216 end if
217
218 if (this%flow_dir .eq. 2) then
219 call cdtp(p_res%x, c_xh%h1, c_xh%drdy, c_xh%dsdy, c_xh%dtdy, c_xh)
220 end if
221
222 if (this%flow_dir .eq. 3) then
223 call cdtp(p_res%x, c_xh%h1, c_xh%drdz, c_xh%dsdz, c_xh%dtdz, c_xh)
224 end if
225
226 call gs_xh%op(p_res, gs_op_add)
227 call bclst_dp%apply_scalar(p_res%x, n)
228 call pc_prs%update()
229 ksp_results(1) = ksp_prs%solve(ax_prs, p_vol, p_res%x, n, &
230 c_xh, bclst_dp, gs_xh, prs_max_iter)
231
232 ! Compute velocity
233
234 call opgrad(u_res%x, v_res%x, w_res%x, p_vol%x, c_xh)
235
236 if (neko_bcknd_device .eq. 1) then
237 call device_opchsign(u_res%x_d, v_res%x_d, w_res%x_d, msh%gdim, n)
238 call device_copy(ta1%x_d, c_xh%B_d, n)
239 call device_copy(ta2%x_d, c_xh%B_d, n)
240 call device_copy(ta3%x_d, c_xh%B_d, n)
241 else
242 call opchsign(u_res%x, v_res%x, w_res%x, msh%gdim, n)
243 call copy(ta1%x, c_xh%B, n)
244 call copy(ta2%x, c_xh%B, n)
245 call copy(ta3%x, c_xh%B, n)
246 end if
247 call bclst_vel_res%apply_vector(ta1%x, ta2%x, ta3%x, n)
248
249 ! add forcing
250
251 if (neko_bcknd_device .eq. 1) then
252 if (this%flow_dir .eq. 1) then
253 call device_add2(u_res%x_d, ta1%x_d, n)
254 else if (this%flow_dir .eq. 2) then
255 call device_add2(v_res%x_d, ta2%x_d, n)
256 else if (this%flow_dir .eq. 3) then
257 call device_add2(w_res%x_d, ta3%x_d, n)
258 end if
259 else
260 if (this%flow_dir .eq. 1) then
261 call add2(u_res%x, ta1%x, n)
262 else if (this%flow_dir .eq. 2) then
263 call add2(v_res%x, ta2%x, n)
264 else if (this%flow_dir .eq. 3) then
265 call add2(w_res%x, ta3%x, n)
266 end if
267 end if
268
269 if (neko_bcknd_device .eq. 1) then
270 call device_copy(c_xh%h1_d, mu%x_d, n)
271 call device_cfill(c_xh%h2_d, rho * (bd / dt), n)
272 else
273 call copy(c_xh%h1, mu%x, n)
274 c_xh%h2 = rho * (bd / dt)
275 end if
276 c_xh%ifh2 = .true.
277
278 call rotate_cyc(u_res%x, v_res%x, w_res%x, 1, c_xh)
279 call gs_xh%op(u_res, gs_op_add)
280 call gs_xh%op(v_res, gs_op_add)
281 call gs_xh%op(w_res, gs_op_add)
282 call rotate_cyc(u_res%x, v_res%x, w_res%x, 0, c_xh)
283
284 call bclst_vel_res%apply_vector(u_res%x, v_res%x, w_res%x, n)
285 call pc_vel%update()
286
287 ksp_results(2:4) = ksp_vel%solve_coupled(ax_vel, &
288 u_vol, v_vol, w_vol, &
289 u_res%x, v_res%x, w_res%x, &
290 n, c_xh, &
291 bclst_du, bclst_dv, bclst_dw, &
292 gs_xh, vel_max_iter)
293
294 if (neko_bcknd_device .eq. 1) then
295 if (this%flow_dir .eq. 1) then
296 this%base_flow = &
297 device_glsc2(u_vol%x_d, c_xh%B_d, n) / this%domain_length
298 end if
299
300 if (this%flow_dir .eq. 2) then
301 this%base_flow = &
302 device_glsc2(v_vol%x_d, c_xh%B_d, n) / this%domain_length
303 end if
304
305 if (this%flow_dir .eq. 3) then
306 this%base_flow = &
307 device_glsc2(w_vol%x_d, c_xh%B_d, n) / this%domain_length
308 end if
309 else
310 if (this%flow_dir .eq. 1) then
311 this%base_flow = glsc2(u_vol%x, c_xh%B, n) / this%domain_length
312 end if
313
314 if (this%flow_dir .eq. 2) then
315 this%base_flow = glsc2(v_vol%x, c_xh%B, n) / this%domain_length
316 end if
317
318 if (this%flow_dir .eq. 3) then
319 this%base_flow = glsc2(w_vol%x, c_xh%B, n) / this%domain_length
320 end if
321 end if
322 end associate
323
324 call neko_scratch_registry%relinquish_field(temp_indices)
325 end subroutine fluid_vol_flow_compute
326
336 subroutine fluid_vol_flow(this, u, v, w, p, u_res, v_res, w_res, p_res, &
337 c_Xh, gs_Xh, ext_bdf, rho, mu, dt, &
338 bclst_dp, bclst_du, bclst_dv, bclst_dw, bclst_vel_res, &
339 Ax_vel, Ax_prs, ksp_prs, ksp_vel, pc_prs, pc_vel, prs_max_iter, &
340 vel_max_iter)
341
342 class(fluid_volflow_t), intent(inout) :: this
343 type(field_t), intent(inout) :: u, v, w, p
344 type(field_t), intent(inout) :: u_res, v_res, w_res, p_res
345 type(coef_t), intent(inout) :: c_Xh
346 type(gs_t), intent(inout) :: gs_Xh
347 type(time_scheme_controller_t), intent(in) :: ext_bdf
348 real(kind=rp), intent(in) :: rho, dt
349 type(field_t) :: mu
350 type(bc_list_t), intent(inout) :: bclst_dp, bclst_du, bclst_dv, bclst_dw
351 type(bc_list_t), intent(inout) :: bclst_vel_res
352 class(ax_t), intent(in) :: Ax_vel
353 class(ax_t), intent(in) :: Ax_prs
354 class(ksp_t), intent(inout) :: ksp_prs, ksp_vel
355 class(pc_t), intent(inout) :: pc_prs, pc_vel
356 integer, intent(in) :: prs_max_iter, vel_max_iter
357 real(kind=rp) :: ifcomp, flow_rate, xsec
358 real(kind=rp) :: current_flow, delta_flow, scale
359 integer :: n, ierr, i
360
361 associate(u_vol => this%u_vol, v_vol => this%v_vol, &
362 w_vol => this%w_vol, p_vol => this%p_vol)
363
364 n = c_xh%dof%size()
365
366 ! If either dt or the backwards difference coefficient change,
367 ! then recompute base flow solution corresponding to unit forcing:
368
369 ifcomp = 0.0_rp
370
371 if ((.not. abscmp(dt, this%dtlag)) .or. &
372 (.not. abscmp(ext_bdf%diffusion_coeffs(1), this%bdlag))) then
373 ifcomp = 1.0_rp
374 end if
375
376 this%dtlag = dt
377 this%bdlag = ext_bdf%diffusion_coeffs(1)
378
379 call mpi_allreduce(mpi_in_place, ifcomp, 1, &
380 mpi_real_precision, mpi_sum, neko_comm, ierr)
381
382 if (ifcomp .gt. 0d0) then
383 call this%compute(u_res, v_res, w_res, p_res, &
384 ext_bdf, gs_xh, c_xh, rho, mu, ext_bdf%diffusion_coeffs(1), dt, &
385 bclst_dp, bclst_du, bclst_dv, bclst_dw, bclst_vel_res, &
386 ax_vel, ax_prs, ksp_prs, ksp_vel, pc_prs, pc_vel, prs_max_iter, &
387 vel_max_iter)
388 end if
389
390 if (neko_bcknd_device .eq. 1) then
391 if (this%flow_dir .eq. 1) then
392 current_flow = &
393 device_glsc2(u%x_d, c_xh%B_d, n) / this%domain_length ! for X
394 else if (this%flow_dir .eq. 2) then
395 current_flow = &
396 device_glsc2(v%x_d, c_xh%B_d, n) / this%domain_length ! for Y
397 else if (this%flow_dir .eq. 3) then
398 current_flow = &
399 device_glsc2(w%x_d, c_xh%B_d, n) / this%domain_length ! for Z
400 end if
401 else
402 if (this%flow_dir .eq. 1) then
403 current_flow = glsc2(u%x, c_xh%B, n) / this%domain_length ! for X
404 else if (this%flow_dir .eq. 2) then
405 current_flow = glsc2(v%x, c_xh%B, n) / this%domain_length ! for Y
406 else if (this%flow_dir .eq. 3) then
407 current_flow = glsc2(w%x, c_xh%B, n) / this%domain_length ! for Z
408 end if
409 end if
410
411 if (this%avflow) then
412 xsec = c_xh%volume / this%domain_length
413 flow_rate = this%flow_rate*xsec
414 else
415 flow_rate = this%flow_rate
416 end if
417
418 delta_flow = flow_rate - current_flow
419 scale = delta_flow / this%base_flow
420
421 if (neko_bcknd_device .eq. 1) then
422 call device_add2s2(u%x_d, u_vol%x_d, scale, n)
423 call device_add2s2(v%x_d, v_vol%x_d, scale, n)
424 call device_add2s2(w%x_d, w_vol%x_d, scale, n)
425 call device_add2s2(p%x_d, p_vol%x_d, scale, n)
426 else
427 do concurrent(i = 1: n)
428 u%x(i,1,1,1) = u%x(i,1,1,1) + scale * u_vol%x(i,1,1,1)
429 v%x(i,1,1,1) = v%x(i,1,1,1) + scale * v_vol%x(i,1,1,1)
430 w%x(i,1,1,1) = w%x(i,1,1,1) + scale * w_vol%x(i,1,1,1)
431 p%x(i,1,1,1) = p%x(i,1,1,1) + scale * p_vol%x(i,1,1,1)
432 end do
433 end if
434 end associate
435
436 end subroutine fluid_vol_flow
437
438
439end module fluid_volflow
Retrieves a parameter by name or throws an error.
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_datatype), public mpi_real_precision
MPI type for working precision of REAL types.
Definition comm.F90:51
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:43
subroutine, public device_add2s2(a_d, b_d, c1, n, strm)
Vector addition with scalar multiplication (multiplication on first argument)
subroutine, public device_add2(a_d, b_d, n, strm)
Vector addition .
subroutine, public device_rzero(a_d, n, strm)
Zero a real vector.
subroutine, public device_copy(a_d, b_d, n, strm)
Copy a vector .
real(kind=rp) function, public device_glsc2(a_d, b_d, n, strm)
Weighted inner product .
subroutine, public device_cfill(a_d, c, n, strm)
Set all elements to a constant c .
subroutine, public device_opchsign(a1_d, a2_d, a3_d, gdim, n)
Defines a mapping of the degrees of freedom.
Definition dofmap.f90:35
Defines a field.
Definition field.f90:34
subroutine fluid_vol_flow(this, u, v, w, p, u_res, v_res, w_res, p_res, c_xh, gs_xh, ext_bdf, rho, mu, dt, bclst_dp, bclst_du, bclst_dv, bclst_dw, bclst_vel_res, ax_vel, ax_prs, ksp_prs, ksp_vel, pc_prs, pc_vel, prs_max_iter, vel_max_iter)
Adjust flow volume.
subroutine fluid_vol_flow_init(this, dm_xh, params)
subroutine fluid_vol_flow_compute(this, u_res, v_res, w_res, p_res, ext_bdf, gs_xh, c_xh, rho, mu, bd, dt, bclst_dp, bclst_du, bclst_dv, bclst_dw, bclst_vel_res, ax_vel, ax_prs, ksp_prs, ksp_vel, pc_prs, pc_vel, prs_max_iter, vel_max_iter)
Compute flow adjustment.
subroutine fluid_vol_flow_free(this)
Gather-scatter.
Utilities for retrieving parameters from the case files.
Implements the base abstract type for Krylov solvers plus helper types.
Definition krylov.f90:34
Definition math.f90:60
real(kind=rp) function, public glsc2(a, b, n)
Weighted inner product .
Definition math.f90:1048
subroutine, public add2(a, b, n)
Vector addition .
Definition math.f90:726
real(kind=rp) function, public glmax(a, n)
Max of a vector of length n.
Definition math.f90:516
subroutine, public copy(a, b, n)
Copy a vector .
Definition math.f90:249
real(kind=rp) function, public glmin(a, n)
Min of a vector of length n.
Definition math.f90:546
Collection of vector field operations operating on and . Note that in general the indices and ....
Definition mathops.f90:65
subroutine, public opchsign(a1, a2, a3, gdim, n)
for and .
Definition mathops.f90:76
Build configurations.
integer, parameter neko_bcknd_device
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public opgrad(ux, uy, uz, u, coef, es, ee)
Compute the weak gradient of a scalar field, i.e. the gradient multiplied by the mass matrix.
subroutine, public cdtp(dtx, x, dr, ds, dt, coef, es, ee)
Apply D^T to a scalar field, where D is the derivative matrix.
Krylov preconditioner.
Definition precon.f90:34
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Compound scheme for the advection and diffusion operators in a transport equation.
Base type for a matrix-vector product providing .
Definition ax.f90:43
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:48
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
Type for storing initial and final residuals in a Krylov solver.
Definition krylov.f90:56
Base abstract type for a canonical Krylov method, solving .
Definition krylov.f90:73
Defines a canonical Krylov preconditioner.
Definition precon.f90:40
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...