Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
scalar_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!
35 use gather_scatter, only : gs_t, gs_op_add
37 use num_types, only : rp
38 use device_math, only : device_col2
40 use field, only : field_t
43 use coefs, only : coef_t
44 use math, only : col2, cfill, cfill_mask
46 use json_module, only : json_file
48 use point_zone, only: point_zone_t
50 use logger, only: neko_log, log_size
52 use fld_file, only: fld_file_t
53 use checkpoint, only: chkp_t
54 use file, only: file_t
57 use space, only: space_t, gll
58 use field_list, only : field_list_t
59 implicit none
60 private
61
62 interface set_scalar_ic
63 module procedure set_scalar_ic_int, set_scalar_ic_usr
64 end interface set_scalar_ic
65
66 public :: set_scalar_ic
67
68contains
69
82 subroutine set_scalar_ic_int(s, coef, gs, type, params, i)
83 type(field_t), intent(inout) :: s
84 type(coef_t), intent(in) :: coef
85 type(gs_t), intent(inout) :: gs
86 character(len=*) :: type
87 type(json_file), intent(inout) :: params
88 integer, intent(in) :: i
89
90 ! Variables for retrieving JSON parameters
91 real(kind=rp) :: ic_value
92 character(len=:), allocatable :: read_str
93 character(len=NEKO_FNAME_LEN) :: fname, mesh_fname
94 real(kind=rp) :: zone_value, tol
95 logical :: interpolate
96
97 if (trim(type) .eq. 'uniform') then
98
99 call json_get(params, 'value', ic_value)
100 call set_scalar_ic_uniform(s, ic_value)
101
102 else if (trim(type) .eq. 'point_zone') then
103
104 call json_get(params, 'base_value', ic_value)
105 call json_get(params, 'zone_name', read_str)
106 call json_get(params, 'zone_value', zone_value)
107
108 call set_scalar_ic_point_zone(s, ic_value, read_str, zone_value)
109
110 else if (trim(type) .eq. 'field') then
111
112 call json_get(params, 'file_name', read_str)
113 fname = trim(read_str)
114 call json_get_or_default(params, 'interpolate', interpolate, &
115 .false.)
116 call json_get_or_default(params, 'tolerance', tol, 0.000001_rp)
117 call json_get_or_default(params, 'mesh_file_name', read_str, &
118 "none")
119 mesh_fname = trim(read_str)
120
121 call set_scalar_ic_fld(s, fname, interpolate, tol, mesh_fname, i)
122
123 else
124 call neko_error('Invalid initial condition')
125 end if
126
127 call set_scalar_ic_common(s, coef, gs)
128
129 end subroutine set_scalar_ic_int
130
138 subroutine set_scalar_ic_usr(scheme_name, s, coef, gs, user_proc)
139 character(len=*), intent(in) :: scheme_name
140 type(field_t), target, intent(inout) :: s
141 type(coef_t), intent(in) :: coef
142 type(gs_t), intent(inout) :: gs
143 procedure(user_initial_conditions_intf) :: user_proc
144 type(field_list_t) :: fields
145
146 call neko_log%message("Type: user")
147
148 call fields%init(1)
149 call fields%assign_to_field(1, s)
150
151 call user_proc(scheme_name, fields)
152 call set_scalar_ic_common(s, coef, gs)
153
154 end subroutine set_scalar_ic_usr
155
162 subroutine set_scalar_ic_common(s, coef, gs)
163 type(field_t), intent(inout) :: s
164 type(coef_t), intent(in) :: coef
165 type(gs_t), intent(inout) :: gs
166 integer :: n
167
168 n = s%dof%size()
169 if (neko_bcknd_device .eq. 1) then
170 call device_memcpy(s%x, s%x_d, n, host_to_device, sync = .false.)
171 end if
172
173 ! Ensure continuity across elements for initial conditions
174 call gs%op(s%x, n, gs_op_add)
175
176 if (neko_bcknd_device .eq. 1) then
177 call device_col2(s%x_d, coef%mult_d, n)
178 else
179 call col2(s%x, coef%mult, n)
180 end if
181
182 end subroutine set_scalar_ic_common
183
188 subroutine set_scalar_ic_uniform(s, ic_value)
189 type(field_t), intent(inout) :: s
190 real(kind=rp), intent(in) :: ic_value
191 integer :: n
192 character(len=LOG_SIZE) :: log_buf
193
194 call neko_log%message("Type : uniform")
195 write (log_buf, '(A,ES12.6)') "Value: ", ic_value
196 call neko_log%message(log_buf)
197
198 s = ic_value
199 n = s%dof%size()
200 if (neko_bcknd_device .eq. 1) then
201 call cfill(s%x, ic_value, n)
202 end if
203
204 end subroutine set_scalar_ic_uniform
205
213 subroutine set_scalar_ic_point_zone(s, base_value, zone_name, zone_value)
214 type(field_t), intent(inout) :: s
215 real(kind=rp), intent(in) :: base_value
216 character(len=*), intent(in) :: zone_name
217 real(kind=rp), intent(in) :: zone_value
218
219 ! Internal variables
220 character(len=LOG_SIZE) :: log_buf
221 class(point_zone_t), pointer :: zone
222 integer :: size
223
224 call neko_log%message("Type : point_zone")
225 write (log_buf, '(A,ES12.6)') "Base value: ", base_value
226 call neko_log%message(log_buf)
227 call neko_log%message("Zone name : " // trim(zone_name))
228 write (log_buf, '(A,ES12.6)') "Zone value: ", zone_value
229 call neko_log%message(log_buf)
230
231 size = s%dof%size()
232 zone => neko_point_zone_registry%get_point_zone(trim(zone_name))
233
234 call set_scalar_ic_uniform(s, base_value)
235 call cfill_mask(s%x, zone_value, size, zone%mask%get(), zone%size)
236
237 end subroutine set_scalar_ic_point_zone
238
254 subroutine set_scalar_ic_fld(s, file_name, &
255 interpolate, tolerance, mesh_file_name, i)
256 type(field_t), intent(inout) :: s
257 character(len=*), intent(in) :: file_name
258 logical, intent(in) :: interpolate
259 real(kind=rp), intent(in) :: tolerance
260 character(len=*), intent(inout) :: mesh_file_name
261 integer, intent(in) :: i
262
263 character(len=LOG_SIZE) :: log_buf
264 integer :: sample_idx, sample_mesh_idx
265 integer :: last_index
266 type(fld_file_data_t) :: fld_data
267 type(file_t) :: f
268 logical :: mesh_mismatch
269
270 ! ---- For the mesh to mesh interpolation
271 type(global_interpolation_t) :: global_interp
272 ! -----
273
274 ! ---- For space to space interpolation
275 type(space_t) :: prev_Xh
276 type(interpolator_t) :: space_interp
277 ! ----
278
279 call neko_log%message("Type : field")
280 call neko_log%message("File name : " // trim(file_name))
281 write (log_buf, '(A,L1)') "Interpolation : ", interpolate
282 call neko_log%message(log_buf)
283
284 ! Extract sample index from the file name
285 sample_idx = extract_fld_file_index(file_name, -1)
286
287 if (sample_idx .eq. -1) then
288 call neko_error("Invalid file name for the initial condition. The " // &
289 "file format must be e.g. 'mean0.f00001'")
290 end if
291
292 ! Change from "field0.f000*" to "field0.fld" for the fld reader
293 call filename_chsuffix(file_name, file_name, 'fld')
294
295 call fld_data%init
296 call f%init(trim(file_name))
297
298 if (interpolate) then
299
300 ! If no mesh file is specified, use the default file name
301 if (mesh_file_name .eq. "none") then
302 mesh_file_name = trim(file_name)
303 sample_mesh_idx = sample_idx
304 else
305
306 ! Extract sample index from the mesh file name
307 sample_mesh_idx = extract_fld_file_index(mesh_file_name, -1)
308
309 if (sample_mesh_idx .eq. -1) then
310 call neko_error("Invalid file name for the initial condition." // &
311 " The file format must be e.g. 'mean0.f00001'")
312 end if
313
314 write (log_buf, '(A,ES12.6)') "Tolerance : ", tolerance
315 call neko_log%message(log_buf)
316 write (log_buf, '(A,A)') "Mesh file : ", &
317 trim(mesh_file_name)
318 call neko_log%message(log_buf)
319
320 end if ! if mesh_file_name .eq. none
321
322 ! Read the mesh coordinates if they are not in our fld file
323 if (sample_mesh_idx .ne. sample_idx) then
324 call f%set_counter(sample_mesh_idx)
325 call f%read(fld_data)
326 end if
327
328 end if
329
330 ! Read the field file containing (u,v,w,p)
331 call f%set_counter(sample_idx)
332 call f%read(fld_data)
333
334 !
335 ! Check that the data in the fld file matches the current case.
336 ! Note that this is a safeguard and there are corner cases where
337 ! two different meshes have the same dimension and same # of elements
338 ! but this should be enough to cover obvious cases.
339 !
340 mesh_mismatch = (fld_data%glb_nelv .ne. s%msh%glb_nelv .or. &
341 fld_data%gdim .ne. s%msh%gdim)
342
343 if (mesh_mismatch .and. .not. interpolate) then
344 call neko_error("The fld file must match the current mesh! " // &
345 "Use 'interpolate': 'true' to enable interpolation.")
346 else if (.not. mesh_mismatch .and. interpolate) then
347 call neko_log%warning("You have activated interpolation but you " // &
348 "might still be using the same mesh.")
349 end if
350
351
352 ! Mesh interpolation if specified
353 if (interpolate) then
354 ! Issue a warning if the mesh is in single precision
355 select type (ft => f%file_type)
356 type is (fld_file_t)
357 if (.not. ft%dp_precision) then
358 call neko_warning("The coordinates read from the field file " // &
359 "are in single precision.")
360 call neko_log%message("It is recommended to use a mesh in " // &
361 "double precision for better interpolation results.")
362 call neko_log%message("If the interpolation does not work, " // &
363 "you can try to increase the tolerance.")
364 end if
365 class default
366 end select
367
368 ! Copy all fld data to device since the reader loads everything on the host
369 call fld_data%x%copy_from(host_to_device, .false.)
370 call fld_data%y%copy_from(host_to_device, .false.)
371 call fld_data%z%copy_from(host_to_device, .false.)
372 call fld_data%t%copy_from(host_to_device, .true.)
373
374 ! Generates an interpolator object and performs the point search
375 call fld_data%generate_interpolator(global_interp, s%dof, s%msh, &
376 tolerance)
377
378 ! Evaluate scalar
379
380 ! i == 0 means it's the temperature field
381 if (i .ne. 0) then
382 call global_interp%evaluate(s%x, fld_data%s(i)%x, .false.)
383 else
384 call global_interp%evaluate(s%x, fld_data%t%x, .false.)
385 end if
386
387 call global_interp%free
388
389 ! Copy back to the host for set_scalar_ic_common
390 call fld_data%t%copy_from(device_to_host, .true.)
391
392 else ! No interpolation, just potentially from different spaces
393
394 ! Build a space_t object from the data in the fld file
395 call prev_xh%init(gll, fld_data%lx, fld_data%ly, fld_data%lz)
396 call space_interp%init(s%Xh, prev_xh)
397
398 ! Do the space-to-space interpolation
399 ! i == 0 means it's the temperature field
400 if (i .ne. 0) then
401 call space_interp%map_host(s%x, fld_data%s(i)%x, fld_data%nelv, s%Xh)
402 else
403 call space_interp%map_host(s%x, fld_data%t%x, fld_data%nelv, s%Xh)
404 end if
405
406 call space_interp%free
407
408 end if
409
410 call fld_data%free
411 call prev_xh%free
412
413 end subroutine set_scalar_ic_fld
414
415end module scalar_ic
Copy data between host and device (or device and device)
Definition device.F90:71
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
Defines a checkpoint.
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
Gather-scatter.
Implements global_interpolation given a dofmap.
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:76
integer, parameter, public log_size
Definition log.f90:46
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
type(point_zone_registry_t), target, public neko_point_zone_registry
Global point_zone registry.
Scalar initial condition.
Definition scalar_ic.f90:34
subroutine set_scalar_ic_uniform(s, ic_value)
Uniform initial condition.
subroutine set_scalar_ic_point_zone(s, base_value, zone_name, zone_value)
Point zone initial condition.
subroutine set_scalar_ic_common(s, coef, gs)
Set scalar initial condition (common)
subroutine set_scalar_ic_int(s, coef, gs, type, params, i)
Set scalar initial condition (builtin)
Definition scalar_ic.f90:83
subroutine set_scalar_ic_fld(s, file_name, interpolate, tolerance, mesh_file_name, i)
Set the initial condition of the scalar based on a field. @detail The field is read from an fld file....
subroutine set_scalar_ic_usr(scheme_name, s, coef, gs, user_proc)
Set scalar intial condition (user defined)
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
subroutine, public filename_suffix(fname, suffix)
Extract a filename's suffix.
Definition utils.f90:108
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:56
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