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
34!
36!
37! I promise I'll write this document properly in the future...
38!
39! But the Borval Peterson (I think) paper had an objective function
40! that had 2 terms, dissipation and this term they claimed represented
41! out of plane stresses.
42! I never really understood that extra term, I also don't think it
43! applies to 3D cases, but everyone includes it anyway.
44!
45! It appears to me to be basically a heuristic penality that targets
46! non-binary designs
47!
48! so let's call
49!
50! F = \int |\nabla u|^2 + K \int \chi \u^2
51!
52! | dissipation | |"lube term"|
53!
54! I say "lube term" because they said it came from lubrication theory...
55! Anyway, we can change all this later (especially the names!)
56
57! If the objective function \int |\nabla u|^2,
58! the corresponding adjoint forcing is \int \nabla v \cdot \nabla u
59!
60! for the lube term, the adjoint forcing is \chi u
61!
62! This has always annoyed me...
63! because now I see one objective and one constraint
64!
66 use objective, only: objective_t
67 use design, only: design_t
68 use brinkman_design, only: brinkman_design_t
69 use simulation_m, only: simulation_t
72 use num_types, only: rp
73 use field, only: field_t
74 use scratch_registry, only: neko_scratch_registry, scratch_registry_t
75 use neko_config, only: neko_bcknd_device
77 use utils, only: neko_error
78 use json_module, only: json_file
79 use json_utils, only: json_get_or_default
80 use registry, only: neko_registry
81 use interpolation, only: interpolator_t
82 use space, only: space_t, gl
83 use coefs, only: coef_t
84 use math, only: glsc2, copy, col2, invcol2
85 use device_math, only: device_copy, device_glsc2, device_col2, device_invcol2
86 use math_ext, only: glsc2_mask
87 use field_math, only: field_col3, field_addcol3, field_cmult, field_col2
88 implicit none
89 private
90
93 type, public, extends(objective_t) :: lube_term_objective_t
94 private
95
97 type(field_t), pointer :: u => null()
99 type(field_t), pointer :: v => null()
101 type(field_t), pointer :: w => null()
103 type(field_t), pointer :: brinkman_amplitude => null()
105 logical :: dealias_sensitivity
107 logical :: dealias_forcing
109 real(kind=rp) :: volume
110
111 ! ---- everything GLL ----
113 type(coef_t), pointer :: c_xh_gll
115 type(space_t), pointer :: xh_gll
116
117 ! ---- everything for GL ----
119 type(space_t), pointer :: xh_gl
121 type(coef_t), pointer :: c_xh_gl
123 type(interpolator_t), pointer :: gll_to_gl
125 type(scratch_registry_t), pointer :: scratch_gl
126
127 contains
128
130 procedure, public, pass(this) :: init_json_sim => lube_term_init_json_sim
132 procedure, public, pass(this) :: init_from_attributes => &
133 lube_term_init_attributes
135 procedure, public, pass(this) :: free => lube_term_free
137 procedure, public, pass(this) :: update_value => &
138 lube_term_update_value
140 procedure, public, pass(this) :: update_sensitivity => &
141 lube_term_update_sensitivity
142
143 end type lube_term_objective_t
144
145contains
146
152 subroutine lube_term_init_json_sim(this, json, design, simulation)
153 class(lube_term_objective_t), intent(inout) :: this
154 type(json_file), intent(inout) :: json
155 class(design_t), intent(in) :: design
156 type(simulation_t), target, intent(inout) :: simulation
157
158 character(len=:), allocatable :: mask_name
159 character(len=:), allocatable :: name
160 real(kind=rp) :: weight
161 logical :: dealias_sensitivity, dealias_forcing
162
163 call json_get_or_default(json, "weight", weight, 1.0_rp)
164 call json_get_or_default(json, "mask_name", mask_name, "")
165 call json_get_or_default(json, "name", name, "Out of plane stresses")
166 call json_get_or_default(json, "dealias_sensitivity", &
167 dealias_sensitivity, .true.)
168 call json_get_or_default(json, "dealias_forcing", &
169 dealias_forcing, .true.)
170
171 call this%init_from_attributes(design, simulation, weight, name, &
172 mask_name, dealias_sensitivity, dealias_forcing)
173 end subroutine lube_term_init_json_sim
174
184 subroutine lube_term_init_attributes(this, design, simulation, weight, &
185 name, mask_name, dealias_sensitivity, dealias_forcing)
186 class(lube_term_objective_t), intent(inout) :: this
187 class(design_t), intent(in) :: design
188 type(simulation_t), target, intent(inout) :: simulation
189 real(kind=rp), intent(in) :: weight
190 character(len=*), intent(in) :: mask_name
191 character(len=*), intent(in) :: name
192 logical, intent(in) :: dealias_sensitivity
193 logical, intent(in) :: dealias_forcing
194 type(adjoint_lube_source_term_t) :: lube_term
195
196 ! Call the base initializer
197 call this%init_base(name, design%size(), weight, mask_name)
198
199 this%dealias_forcing = dealias_forcing
200 this%dealias_sensitivity = dealias_sensitivity
201
202 ! Grab the brinkman amplitude for the lube term
203 select type (design)
204 type is (brinkman_design_t)
205 this%brinkman_amplitude => &
206 neko_registry%get_field("brinkman_amplitude")
207
208
209 class default
210 call neko_error('Minimum dissipation only works with brinkman_design')
211 end select
212
213 this%u => neko_registry%get_field('u')
214 this%v => neko_registry%get_field('v')
215 this%w => neko_registry%get_field('w')
216
217 ! GLL
218 this%c_Xh_GLL => simulation%neko_case%fluid%c_Xh
219 this%Xh_GLL => simulation%neko_case%fluid%c_Xh%Xh
220
221 ! GL
222 this%c_Xh_GL => simulation%adjoint_case%fluid_adj%c_Xh_GL
223 this%Xh_GL => this%c_Xh_GL%Xh
224
225 ! GLL to GL
226 this%GLL_to_GL => simulation%adjoint_case%fluid_adj%GLL_to_GL
227
228 this%scratch_GL => simulation%adjoint_case%fluid_adj%scratch_GL
229
230 ! compute the volume of the objective domain
231 if (this%has_mask) then
232 this%volume = compute_masked_volume(this%mask, this%c_Xh_GLL)
233 else
234 this%volume = this%c_Xh_GLL%volume
235 end if
236
237 ! if we have the lube term we need to initialize and append that too
238
239 associate(f_adj_x => simulation%adjoint_fluid%f_adj_x, &
240 f_adj_y => simulation%adjoint_fluid%f_adj_y, &
241 f_adj_z => simulation%adjoint_fluid%f_adj_z)
242
243 call lube_term%init_from_components(f_adj_x, f_adj_y, f_adj_z, design, &
244 this%weight, this%u, this%v, this%w, this%mask, &
245 this%has_mask, this%c_Xh_GLL, this%c_Xh_GL, this%GLL_to_GL, &
246 this%dealias_forcing, this%volume, &
247 this%scratch_GL)
248
249 end associate
250
251 ! append adjoint forcing term based on objective function
252 select type (f => simulation%adjoint_fluid)
253 type is (adjoint_fluid_pnpn_t)
254 call f%source_term%add_source_term(lube_term)
255 end select
256
257 end subroutine lube_term_init_attributes
258
260 subroutine lube_term_free(this)
261 class(lube_term_objective_t), intent(inout) :: this
262 call this%free_base()
263
264 this%u => null()
265 this%v => null()
266 this%w => null()
267 this%c_Xh_GLL => null()
268 this%brinkman_amplitude => null()
269 nullify(this%c_Xh_GL)
270 nullify(this%Xh_GL)
271 nullify(this%Xh_GLL)
272 nullify(this%GLL_to_GL)
273 nullify(this%scratch_GL)
274
275 end subroutine lube_term_free
276
280 subroutine lube_term_update_value(this, design)
281 class(lube_term_objective_t), intent(inout) :: this
282 class(design_t), intent(in) :: design
283 type(field_t), pointer :: work
284 integer :: temp_indices(1)
285
286 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
287
288 call field_col3(work, this%u, this%u)
289 call field_addcol3(work, this%v, this%v)
290 call field_addcol3(work, this%w, this%w)
291 call field_col2(work, this%brinkman_amplitude)
292
293 if (this%has_mask) then
294 if (neko_bcknd_device .eq. 1) then
295 ! note, this could be done more elagantly by writing
296 ! device_glsc2_mask
297 call mask_exterior_const(work, this%mask, 0.0_rp)
298 this%value = device_glsc2(work%x_d, this%c_Xh_GLL%B_d, design%size())
299 else
300 this%value = glsc2_mask(work%x, this%c_Xh_GLL%B, design%size(), &
301 this%mask%mask%get(), this%mask%size)
302 end if
303 else
304 if (neko_bcknd_device .eq. 1) then
305 this%value = device_glsc2(work%x_d, this%c_Xh_GLL%B_d, design%size())
306 else
307 this%value = glsc2(work%x, this%c_Xh_GLL%B, design%size())
308 end if
309 end if
310 this%value = 0.5_rp * this%value / this%volume
311
312 call neko_scratch_registry%relinquish_field(temp_indices)
313
314 end subroutine lube_term_update_value
315
320 subroutine lube_term_update_sensitivity(this, design)
321 class(lube_term_objective_t), intent(inout) :: this
322 class(design_t), intent(in) :: design
323 type(field_t), pointer :: work
324 integer :: temp_indices(1)
325 integer :: n_GL, nel
326 type(field_t), pointer :: accumulate, fld_GL
327 integer :: temp_indices_GL(2)
328
329 ! if we have the lube term we also get an extra term in the sensitivity
330
331 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
332
333 if(this%dealias_sensitivity) then
334 nel = this%c_Xh_GLL%msh%nelv
335 n_gl = nel * this%Xh_GL%lxyz
336 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1), &
337 .false.)
338 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2), .false.)
339
340 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
341 call field_col3(accumulate, fld_gl, fld_gl)
342 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
343 call field_addcol3(accumulate, fld_gl, fld_gl)
344 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
345 call field_addcol3(accumulate, fld_gl, fld_gl)
346 ! scale
347 call field_cmult(accumulate, this%weight * 0.5_rp / this%volume)
348
349 ! Evaluate term on GL and preempt the GLL premultiplication
350 if (neko_bcknd_device .eq. 1) then
351 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
352 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
353 call device_invcol2(work%x_d, this%c_Xh_GLL%B_d, work%size())
354 else
355 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
356 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
357 call invcol2(work%x, this%c_Xh_GLL%B, work%size())
358 end if
359
360 call this%scratch_GL%relinquish_field(temp_indices_gl)
361
362 else
363 call field_col3(work, this%u, this%u)
364 call field_addcol3(work, this%v, this%v)
365 call field_addcol3(work, this%w, this%w)
366 ! scale
367 call field_cmult(work, this%weight * 0.5_rp / this%volume)
368 end if
369
370 if (neko_bcknd_device .eq. 1) then
371 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
372 else
373 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
374 end if
375
376 call neko_scratch_registry%relinquish_field(temp_indices)
377
378 end subroutine lube_term_update_sensitivity
379
380end module lube_term_objective
Adjoint Pn/Pn formulation.
Implements the adjoint_lube_source_term_t type.
Implements the design_t.
Definition design.f90:36
Implements the lube_term_objective_t type.
subroutine lube_term_init_json_sim(this, json, design, simulation)
The common constructor using a JSON object.
Some common Masking operations we may need.
Definition mask_ops.f90:36
real(kind=rp) function, public compute_masked_volume(mask, coef)
Compute the volume of the domain contained within the mask.
Definition mask_ops.f90:182
Implements the objective_t type.
Definition objective.f90:36
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:54
An objective function corresponding to out of plane stresses .
The abstract objective type.
Definition objective.f90:52