Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
import_field_utils.f90
Go to the documentation of this file.
1! Copyright (c) 2019-2026, 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!
38 use file, only : file_t
39 use num_types, only : rp, dp
40 use field, only : field_t
41 use field_list, only : field_list_t
45 use logger, only : log_size, neko_log
46 use device, only : host_to_device
47 use json_module, only : json_file
49 implicit none
50 private
51
52 public :: import_fields
53
56 end interface import_fields
57
58contains
59
91 subroutine import_fields_from_json(fname, global_interp_subdict, mesh_fname, &
92 u, v, w, p, t, s_target_list, s_index_list, interpolate)
93 character(len=*), intent(in) :: fname
94 type(json_file), intent(inout) :: global_interp_subdict
95 character(len=*), intent(in), optional :: mesh_fname
96 type(field_t), pointer, intent(inout), optional :: u,v,w,p,t
97 type(field_list_t), intent(inout), optional :: s_target_list
98 integer, intent(in), optional :: s_index_list(:)
99 logical, intent(in), optional :: interpolate
100
101 real(kind=dp) :: tolerance, padding
102
103 call json_get_or_default(global_interp_subdict, "tolerance", &
104 tolerance, glob_interp_tol)
105 call json_get_or_default(global_interp_subdict, "padding", &
106 padding, glob_interp_pad)
107
108 call import_fields_from_params(fname, mesh_fname, &
109 u, v, w, p, t, s_target_list, s_index_list, &
110 interpolate, tolerance = tolerance, padding = padding)
111
112 end subroutine import_fields_from_json
113
147 subroutine import_fields_from_params(fname, mesh_fname, u, v, w, p, t, &
148 s_target_list, s_index_list, interpolate, tolerance, padding)
149 character(len=*), intent(in) :: fname
150 character(len=*), intent(in), optional :: mesh_fname
151 type(field_t), pointer, intent(inout), optional :: u,v,w,p,t
152 type(field_list_t), intent(inout), optional :: s_target_list
153 integer, intent(in), optional :: s_index_list(:)
154 logical, intent(in), optional :: interpolate
155 real(kind=dp), intent(in), optional :: tolerance
156 real(kind=dp), intent(in), optional :: padding
157
158 character(len=LOG_SIZE) :: log_buf
159 integer :: sample_idx, sample_mesh_idx, i
160 character(len=NEKO_FNAME_LEN) :: fname_, mesh_fname_
161
162 logical :: interpolate_
163
164 type(file_t) :: f
165 type(fld_file_data_t) :: fld_data
166
167 ! ---- Default values
168 interpolate_ = .false.
169 if (present(interpolate)) interpolate_ = interpolate
170 mesh_fname_ = "none"
171 if (present(mesh_fname)) mesh_fname_ = trim(mesh_fname)
172 ! ----
173
174 call neko_log%section("Import fields")
175 call neko_log%message("File name : " // trim(fname))
176 write (log_buf, '(A,L1)') "Interpolation : ", interpolate_
177 call neko_log%message(log_buf)
178
179 !
180 ! Handling of file names and reading of data
181 !
182
183 ! Extract sample index from the file name
184 sample_idx = extract_fld_file_index(fname, -1)
185
186 if (sample_idx .eq. -1) &
187 call neko_error("Invalid file name for the initial condition. The&
188 & file format must be e.g. 'mean0.f00001'")
189
190 ! Change from "field0.f000*" to "field0.fld" for the fld reader
191 call filename_chsuffix(fname, fname_, 'fld')
192
193 ! Initialize file object
194 call f%init(trim(fname_))
195
196 ! If interpolate, check if we need to read the mesh file
197 if (interpolate_) then
198
199 ! If no mesh file is specified, use the default file name
200 if (mesh_fname_ .eq. "none") then
201 mesh_fname_ = trim(fname_)
202 sample_mesh_idx = sample_idx
203 else
204 mesh_fname_ = trim(mesh_fname)
205
206 ! Extract sample index from the mesh file name
207 sample_mesh_idx = extract_fld_file_index(mesh_fname_, -1)
208
209 if (sample_mesh_idx .eq. -1) then
210 call neko_error("Invalid file name for the initial condition." // &
211 "The file format must be e.g. 'mean0.f00001'")
212 end if
213
214 write (log_buf, '(A,A)') "Mesh file : ", &
215 trim(mesh_fname_)
216 call neko_log%message(log_buf)
217
218 end if ! if mesh_file_name .eq. none
219
220 ! Read the mesh coordinates if they are not in our fld file
221 if (sample_mesh_idx .ne. sample_idx) then
222 call f%set_counter(sample_mesh_idx)
223 call f%read(fld_data)
224 end if
225
226 end if
227
228 ! Read the field file containing (u,v,w,p)
229 call f%set_counter(sample_idx)
230 call f%read(fld_data)
231
232 !
233 ! Copy all vectors to device (GPU) since everything is read on the CPU
234 !
235 if (present(u)) call fld_data%u%copy_from(host_to_device, .true.)
236 if (present(v)) call fld_data%v%copy_from(host_to_device, .true.)
237 if (present(w)) call fld_data%w%copy_from(host_to_device, .true.)
238 if (present(p)) call fld_data%p%copy_from(host_to_device, .true.)
239 if (present(t)) call fld_data%t%copy_from(host_to_device, .true.)
240 if (present(s_target_list)) then
241
242 if (present(s_index_list)) then
243 if (size(s_index_list) .ne. s_target_list%size()) then
244 call neko_error("Scalar lists must have same size!")
245 end if
246
247 do i = 1, size(s_index_list)
248 ! Take care that if we set i=0 we want temperature
249 if (s_index_list(i) .eq. 0) then
250 call fld_data%t%copy_from(host_to_device, .true.)
251 else
252 ! For scalar fields, require indices in 1:this%n_scalars
253 if (s_index_list(i) < 1 .or. &
254 s_index_list(i) > fld_data%n_scalars) then
255 call neko_error("s_index_list entry out of bounds")
256 end if
257 call fld_data%s(s_index_list(i))%copy_from(host_to_device, &
258 .true.)
259 end if
260 end do
261 else
262 do i = 1, s_target_list%size()
263 call fld_data%s(i)%copy_from(host_to_device, .true.)
264 end do
265 end if ! if present s_index_list
266
267 end if ! present s_tgt
268
269
270 if (interpolate_) then
271 ! Sync coordinates to device for the interpolation
272 call fld_data%x%copy_from(host_to_device, .false.)
273 call fld_data%y%copy_from(host_to_device, .false.)
274 call fld_data%z%copy_from(host_to_device, .true.)
275 end if
276
277 ! Call the import of fields
278 call fld_data%import_fields(u, v, w, p, t, s_target_list, s_index_list, &
279 interpolate_, tolerance = tolerance, padding = padding)
280
281 call neko_log%end_section()
282 call fld_data%free()
283
284 end subroutine import_fields_from_params
285
286end module import_field_utils
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Device abstraction, common interface for various accelerators.
Definition device.F90:34
integer, parameter, public host_to_device
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...
Implements global_interpolation given a dofmap.
real(kind=dp), parameter, public glob_interp_tol
real(kind=dp), parameter, public glob_interp_pad
Importation of fields from fld files.
subroutine import_fields_from_json(fname, global_interp_subdict, mesh_fname, u, v, w, p, t, s_target_list, s_index_list, interpolate)
Imports fields from an fld file, potentially with interpolation, with parameters provided in a JSON s...
subroutine import_fields_from_params(fname, mesh_fname, u, v, w, p, t, s_target_list, s_index_list, interpolate, tolerance, padding)
Imports fields from an fld file, potentially with interpolation.
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
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
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 filename_chsuffix(fname, new_fname, new_suffix)
Change a filename's suffix.
Definition utils.f90:140
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