Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
field_writer.f90
Go to the documentation of this file.
1! Copyright (c) 2024-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!
33!
35
37 use num_types, only : rp, dp, sp
38 use json_module, only : json_file
40 use time_state, only : time_state_t
41 use registry, only : neko_registry
42 use case, only : case_t
44 use fld_file, only : fld_file_t
45 use point_zone, only : point_zone_t
50 use logger, only : neko_log
51 implicit none
52 private
53
55 type, public, extends(simulation_component_t) :: field_writer_t
57 type(field_output_t), private :: output
58
59 ! Default values for optional parameters in the constructor.
60 character(len=20), private :: default_name = "field_writer"
61 character(len=20), private :: default_precision = "single"
62 character(len=20), private :: default_filename = ""
63 character(len=20), private :: default_format = "nek5000"
64 logical, private :: default_subdivide = .false.
65
66 contains
68 procedure, pass(this) :: init => field_writer_init_from_json
70 generic :: init_from_components => &
71 init_from_controllers, init_from_controllers_properties
73 procedure, pass(this) :: init_from_controllers => &
77 procedure, pass(this) :: init_from_controllers_properties => &
80 procedure, private, pass(this) :: init_common => field_writer_init_common
82 procedure, pass(this) :: free => field_writer_free
84 procedure, pass(this) :: compute_ => field_writer_compute
85 end type field_writer_t
86
87contains
88
92 subroutine field_writer_init_from_json(this, json, case)
93 class(field_writer_t), intent(inout), target :: this
94 type(json_file), intent(inout) :: json
95 class(case_t), intent(inout), target :: case
96 character(len=:), allocatable :: filename
97 character(len=:), allocatable :: name
98 character(len=:), allocatable :: precision
99 character(len=:), allocatable :: format
100 character(len=:), allocatable :: point_zone_name
101 character(len=NEKO_VARNAME_LEN), allocatable :: fields(:)
102 integer :: precision_value
103 logical :: subdivide
104 class(point_zone_t), pointer :: point_zone
105
106 call this%init_base(json, case)
107 call json_get(json, "fields", fields)
108 call json_get_or_default(json, "name", name, this%default_name)
109 call json_get_or_default(json, "output_filename", filename, &
110 this%default_filename)
111 call json_get_or_default(json, "output_precision", precision, &
112 this%default_precision)
113 call json_get_or_default(json, "output_format", format, this%default_format)
114
115 if (precision .eq. "single") then
116 precision_value = sp
117 else if (precision .eq. "double") then
118 precision_value = dp
119 else
120 call neko_error("Invalid precision specified for field_writer: " &
121 // trim(precision))
122 end if
123
124 call json_get_or_default(json, "output_subdivide", subdivide, &
125 this%default_subdivide)
126
127 if (json%valid_path('point_zone')) then
128 call json_get(json, 'point_zone', point_zone_name)
129 point_zone => neko_point_zone_registry%get_point_zone(point_zone_name)
130 call this%init_common(name, fields, filename, precision_value, format, &
131 subdivide, point_zone)
132 else
133 call this%init_common(name, fields, filename, precision_value, format, &
134 subdivide)
135 end if
136
137 end subroutine field_writer_init_from_json
138
155 subroutine field_writer_init_from_controllers(this, name, case, order, &
156 preprocess_controller, compute_controller, output_controller, &
157 fields, filename, precision, format, subdivide, point_zone)
158 class(field_writer_t), intent(inout) :: this
159 character(len=*), intent(in) :: name
160 class(case_t), intent(inout), target :: case
161 integer :: order
162 type(time_based_controller_t), intent(in) :: preprocess_controller
163 type(time_based_controller_t), intent(in) :: compute_controller
164 type(time_based_controller_t), intent(in) :: output_controller
165 character(len=*), intent(in) :: fields(:)
166 character(len=*), intent(in), optional :: filename
167 integer, intent(in), optional :: precision
168 character(len=20), intent(in), optional :: format
169 logical, intent(in), optional :: subdivide
170 class(point_zone_t), intent(inout), optional :: point_zone
171
172 call this%init_base_from_components(case, order, preprocess_controller, &
173 compute_controller, output_controller)
174 call this%init_common(name, fields, filename, precision, format, &
175 subdivide, point_zone)
176
178
200 case, order, preprocess_control, preprocess_value, compute_control, &
201 compute_value, output_control, output_value, fields, filename, &
202 precision, format, subdivide, point_zone)
203 class(field_writer_t), intent(inout) :: this
204 character(len=*), intent(in) :: name
205 class(case_t), intent(inout), target :: case
206 integer :: order
207 character(len=*), intent(in) :: preprocess_control
208 real(kind=rp), intent(in) :: preprocess_value
209 character(len=*), intent(in) :: compute_control
210 real(kind=rp), intent(in) :: compute_value
211 character(len=*), intent(in) :: output_control
212 real(kind=rp), intent(in) :: output_value
213 character(len=*), intent(in) :: fields(:)
214 character(len=*), intent(in), optional :: filename
215 integer, intent(in), optional :: precision
216 character(len=*), intent(in), optional :: format
217 logical, intent(in), optional :: subdivide
218 class(point_zone_t), intent(inout), optional :: point_zone
219
220 call this%init_base_from_components(case, order, preprocess_control, &
221 preprocess_value, compute_control, compute_value, output_control, &
222 output_value)
223 call this%init_common(name, fields, filename, precision, format, &
224 subdivide, point_zone)
225
227
239 subroutine field_writer_init_common(this, name, fields, filename, precision, &
240 format, subdivide, point_zone)
241 class(field_writer_t), intent(inout) :: this
242 character(len=*), intent(in) :: name
243 character(len=*), intent(in) :: fields(:)
244 character(len=*), intent(in), optional :: filename
245 integer, intent(in), optional :: precision
246 character(len=*), intent(in), optional :: format
247 logical, intent(in), optional :: subdivide
248 class(point_zone_t), intent(inout), optional :: point_zone
249
250 logical :: filename_provided
251 integer :: i
252
253 this%name = name
254 ! Register fields if they don't exist.
255 do i = 1, size(fields)
256 call neko_registry%add_field(this%case%fluid%dm_Xh, trim(fields(i)), &
257 ignore_existing = .true.)
258 end do
259
260 ! Throw an error if no output filename is provided to avoid dumping masked
261 ! fields in the default fluid output.
262 if (present(point_zone) .and. &
263 (.not. present(filename) .or. len_trim(filename) .eq. 0 ) ) then
264 call neko_error("Please provide a filename for the " // &
265 "field_writer when using a point_zone with the key " // &
266 "'output_filename'.")
267 end if
268
269 filename_provided = .false.
270 if (present(filename)) then
271 if (len_trim(filename) .ne. 0) then
272 filename_provided = .true.
273 call this%output%init(trim(filename), size(fields), &
274 precision = precision, format = format, &
275 path = this%case%output_directory)
276
277 if (present(subdivide)) then
278 call this%output%file_%set_subdivide(subdivide)
279 end if
280
281 do i = 1, size(fields)
282 call this%output%fields%assign(i, &
283 neko_registry%get_field(trim(fields(i))))
284 end do
285
286 call this%case%output_controller%add(this%output, &
287 this%output_controller%control_value, &
288 this%output_controller%control_mode)
289
290 end if
291 end if
292
293 if (.not. filename_provided) then
294 do i = 1, size(fields)
295 call this%case%f_out%fluid%append( &
296 neko_registry%get_field(trim(fields(i))))
297 end do
298 end if
299
300 !
301 ! Set the mask for subsampling is a point zone is provided. Only supported
302 ! for fld/nek5000 files
303 !
304 if (present(point_zone)) then
305 select type (ft => this%output%file_%file_type)
306 class is (fld_file_t)
307 call ft%set_mask(point_zone%mask)
308 class default
309 call neko_error("point_zone can only be used with nek5000/fld files")
310 end select
311 end if
312
313 end subroutine field_writer_init_common
314
316 subroutine field_writer_free(this)
317 class(field_writer_t), intent(inout) :: this
318 call this%free_base()
319 call this%output%free()
320 end subroutine field_writer_free
321
323 subroutine field_writer_compute(this, time)
324 class(field_writer_t), intent(inout) :: this
325 type(time_state_t), intent(in) :: time
326
327 end subroutine field_writer_compute
328
329end module field_writer
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.
Defines a simulation case.
Definition case.f90:34
Implements field_output_t.
Implements the field_writer_t type.
subroutine field_writer_compute(this, time)
Here to comply with the interface, does nothing.
subroutine field_writer_free(this)
Destructor.
subroutine field_writer_init_from_json(this, json, case)
Constructor from json.
subroutine field_writer_init_from_controllers_properties(this, name, case, order, preprocess_control, preprocess_value, compute_control, compute_value, output_control, output_value, fields, filename, precision, format, subdivide, point_zone)
Constructor from components, passing properties to the time_based_controllercomponents in the base ty...
subroutine field_writer_init_common(this, name, fields, filename, precision, format, subdivide, point_zone)
Common part of both constructors.
subroutine field_writer_init_from_controllers(this, name, case, order, preprocess_controller, compute_controller, output_controller, fields, filename, precision, format, subdivide, point_zone)
Constructor from components, passing controllers.
NEKTON fld file format.
Definition fld_file.f90:35
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 dp
Definition num_types.f90:9
integer, parameter, public sp
Definition num_types.f90:8
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Implements output_controller_t
Defines an output.
Definition output.f90:34
type(point_zone_registry_t), target, public neko_point_zone_registry
Global point_zone registry.
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
Simulation components are objects that encapsulate functionality that can be fit to a particular comp...
subroutine compute_(this, time)
Dummy compute function.
Contains the time_based_controller_t type.
Module with things related to the simulation time.
Utilities.
Definition utils.f90:35
integer, parameter, public neko_varname_len
Definition utils.f90:42
A simple output saving a list of fields to a file.
A simulation component that writes a 3d field to a file.
Interface for NEKTON fld files.
Definition fld_file.f90:66
Base abstract type for point zones.
Base abstract class for simulation components.
A utility type for determining whether an action should be executed based on the current time value....
A struct that contains all info about the time, expand as needed.