Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
base_functional.f90
Go to the documentation of this file.
1
34
37 use design, only: design_t
38 use json_module, only: json_file
39 use json_utils, only: json_get
40 use num_types, only: rp
41 use point_zone, only: point_zone_t
42 use simulation_m, only: simulation_t
43 use vector, only: vector_t
44 use utils, only: neko_error
45 use vector_math, only: vector_copy, vector_add2s1
46 implicit none
47 private
48
59 type, abstract, public :: base_functional_t
60
62 real(kind=rp) :: value = 0.0_rp
64 real(kind=rp) :: value_old = 0.0_rp
66 type(vector_t) :: sensitivity
68 type(vector_t) :: sensitivity_old
70 character(len=25) :: name = ""
72 logical :: has_mask = .false.
74 class(point_zone_t), pointer :: mask => null()
75
76 contains
77
78 ! ----------------------------------------------------------------------- !
79 ! Derived class interfaces
80
82 generic :: init => init_json, init_json_sim
83
85 procedure, pass(this) :: init_json => functional_init_json
87 procedure, pass(this) :: init_json_sim => functional_init_json_sim
89 procedure(functional_free), pass(this), deferred :: free
90
92 procedure(functional_update_value), pass(this), deferred :: update_value
94 procedure(functional_update_sensitivity), pass(this), deferred :: &
95 update_sensitivity
96
98 procedure, pass(this) :: get_value => functional_get_value
100 procedure, pass(this) :: get_sensitivity => functional_get_sensitivity
102 procedure, pass(this) :: get_log_size => functional_get_log_size
104 procedure, pass(this) :: get_log_headers => functional_get_log_headers
106 procedure, pass(this) :: get_log_values => functional_get_log_values
108 procedure, pass(this) :: reset_value => functional_reset_value
110 procedure, pass(this) :: reset_sensitivity => functional_reset_sensitivity
112 procedure, pass(this) :: accumulate_value => functional_accumulate_value
114 procedure, pass(this) :: accumulate_sensitivity => &
115 functional_accumulate_sensitivity
116
117 end type base_functional_t
118
119 ! -------------------------------------------------------------------------- !
120 ! Interface specifications for the derived types, these are the constructors
121 ! for the different types of objective functions.
122
123 abstract interface
124
126 subroutine functional_free(this)
127 import base_functional_t
128 class(base_functional_t), intent(inout) :: this
129 end subroutine functional_free
130
132 subroutine functional_update_value(this, design)
134 class(base_functional_t), intent(inout) :: this
135 class(design_t), intent(in) :: design
136 end subroutine functional_update_value
137
139 subroutine functional_update_sensitivity(this, design)
141 class(base_functional_t), intent(inout) :: this
142 class(design_t), intent(in) :: design
143 end subroutine functional_update_sensitivity
144
145 end interface
146
147contains
148
150 subroutine functional_init_json(this, json, design)
151 class(base_functional_t), intent(inout) :: this
152 type(json_file), intent(inout) :: json
153 class(design_t), intent(in) :: design
154 character(len=:), allocatable :: type
155
156 call json_get(json, 'type', type)
157 call neko_error("Functional type: '" // type // &
158 "' does not support initialization without simulation")
159 end subroutine functional_init_json
160
162 subroutine functional_init_json_sim(this, json, design, simulation)
163 class(base_functional_t), intent(inout) :: this
164 type(json_file), intent(inout) :: json
165 class(design_t), intent(in) :: design
166 type(simulation_t), target, intent(inout) :: simulation
167 character(len=:), allocatable :: type
168
169 call json_get(json, 'type', type)
170 call neko_error("Functional type: '" // type // &
171 "' does not support initialization with simulation")
172 end subroutine functional_init_json_sim
173
174
176 function functional_get_value(this) result(v)
177 class(base_functional_t), intent(in) :: this
178 real(kind=rp) :: v
179
180 v = this%value
181 end function functional_get_value
182
184 subroutine functional_get_sensitivity(this, sensitivity)
185 class(base_functional_t), intent(in) :: this
186 type(vector_t), intent(inout) :: sensitivity
187
188 sensitivity = this%sensitivity
189 end subroutine functional_get_sensitivity
190
194 function functional_get_log_size(this) result(n)
195 class(base_functional_t), intent(in) :: this
196 integer :: n
197
198 n = 1
199 end function functional_get_log_size
200
204 subroutine functional_get_log_headers(this, headers)
205 class(base_functional_t), intent(in) :: this
206 character(len=*), intent(out) :: headers(:)
207
208 if (size(headers) .eq. 0) return
209 headers(1) = trim(this%name)
210 end subroutine functional_get_log_headers
211
215 subroutine functional_get_log_values(this, values)
216 class(base_functional_t), intent(in) :: this
217 real(kind=rp), intent(out) :: values(:)
218
219 if (size(values) .eq. 0) return
220 values(1) = this%value
221 end subroutine functional_get_log_values
222
224 subroutine functional_reset_value(this)
225 class(base_functional_t), intent(inout) :: this
226
227 this%value = 0.0_rp
228 end subroutine functional_reset_value
229
231 subroutine functional_reset_sensitivity(this)
232 class(base_functional_t), intent(inout) :: this
233
234 this%sensitivity = 0.0_rp
235 end subroutine functional_reset_sensitivity
236
238 subroutine functional_accumulate_value(this, design, dt)
239 class(base_functional_t), intent(inout) :: this
240 class(design_t), intent(in) :: design
241 real(kind=rp), intent(in) :: dt
242
243 this%value_old = this%value
244 call this%update_value(design)
245
246 ! could potentially use higher order trapezoidal/Simpson etc, but this
247 ! should suffice
248 this%value = this%value_old + this%value * dt
249 end subroutine functional_accumulate_value
250
252 subroutine functional_accumulate_sensitivity(this, design, dt)
253 class(base_functional_t), intent(inout) :: this
254 class(design_t), intent(in) :: design
255 real(kind=rp), intent(in) :: dt
256
257 call vector_copy(this%sensitivity_old, this%sensitivity)
258 call this%update_sensitivity(design)
259
260 ! could potentially use higher order trapezoidal/Simpson etc, but this
261 ! should suffice
262 call vector_add2s1(this%sensitivity, this%sensitivity_old, dt)
263 end subroutine functional_accumulate_sensitivity
264end module base_functional
Defines the abstract the base_functional_t type.
Implements the design_t.
Definition design.f90:36
Implements the steady_problem_t type.
An abstract design type.
Definition design.f90:54