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
54 use point_zone, only : point_zone_t
57 use fld_file, only : fld_file_t
58 use file, only : file_t
61 use space, only : space_t, gll
62 use field_list, only : field_list_t
63 use operators, only : rotate_cyc
64 implicit none
65 private
66
70 end interface set_flow_ic
71
73
74contains
75
77 subroutine set_flow_ic_int(u, v, w, p, coef, gs, type, params)
78 type(field_t), intent(inout) :: u
79 type(field_t), intent(inout) :: v
80 type(field_t), intent(inout) :: w
81 type(field_t), intent(inout) :: p
82 type(coef_t), intent(in) :: coef
83 type(gs_t), intent(inout) :: gs
84 character(len=*) :: type
85 type(json_file), intent(inout) :: params
86 real(kind=rp) :: delta
87 real(kind=rp), allocatable :: uinf(:)
88 real(kind=rp), allocatable :: zone_value(:)
89 character(len=:), allocatable :: read_str
90
91
92 !
93 ! Uniform (Uinf, Vinf, Winf)
94 !
95 if (trim(type) .eq. 'uniform') then
96
97 call json_get_or_lookup(params, 'value', uinf)
98 call set_flow_ic_uniform(u, v, w, uinf)
99
100 !
101 ! Blasius boundary layer
102 !
103 else if (trim(type) .eq. 'blasius') then
104
105 call json_get_or_lookup(params, 'delta', delta)
106 call json_get(params, 'approximation', read_str)
107 call json_get_or_lookup(params, 'freestream_velocity', uinf)
108
109 call set_flow_ic_blasius(u, v, w, delta, uinf, read_str)
110
111 !
112 ! Point zone initial condition
113 !
114 else if (trim(type) .eq. 'point_zone') then
115
116 call json_get_or_lookup(params, 'base_value', uinf)
117 call json_get(params, 'zone_name', read_str)
118 call json_get_or_lookup(params, 'zone_value', zone_value)
119
120 call set_flow_ic_point_zone(u, v, w, uinf, read_str, zone_value)
121
122 !
123 ! Field initial condition (from fld file)
124 !
125 else if (trim(type) .eq. 'field') then
126
127 block
128 character(len=NEKO_FNAME_LEN) :: fname, mesh_fname
129 logical :: interpolate
130 type(json_file) :: interp_subdict
131
132 call json_get(params, 'file_name', read_str)
133 fname = trim(read_str)
134
135 call json_get_or_default(params, 'interpolate', interpolate, &
136 .false.)
137
138 call json_get_or_default(params, 'mesh_file_name', read_str, "none")
139 mesh_fname = trim(read_str)
140
141 call json_get_subdict_or_empty(params, "interpolation", &
142 interp_subdict)
143 call set_flow_ic_fld(u, v, w, p, fname, interpolate, &
144 mesh_fname, interp_subdict)
145 end block
146
147 else
148 call neko_error('Invalid initial condition')
149 end if
150
151 call set_flow_ic_common(u, v, w, p, coef, gs)
152
153 end subroutine set_flow_ic_int
154
156 subroutine set_flow_ic_usr(u, v, w, p, coef, gs, user_proc, scheme_name)
157 type(field_t), target, intent(inout) :: u
158 type(field_t), target, intent(inout) :: v
159 type(field_t), target, intent(inout) :: w
160 type(field_t), target, intent(inout) :: p
161 type(coef_t), intent(in) :: coef
162 type(gs_t), intent(inout) :: gs
163 procedure(user_initial_conditions_intf) :: user_proc
164 character(len=*), intent(in) :: scheme_name
165
166 type(field_list_t) :: fields
167
168
169 call neko_log%message("Type: user")
170
171 call fields%init(4)
172 call fields%assign_to_field(1, u)
173 call fields%assign_to_field(2, v)
174 call fields%assign_to_field(3, w)
175 call fields%assign_to_field(4, p)
176
177 call user_proc(scheme_name, fields)
178
179 call set_flow_ic_common(u, v, w, p, coef, gs)
180
181 end subroutine set_flow_ic_usr
182
185 subroutine set_compressible_flow_ic_usr(rho, u, v, w, p, coef, gs, &
186 user_proc, scheme_name)
187 type(field_t), target, intent(inout) :: rho
188 type(field_t), target, intent(inout) :: u
189 type(field_t), target, intent(inout) :: v
190 type(field_t), target, intent(inout) :: w
191 type(field_t), target, intent(inout) :: p
192 type(coef_t), intent(in) :: coef
193 type(gs_t), intent(inout) :: gs
194 procedure(user_initial_conditions_intf) :: user_proc
195 character(len=*), intent(in) :: scheme_name
196 integer :: n
197 type(field_list_t) :: fields
198
199
200 call neko_log%message("Type: user (compressible flows)")
201
202 call fields%init(5)
203 call fields%assign_to_field(1, rho)
204 call fields%assign_to_field(2, u)
205 call fields%assign_to_field(3, v)
206 call fields%assign_to_field(4, w)
207 call fields%assign_to_field(5, p)
208 call user_proc(scheme_name, fields)
209
210 call set_flow_ic_common(u, v, w, p, coef, gs)
211
212 n = u%dof%size()
213
214 if (neko_bcknd_device .eq. 1) then
215 call device_memcpy(p%x, p%x_d, n, host_to_device, sync = .false.)
216 call device_memcpy(rho%x, rho%x_d, n, host_to_device, sync = .false.)
217 end if
218
219 ! Ensure continuity across elements for initial conditions
220 ! These variables are not treated in the common constructor
221 call gs%op(p%x, p%dof%size(), gs_op_add)
222 call gs%op(rho%x, rho%dof%size(), gs_op_add)
223
224 if (neko_bcknd_device .eq. 1) then
225 call device_col2(rho%x_d, coef%mult_d, rho%dof%size())
226 call device_col2(p%x_d, coef%mult_d, p%dof%size())
227 else
228 call col2(rho%x, coef%mult, rho%dof%size())
229 call col2(p%x, coef%mult, p%dof%size())
230 end if
231
232 end subroutine set_compressible_flow_ic_usr
233
234 subroutine set_flow_ic_common(u, v, w, p, coef, gs)
235 type(field_t), intent(inout) :: u
236 type(field_t), intent(inout) :: v
237 type(field_t), intent(inout) :: w
238 type(field_t), intent(inout) :: p
239 type(coef_t), intent(in) :: coef
240 type(gs_t), intent(inout) :: gs
241 integer :: n
242
243 n = u%dof%size()
244
245 if (neko_bcknd_device .eq. 1) then
246 call u%copy_from(host_to_device, sync = .false.)
247 call v%copy_from(host_to_device, sync = .false.)
248 call w%copy_from(host_to_device, sync = .false.)
249
250 ! also copy pressure for consistency
251 call p%copy_from(host_to_device, sync = .true.)
252 end if
253
254 ! Ensure continuity across elements for initial conditions
255 call rotate_cyc(u%x, v%x, w%x, 1, coef)
256 call gs%op(u%x, u%dof%size(), gs_op_add)
257 call gs%op(v%x, v%dof%size(), gs_op_add)
258 call gs%op(w%x, w%dof%size(), gs_op_add)
259 call rotate_cyc(u%x, v%x, w%x, 0, coef)
260
261 if (neko_bcknd_device .eq. 1) then
262 call device_col2(u%x_d, coef%mult_d, u%dof%size())
263 call device_col2(v%x_d, coef%mult_d, v%dof%size())
264 call device_col2(w%x_d, coef%mult_d, w%dof%size())
265 else
266 call col2(u%x, coef%mult, u%dof%size())
267 call col2(v%x, coef%mult, v%dof%size())
268 call col2(w%x, coef%mult, w%dof%size())
269 end if
270
271 end subroutine set_flow_ic_common
272
274 subroutine set_flow_ic_uniform(u, v, w, uinf)
275 type(field_t), intent(inout) :: u
276 type(field_t), intent(inout) :: v
277 type(field_t), intent(inout) :: w
278 real(kind=rp), intent(in) :: uinf(3)
279 integer :: n, i
280 character(len=LOG_SIZE) :: log_buf
281
282 call neko_log%message("Type : uniform")
283 write (log_buf, '(A, 3(ES12.6, A))') "Value: [", &
284 (uinf(i), ", ", i = 1, 2), uinf(3), "]"
285 call neko_log%message(log_buf)
286
287 u = uinf(1)
288 v = uinf(2)
289 w = uinf(3)
290 n = u%dof%size()
291 if (neko_bcknd_device .eq. 1) then
292 call cfill(u%x, uinf(1), n)
293 call cfill(v%x, uinf(2), n)
294 call cfill(w%x, uinf(3), n)
295 end if
296
297 end subroutine set_flow_ic_uniform
298
301 subroutine set_flow_ic_blasius(u, v, w, delta, uinf, type)
302 type(field_t), intent(inout) :: u
303 type(field_t), intent(inout) :: v
304 type(field_t), intent(inout) :: w
305 real(kind=rp), intent(in) :: delta
306 real(kind=rp), intent(in) :: uinf(3)
307 character(len=*), intent(in) :: type
308 procedure(blasius_profile), pointer :: bla => null()
309 integer :: i
310 character(len=LOG_SIZE) :: log_buf
311
312 call neko_log%message("Type : blasius")
313 write (log_buf, '(A,ES12.6)') "delta : ", delta
314 call neko_log%message(log_buf)
315 call neko_log%message("Approximation : " // trim(type))
316 write (log_buf, '(A,"[",2(ES12.6,","),ES12.6,"]")') "Value : ", &
317 uinf(1), uinf(2), uinf(3)
318 call neko_log%message(log_buf)
319
320 select case (trim(type))
321 case ('linear')
322 bla => blasius_linear
323 case ('quadratic')
324 bla => blasius_quadratic
325 case ('cubic')
326 bla => blasius_cubic
327 case ('quartic')
328 bla => blasius_quartic
329 case ('sin')
330 bla => blasius_sin
331 case ('tanh')
332 bla => blasius_tanh
333 case default
334 call neko_error('Invalid Blasius approximation')
335 end select
336
337 if ((uinf(1) .gt. 0.0_rp) .and. abscmp(uinf(2), 0.0_rp) &
338 .and. abscmp(uinf(3), 0.0_rp)) then
339 do i = 1, u%dof%size()
340 u%x(i,1,1,1) = bla(u%dof%z(i,1,1,1), delta, uinf(1))
341 v%x(i,1,1,1) = 0.0_rp
342 w%x(i,1,1,1) = 0.0_rp
343 end do
344 else if (abscmp(uinf(1), 0.0_rp) .and. (uinf(2) .gt. 0.0_rp) &
345 .and. abscmp(uinf(3), 0.0_rp)) then
346 do i = 1, u%dof%size()
347 u%x(i,1,1,1) = 0.0_rp
348 v%x(i,1,1,1) = bla(u%dof%x(i,1,1,1), delta, uinf(2))
349 w%x(i,1,1,1) = 0.0_rp
350 end do
351 else if (abscmp(uinf(1), 0.0_rp) .and. abscmp(uinf(2), 0.0_rp) &
352 .and. (uinf(3) .gt. 0.0_rp)) then
353 do i = 1, u%dof%size()
354 u%x(i,1,1,1) = 0.0_rp
355 v%x(i,1,1,1) = 0.0_rp
356 w%x(i,1,1,1) = bla(u%dof%y(i,1,1,1), delta, uinf(3))
357 end do
358 end if
359
360 end subroutine set_flow_ic_blasius
361
371 subroutine set_flow_ic_point_zone(u, v, w, base_value, zone_name, zone_value)
372 type(field_t), intent(inout) :: u
373 type(field_t), intent(inout) :: v
374 type(field_t), intent(inout) :: w
375 real(kind=rp), intent(in), dimension(3) :: base_value
376 character(len=*), intent(in) :: zone_name
377 real(kind=rp), intent(in) :: zone_value(:)
378 character(len=LOG_SIZE) :: log_buf
379
380 ! Internal variables
381 class(point_zone_t), pointer :: zone
382 integer :: size
383
384 call neko_log%message("Type : point_zone")
385 write (log_buf, '(A,ES12.6)') "Base value : ", base_value
386 call neko_log%message(log_buf)
387 call neko_log%message("Zone name : " // trim(zone_name))
388 write (log_buf, '(A,"[",2(ES12.6,","),ES12.6," ]")') "Value : ", &
389 zone_value(1), zone_value(2), zone_value(3)
390 call neko_log%message(log_buf)
391
392 call set_flow_ic_uniform(u, v, w, base_value)
393 size = u%dof%size()
394
395 zone => neko_point_zone_registry%get_point_zone(trim(zone_name))
396
397 call cfill_mask(u%x, zone_value(1), size, zone%mask%get(), zone%size)
398 call cfill_mask(v%x, zone_value(2), size, zone%mask%get(), zone%size)
399 call cfill_mask(w%x, zone_value(3), size, zone%mask%get(), zone%size)
400
401 end subroutine set_flow_ic_point_zone
402
418 subroutine set_flow_ic_fld(u, v, w, p, file_name, &
419 interpolate, mesh_file_name, global_interp_subdict)
420 type(field_t), target, intent(inout) :: u
421 type(field_t), target, intent(inout) :: v
422 type(field_t), target, intent(inout) :: w
423 type(field_t), target, intent(inout) :: p
424 character(len=*), intent(inout) :: file_name
425 logical, intent(in) :: interpolate
426 character(len=*), intent(inout) :: mesh_file_name
427 type(json_file), intent(inout) :: global_interp_subdict
428
429 type(field_t), pointer :: us, vs, ws, ps
430
431 us => u
432 vs => v
433 ws => w
434 ps => p
435
436 call import_fields(file_name, global_interp_subdict, mesh_file_name, &
437 u = us, v = vs, w = ws, p = ps, &
438 interpolate = interpolate)
439
440 nullify(us, vs, ws, ps)
441
442 ! If we are on GPU we need to move (u,v,w) and p back to the host
443 ! since set_flow_ic_common copies it again to the device.
444 call u%copy_from(device_to_host, sync = .false.)
445 call v%copy_from(device_to_host, sync = .false.)
446 call w%copy_from(device_to_host, sync = .false.)
447 call p%copy_from(device_to_host, sync = .true.)
448
449 end subroutine set_flow_ic_fld
450
451end 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:113
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:67
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:157
subroutine, public set_flow_ic_fld(u, v, w, p, file_name, interpolate, mesh_file_name, global_interp_subdict)
Set the initial condition of the flow based on a field. @detail The fields are read from an fld file....
Definition flow_ic.f90:420
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:372
subroutine set_flow_ic_int(u, v, w, p, coef, gs, type, params)
Set initial flow condition (builtin)
Definition flow_ic.f90:78
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:187
subroutine set_flow_ic_uniform(u, v, w, uinf)
Uniform initial condition.
Definition flow_ic.f90:275
subroutine set_flow_ic_common(u, v, w, p, coef, gs)
Definition flow_ic.f90:235
subroutine set_flow_ic_blasius(u, v, w, delta, uinf, type)
Set a Blasius profile as initial condition.
Definition flow_ic.f90:302
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.
Routines to interpolate between different spaces.
Utilities for retrieving parameters from the case files.
subroutine, public json_get_subdict_or_empty(json, key, output)
Extract a sub-object from a json object and returns an empty object if the key is missing.
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:488
subroutine, public col2(a, b, n)
Vector multiplication .
Definition math.f90:855
subroutine, public cfill_mask(a, c, n, mask, n_mask)
Fill a constant to a masked vector. .
Definition math.f90:398
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:62
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:56
Interface for NEKTON fld files.
Definition fld_file.f90:65
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