Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
flow_ic.f90
Go to the documentation of this file.
1! Copyright (c) 2021-2025, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
34module flow_ic
35 use num_types, only : rp
36 use logger, only : neko_log, log_size
37 use gather_scatter, only : gs_t, gs_op_add
43 use field, only : field_t
44 use utils, only : neko_error, filename_chsuffix, &
46 use coefs, only : coef_t
47 use math, only : col2, cfill, cfill_mask, abscmp
48 use device_math, only : device_col2
50 use json_module, only : json_file
53 use point_zone, only : point_zone_t
56 use fld_file, only : fld_file_t
57 use file, only : file_t
60 use space, only : space_t, gll
61 use field_list, only : field_list_t
62 use operators, only : rotate_cyc
63 implicit none
64 private
65
69 end interface set_flow_ic
70
72
73contains
74
76 subroutine set_flow_ic_int(u, v, w, p, coef, gs, type, params)
77 type(field_t), intent(inout) :: u
78 type(field_t), intent(inout) :: v
79 type(field_t), intent(inout) :: w
80 type(field_t), intent(inout) :: p
81 type(coef_t), intent(in) :: coef
82 type(gs_t), intent(inout) :: gs
83 character(len=*) :: type
84 type(json_file), intent(inout) :: params
85 real(kind=rp) :: delta, tol
86 real(kind=rp), allocatable :: uinf(:)
87 real(kind=rp), allocatable :: zone_value(:)
88 character(len=:), allocatable :: read_str
89 character(len=NEKO_FNAME_LEN) :: fname, mesh_fname
90 logical :: interpolate
91
92
93 !
94 ! Uniform (Uinf, Vinf, Winf)
95 !
96 if (trim(type) .eq. 'uniform') then
97
98 call json_get_or_lookup(params, 'value', uinf)
99 call set_flow_ic_uniform(u, v, w, uinf)
100
101 !
102 ! Blasius boundary layer
103 !
104 else if (trim(type) .eq. 'blasius') then
105
106 call json_get_or_lookup(params, 'delta', delta)
107 call json_get(params, 'approximation', read_str)
108 call json_get_or_lookup(params, 'freestream_velocity', uinf)
109
110 call set_flow_ic_blasius(u, v, w, delta, uinf, read_str)
111
112 !
113 ! Point zone initial condition
114 !
115 else if (trim(type) .eq. 'point_zone') then
116
117 call json_get_or_lookup(params, 'base_value', uinf)
118 call json_get(params, 'zone_name', read_str)
119 call json_get_or_lookup(params, 'zone_value', zone_value)
120
121 call set_flow_ic_point_zone(u, v, w, uinf, read_str, zone_value)
122
123 !
124 ! Field initial condition (from fld file)
125 !
126 else if (trim(type) .eq. 'field') then
127
128 call json_get(params, 'file_name', read_str)
129 fname = trim(read_str)
130 call json_get_or_default(params, 'interpolate', interpolate, &
131 .false.)
132 call json_get_or_lookup_or_default(params, 'tolerance', tol, &
133 0.000001_rp)
134 call json_get_or_default(params, 'mesh_file_name', read_str, "none")
135 mesh_fname = trim(read_str)
136
137 call set_flow_ic_fld(u, v, w, p, fname, interpolate, tol, mesh_fname)
138
139 else
140 call neko_error('Invalid initial condition')
141 end if
142
143 call set_flow_ic_common(u, v, w, p, coef, gs)
144
145 end subroutine set_flow_ic_int
146
148 subroutine set_flow_ic_usr(u, v, w, p, coef, gs, user_proc, scheme_name)
149 type(field_t), target, intent(inout) :: u
150 type(field_t), target, intent(inout) :: v
151 type(field_t), target, intent(inout) :: w
152 type(field_t), target, intent(inout) :: p
153 type(coef_t), intent(in) :: coef
154 type(gs_t), intent(inout) :: gs
155 procedure(user_initial_conditions_intf) :: user_proc
156 character(len=*), intent(in) :: scheme_name
157
158 type(field_list_t) :: fields
159
160
161 call neko_log%message("Type: user")
162
163 call fields%init(4)
164 call fields%assign_to_field(1, u)
165 call fields%assign_to_field(2, v)
166 call fields%assign_to_field(3, w)
167 call fields%assign_to_field(4, p)
168
169 call user_proc(scheme_name, fields)
170
171 call set_flow_ic_common(u, v, w, p, coef, gs)
172
173 end subroutine set_flow_ic_usr
174
177 subroutine set_compressible_flow_ic_usr(rho, u, v, w, p, coef, gs, &
178 user_proc, scheme_name)
179 type(field_t), target, intent(inout) :: rho
180 type(field_t), target, intent(inout) :: u
181 type(field_t), target, intent(inout) :: v
182 type(field_t), target, intent(inout) :: w
183 type(field_t), target, intent(inout) :: p
184 type(coef_t), intent(in) :: coef
185 type(gs_t), intent(inout) :: gs
186 procedure(user_initial_conditions_intf) :: user_proc
187 character(len=*), intent(in) :: scheme_name
188 integer :: n
189 type(field_list_t) :: fields
190
191
192 call neko_log%message("Type: user (compressible flows)")
193
194 call fields%init(5)
195 call fields%assign_to_field(1, rho)
196 call fields%assign_to_field(2, u)
197 call fields%assign_to_field(3, v)
198 call fields%assign_to_field(4, w)
199 call fields%assign_to_field(5, p)
200 call user_proc(scheme_name, fields)
201
202 call set_flow_ic_common(u, v, w, p, coef, gs)
203
204 n = u%dof%size()
205
206 if (neko_bcknd_device .eq. 1) then
207 call device_memcpy(p%x, p%x_d, n, host_to_device, sync = .false.)
208 call device_memcpy(rho%x, rho%x_d, n, host_to_device, sync = .false.)
209 end if
210
211 ! Ensure continuity across elements for initial conditions
212 ! These variables are not treated in the common constructor
213 call gs%op(p%x, p%dof%size(), gs_op_add)
214 call gs%op(rho%x, rho%dof%size(), gs_op_add)
215
216 if (neko_bcknd_device .eq. 1) then
217 call device_col2(rho%x_d, coef%mult_d, rho%dof%size())
218 call device_col2(p%x_d, coef%mult_d, p%dof%size())
219 else
220 call col2(rho%x, coef%mult, rho%dof%size())
221 call col2(p%x, coef%mult, p%dof%size())
222 end if
223
224 end subroutine set_compressible_flow_ic_usr
225
226 subroutine set_flow_ic_common(u, v, w, p, coef, gs)
227 type(field_t), intent(inout) :: u
228 type(field_t), intent(inout) :: v
229 type(field_t), intent(inout) :: w
230 type(field_t), intent(inout) :: p
231 type(coef_t), intent(in) :: coef
232 type(gs_t), intent(inout) :: gs
233 integer :: n
234
235 n = u%dof%size()
236
237 if (neko_bcknd_device .eq. 1) then
238 call device_memcpy(u%x, u%x_d, n, &
239 host_to_device, sync = .false.)
240 call device_memcpy(v%x, v%x_d, n, &
241 host_to_device, sync = .false.)
242 call device_memcpy(w%x, w%x_d, n, &
243 host_to_device, sync = .false.)
244 end if
245
246 ! Ensure continuity across elements for initial conditions
247 call rotate_cyc(u%x, v%x, w%x, 1, coef)
248 call gs%op(u%x, u%dof%size(), gs_op_add)
249 call gs%op(v%x, v%dof%size(), gs_op_add)
250 call gs%op(w%x, w%dof%size(), gs_op_add)
251 call rotate_cyc(u%x, v%x, w%x, 0, coef)
252
253 if (neko_bcknd_device .eq. 1) then
254 call device_col2(u%x_d, coef%mult_d, u%dof%size())
255 call device_col2(v%x_d, coef%mult_d, v%dof%size())
256 call device_col2(w%x_d, coef%mult_d, w%dof%size())
257 else
258 call col2(u%x, coef%mult, u%dof%size())
259 call col2(v%x, coef%mult, v%dof%size())
260 call col2(w%x, coef%mult, w%dof%size())
261 end if
262
263 end subroutine set_flow_ic_common
264
266 subroutine set_flow_ic_uniform(u, v, w, uinf)
267 type(field_t), intent(inout) :: u
268 type(field_t), intent(inout) :: v
269 type(field_t), intent(inout) :: w
270 real(kind=rp), intent(in) :: uinf(3)
271 integer :: n, i
272 character(len=LOG_SIZE) :: log_buf
273
274 call neko_log%message("Type : uniform")
275 write (log_buf, '(A, 3(ES12.6, A))') "Value: [", &
276 (uinf(i), ", ", i = 1, 2), uinf(3), "]"
277 call neko_log%message(log_buf)
278
279 u = uinf(1)
280 v = uinf(2)
281 w = uinf(3)
282 n = u%dof%size()
283 if (neko_bcknd_device .eq. 1) then
284 call cfill(u%x, uinf(1), n)
285 call cfill(v%x, uinf(2), n)
286 call cfill(w%x, uinf(3), n)
287 end if
288
289 end subroutine set_flow_ic_uniform
290
293 subroutine set_flow_ic_blasius(u, v, w, delta, uinf, type)
294 type(field_t), intent(inout) :: u
295 type(field_t), intent(inout) :: v
296 type(field_t), intent(inout) :: w
297 real(kind=rp), intent(in) :: delta
298 real(kind=rp), intent(in) :: uinf(3)
299 character(len=*), intent(in) :: type
300 procedure(blasius_profile), pointer :: bla => null()
301 integer :: i
302 character(len=LOG_SIZE) :: log_buf
303
304 call neko_log%message("Type : blasius")
305 write (log_buf, '(A,ES12.6)') "delta : ", delta
306 call neko_log%message(log_buf)
307 call neko_log%message("Approximation : " // trim(type))
308 write (log_buf, '(A,"[",2(ES12.6,","),ES12.6,"]")') "Value : ", &
309 uinf(1), uinf(2), uinf(3)
310 call neko_log%message(log_buf)
311
312 select case (trim(type))
313 case ('linear')
314 bla => blasius_linear
315 case ('quadratic')
316 bla => blasius_quadratic
317 case ('cubic')
318 bla => blasius_cubic
319 case ('quartic')
320 bla => blasius_quartic
321 case ('sin')
322 bla => blasius_sin
323 case ('tanh')
324 bla => blasius_tanh
325 case default
326 call neko_error('Invalid Blasius approximation')
327 end select
328
329 if ((uinf(1) .gt. 0.0_rp) .and. abscmp(uinf(2), 0.0_rp) &
330 .and. abscmp(uinf(3), 0.0_rp)) then
331 do i = 1, u%dof%size()
332 u%x(i,1,1,1) = bla(u%dof%z(i,1,1,1), delta, uinf(1))
333 v%x(i,1,1,1) = 0.0_rp
334 w%x(i,1,1,1) = 0.0_rp
335 end do
336 else if (abscmp(uinf(1), 0.0_rp) .and. (uinf(2) .gt. 0.0_rp) &
337 .and. abscmp(uinf(3), 0.0_rp)) then
338 do i = 1, u%dof%size()
339 u%x(i,1,1,1) = 0.0_rp
340 v%x(i,1,1,1) = bla(u%dof%x(i,1,1,1), delta, uinf(2))
341 w%x(i,1,1,1) = 0.0_rp
342 end do
343 else if (abscmp(uinf(1), 0.0_rp) .and. abscmp(uinf(2), 0.0_rp) &
344 .and. (uinf(3) .gt. 0.0_rp)) then
345 do i = 1, u%dof%size()
346 u%x(i,1,1,1) = 0.0_rp
347 v%x(i,1,1,1) = 0.0_rp
348 w%x(i,1,1,1) = bla(u%dof%y(i,1,1,1), delta, uinf(3))
349 end do
350 end if
351
352 end subroutine set_flow_ic_blasius
353
363 subroutine set_flow_ic_point_zone(u, v, w, base_value, zone_name, zone_value)
364 type(field_t), intent(inout) :: u
365 type(field_t), intent(inout) :: v
366 type(field_t), intent(inout) :: w
367 real(kind=rp), intent(in), dimension(3) :: base_value
368 character(len=*), intent(in) :: zone_name
369 real(kind=rp), intent(in) :: zone_value(:)
370 character(len=LOG_SIZE) :: log_buf
371
372 ! Internal variables
373 class(point_zone_t), pointer :: zone
374 integer :: size
375
376 call neko_log%message("Type : point_zone")
377 write (log_buf, '(A,ES12.6)') "Base value : ", base_value
378 call neko_log%message(log_buf)
379 call neko_log%message("Zone name : " // trim(zone_name))
380 write (log_buf, '(A,"[",2(ES12.6,","),ES12.6," ]")') "Value : ", &
381 zone_value(1), zone_value(2), zone_value(3)
382 call neko_log%message(log_buf)
383
384 call set_flow_ic_uniform(u, v, w, base_value)
385 size = u%dof%size()
386
387 zone => neko_point_zone_registry%get_point_zone(trim(zone_name))
388
389 call cfill_mask(u%x, zone_value(1), size, zone%mask%get(), zone%size)
390 call cfill_mask(v%x, zone_value(2), size, zone%mask%get(), zone%size)
391 call cfill_mask(w%x, zone_value(3), size, zone%mask%get(), zone%size)
392
393 end subroutine set_flow_ic_point_zone
394
412 subroutine set_flow_ic_fld(u, v, w, p, file_name, &
413 interpolate, tolerance, mesh_file_name)
414 type(field_t), target, intent(inout) :: u
415 type(field_t), target, intent(inout) :: v
416 type(field_t), target, intent(inout) :: w
417 type(field_t), target, intent(inout) :: p
418 character(len=*), intent(inout) :: file_name
419 logical, intent(in) :: interpolate
420 real(kind=rp), intent(in) :: tolerance
421 character(len=*), intent(inout) :: mesh_file_name
422
423 type(field_t), pointer :: us, vs, ws, ps
424
425 us => u
426 vs => v
427 ws => w
428 ps => p
429
430 call import_fields(file_name, mesh_file_name, &
431 u = us, v = vs, w = ws, p = ps, &
432 interpolate = interpolate, tolerance = tolerance)
433
434 nullify(us, vs, ws, ps)
435
436 ! If we are on GPU we need to move (u,v,w) back to the host
437 ! since set_flow_ic_common copies it again to the device.
438 ! NOTE: p can stay on the device since it is not copied by
439 ! set_flow_ic_common
440 call u%copy_from(device_to_host, .false.)
441 call v%copy_from(device_to_host, .false.)
442 call w%copy_from(device_to_host, .true.)
443
444 end subroutine set_flow_ic_fld
445
446end module flow_ic
Copy data between host and device (or device and device)
Definition device.F90:71
Synchronize a device or stream.
Definition device.F90:107
Abstract interface for computing a Blasius flow profile.
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Abstract interface for user defined initial conditions.
Definition user_intf.f90:65
Coefficients.
Definition coef.f90:34
subroutine, public device_col2(a_d, b_d, n, strm)
Vector multiplication .
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
Definition device.F90:47
integer, parameter, public device_to_host
Definition device.F90:47
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
Simple module to handle fld file series. Provides an interface to the different fields sotred in a fl...
NEKTON fld file format.
Definition fld_file.f90:35
Initial flow condition.
Definition flow_ic.f90:34
subroutine set_flow_ic_usr(u, v, w, p, coef, gs, user_proc, scheme_name)
Set intial flow condition (user defined)
Definition flow_ic.f90:149
subroutine set_flow_ic_point_zone(u, v, w, base_value, zone_name, zone_value)
Set the initial condition of the flow based on a point zone.
Definition flow_ic.f90:364
subroutine set_flow_ic_int(u, v, w, p, coef, gs, type, params)
Set initial flow condition (builtin)
Definition flow_ic.f90:77
subroutine set_compressible_flow_ic_usr(rho, u, v, w, p, coef, gs, user_proc, scheme_name)
Set intial flow condition (user defined) for compressible flows.
Definition flow_ic.f90:179
subroutine set_flow_ic_uniform(u, v, w, uinf)
Uniform initial condition.
Definition flow_ic.f90:267
subroutine, public set_flow_ic_fld(u, v, w, p, file_name, interpolate, tolerance, mesh_file_name)
Set the initial condition of the flow based on a field. @detail The fields are read from an fld file....
Definition flow_ic.f90:414
subroutine set_flow_ic_common(u, v, w, p, coef, gs)
Definition flow_ic.f90:227
subroutine set_flow_ic_blasius(u, v, w, delta, uinf, type)
Set a Blasius profile as initial condition.
Definition flow_ic.f90:294
Defines a flow profile.
real(kind=rp) function, public blasius_quadratic(y, delta, u)
Quadratic approximate Blasius Profile .
real(kind=rp) function, public blasius_quartic(y, delta, u)
Quartic approximate Blasius Profile .
real(kind=rp) function, public blasius_sin(y, delta, u)
Sinusoidal approximate Blasius Profile .
real(kind=rp) function, public blasius_cubic(y, delta, u)
Cubic approximate Blasius Profile .
real(kind=rp) function, public blasius_tanh(y, delta, u)
Hyperbolic tangent approximate Blasius Profile from O. Savas (2012) where is the 99 percent thickne...
real(kind=rp) function, public blasius_linear(y, delta, u)
Linear approximate Blasius profile .
Gather-scatter.
Implements global_interpolation given a dofmap.
Importation of fields from fld files.
subroutine, public import_fields(fname, mesh_fname, u, v, w, p, t, s_target_list, s_index_list, interpolate, tolerance)
Imports fields from an fld file, potentially with interpolation.
Routines to interpolate between different spaces.
Utilities for retrieving parameters from the case files.
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:79
integer, parameter, public log_size
Definition log.f90:45
Definition math.f90:60
subroutine, public cfill(a, c, n)
Set all elements to a constant c .
Definition math.f90:487
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:854
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
Definition math.f90:397
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
type(point_zone_registry_t), target, public neko_point_zone_registry
Global point_zone registry.
Defines a function space.
Definition space.f90:34
integer, parameter, public gll
Definition space.f90:49
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
Utilities.
Definition utils.f90:35
integer function, public extract_fld_file_index(fld_filename, default_index)
Extracts the index of a field file. For example, "myfield.f00045" will return 45. If the suffix of th...
Definition utils.f90:157
integer, parameter, public neko_fname_len
Definition utils.f90:40
subroutine, public neko_warning(warning_msg)
Reports a warning to standard output.
Definition utils.f90:346
subroutine, public filename_chsuffix(fname, new_fname, new_suffix)
Change a filename's suffix.
Definition utils.f90:140
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:58
field_list_t, To be able to group fields together
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition file.f90:55
Interface for NEKTON fld files.
Definition fld_file.f90:64
Implements global interpolation for arbitrary points in the domain.
Interpolation between two space::space_t.
Base abstract type for point zones.
The function space for the SEM solution fields.
Definition space.f90:63