Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
lube_term_objective.f90
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
66 use brinkman_design, only: brinkman_design_t
67 use simulation_m, 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_sim => lube_term_init_json_sim
113 procedure, public, pass(this) :: init_from_attributes => &
114 lube_term_init_attributes
116 procedure, public, pass(this) :: free => lube_term_free
118 procedure, public, pass(this) :: update_value => &
119 lube_term_update_value
121 procedure, public, pass(this) :: update_sensitivity => &
122 lube_term_update_sensitivity
123
124 end type lube_term_objective_t
125
126contains
127
133 subroutine lube_term_init_json_sim(this, json, design, simulation)
134 class(lube_term_objective_t), intent(inout) :: this
135 type(json_file), intent(inout) :: json
136 class(design_t), intent(in) :: design
137 type(simulation_t), target, intent(inout) :: simulation
138
139 character(len=:), allocatable :: mask_name
140 character(len=:), allocatable :: name
141 real(kind=rp) :: weight, k
142
143 call json_get_or_default(json, "weight", weight, 1.0_rp)
144 call json_get_or_default(json, "mask_name", mask_name, "")
145 call json_get_or_default(json, "name", name, "Out of plane stresses")
146 call json_get_or_default(json, "K", k, 1.0_rp)
147
148 call this%init_from_attributes(design, simulation, weight, name, &
149 mask_name, k)
150 end subroutine lube_term_init_json_sim
151
160 subroutine lube_term_init_attributes(this, design, simulation, weight, &
161 name, mask_name, K)
162 class(lube_term_objective_t), intent(inout) :: this
163 class(design_t), intent(in) :: design
164 type(simulation_t), target, intent(inout) :: simulation
165 real(kind=rp), intent(in) :: weight
166 character(len=*), intent(in) :: mask_name
167 character(len=*), intent(in) :: name
168 real(kind=rp), intent(in) :: k
169 type(adjoint_lube_source_term_t) :: lube_term
170
171 ! Call the base initializer
172 call this%init_base(name, design%size(), weight, mask_name)
173
174 ! Set the coefficient for the lube term
175 this%K = k
176
177 ! Grab the brinkman amplitude for the lube term
178 select type (design)
179 type is (brinkman_design_t)
180 this%brinkman_amplitude => &
181 neko_field_registry%get_field("brinkman_amplitude")
182
183
184 class default
185 call neko_error('Minimum dissipation only works with brinkman_design')
186 end select
187
188 this%u => neko_field_registry%get_field('u')
189 this%v => neko_field_registry%get_field('v')
190 this%w => neko_field_registry%get_field('w')
191 this%c_Xh => simulation%neko_case%fluid%c_Xh
192
193 ! if we have the lube term we need to initialize and append that too
194
195 associate(f_adj_x => simulation%adjoint_fluid%f_adj_x, &
196 f_adj_y => simulation%adjoint_fluid%f_adj_y, &
197 f_adj_z => simulation%adjoint_fluid%f_adj_z, &
198 c_xh => simulation%adjoint_fluid%c_Xh)
199
200 call lube_term%init_from_components(f_adj_x, f_adj_y, f_adj_z, design, &
201 this%k * this%weight, this%u, this%v, this%w, this%mask, &
202 this%has_mask, c_xh)
203
204 end associate
205
206 ! append adjoint forcing term based on objective function
207 call simulation%adjoint_fluid%source_term%add_source_term(lube_term)
208
209 end subroutine lube_term_init_attributes
210
212 subroutine lube_term_free(this)
213 class(lube_term_objective_t), intent(inout) :: this
214 call this%free_base()
215
216 this%u => null()
217 this%v => null()
218 this%w => null()
219 this%c_Xh => null()
220 this%brinkman_amplitude => null()
221
222 end subroutine lube_term_free
223
227 subroutine lube_term_update_value(this, design)
228 class(lube_term_objective_t), intent(inout) :: this
229 class(design_t), intent(in) :: design
230 type(field_t), pointer :: work
231 integer :: temp_indices(1)
232
233 call neko_scratch_registry%request_field(work, temp_indices(1))
234
235 ! it's becoming so stupid to pass the whole fluid and adjoint and
236 ! design through
237 ! I feel like every objective function should have internal pointers to
238 ! u,v,w and u_adj, v_adj, w_adj and perhaps the design
239 ! (the whole design, so we get all the coeffients)
240 call field_col3(work, this%u, this%brinkman_amplitude)
241 call field_addcol3(work, this%v, this%brinkman_amplitude)
242 call field_addcol3(work, this%w, this%brinkman_amplitude)
243
244 if (this%has_mask) then
245 if (neko_bcknd_device .eq. 1) then
246 ! note, this could be done more elagantly by writing
247 ! device_glsc2_mask
248 call mask_exterior_const(work, this%mask, 0.0_rp)
249 this%value = device_glsc2(work%x_d, this%c_Xh%B_d, design%size())
250 else
251 this%value = glsc2_mask(work%x, this%c_Xh%B, design%size(), &
252 this%mask%mask, this%mask%size)
253 end if
254 else
255 if (neko_bcknd_device .eq. 1) then
256 this%value = device_glsc2(work%x_d, this%c_Xh%B_d, design%size())
257 else
258 this%value = glsc2(work%x, this%c_Xh%B, design%size())
259 end if
260 end if
261 this%value = 0.5 * this%K * this%value
262
263 call neko_scratch_registry%relinquish_field(temp_indices)
264
265 end subroutine lube_term_update_value
266
270 subroutine lube_term_update_sensitivity(this, design)
271 class(lube_term_objective_t), intent(inout) :: this
272 class(design_t), intent(in) :: design
273 type(field_t), pointer :: work
274 integer :: temp_indices(1)
275
276 ! if we have the lube term we also get an extra term in the sensitivity
277 ! K * u^2
278 ! TODO
279 ! omfg be so careful with non-dimensionalization etc
280 ! I bet this is scaled a smidge wrong (ie, track if it's 1/2 or not etc)
281 ! do this later
282
283 call neko_scratch_registry%request_field(work, temp_indices(1))
284
285 call field_col3(work, this%u, this%u)
286 call field_addcol3(work, this%v, this%v)
287 call field_addcol3(work, this%w, this%w)
288 call field_cmult(work, this%K)
289
290 if (neko_bcknd_device .eq. 1) then
291 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
292 else
293 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
294 end if
295
296 call neko_scratch_registry%relinquish_field(temp_indices)
297
298 end subroutine lube_term_update_sensitivity
299
300end 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.
Some common Masking operations we may need.
Definition mask_ops.f90:34
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:50
An objective function corresponding to minimum dissipation .
The abstract objective type.
Definition objective.f90:51