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
40 use field, only : field_t
41 use field_list, only : field_list_t
44 use logger, only : log_size, neko_log
45 use device, only : host_to_device
46 implicit none
47 private
48
49 public :: import_fields
50
51
52contains
53
86 subroutine import_fields(fname, mesh_fname, u, v, w, p, t, s_target_list, &
87 s_index_list, interpolate, tolerance)
88 character(len=*), intent(in) :: fname
89 character(len=*), intent(in), optional :: mesh_fname
90 type(field_t), pointer, intent(inout), optional :: u,v,w,p,t
91 type(field_list_t), intent(inout), optional :: s_target_list
92 integer, intent(in), optional :: s_index_list(:)
93 logical, intent(in), optional :: interpolate
94 real(kind=rp), intent(in), optional :: tolerance
95
96 character(len=LOG_SIZE) :: log_buf
97 integer :: sample_idx, sample_mesh_idx, i
98 character(len=NEKO_FNAME_LEN) :: fname_, mesh_fname_
99
100 logical :: interpolate_
101
102 type(file_t) :: f
103 type(fld_file_data_t) :: fld_data
104
105 ! ---- Default values
106 interpolate_ = .false.
107 if (present(interpolate)) interpolate_ = interpolate
108 mesh_fname_ = "none"
109 if (present(mesh_fname)) mesh_fname_ = trim(mesh_fname)
110 ! ----
111
112 call neko_log%section("Import fields")
113 call neko_log%message("File name : " // trim(fname))
114 write (log_buf, '(A,L1)') "Interpolation : ", interpolate_
115 call neko_log%message(log_buf)
116
117 !
118 ! Handling of file names and reading of data
119 !
120
121 ! Extract sample index from the file name
122 sample_idx = extract_fld_file_index(fname, -1)
123
124 if (sample_idx .eq. -1) &
125 call neko_error("Invalid file name for the initial condition. The&
126 & file format must be e.g. 'mean0.f00001'")
127
128 ! Change from "field0.f000*" to "field0.fld" for the fld reader
129 call filename_chsuffix(fname, fname_, 'fld')
130
131 ! Initialize file object
132 call f%init(trim(fname_))
133
134 ! If interpolate, check if we need to read the mesh file
135 if (interpolate_) then
136
137 if (present(tolerance)) then
138 write (log_buf, '(A,ES12.6)') "Tolerance : ", tolerance
139 call neko_log%message(log_buf)
140 end if
141
142 ! If no mesh file is specified, use the default file name
143 if (mesh_fname_ .eq. "none") then
144 mesh_fname_ = trim(fname_)
145 sample_mesh_idx = sample_idx
146 else
147 mesh_fname_ = trim(mesh_fname)
148
149 ! Extract sample index from the mesh file name
150 sample_mesh_idx = extract_fld_file_index(mesh_fname_, -1)
151
152 if (sample_mesh_idx .eq. -1) then
153 call neko_error("Invalid file name for the initial condition." // &
154 "The file format must be e.g. 'mean0.f00001'")
155 end if
156
157 write (log_buf, '(A,A)') "Mesh file : ", &
158 trim(mesh_fname_)
159 call neko_log%message(log_buf)
160
161 end if ! if mesh_file_name .eq. none
162
163 ! Read the mesh coordinates if they are not in our fld file
164 if (sample_mesh_idx .ne. sample_idx) then
165 call f%set_counter(sample_mesh_idx)
166 call f%read(fld_data)
167 end if
168
169 end if
170
171 ! Read the field file containing (u,v,w,p)
172 call f%set_counter(sample_idx)
173 call f%read(fld_data)
174
175 !
176 ! Copy all vectors to device (GPU) since everything is read on the CPU
177 !
178 if (present(u)) call fld_data%u%copy_from(host_to_device, .true.)
179 if (present(v)) call fld_data%v%copy_from(host_to_device, .true.)
180 if (present(w)) call fld_data%w%copy_from(host_to_device, .true.)
181 if (present(p)) call fld_data%p%copy_from(host_to_device, .true.)
182 if (present(t)) call fld_data%t%copy_from(host_to_device, .true.)
183 if (present(s_target_list)) then
184
185 if (present(s_index_list)) then
186 if (size(s_index_list) .ne. s_target_list%size()) then
187 call neko_error("Scalar lists must have same size!")
188 end if
189
190 do i = 1, size(s_index_list)
191 ! Take care that if we set i=0 we want temperature
192 if (s_index_list(i) .eq. 0) then
193 call fld_data%t%copy_from(host_to_device, .true.)
194 else
195 ! For scalar fields, require indices in 1:this%n_scalars
196 if (s_index_list(i) < 1 .or. &
197 s_index_list(i) > fld_data%n_scalars) then
198 call neko_error("s_index_list entry out of bounds")
199 end if
200 call fld_data%s(s_index_list(i))%copy_from(host_to_device, &
201 .true.)
202 end if
203 end do
204 else
205 do i = 1, s_target_list%size()
206 call fld_data%s(i)%copy_from(host_to_device, .true.)
207 end do
208 end if ! if present s_index_list
209
210 end if ! present s_tgt
211
212
213 if (interpolate_) then
214 ! Sync coordinates to device for the interpolation
215 call fld_data%x%copy_from(host_to_device, .false.)
216 call fld_data%y%copy_from(host_to_device, .false.)
217 call fld_data%z%copy_from(host_to_device, .true.)
218 end if
219
220 ! Call the import of fields
221 call fld_data%import_fields(u, v, w, p, t, s_target_list, s_index_list, &
222 interpolate_, tolerance = tolerance)
223
224 call neko_log%end_section()
225 call fld_data%free()
226
227 end subroutine import_fields
228
229end module import_field_utils
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...
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.
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 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:55