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 time_state, only: time_state_t
44 use vector, only: vector_t
45 use utils, only: neko_error
46 use vector_math, only: vector_copy, vector_add2s1
47 implicit none
48 private
49
60 type, abstract, public :: base_functional_t
61
63 real(kind=rp) :: value = 0.0_rp
65 real(kind=rp) :: value_old = 0.0_rp
67 type(vector_t) :: sensitivity
69 type(vector_t) :: sensitivity_old
71 character(len=25) :: name = ""
73 logical :: has_mask = .false.
75 class(point_zone_t), pointer :: mask => null()
77 real(kind=rp) :: start_time = 0.0_rp
79 real(kind=rp) :: end_time = huge(0.0_rp)
80
81 contains
82
83 ! ----------------------------------------------------------------------- !
84 ! Derived class interfaces
85
87 generic :: init => init_json, init_json_sim
88
90 procedure, pass(this) :: init_json => functional_init_json
92 procedure, pass(this) :: init_json_sim => functional_init_json_sim
94 procedure(functional_free), pass(this), deferred :: free
95
97 procedure(functional_update_value), pass(this), deferred :: update_value
99 procedure(functional_update_sensitivity), pass(this), deferred :: &
100 update_sensitivity
101
103 procedure, pass(this) :: get_value => functional_get_value
105 procedure, pass(this) :: get_sensitivity => functional_get_sensitivity
107 procedure, pass(this) :: get_log_size => functional_get_log_size
109 procedure, pass(this) :: get_log_headers => functional_get_log_headers
111 procedure, pass(this) :: get_log_values => functional_get_log_values
113 procedure, pass(this) :: reset_value => functional_reset_value
115 procedure, pass(this) :: reset_sensitivity => functional_reset_sensitivity
117 procedure, pass(this) :: accumulate_value => functional_accumulate_value
119 procedure, pass(this) :: accumulate_sensitivity => &
120 functional_accumulate_sensitivity
121
122 end type base_functional_t
123
124 ! -------------------------------------------------------------------------- !
125 ! Interface specifications for the derived types, these are the constructors
126 ! for the different types of objective functions.
127
128 abstract interface
129
131 subroutine functional_free(this)
132 import base_functional_t
133 class(base_functional_t), intent(inout) :: this
134 end subroutine functional_free
135
137 subroutine functional_update_value(this, design)
139 class(base_functional_t), intent(inout) :: this
140 class(design_t), intent(in) :: design
141 end subroutine functional_update_value
142
144 subroutine functional_update_sensitivity(this, design)
146 class(base_functional_t), intent(inout) :: this
147 class(design_t), intent(in) :: design
148 end subroutine functional_update_sensitivity
149
150 end interface
151
152contains
153
155 subroutine functional_init_json(this, json, design)
156 class(base_functional_t), intent(inout) :: this
157 type(json_file), intent(inout) :: json
158 class(design_t), intent(in) :: design
159 character(len=:), allocatable :: type
160
161 call json_get(json, 'type', type)
162 call neko_error("Functional type: '" // type // &
163 "' does not support initialization without simulation")
164 end subroutine functional_init_json
165
167 subroutine functional_init_json_sim(this, json, design, simulation)
168 class(base_functional_t), intent(inout) :: this
169 type(json_file), intent(inout) :: json
170 class(design_t), intent(in) :: design
171 type(simulation_t), target, intent(inout) :: simulation
172 character(len=:), allocatable :: type
173
174 call json_get(json, 'type', type)
175 call neko_error("Functional type: '" // type // &
176 "' does not support initialization with simulation")
177 end subroutine functional_init_json_sim
178
180 function functional_get_value(this) result(v)
181 class(base_functional_t), intent(in) :: this
182 real(kind=rp) :: v
183
184 v = this%value
185 end function functional_get_value
186
188 subroutine functional_get_sensitivity(this, sensitivity)
189 class(base_functional_t), intent(in) :: this
190 type(vector_t), intent(inout) :: sensitivity
191
192 sensitivity = this%sensitivity
193 end subroutine functional_get_sensitivity
194
198 function functional_get_log_size(this) result(n)
199 class(base_functional_t), intent(in) :: this
200 integer :: n
201
202 n = 1
203 end function functional_get_log_size
204
208 subroutine functional_get_log_headers(this, headers)
209 class(base_functional_t), intent(in) :: this
210 character(len=*), intent(out) :: headers(:)
211
212 if (size(headers) .eq. 0) return
213 headers(1) = trim(this%name)
214 end subroutine functional_get_log_headers
215
219 subroutine functional_get_log_values(this, values)
220 class(base_functional_t), intent(in) :: this
221 real(kind=rp), intent(out) :: values(:)
222
223 if (size(values) .eq. 0) return
224 values(1) = this%value
225 end subroutine functional_get_log_values
226
228 subroutine functional_reset_value(this)
229 class(base_functional_t), intent(inout) :: this
230
231 this%value = 0.0_rp
232 end subroutine functional_reset_value
233
235 subroutine functional_reset_sensitivity(this)
236 class(base_functional_t), intent(inout) :: this
237
238 this%sensitivity = 0.0_rp
239 end subroutine functional_reset_sensitivity
240
242 subroutine functional_accumulate_value(this, design, time)
243 class(base_functional_t), intent(inout) :: this
244 class(design_t), intent(in) :: design
245 type(time_state_t), intent(in) :: time
246
247 if (time%t .lt. this%start_time .or. time%t .gt. this%end_time) return
248
249 this%value_old = this%value
250 call this%update_value(design)
251
252 ! could potentially use higher order trapezoidal/Simpson etc, but this
253 ! should suffice
254 this%value = this%value_old + this%value * time%dt
255 end subroutine functional_accumulate_value
256
258 subroutine functional_accumulate_sensitivity(this, design, time)
259 class(base_functional_t), intent(inout) :: this
260 class(design_t), intent(in) :: design
261 type(time_state_t), intent(in) :: time
262
263 if (time%t .lt. this%start_time .or. time%t .gt. this%end_time) return
264
265 call vector_copy(this%sensitivity_old, this%sensitivity)
266 call this%update_sensitivity(design)
267
268 ! could potentially use higher order trapezoidal/Simpson etc, but this
269 ! should suffice
270 call vector_add2s1(this%sensitivity, this%sensitivity_old, time%dt)
271 end subroutine functional_accumulate_sensitivity
272end 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:53