Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
lube_term_objective.f90
Go to the documentation of this file.
1! Copyright (c) 2023, 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!
34!
35! I promise I'll write this document properly in the future...
36!
37! But the Borval Peterson (I think) paper had an objective function
38! that had 2 terms, dissipation and this term they claimed represented
39! out of plane stresses.
40! I never really understood that extra term, I also don't think it
41! applies to 3D cases, but everyone includes it anyway.
42!
43! It appears to me to be basically a heuristic penality that targets
44! non-binary designs
45!
46! so let's call
47!
48! F = \int |\nabla u|^2 + K \int \chi \u^2
49!
50! | dissipation | |"lube term"|
51!
52! I say "lube term" because they said it came from lubrication theory...
53! Anyway, we can change all this later (especially the names!)
54
55! If the objective function \int |\nabla u|^2,
56! the corresponding adjoint forcing is \int \nabla v \cdot \nabla u
57!
58! for the lube term, the adjoint forcing is \chi u
59!
60! This has always annoyed me...
61! because now I see one objective and one constraint
62!
64 use objective, only: objective_t
65 use design, only: design_t
67 use simulation, only: simulation_t
69
70 use num_types, only: rp
71 use field, only: field_t
72 use scratch_registry, only: neko_scratch_registry
73 use neko_config, only: neko_bcknd_device
75 use utils, only: neko_error
76 use json_module, only: json_file
77 use json_utils, only: json_get_or_default
78 use field_registry, only: neko_field_registry
79 use coefs, only: coef_t
80 use math, only: glsc2, copy
81 use device_math, only: device_copy, device_glsc2
82 use math_ext, only: glsc2_mask
83 use field_math, only: field_col3, field_addcol3, field_cmult
84 use, intrinsic :: iso_c_binding, only: c_ptr, c_null_ptr
85 implicit none
86 private
87
91 type, public, extends(objective_t) :: lube_term_objective_t
92 private
93
95 real(kind=rp) :: k
96
98 type(field_t), pointer :: u => null()
100 type(field_t), pointer :: v => null()
102 type(field_t), pointer :: w => null()
104 type(coef_t), pointer :: c_xh => null()
106 type(field_t), pointer :: brinkman_amplitude => null()
107
108 contains
109
111 procedure, public, pass(this) :: init_json => lube_term_init_json
113 procedure, public, pass(this) :: init_from_attributes => &
116 procedure, public, pass(this) :: free => lube_term_free
118 procedure, public, pass(this) :: update_value => &
121 procedure, public, pass(this) :: update_sensitivity => &
123
124 end type lube_term_objective_t
125
126contains
127
132 subroutine lube_term_init_json(this, json, design, simulation)
133 class(lube_term_objective_t), intent(inout) :: this
134 type(json_file), intent(inout) :: json
135 class(design_t), intent(in) :: design
136 type(simulation_t), target, intent(inout) :: simulation
137 character(len=:), allocatable :: mask_name
138 character(len=:), allocatable :: name
139 real(kind=rp) :: weight, k
140
141 call json_get_or_default(json, "weight", weight, 1.0_rp)
142 call json_get_or_default(json, "mask_name", mask_name, "")
143 call json_get_or_default(json, "name", name, "Out of plane stresses")
144 call json_get_or_default(json, "K", k, 1.0_rp)
145
146 call this%init_from_attributes(design, simulation, weight, name, &
147 mask_name, k)
148 end subroutine lube_term_init_json
149
156 subroutine lube_term_init_attributes(this, design, simulation, weight, &
157 name, mask_name, K)
158 class(lube_term_objective_t), intent(inout) :: this
159 class(design_t), intent(in) :: design
160 type(simulation_t), target, intent(inout) :: simulation
161 real(kind=rp), intent(in) :: weight
162 character(len=*), intent(in) :: mask_name
163 character(len=*), intent(in) :: name
164 real(kind=rp), intent(in) :: k
165 type(adjoint_lube_source_term_t) :: lube_term
166
167 ! Call the base initializer
168 call this%init_base(name, design%size(), weight, mask_name)
169
170 ! Set the coefficient for the lube term
171 this%K = k
172
173 ! Grab the brinkman amplitude for the lube term
174 select type (design)
175 type is (brinkman_design_t)
176 this%brinkman_amplitude => &
177 neko_field_registry%get_field("brinkman_amplitude")
178
179
180 class default
181 call neko_error('Minimum dissipation only works with brinkman_design')
182 end select
183
184 this%u => neko_field_registry%get_field('u')
185 this%v => neko_field_registry%get_field('v')
186 this%w => neko_field_registry%get_field('w')
187 this%c_Xh => simulation%neko_case%fluid%c_Xh
188
189 ! if we have the lube term we need to initialize and append that too
190
191 associate(f_adj_x => simulation%adjoint_case%scheme%f_adj_x, &
192 f_adj_y => simulation%adjoint_case%scheme%f_adj_y, &
193 f_adj_z => simulation%adjoint_case%scheme%f_adj_z, &
194 c_xh => simulation%adjoint_case%scheme%c_Xh)
195
196 call lube_term%init_from_components(f_adj_x, f_adj_y, f_adj_z, design, &
197 this%k * this%weight, &
198 this%u, this%v, this%w, &
199 this%mask, this%has_mask, c_xh)
200
201 end associate
202
203 ! append adjoint forcing term based on objective function
204 call simulation%adjoint_case%scheme%source_term%add_source_term(lube_term)
205
206 end subroutine lube_term_init_attributes
207
209 subroutine lube_term_free(this)
210 class(lube_term_objective_t), intent(inout) :: this
211 call this%free_base()
212
213 this%u => null()
214 this%v => null()
215 this%w => null()
216 this%c_Xh => null()
217 this%brinkman_amplitude => null()
218
219 end subroutine lube_term_free
220
225 subroutine lube_term_update_value(this, design)
226 class(lube_term_objective_t), intent(inout) :: this
227 class(design_t), intent(in) :: design
228 type(field_t), pointer :: work
229 integer :: temp_indices(1)
230
231 call neko_scratch_registry%request_field(work, temp_indices(1))
232
233 ! it's becoming so stupid to pass the whole fluid and adjoint and
234 ! design through
235 ! I feel like every objective function should have internal pointers to
236 ! u,v,w and u_adj, v_adj, w_adj and perhaps the design
237 ! (the whole design, so we get all the coeffients)
238 call field_col3(work, this%u, this%brinkman_amplitude)
239 call field_addcol3(work, this%v, this%brinkman_amplitude)
240 call field_addcol3(work, this%w, this%brinkman_amplitude)
241
242 if (this%has_mask) then
243 if (neko_bcknd_device .eq. 1) then
244 ! note, this could be done more elagantly by writing
245 ! device_glsc2_mask
246 call mask_exterior_const(work, this%mask, 0.0_rp)
247 this%value = device_glsc2(work%x_d, this%c_Xh%B_d, design%size())
248 else
249 this%value = glsc2_mask(work%x, this%c_Xh%B, design%size(), &
250 this%mask%mask, this%mask%size)
251 end if
252 else
253 if (neko_bcknd_device .eq. 1) then
254 this%value = device_glsc2(work%x_d, this%c_Xh%B_d, design%size())
255 else
256 this%value = glsc2(work%x, this%c_Xh%B, design%size())
257 end if
258 end if
259 this%value = 0.5 * this%K * this%value
260
261 call neko_scratch_registry%relinquish_field(temp_indices)
262
263 end subroutine lube_term_update_value
264
267 subroutine lube_term_update_sensitivity(this, design)
268 class(lube_term_objective_t), intent(inout) :: this
269 class(design_t), intent(in) :: design
270 type(field_t), pointer :: work
271 integer :: temp_indices(1)
272
273 ! if we have the lube term we also get an extra term in the sensitivity
274 ! K * u^2
275 ! TODO
276 ! omfg be so careful with non-dimensionalization etc
277 ! I bet this is scaled a smidge wrong (ie, track if it's 1/2 or not etc)
278 ! do this later
279
280 call neko_scratch_registry%request_field(work, temp_indices(1))
281
282 call field_col3(work, this%u, this%u)
283 call field_addcol3(work, this%v, this%v)
284 call field_addcol3(work, this%w, this%w)
285 call field_cmult(work, this%K)
286
287 if (neko_bcknd_device .eq. 1) then
288 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
289 else
290 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
291 end if
292
293 call neko_scratch_registry%relinquish_field(temp_indices)
294
295 end subroutine lube_term_update_sensitivity
296
297end module lube_term_objective
Implements the adjoint_lube_source_term_t type.
Implements the design_t.
Definition design.f90:34
Implements the lube_term_objective_t type.
subroutine lube_term_init_json(this, json, design, simulation)
The common constructor using a JSON object.
subroutine lube_term_update_value(this, design)
Compute the objective function.
subroutine lube_term_init_attributes(this, design, simulation, weight, name, mask_name, k)
The actual constructor.
subroutine lube_term_free(this)
Destructor.
subroutine lube_term_update_sensitivity(this, design)
update_value the sensitivity of the objective function with respect to $\chi$
Some common Masking operations we may need.
Definition mask_ops.f90:34
real(kind=rp) function glsc2_mask(a, b, size, mask, mask_size)
Weighted inner product for indices in the mask.
Definition math_ext.f90:117
Implements the objective_t type.
Definition objective.f90:35
Implements the steady_problem_t type.
A adjoint source term corresponding to an objective of.
A topology optimization design variable.
An abstract design type.
Definition design.f90:48
An objective function corresponding to minimum dissipation $ F = \int_\Omega |\nabla u|^2 d \Omega + ...
The abstract objective type.
Definition objective.f90:51