Neko 1.99.1
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
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
101
102 call json_get_or_default(json, "computed_field", computed_field, "divergence")
103 call json_get(json, "fields", field_names)
104
105 if (size(field_names) .ne. 3) then
106 call neko_error("The divergence simcomp requires exactly 3 entries in " // &
107 "fieldes.")
108 end if
109
110 fields(1) = trim(computed_field)
111
112 ! This is needed for the field writer to pick up the fields.
113 call json%add("fields", fields)
114
115 call this%init_base(json, case)
116 call this%writer%init(json, case)
117
118 call divergence_init_common(this, field_names, computed_field)
119 end subroutine divergence_init_from_json
120
124 subroutine divergence_init_common(this, field_names, computed_field)
125 class(divergence_t), intent(inout) :: this
126 character(len=*) :: field_names(3)
127 character(len=*) :: computed_field
128
129 this%u => neko_field_registry%get_field_by_name(field_names(1))
130 this%v => neko_field_registry%get_field_by_name(field_names(2))
131 this%w => neko_field_registry%get_field_by_name(field_names(3))
132
133 this%divergence => neko_field_registry%get_field_by_name(computed_field)
134
135 end subroutine divergence_init_common
136
148 subroutine divergence_init_from_controllers(this, case, order, &
149 preprocess_controller, compute_controller, output_controller, &
150 field_names, computed_field, filename, precision)
151 class(divergence_t), intent(inout) :: this
152 class(case_t), intent(inout), target :: case
153 integer :: order
154 type(time_based_controller_t), intent(in) :: preprocess_controller
155 type(time_based_controller_t), intent(in) :: compute_controller
156 type(time_based_controller_t), intent(in) :: output_controller
157 character(len=*) :: field_names(3)
158 character(len=*) :: computed_field
159 character(len=*), intent(in), optional :: filename
160 integer, intent(in), optional :: precision
161
162 character(len=20) :: fields(1)
163
164 fields(1) = trim(computed_field)
165
166 call this%init_base_from_components(case, order, preprocess_controller, &
167 compute_controller, output_controller)
168 call this%writer%init_from_components(case, order, preprocess_controller, &
169 compute_controller, output_controller, fields, filename, precision)
170 call this%init_common(field_names, computed_field)
171
173
191 case, order, preprocess_control, preprocess_value, compute_control, &
192 compute_value, output_control, output_value, field_names, &
193 computed_field, filename, precision)
194 class(divergence_t), intent(inout) :: this
195 class(case_t), intent(inout), target :: case
196 integer :: order
197 character(len=*), intent(in) :: preprocess_control
198 real(kind=rp), intent(in) :: preprocess_value
199 character(len=*), intent(in) :: compute_control
200 real(kind=rp), intent(in) :: compute_value
201 character(len=*), intent(in) :: output_control
202 real(kind=rp), intent(in) :: output_value
203 character(len=*) :: field_names(3)
204 character(len=*) :: computed_field
205 character(len=*), intent(in), optional :: filename
206 integer, intent(in), optional :: precision
207
208 character(len=20) :: fields(1)
209
210 fields(1) = trim(computed_field) // "_x"
211
212 call this%init_base_from_components(case, order, preprocess_control, &
213 preprocess_value, compute_control, compute_value, output_control, &
214 output_value)
215 call this%writer%init_from_components(case, order, preprocess_control, &
216 preprocess_value, compute_control, compute_value, output_control, &
217 output_value, fields, filename, precision)
218 call this%init_common(field_names, computed_field)
219
221
223 subroutine divergence_free(this)
224 class(divergence_t), intent(inout) :: this
225 call this%free_base()
226 call this%writer%free()
227
228 nullify(this%u)
229 nullify(this%v)
230 nullify(this%w)
231 nullify(this%divergence)
232 end subroutine divergence_free
233
235 subroutine divergence_compute(this, time)
236 class(divergence_t), intent(inout) :: this
237 type(time_state_t), intent(in) :: time
238
239 call div(this%divergence%x, this%u%x, this%v%x, this%w%x, this%case%fluid%c_Xh)
240
241 end subroutine divergence_compute
242
243end 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_json(this, json, case)
Constructor from json.
subroutine divergence_init_from_controllers_properties(this, 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_init_from_controllers(this, case, order, preprocess_controller, compute_controller, output_controller, field_names, computed_field, filename, precision)
Constructor from components, passing controllers.
subroutine divergence_free(this)
Destructor.
subroutine divergence_init_common(this, field_names, computed_field)
Actual constructor.
subroutine divergence_compute(this, time)
Compute the divergence field.
Defines a registry for storing solution fields.
type(field_registry_t), target, public neko_field_registry
Global field registry.
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
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.