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
50 use point_zone, only : point_zone_t
52 use logger, only : neko_log, log_size
54 use fld_file, only : fld_file_t
55 use checkpoint, only : chkp_t
56 use file, only : file_t
59 use space, only : space_t, gll
60 use field_list, only : field_list_t
62 implicit none
63 private
64
65 interface set_scalar_ic
66 module procedure set_scalar_ic_int, set_scalar_ic_usr
67 end interface set_scalar_ic
68
69 public :: set_scalar_ic
70
71contains
72
85 subroutine set_scalar_ic_int(s, coef, gs, type, params, i)
86 type(field_t), intent(inout) :: s
87 type(coef_t), intent(in) :: coef
88 type(gs_t), intent(inout) :: gs
89 character(len=*) :: type
90 type(json_file), intent(inout) :: params
91 integer, intent(in) :: i
92
93 ! Variables for retrieving JSON parameters
94 real(kind=rp) :: ic_value
95 character(len=:), allocatable :: read_str
96 real(kind=rp) :: zone_value
97
98 if (trim(type) .eq. 'uniform') then
99
100 call json_get_or_lookup(params, 'value', ic_value)
101 call set_scalar_ic_uniform(s, ic_value)
102
103 else if (trim(type) .eq. 'point_zone') then
104
105 call json_get_or_lookup(params, 'base_value', ic_value)
106 call json_get(params, 'zone_name', read_str)
107 call json_get_or_lookup(params, 'zone_value', zone_value)
108
109 call set_scalar_ic_point_zone(s, ic_value, read_str, zone_value)
110
111 else if (trim(type) .eq. 'field') then
112
113 block
114 character(len=NEKO_FNAME_LEN) :: fname, mesh_fname
115 logical :: interpolate
116 type(json_file) :: interp_subdict
117 integer :: tgt_scal_idx
118
119 call json_get(params, 'file_name', read_str)
120 fname = trim(read_str)
121
122 call json_get_or_default(params, 'interpolate', interpolate, &
123 .false.)
124
125 call json_get_or_default(params, 'mesh_file_name', read_str, &
126 "none")
127 mesh_fname = trim(read_str)
128
129 ! Give the user the option to select which scalar they want to import
130 ! the values from, in the fld file. 0 corresponds to temperature.
131 call json_get_or_default(params, 'target_index', tgt_scal_idx, i)
132
133 call json_get_subdict_or_empty(params, "interpolation", &
134 interp_subdict)
135
136 call set_scalar_ic_fld(s, fname, interpolate, mesh_fname, i, &
137 tgt_scal_idx, interp_subdict)
138
139 end block
140
141 else
142 call neko_error('Invalid initial condition')
143 end if
144
145 call set_scalar_ic_common(s, coef, gs)
146
147 end subroutine set_scalar_ic_int
148
156 subroutine set_scalar_ic_usr(scheme_name, s, coef, gs, user_proc)
157 character(len=*), intent(in) :: scheme_name
158 type(field_t), target, intent(inout) :: s
159 type(coef_t), intent(in) :: coef
160 type(gs_t), intent(inout) :: gs
161 procedure(user_initial_conditions_intf) :: user_proc
162 type(field_list_t) :: fields
163
164 call neko_log%message("Type: user")
165
166 call fields%init(1)
167 call fields%assign_to_field(1, s)
168
169 call user_proc(scheme_name, fields)
170 call set_scalar_ic_common(s, coef, gs)
171
172 end subroutine set_scalar_ic_usr
173
180 subroutine set_scalar_ic_common(s, coef, gs)
181 type(field_t), intent(inout) :: s
182 type(coef_t), intent(in) :: coef
183 type(gs_t), intent(inout) :: gs
184 integer :: n
185
186 n = s%dof%size()
187 if (neko_bcknd_device .eq. 1) then
188 call device_memcpy(s%x, s%x_d, n, host_to_device, sync = .false.)
189 end if
190
191 ! Ensure continuity across elements for initial conditions
192 call gs%op(s%x, n, gs_op_add)
193
194 if (neko_bcknd_device .eq. 1) then
195 call device_col2(s%x_d, coef%mult_d, n)
196 else
197 call col2(s%x, coef%mult, n)
198 end if
199
200 end subroutine set_scalar_ic_common
201
206 subroutine set_scalar_ic_uniform(s, ic_value)
207 type(field_t), intent(inout) :: s
208 real(kind=rp), intent(in) :: ic_value
209 integer :: n
210 character(len=LOG_SIZE) :: log_buf
211
212 call neko_log%message("Type : uniform")
213 write (log_buf, '(A,ES12.6)') "Value: ", ic_value
214 call neko_log%message(log_buf)
215
216 s = ic_value
217 n = s%dof%size()
218 if (neko_bcknd_device .eq. 1) then
219 call cfill(s%x, ic_value, n)
220 end if
221
222 end subroutine set_scalar_ic_uniform
223
231 subroutine set_scalar_ic_point_zone(s, base_value, zone_name, zone_value)
232 type(field_t), intent(inout) :: s
233 real(kind=rp), intent(in) :: base_value
234 character(len=*), intent(in) :: zone_name
235 real(kind=rp), intent(in) :: zone_value
236
237 ! Internal variables
238 character(len=LOG_SIZE) :: log_buf
239 class(point_zone_t), pointer :: zone
240 integer :: size
241
242 call neko_log%message("Type : point_zone")
243 write (log_buf, '(A,ES12.6)') "Base value: ", base_value
244 call neko_log%message(log_buf)
245 call neko_log%message("Zone name : " // trim(zone_name))
246 write (log_buf, '(A,ES12.6)') "Zone value: ", zone_value
247 call neko_log%message(log_buf)
248
249 size = s%dof%size()
250 zone => neko_point_zone_registry%get_point_zone(trim(zone_name))
251
252 call set_scalar_ic_uniform(s, base_value)
253 call cfill_mask(s%x, zone_value, size, zone%mask%get(), zone%size)
254
255 end subroutine set_scalar_ic_point_zone
256
272 subroutine set_scalar_ic_fld(s, file_name, &
273 interpolate, mesh_file_name, i, target_idx, global_interp_subdict)
274 type(field_t), target, intent(inout) :: s
275 character(len=*), intent(in) :: file_name
276 logical, intent(in) :: interpolate
277 character(len=*), intent(inout) :: mesh_file_name
278 integer, intent(in) :: i
279 integer, intent(in) :: target_idx
280 type(json_file), intent(inout) :: global_interp_subdict
281
282 character(len=LOG_SIZE) :: log_buf
283 type(field_t), pointer :: ss
284 type(field_list_t) :: s_tgt_list
285
286 if (i .ne. target_idx) then
287 write (log_buf, '(A,I0,A,I0)') "Loading scalar #", target_idx, &
288 " into scalar #", i
289 call neko_log%message(log_buf)
290 end if
291
292 ! use a pointer since import_fields needs a pointer as input
293 ss => s
294
295 ! Put ss in a field list of 1 element
296 call s_tgt_list%init(1)
297 call s_tgt_list%assign(1, ss)
298
299 call import_fields(file_name, global_interp_subdict, mesh_file_name, &
300 s_target_list = s_tgt_list, & ! The target field
301 s_index_list = [target_idx], & ! Take values from target scalar
302 interpolate = interpolate)
303
304 call s_tgt_list%free()
305
306 nullify(ss)
307
308 ! If we are on GPU we need to move s back to the host
309 ! since set_scalar_ic_common copies it again to the device.
310 call s%copy_from(device_to_host, .true.)
311
312 end subroutine set_scalar_ic_fld
313
314end 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:67
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.
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
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_fld(s, file_name, interpolate, mesh_file_name, i, target_idx, global_interp_subdict)
Set the initial condition of the scalar based on a field. @detail The field is read from an fld file....
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:86
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: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