Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
weak_gradient_simcomp.f90
Go to the documentation of this file.
1! Copyright (c) 2024, 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 operators, only : opgrad
43 use time_state, only : time_state_t
44 use case, only : case_t
48 use utils, only : neko_varname_len
50 implicit none
51 private
52
55 type, public, extends(simulation_component_t) :: weak_gradient_t
57 type(field_t), pointer :: u
59 type(field_t), pointer :: gradient_x
61 type(field_t), pointer :: gradient_y
63 type(field_t), pointer :: gradient_z
65 type(field_writer_t) :: writer
66
67 contains
69 procedure, pass(this) :: init => weak_gradient_init_from_json
71 generic :: init_from_components => &
72 init_from_controllers, init_from_controllers_properties
74 procedure, pass(this) :: init_from_controllers => &
78 procedure, pass(this) :: init_from_controllers_properties => &
81 procedure, private, pass(this) :: init_common => weak_gradient_init_common
83 procedure, pass(this) :: free => weak_gradient_free
85 procedure, pass(this) :: compute_ => weak_gradient_compute
86 end type weak_gradient_t
87
88contains
89
91 subroutine weak_gradient_init_from_json(this, json, case)
92 class(weak_gradient_t), intent(inout), target :: this
93 type(json_file), intent(inout) :: json
94 class(case_t), intent(inout), target :: case
95 character(len=:), allocatable :: field_name
96 character(len=:), allocatable :: name
97 character(len=NEKO_VARNAME_LEN) :: fields(3)
98 character(len=:), allocatable :: computed_field
99
100 ! Add fields keyword to the json so that the field_writer picks it up.
101 ! Will also add fields to the registry.
102 call json_get_or_default(json, "name", name, "weak_gradient")
103 call json_get(json, "field", field_name)
104 call json_get_or_default(json, "computed_field", computed_field, &
105 "weak_gradient" // trim(field_name))
106
107 fields(1) = computed_field // "_x"
108 fields(2) = computed_field // "_y"
109 fields(3) = computed_field // "_z"
110
111 call json%add("fields", fields)
112
113 call this%init_base(json, case)
114 call this%writer%init(json, case)
115
116 call this%init_common(name, field_name, computed_field)
117 end subroutine weak_gradient_init_from_json
118
121 subroutine weak_gradient_init_common(this, name, field_name, computed_field)
122 class(weak_gradient_t), intent(inout) :: this
123 character(len=*) :: name
124 character(len=*) :: field_name
125 character(len=*) :: computed_field
126
127 this%name = name
128 this%u => neko_registry%get_field_by_name(trim(field_name))
129
130 this%gradient_x => neko_registry%get_field_by_name( &
131 computed_field // "_x")
132 this%gradient_y => neko_registry%get_field_by_name( &
133 computed_field // "_y")
134 this%gradient_z => neko_registry%get_field_by_name( &
135 computed_field // "_z")
136
137
138 end subroutine weak_gradient_init_common
139
152 subroutine weak_gradient_init_from_controllers(this, name, case, order, &
153 preprocess_controller, compute_controller, output_controller, &
154 field_name, computed_field, filename, precision)
155 class(weak_gradient_t), intent(inout) :: this
156 character(len=*), intent(in) :: name
157 class(case_t), intent(inout), target :: case
158 integer :: order
159 type(time_based_controller_t), intent(in) :: preprocess_controller
160 type(time_based_controller_t), intent(in) :: compute_controller
161 type(time_based_controller_t), intent(in) :: output_controller
162 character(len=*) :: field_name
163 character(len=*) :: computed_field
164 character(len=*), intent(in), optional :: filename
165 integer, intent(in), optional :: precision
166
167 character(len=NEKO_VARNAME_LEN) :: fields(3)
168
169 fields(1) = trim(computed_field) // "_x"
170 fields(2) = trim(computed_field) // "_y"
171 fields(3) = trim(computed_field) // "_z"
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_name, computed_field)
179
181
200 case, order, preprocess_control, preprocess_value, compute_control, &
201 compute_value, output_control, output_value, field_name, &
202 computed_field, filename, precision)
203 class(weak_gradient_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_name
214 character(len=*) :: computed_field
215 character(len=*), intent(in), optional :: filename
216 integer, intent(in), optional :: precision
217
218 character(len=NEKO_VARNAME_LEN) :: fields(3)
219
220 fields(1) = trim(computed_field) // "_x"
221 fields(2) = trim(computed_field) // "_y"
222 fields(3) = trim(computed_field) // "_z"
223
224 call this%init_base_from_components(case, order, preprocess_control, &
225 preprocess_value, compute_control, compute_value, output_control, &
226 output_value)
227 call this%writer%init_from_components("field_writer", case, order, &
228 preprocess_control, preprocess_value, compute_control, compute_value, &
229 output_control, output_value, fields, filename, precision)
230 call this%init_common(name, field_name, computed_field)
231
233
235 subroutine weak_gradient_free(this)
236 class(weak_gradient_t), intent(inout) :: this
237 call this%free_base()
238 call this%writer%free()
239 nullify(this%gradient_x)
240 nullify(this%gradient_y)
241 nullify(this%gradient_z)
242 nullify(this%u)
243 end subroutine weak_gradient_free
244
247 subroutine weak_gradient_compute(this, time)
248 class(weak_gradient_t), intent(inout) :: this
249 type(time_state_t), intent(in) :: time
250
251 call opgrad(this%gradient_x%x, this%gradient_y%x, this%gradient_z%x, &
252 this%u%x, this%case%fluid%c_Xh)
253 end subroutine weak_gradient_compute
254
255end module weak_gradient_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 field_writer_t type.
Defines a field.
Definition field.f90:34
Implements fld_file_output_t.
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 opgrad(ux, uy, uz, u, coef, es, ee)
Compute the weak gradient of a scalar field, i.e. the gradient multiplied by the mass matrix.
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: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
Implements the weak_gradient_t type.
subroutine weak_gradient_init_from_controllers(this, name, case, order, preprocess_controller, compute_controller, output_controller, field_name, computed_field, filename, precision)
Constructor from components, passing controllers.
subroutine weak_gradient_free(this)
Destructor.
subroutine weak_gradient_init_from_controllers_properties(this, name, case, order, preprocess_control, preprocess_value, compute_control, compute_value, output_control, output_value, field_name, computed_field, filename, precision)
Constructor from components, passing properties to the time_based_controller` components in the base ...
subroutine weak_gradient_init_common(this, name, field_name, computed_field)
Common part of the constructors.
subroutine weak_gradient_init_from_json(this, json, case)
Constructor from json.
subroutine weak_gradient_compute(this, time)
Compute the weak_gradient field.
A simulation component that writes a 3d field to a file.
A simple output saving a list of fields to a .fld 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.
A simulation component that computes the weak gradient of a field. Wraps the opgradient operator.