Neko 1.99.2
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
49 implicit none
50 private
51
54 type, public, extends(simulation_component_t) :: weak_gradient_t
56 type(field_t), pointer :: u
58 type(field_t), pointer :: gradient_x
60 type(field_t), pointer :: gradient_y
62 type(field_t), pointer :: gradient_z
64 type(field_writer_t) :: writer
65
66 contains
68 procedure, pass(this) :: init => weak_gradient_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 => weak_gradient_init_common
82 procedure, pass(this) :: free => weak_gradient_free
84 procedure, pass(this) :: compute_ => weak_gradient_compute
85 end type weak_gradient_t
86
87contains
88
90 subroutine weak_gradient_init_from_json(this, json, case)
91 class(weak_gradient_t), intent(inout), target :: this
92 type(json_file), intent(inout) :: json
93 class(case_t), intent(inout), target :: case
94 character(len=:), allocatable :: field_name
95 character(len=:), allocatable :: name
96 character(len=20) :: fields(3)
97 character(len=:), allocatable :: computed_field
98
99 ! Add fields keyword to the json so that the field_writer picks it up.
100 ! Will also add fields to the registry.
101 call json_get_or_default(json, "name", name, "weak_gradient")
102 call json_get(json, "field", field_name)
103 call json_get_or_default(json, "computed_field", computed_field, &
104 "weak_gradient" // trim(field_name))
105
106 fields(1) = computed_field // "_x"
107 fields(2) = computed_field // "_y"
108 fields(3) = computed_field // "_z"
109
110 call json%add("fields", fields)
111
112 call this%init_base(json, case)
113 call this%writer%init(json, case)
114
115 call this%init_common(name, field_name, computed_field)
116 end subroutine weak_gradient_init_from_json
117
120 subroutine weak_gradient_init_common(this, name, field_name, computed_field)
121 class(weak_gradient_t), intent(inout) :: this
122 character(len=*) :: name
123 character(len=*) :: field_name
124 character(len=*) :: computed_field
125
126 this%name = name
127 this%u => neko_registry%get_field_by_name(trim(field_name))
128
129 this%gradient_x => neko_registry%get_field_by_name( &
130 computed_field // "_x")
131 this%gradient_y => neko_registry%get_field_by_name( &
132 computed_field // "_y")
133 this%gradient_z => neko_registry%get_field_by_name( &
134 computed_field // "_z")
135
136
137 end subroutine weak_gradient_init_common
138
151 subroutine weak_gradient_init_from_controllers(this, name, case, order, &
152 preprocess_controller, compute_controller, output_controller, &
153 field_name, computed_field, filename, precision)
154 class(weak_gradient_t), intent(inout) :: this
155 character(len=*), intent(in) :: name
156 class(case_t), intent(inout), target :: case
157 integer :: order
158 type(time_based_controller_t), intent(in) :: preprocess_controller
159 type(time_based_controller_t), intent(in) :: compute_controller
160 type(time_based_controller_t), intent(in) :: output_controller
161 character(len=*) :: field_name
162 character(len=*) :: computed_field
163 character(len=*), intent(in), optional :: filename
164 integer, intent(in), optional :: precision
165
166 character(len=20) :: fields(3)
167
168 fields(1) = trim(computed_field) // "_x"
169 fields(2) = trim(computed_field) // "_y"
170 fields(3) = trim(computed_field) // "_z"
171
172 call this%init_base_from_components(case, order, preprocess_controller, &
173 compute_controller, output_controller)
174 call this%writer%init_from_components("field_writer", case, order, &
175 preprocess_controller, compute_controller, output_controller, fields, &
176 filename, precision)
177 call this%init_common(name, field_name, computed_field)
178
180
199 case, order, preprocess_control, preprocess_value, compute_control, &
200 compute_value, output_control, output_value, field_name, &
201 computed_field, filename, precision)
202 class(weak_gradient_t), intent(inout) :: this
203 character(len=*), intent(in) :: name
204 class(case_t), intent(inout), target :: case
205 integer :: order
206 character(len=*), intent(in) :: preprocess_control
207 real(kind=rp), intent(in) :: preprocess_value
208 character(len=*), intent(in) :: compute_control
209 real(kind=rp), intent(in) :: compute_value
210 character(len=*), intent(in) :: output_control
211 real(kind=rp), intent(in) :: output_value
212 character(len=*) :: field_name
213 character(len=*) :: computed_field
214 character(len=*), intent(in), optional :: filename
215 integer, intent(in), optional :: precision
216
217 character(len=20) :: fields(3)
218
219 fields(1) = trim(computed_field) // "_x"
220 fields(2) = trim(computed_field) // "_y"
221 fields(3) = trim(computed_field) // "_z"
222
223 call this%init_base_from_components(case, order, preprocess_control, &
224 preprocess_value, compute_control, compute_value, output_control, &
225 output_value)
226 call this%writer%init_from_components("field_writer", case, order, &
227 preprocess_control, preprocess_value, compute_control, compute_value, &
228 output_control, output_value, fields, filename, precision)
229 call this%init_common(name, field_name, computed_field)
230
232
234 subroutine weak_gradient_free(this)
235 class(weak_gradient_t), intent(inout) :: this
236 call this%free_base()
237 call this%writer%free()
238 nullify(this%gradient_x)
239 nullify(this%gradient_y)
240 nullify(this%gradient_z)
241 nullify(this%u)
242 end subroutine weak_gradient_free
243
246 subroutine weak_gradient_compute(this, time)
247 class(weak_gradient_t), intent(inout) :: this
248 type(time_state_t), intent(in) :: time
249
250 call opgrad(this%gradient_x%x, this%gradient_y%x, this%gradient_z%x, &
251 this%u%x, this%case%fluid%c_Xh)
252 end subroutine weak_gradient_compute
253
254end 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: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.
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.