Neko 1.99.3
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
49 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
61 implicit none
62 private
63
64 interface set_scalar_ic
65 module procedure set_scalar_ic_int, set_scalar_ic_usr
66 end interface set_scalar_ic
67
68 public :: set_scalar_ic
69
70contains
71
84 subroutine set_scalar_ic_int(s, coef, gs, type, params, i)
85 type(field_t), intent(inout) :: s
86 type(coef_t), intent(in) :: coef
87 type(gs_t), intent(inout) :: gs
88 character(len=*) :: type
89 type(json_file), intent(inout) :: params
90 integer, intent(in) :: i
91
92 ! Variables for retrieving JSON parameters
93 real(kind=rp) :: ic_value
94 character(len=:), allocatable :: read_str
95 character(len=NEKO_FNAME_LEN) :: fname, mesh_fname
96 real(kind=rp) :: zone_value, tol
97 logical :: interpolate
98 integer :: tgt_scal_idx
99
100 if (trim(type) .eq. 'uniform') then
101
102 call json_get_or_lookup(params, 'value', ic_value)
103 call set_scalar_ic_uniform(s, ic_value)
104
105 else if (trim(type) .eq. 'point_zone') then
106
107 call json_get_or_lookup(params, 'base_value', ic_value)
108 call json_get(params, 'zone_name', read_str)
109 call json_get_or_lookup(params, 'zone_value', zone_value)
110
111 call set_scalar_ic_point_zone(s, ic_value, read_str, zone_value)
112
113 else if (trim(type) .eq. 'field') then
114
115 call json_get(params, 'file_name', read_str)
116 fname = trim(read_str)
117 call json_get_or_default(params, 'interpolate', interpolate, &
118 .false.)
119 call json_get_or_lookup_or_default(params, 'tolerance', tol, 0.000001_rp)
120 call json_get_or_default(params, 'mesh_file_name', read_str, &
121 "none")
122 mesh_fname = trim(read_str)
123
124 ! Give the user the option to select which scalar they want to import
125 ! the values from, in the fld file. 0 corresponds to temperature.
126 call json_get_or_default(params, 'target_index', tgt_scal_idx, i)
127
128 call set_scalar_ic_fld(s, fname, interpolate, tol, mesh_fname, i, &
129 tgt_scal_idx)
130
131 else
132 call neko_error('Invalid initial condition')
133 end if
134
135 call set_scalar_ic_common(s, coef, gs)
136
137 end subroutine set_scalar_ic_int
138
146 subroutine set_scalar_ic_usr(scheme_name, s, coef, gs, user_proc)
147 character(len=*), intent(in) :: scheme_name
148 type(field_t), target, intent(inout) :: s
149 type(coef_t), intent(in) :: coef
150 type(gs_t), intent(inout) :: gs
151 procedure(user_initial_conditions_intf) :: user_proc
152 type(field_list_t) :: fields
153
154 call neko_log%message("Type: user")
155
156 call fields%init(1)
157 call fields%assign_to_field(1, s)
158
159 call user_proc(scheme_name, fields)
160 call set_scalar_ic_common(s, coef, gs)
161
162 end subroutine set_scalar_ic_usr
163
170 subroutine set_scalar_ic_common(s, coef, gs)
171 type(field_t), intent(inout) :: s
172 type(coef_t), intent(in) :: coef
173 type(gs_t), intent(inout) :: gs
174 integer :: n
175
176 n = s%dof%size()
177 if (neko_bcknd_device .eq. 1) then
178 call device_memcpy(s%x, s%x_d, n, host_to_device, sync = .false.)
179 end if
180
181 ! Ensure continuity across elements for initial conditions
182 call gs%op(s%x, n, gs_op_add)
183
184 if (neko_bcknd_device .eq. 1) then
185 call device_col2(s%x_d, coef%mult_d, n)
186 else
187 call col2(s%x, coef%mult, n)
188 end if
189
190 end subroutine set_scalar_ic_common
191
196 subroutine set_scalar_ic_uniform(s, ic_value)
197 type(field_t), intent(inout) :: s
198 real(kind=rp), intent(in) :: ic_value
199 integer :: n
200 character(len=LOG_SIZE) :: log_buf
201
202 call neko_log%message("Type : uniform")
203 write (log_buf, '(A,ES12.6)') "Value: ", ic_value
204 call neko_log%message(log_buf)
205
206 s = ic_value
207 n = s%dof%size()
208 if (neko_bcknd_device .eq. 1) then
209 call cfill(s%x, ic_value, n)
210 end if
211
212 end subroutine set_scalar_ic_uniform
213
221 subroutine set_scalar_ic_point_zone(s, base_value, zone_name, zone_value)
222 type(field_t), intent(inout) :: s
223 real(kind=rp), intent(in) :: base_value
224 character(len=*), intent(in) :: zone_name
225 real(kind=rp), intent(in) :: zone_value
226
227 ! Internal variables
228 character(len=LOG_SIZE) :: log_buf
229 class(point_zone_t), pointer :: zone
230 integer :: size
231
232 call neko_log%message("Type : point_zone")
233 write (log_buf, '(A,ES12.6)') "Base value: ", base_value
234 call neko_log%message(log_buf)
235 call neko_log%message("Zone name : " // trim(zone_name))
236 write (log_buf, '(A,ES12.6)') "Zone value: ", zone_value
237 call neko_log%message(log_buf)
238
239 size = s%dof%size()
240 zone => neko_point_zone_registry%get_point_zone(trim(zone_name))
241
242 call set_scalar_ic_uniform(s, base_value)
243 call cfill_mask(s%x, zone_value, size, zone%mask%get(), zone%size)
244
245 end subroutine set_scalar_ic_point_zone
246
264 subroutine set_scalar_ic_fld(s, file_name, &
265 interpolate, tolerance, mesh_file_name, i, target_idx)
266 type(field_t), target, intent(inout) :: s
267 character(len=*), intent(in) :: file_name
268 logical, intent(in) :: interpolate
269 real(kind=rp), intent(in) :: tolerance
270 character(len=*), intent(inout) :: mesh_file_name
271 integer, intent(in) :: i
272 integer, intent(in) :: target_idx
273
274 character(len=LOG_SIZE) :: log_buf
275 type(field_t), pointer :: ss
276 type(field_list_t) :: s_tgt_list
277
278 if (i .ne. target_idx) then
279 write (log_buf, '(A,I0,A,I0)') "Loading scalar #", target_idx, &
280 " into scalar #", i
281 call neko_log%message(log_buf)
282 end if
283
284 ! use a pointer since import_fields needs a pointer as input
285 ss => s
286
287 ! Put ss in a field list of 1 element
288 call s_tgt_list%init(1)
289 call s_tgt_list%assign(1, ss)
290
291 call import_fields(file_name, mesh_file_name, &
292 s_target_list = s_tgt_list, & ! The target field
293 s_index_list = [target_idx], & ! Take values from target scalar
294 interpolate = interpolate, tolerance = tolerance)
295
296 call s_tgt_list%free()
297
298 nullify(ss)
299
300 ! If we are on GPU we need to move s back to the host
301 ! since set_scalar_ic_common copies it again to the device.
302 call s%copy_from(device_to_host, .true.)
303
304 end subroutine set_scalar_ic_fld
305
306end 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.
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
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:85
subroutine set_scalar_ic_fld(s, file_name, interpolate, tolerance, mesh_file_name, i, target_idx)
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: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