Neko 1.99.2
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
divergence_simcomp.f90
Go to the documentation of this file.
1! Copyright (c) 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!
33!
35
37 use num_types, only : rp, dp, sp
38 use json_module, only : json_file
40 use registry, only : neko_registry
41 use field, only : field_t
42 use time_state, only : time_state_t
43 use operators, only : div
44 use case, only : case_t
48 use utils, only : neko_error
49 implicit none
50 private
51
55 type, public, extends(simulation_component_t) :: divergence_t
57 type(field_t), pointer :: u
59 type(field_t), pointer :: v
61 type(field_t), pointer :: w
62
64 type(field_t), pointer :: divergence
65
67 type(field_writer_t) :: writer
68
69 contains
71 procedure, pass(this) :: init => divergence_init_from_json
73 generic :: init_from_components => &
74 init_from_controllers, init_from_controllers_properties
76 procedure, pass(this) :: init_from_controllers => &
80 procedure, pass(this) :: init_from_controllers_properties => &
83 procedure, private, pass(this) :: init_common => divergence_init_common
85 procedure, pass(this) :: free => divergence_free
87 procedure, pass(this) :: compute_ => divergence_compute
88 end type divergence_t
89
90contains
91
93 subroutine divergence_init_from_json(this, json, case)
94 class(divergence_t), intent(inout), target :: this
95 type(json_file), intent(inout) :: json
96 class(case_t), intent(inout), target :: case
97 character(len=20) :: fields(1)
98 character(len=20), allocatable :: field_names(:)
99 character(len=:), allocatable :: computed_field
100 character(len=:), allocatable :: name
101
102 call json_get_or_default(json, "name", name, "divergence")
103 call json_get_or_default(json, "computed_field", computed_field, &
104 "divergence")
105 call json_get(json, "fields", field_names)
106
107 if (size(field_names) .ne. 3) then
108 call neko_error("The divergence simcomp requires exactly 3 entries in " &
109 // "fieldes.")
110 end if
111
112 fields(1) = trim(computed_field)
113
114 ! This is needed for the field writer to pick up the fields.
115 call json%add("fields", fields)
116
117 call this%init_base(json, case)
118 call this%writer%init(json, case)
119
120 call divergence_init_common(this, name, field_names, computed_field)
121 end subroutine divergence_init_from_json
122
127 subroutine divergence_init_common(this, name, field_names, computed_field)
128 class(divergence_t), intent(inout) :: this
129 character(len=*), intent(in) :: name
130 character(len=*), intent(in) :: field_names(3)
131 character(len=*), intent(in) :: computed_field
132
133 this%name = name
134 this%u => neko_registry%get_field_by_name(field_names(1))
135 this%v => neko_registry%get_field_by_name(field_names(2))
136 this%w => neko_registry%get_field_by_name(field_names(3))
137
138 this%divergence => neko_registry%get_field_by_name(computed_field)
139
140 end subroutine divergence_init_common
141
154 subroutine divergence_init_from_controllers(this, name, case, order, &
155 preprocess_controller, compute_controller, output_controller, &
156 field_names, computed_field, filename, precision)
157 class(divergence_t), intent(inout) :: this
158 character(len=*), intent(in) :: name
159 class(case_t), intent(inout), target :: case
160 integer :: order
161 type(time_based_controller_t), intent(in) :: preprocess_controller
162 type(time_based_controller_t), intent(in) :: compute_controller
163 type(time_based_controller_t), intent(in) :: output_controller
164 character(len=*), intent(in) :: field_names(3)
165 character(len=*), intent(in) :: computed_field
166 character(len=*), intent(in), optional :: filename
167 integer, intent(in), optional :: precision
168
169 character(len=20) :: fields(1)
170
171 fields(1) = trim(computed_field)
172
173 call this%init_base_from_components(case, order, preprocess_controller, &
174 compute_controller, output_controller)
175 call this%writer%init_from_components("field_writer", case, order, &
176 preprocess_controller, compute_controller, output_controller, fields, &
177 filename, precision)
178 call this%init_common(name, field_names, computed_field)
179
181
200 case, order, preprocess_control, preprocess_value, compute_control, &
201 compute_value, output_control, output_value, field_names, &
202 computed_field, filename, precision)
203 class(divergence_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=*) :: field_names(3)
214 character(len=*) :: computed_field
215 character(len=*), intent(in), optional :: filename
216 integer, intent(in), optional :: precision
217
218 character(len=20) :: fields(1)
219
220 fields(1) = trim(computed_field) // "_x"
221
222 call this%init_base_from_components(case, order, preprocess_control, &
223 preprocess_value, compute_control, compute_value, output_control, &
224 output_value)
225 call this%writer%init_from_components("field_writer", case, order, &
226 preprocess_control, preprocess_value, compute_control, compute_value, &
227 output_control, output_value, fields, filename, precision)
228 call this%init_common(name, field_names, computed_field)
229
231
233 subroutine divergence_free(this)
234 class(divergence_t), intent(inout) :: this
235 call this%free_base()
236 call this%writer%free()
237
238 nullify(this%u)
239 nullify(this%v)
240 nullify(this%w)
241 nullify(this%divergence)
242 end subroutine divergence_free
243
245 subroutine divergence_compute(this, time)
246 class(divergence_t), intent(inout) :: this
247 type(time_state_t), intent(in) :: time
248
249 call div(this%divergence%x, this%u%x, this%v%x, this%w%x, &
250 this%case%fluid%c_Xh)
251
252 end subroutine divergence_compute
253
254end module divergence_simcomp
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 the divergence_t type.
subroutine divergence_init_from_controllers(this, name, case, order, preprocess_controller, compute_controller, output_controller, field_names, computed_field, filename, precision)
Constructor from components, passing controllers.
subroutine divergence_init_from_json(this, json, case)
Constructor from json.
subroutine divergence_init_common(this, name, field_names, computed_field)
Actual constructor.
subroutine divergence_init_from_controllers_properties(this, name, case, order, preprocess_control, preprocess_value, compute_control, compute_value, output_control, output_value, field_names, computed_field, filename, precision)
Constructor from components, passing properties to the time_based_controller` components in the base ...
subroutine divergence_free(this)
Destructor.
subroutine divergence_compute(this, time)
Compute the divergence field.
Implements the field_writer_t type.
Defines a field.
Definition field.f90:34
Utilities for retrieving parameters from the case files.
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
Operators.
Definition operators.f90:34
subroutine, public div(res, ux, uy, uz, coef)
Compute the divergence of a vector field.
Implements output_controller_t
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:149
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
A simulation component that computes the divergence of a vector field. Added to the field registry as...
A simulation component that writes a 3d field to a file.
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.