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