Neko 1.99.1
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
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=20) :: fields(3)
96 character(len=:), allocatable :: computed_field
97
98 ! Add fields keyword to the json so that the field_writer picks it up.
99 ! Will also add fields to the registry.
100 call json_get(json, "field", field_name)
101 call json_get_or_default(json, "computed_field", computed_field, &
102 "weak_gradient" // trim(field_name))
103
104 fields(1) = computed_field // "_x"
105 fields(2) = computed_field // "_y"
106 fields(3) = computed_field // "_z"
107
108 call json%add("fields", fields)
109
110 call this%init_base(json, case)
111 call this%writer%init(json, case)
112
113 call this%init_common(field_name, computed_field)
114 end subroutine weak_gradient_init_from_json
115
117 subroutine weak_gradient_init_common(this, field_name, computed_field)
118 class(weak_gradient_t), intent(inout) :: this
119 character(len=*) :: field_name
120 character(len=*) :: computed_field
121
122 this%u => neko_field_registry%get_field_by_name(trim(field_name))
123
124 this%gradient_x => neko_field_registry%get_field_by_name( &
125 computed_field // "_x")
126 this%gradient_y => neko_field_registry%get_field_by_name( &
127 computed_field // "_y")
128 this%gradient_z => neko_field_registry%get_field_by_name( &
129 computed_field // "_z")
130
131
132 end subroutine weak_gradient_init_common
133
145 subroutine weak_gradient_init_from_controllers(this, case, order, &
146 preprocess_controller, compute_controller, output_controller, &
147 field_name, computed_field, filename, precision)
148 class(weak_gradient_t), intent(inout) :: this
149 class(case_t), intent(inout), target :: case
150 integer :: order
151 type(time_based_controller_t), intent(in) :: preprocess_controller
152 type(time_based_controller_t), intent(in) :: compute_controller
153 type(time_based_controller_t), intent(in) :: output_controller
154 character(len=*) :: field_name
155 character(len=*) :: computed_field
156 character(len=*), intent(in), optional :: filename
157 integer, intent(in), optional :: precision
158
159 character(len=20) :: fields(3)
160
161 fields(1) = trim(computed_field) // "_x"
162 fields(2) = trim(computed_field) // "_y"
163 fields(3) = trim(computed_field) // "_z"
164
165 call this%init_base_from_components(case, order, preprocess_controller, &
166 compute_controller, output_controller)
167 call this%writer%init_from_components(case, order, preprocess_controller, &
168 compute_controller, output_controller, fields, filename, precision)
169 call this%init_common(field_name, computed_field)
170
172
190 case, order, preprocess_control, preprocess_value, compute_control, &
191 compute_value, output_control, output_value, field_name, &
192 computed_field, filename, precision)
193 class(weak_gradient_t), intent(inout) :: this
194 class(case_t), intent(inout), target :: case
195 integer :: order
196 character(len=*), intent(in) :: preprocess_control
197 real(kind=rp), intent(in) :: preprocess_value
198 character(len=*), intent(in) :: compute_control
199 real(kind=rp), intent(in) :: compute_value
200 character(len=*), intent(in) :: output_control
201 real(kind=rp), intent(in) :: output_value
202 character(len=*) :: field_name
203 character(len=*) :: computed_field
204 character(len=*), intent(in), optional :: filename
205 integer, intent(in), optional :: precision
206
207 character(len=20) :: fields(3)
208
209 fields(1) = trim(computed_field) // "_x"
210 fields(2) = trim(computed_field) // "_y"
211 fields(3) = trim(computed_field) // "_z"
212
213 call this%init_base_from_components(case, order, preprocess_control, &
214 preprocess_value, compute_control, compute_value, output_control, &
215 output_value)
216 call this%writer%init_from_components(case, order, preprocess_control, &
217 preprocess_value, compute_control, compute_value, output_control, &
218 output_value, fields, filename, precision)
219 call this%init_common(field_name, computed_field)
220
222
224 subroutine weak_gradient_free(this)
225 class(weak_gradient_t), intent(inout) :: this
226 call this%free_base()
227 call this%writer%free()
228 nullify(this%gradient_x)
229 nullify(this%gradient_y)
230 nullify(this%gradient_z)
231 nullify(this%u)
232 end subroutine weak_gradient_free
233
236 subroutine weak_gradient_compute(this, time)
237 class(weak_gradient_t), intent(inout) :: this
238 type(time_state_t), intent(in) :: time
239
240 call opgrad(this%gradient_x%x, this%gradient_y%x, this%gradient_z%x, this%u%x,&
241 this%case%fluid%c_Xh)
242 end subroutine weak_gradient_compute
243
244end 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
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
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
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_common(this, field_name, computed_field)
Common part of the constructors.
subroutine weak_gradient_free(this)
Destructor.
subroutine weak_gradient_init_from_controllers_properties(this, 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_from_json(this, json, case)
Constructor from json.
subroutine weak_gradient_compute(this, time)
Compute the weak_gradient field.
subroutine weak_gradient_init_from_controllers(this, case, order, preprocess_controller, compute_controller, output_controller, field_name, computed_field, filename, precision)
Constructor from components, passing controllers.
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.