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
70 use num_types, only: rp
71 use field, only: field_t
72 use scratch_registry, only: scratch_registry_t, neko_scratch_registry
73 use neko_config, only: neko_bcknd_device
74 use mask_ops, only: mask_exterior_const, compute_masked_volume
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 interpolation, only: interpolator_t
80 use space, only: space_t, gl
81 use coefs, only: coef_t
82 use math, only: glsc2, copy, col2, invcol2
83 use device_math, only: device_copy, device_glsc2, device_col2, device_invcol2
84 use math_ext, only: glsc2_mask
85 use field_math, only: field_col3, field_addcol3, field_cmult, field_col2
86 implicit none
87 private
88
91 type, public, extends(objective_t) :: lube_term_objective_t
92 private
93
95 type(field_t), pointer :: u => null()
97 type(field_t), pointer :: v => null()
99 type(field_t), pointer :: w => null()
101 type(field_t), pointer :: brinkman_amplitude => null()
103 logical :: dealias_sensitivity
105 logical :: dealias_forcing
107 real(kind=rp) :: volume
108
109 ! ---- everything GLL ----
111 type(coef_t), pointer :: c_xh_gll
113 type(space_t), pointer :: xh_gll
114
115 ! ---- everything for GL ----
117 type(space_t), pointer :: xh_gl
119 type(coef_t), pointer :: c_xh_gl
121 type(interpolator_t), pointer :: gll_to_gl
123 type(scratch_registry_t), pointer :: scratch_gl
124
125 contains
126
128 procedure, public, pass(this) :: init_json_sim => lube_term_init_json_sim
130 procedure, public, pass(this) :: init_from_attributes => &
131 lube_term_init_attributes
133 procedure, public, pass(this) :: free => lube_term_free
135 procedure, public, pass(this) :: update_value => &
136 lube_term_update_value
138 procedure, public, pass(this) :: update_sensitivity => &
139 lube_term_update_sensitivity
140
141 end type lube_term_objective_t
142
143contains
144
150 subroutine lube_term_init_json_sim(this, json, design, simulation)
151 class(lube_term_objective_t), intent(inout) :: this
152 type(json_file), intent(inout) :: json
153 class(design_t), intent(in) :: design
154 type(simulation_t), target, intent(inout) :: simulation
155
156 character(len=:), allocatable :: mask_name
157 character(len=:), allocatable :: name
158 real(kind=rp) :: weight
159 logical :: dealias_sensitivity, dealias_forcing
160
161 call json_get_or_default(json, "weight", weight, 1.0_rp)
162 call json_get_or_default(json, "mask_name", mask_name, "")
163 call json_get_or_default(json, "name", name, "Out of plane stresses")
164 call json_get_or_default(json, "dealias_sensitivity", &
165 dealias_sensitivity, .true.)
166 call json_get_or_default(json, "dealias_forcing", &
167 dealias_forcing, .true.)
168
169 call this%init_from_attributes(design, simulation, weight, name, &
170 mask_name, dealias_sensitivity, dealias_forcing)
171 end subroutine lube_term_init_json_sim
172
182 subroutine lube_term_init_attributes(this, design, simulation, weight, &
183 name, mask_name, dealias_sensitivity, dealias_forcing)
184 class(lube_term_objective_t), intent(inout) :: this
185 class(design_t), intent(in) :: design
186 type(simulation_t), target, intent(inout) :: simulation
187 real(kind=rp), intent(in) :: weight
188 character(len=*), intent(in) :: mask_name
189 character(len=*), intent(in) :: name
190 logical, intent(in) :: dealias_sensitivity
191 logical, intent(in) :: dealias_forcing
192 type(adjoint_lube_source_term_t) :: lube_term
193
194 ! Call the base initializer
195 call this%init_base(name, design%size(), weight, mask_name)
196
197 this%dealias_forcing = dealias_forcing
198 this%dealias_sensitivity = dealias_sensitivity
199
200 ! Grab the brinkman amplitude for the lube term
201 select type (design)
202 type is (brinkman_design_t)
203 this%brinkman_amplitude => &
204 neko_field_registry%get_field("brinkman_amplitude")
205
206
207 class default
208 call neko_error('Minimum dissipation only works with brinkman_design')
209 end select
210
211 this%u => neko_field_registry%get_field('u')
212 this%v => neko_field_registry%get_field('v')
213 this%w => neko_field_registry%get_field('w')
214
215 ! GLL
216 this%c_Xh_GLL => simulation%neko_case%fluid%c_Xh
217 this%Xh_GLL => simulation%neko_case%fluid%c_Xh%Xh
218
219 ! GL
220 this%c_Xh_GL => simulation%adjoint_case%fluid_adj%c_Xh_GL
221 this%Xh_GL => this%c_Xh_GL%Xh
222
223 ! GLL to GL
224 this%GLL_to_GL => simulation%adjoint_case%fluid_adj%GLL_to_GL
225
226 this%scratch_GL => simulation%adjoint_case%fluid_adj%scratch_GL
227
228 ! compute the volume of the objective domain
229 if (this%has_mask) then
230 this%volume = compute_masked_volume(this%mask, this%c_Xh_GLL)
231 else
232 this%volume = this%c_Xh_GLL%volume
233 end if
234
235 ! if we have the lube term we need to initialize and append that too
236
237 associate(f_adj_x => simulation%adjoint_fluid%f_adj_x, &
238 f_adj_y => simulation%adjoint_fluid%f_adj_y, &
239 f_adj_z => simulation%adjoint_fluid%f_adj_z)
240
241 call lube_term%init_from_components(f_adj_x, f_adj_y, f_adj_z, design, &
242 this%weight, this%u, this%v, this%w, this%mask, &
243 this%has_mask, this%c_Xh_GLL, this%c_Xh_GL, this%GLL_to_GL, &
244 this%dealias_forcing, this%volume, &
245 this%scratch_GL)
246
247 end associate
248
249 ! append adjoint forcing term based on objective function
250 select type (f => simulation%adjoint_fluid)
251 type is (adjoint_fluid_pnpn_t)
252 call f%source_term%add_source_term(lube_term)
253 end select
254
255 end subroutine lube_term_init_attributes
256
258 subroutine lube_term_free(this)
259 class(lube_term_objective_t), intent(inout) :: this
260 call this%free_base()
261
262 this%u => null()
263 this%v => null()
264 this%w => null()
265 this%c_Xh_GLL => null()
266 this%brinkman_amplitude => null()
267 nullify(this%c_Xh_GL)
268 nullify(this%Xh_GL)
269 nullify(this%Xh_GLL)
270 nullify(this%GLL_to_GL)
271 nullify(this%scratch_GL)
272
273 end subroutine lube_term_free
274
278 subroutine lube_term_update_value(this, design)
279 class(lube_term_objective_t), intent(inout) :: this
280 class(design_t), intent(in) :: design
281 type(field_t), pointer :: work
282 integer :: temp_indices(1)
283
284 call neko_scratch_registry%request_field(work, temp_indices(1))
285
286 call field_col3(work, this%u, this%u)
287 call field_addcol3(work, this%v, this%v)
288 call field_addcol3(work, this%w, this%w)
289 call field_col2(work, this%brinkman_amplitude)
290
291 if (this%has_mask) then
292 if (neko_bcknd_device .eq. 1) then
293 ! note, this could be done more elagantly by writing
294 ! device_glsc2_mask
295 call mask_exterior_const(work, this%mask, 0.0_rp)
296 this%value = device_glsc2(work%x_d, this%c_Xh_GLL%B_d, design%size())
297 else
298 this%value = glsc2_mask(work%x, this%c_Xh_GLL%B, design%size(), &
299 this%mask%mask%get(), this%mask%size)
300 end if
301 else
302 if (neko_bcknd_device .eq. 1) then
303 this%value = device_glsc2(work%x_d, this%c_Xh_GLL%B_d, design%size())
304 else
305 this%value = glsc2(work%x, this%c_Xh_GLL%B, design%size())
306 end if
307 end if
308 this%value = 0.5_rp * this%value / this%volume
309
310 call neko_scratch_registry%relinquish_field(temp_indices)
311
312 end subroutine lube_term_update_value
313
318 subroutine lube_term_update_sensitivity(this, design)
319 class(lube_term_objective_t), intent(inout) :: this
320 class(design_t), intent(in) :: design
321 type(field_t), pointer :: work
322 integer :: temp_indices(1)
323 integer :: n_GL, nel
324 type(field_t), pointer :: accumulate, fld_GL
325 integer :: temp_indices_GL(2)
326
327 ! if we have the lube term we also get an extra term in the sensitivity
328
329 call neko_scratch_registry%request_field(work, temp_indices(1))
330
331 if(this%dealias_sensitivity) then
332 nel = this%c_Xh_GLL%msh%nelv
333 n_gl = nel * this%Xh_GL%lxyz
334 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1))
335 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2))
336
337 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
338 call field_col3(accumulate, fld_gl, fld_gl)
339 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
340 call field_addcol3(accumulate, fld_gl, fld_gl)
341 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
342 call field_addcol3(accumulate, fld_gl, fld_gl)
343 ! scale
344 call field_cmult(accumulate, this%weight * 0.5_rp / this%volume)
345
346 ! Evaluate term on GL and preempt the GLL premultiplication
347 if (neko_bcknd_device .eq. 1) then
348 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
349 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
350 call device_invcol2(work%x_d, this%c_Xh_GLL%B_d, work%size())
351 else
352 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
353 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
354 call invcol2(work%x, this%c_Xh_GLL%B, work%size())
355 end if
356
357 call this%scratch_GL%relinquish_field(temp_indices_gl)
358
359 else
360 call field_col3(work, this%u, this%u)
361 call field_addcol3(work, this%v, this%v)
362 call field_addcol3(work, this%w, this%w)
363 ! scale
364 call field_cmult(work, this%weight * 0.5_rp / this%volume)
365 end if
366
367 if (neko_bcknd_device .eq. 1) then
368 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
369 else
370 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
371 end if
372
373 call neko_scratch_registry%relinquish_field(temp_indices)
374
375 end subroutine lube_term_update_sensitivity
376
377end module lube_term_objective
Adjoint Pn/Pn formulation.
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:52
An objective function corresponding to out of plane stresses .
The abstract objective type.
Definition objective.f90:51