Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_lube_source_term.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!
36! I know this is a stupid naming convention...
37! The `lube` aspect came from a paper that attributed this term to out of plane
38! stresses based on lubrication theory.
39!
40! I preffer to think of it as a constraint that penalizes non-binary designs
41!
42! The term is $K \int_\Omega \frac{1}{2}\chi|\mathbf{u}|^2$
43!
44! the corresponding adjoint forcing is $K \chi \mathbf{u}$
46 use num_types, only: rp
47 use field_list, only: field_list_t
48 use json_module, only: json_file
49 use source_term, only: source_term_t
50 use coefs, only: coef_t
51 use interpolation, only: interpolator_t
52 use space, only: space_t, gl
53 use field, only: field_t
54 use time_state, only: time_state_t
55 use design, only: design_t
56 use brinkman_design, only: brinkman_design_t
57 use field_math, only: field_addcol3, field_copy, field_cmult, field_add2, &
58 field_col3
60 use point_zone, only: point_zone_t
61 use utils, only: neko_error
62 use field_registry, only : neko_field_registry
63 use scratch_registry, only: scratch_registry_t, neko_scratch_registry
64 use neko_config, only: neko_bcknd_device
65 use math, only: col2, invcol2
66 use device_math, only: device_col2, device_invcol2
67 implicit none
68 private
69
71 ! $K \int_\Omega \frac{1}{2}\chi|\mathbf{u}|^2$.
72 type, public, extends(source_term_t) :: adjoint_lube_source_term_t
73
75 type(field_t), pointer :: u => null()
77 type(field_t), pointer :: v => null()
79 type(field_t), pointer :: w => null()
81 type(field_t), pointer :: chi => null()
83 real(kind=rp) :: k
85 class(point_zone_t), pointer :: mask => null()
87 logical :: if_mask
89 logical :: dealias
91 type(space_t), pointer :: xh_gll
93 type(space_t), pointer :: xh_gl
95 type(coef_t), pointer :: c_xh_gl
97 type(interpolator_t), pointer :: gll_to_gl
99 type(scratch_registry_t), pointer :: scratch_gl
101 real(kind=rp) :: volume
102
103 contains
105 procedure, pass(this) :: init => adjoint_lube_source_term_init_from_json
107 procedure, pass(this) :: init_from_components => &
108 adjoint_lube_source_term_init_from_components
110 procedure, pass(this) :: free => adjoint_lube_source_term_free
112 procedure, pass(this) :: compute_ => adjoint_lube_source_term_compute
114
115contains
122 subroutine adjoint_lube_source_term_init_from_json(this, json, fields, coef, &
123 variable_name)
124 class(adjoint_lube_source_term_t), intent(inout) :: this
125 type(json_file), intent(inout) :: json
126 type(field_list_t), intent(in), target :: fields
127 type(coef_t), intent(in), target :: coef
128 character(len=*), intent(in) :: variable_name
129 ! real(kind=rp), allocatable :: values(:)
130 ! real(kind=rp) :: start_time, end_time
131
132
133 ! we shouldn't be initializing this from JSON
134 ! maybe throw an error?
135
136
137 end subroutine adjoint_lube_source_term_init_from_json
138
153 subroutine adjoint_lube_source_term_init_from_components(this, &
154 f_x, f_y, f_z, design, K, &
155 u, v, w, &
156 mask, if_mask, &
157 coef, c_Xh_GL, GLL_to_GL, dealias, volume, scratch_GL)
158 class(adjoint_lube_source_term_t), intent(inout) :: this
159 type(field_t), pointer, intent(in) :: f_x, f_y, f_z
160 class(design_t), intent(in), target :: design
161 real(kind=rp), intent(in) :: k
162 type(field_t), intent(in), target :: u, v, w
163 class(point_zone_t), intent(in), target :: mask
164 logical :: if_mask
165 type(coef_t), intent(in) :: coef
166 type(coef_t), intent(in), target :: c_Xh_GL
167 type(interpolator_t), intent(in), target :: GLL_to_GL
168 logical, intent(in) :: dealias
169 real(kind=rp), intent(in) :: volume
170 type(scratch_registry_t), intent(in), target :: scratch_gl
171 real(kind=rp) :: start_time
172 real(kind=rp) :: end_time
173 type(field_list_t) :: fields
174
175 ! I wish you didn't need a start time and end time...
176 ! but I'm just going to set a super big number...
177 start_time = 0.0_rp
178 end_time = 100000000.0_rp
179
180 call this%free()
181
182 ! this is copying the fluid source term init
183 ! We package the fields for the source term to operate on in a field list.
184 call fields%init(3)
185 call fields%assign(1, f_x)
186 call fields%assign(2, f_y)
187 call fields%assign(3, f_z)
188
189 call this%init_base(fields, coef, start_time, end_time)
190
191 ! point everything in the correct places
192 this%c_Xh_GL => c_xh_gl
193 this%Xh_GL => this%c_Xh_GL%Xh
194 this%Xh_GLL => this%coef%Xh
195 this%GLL_to_GL => gll_to_gl
196 this%dealias = dealias
197 this%scratch_GL => scratch_gl
198 this%u => u
199 this%v => v
200 this%w => w
201
202 this%volume = volume
203
204 select type (design)
205 type is (brinkman_design_t)
206 this%chi => neko_field_registry%get_field("brinkman_amplitude")
207 class default
208 call neko_error('Unknown design type')
209 end select
210
211 this%K = k
212 this%if_mask = if_mask
213 if (this%if_mask) then
214 this%mask => mask
215 end if
216
217 end subroutine adjoint_lube_source_term_init_from_components
218
220 subroutine adjoint_lube_source_term_free(this)
221 class(adjoint_lube_source_term_t), intent(inout) :: this
222
223 call this%free_base()
224 nullify(this%u)
225 nullify(this%v)
226 nullify(this%w)
227 nullify(this%c_Xh_GL)
228 nullify(this%Xh_GL)
229 nullify(this%Xh_GLL)
230 nullify(this%GLL_to_GL)
231 nullify(this%mask)
232 nullify(this%scratch_GL)
233
234 end subroutine adjoint_lube_source_term_free
235
239 subroutine adjoint_lube_source_term_compute(this, time)
240 class(adjoint_lube_source_term_t), intent(inout) :: this
241 type(time_state_t), intent(in) :: time
242 type(field_t), pointer :: fu, fv, fw
243 type(field_t), pointer :: work
244 integer :: temp_indices(1)
245 type(field_t), pointer :: accumulate, fld_GL, chi_GL
246 integer :: temp_indices_GL(3)
247 integer :: n_GL, nel
248
249 fu => this%fields%get_by_index(1)
250 fv => this%fields%get_by_index(2)
251 fw => this%fields%get_by_index(3)
252
253 ! BE SO CAREFUL HERE
254 ! It make look the same as the Brinkman term, but it's assumed that
255 ! this source term acts on the adjoint, and the u,v,w here come from
256 ! the primal
257 call neko_scratch_registry%request_field(work, temp_indices(1))
258 call field_copy(work, this%chi)
259
260 ! scale by K and volume
261 call field_cmult(work, this%K / this%volume)
262
263 ! mask
264 if (this%if_mask) then
265 call mask_exterior_const(work, this%mask, 0.0_rp)
266 end if
267
268 if (this%dealias) then
269 nel = this%coef%msh%nelv
270 n_gl = nel * this%Xh_GL%lxyz
271 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1))
272 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2))
273 call this%scratch_GL%request_field(chi_gl, temp_indices_gl(3))
274
275 call this%GLL_to_GL%map(chi_gl%x, work%x, nel, this%Xh_GL)
276
277 ! u
278 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
279 call field_col3(accumulate, chi_gl, fld_gl)
280 if (neko_bcknd_device .eq. 1) then
281 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
282 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
283 call device_invcol2(work%x_d, this%coef%B_d, work%size())
284 else
285 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
286 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
287 call invcol2(work%x, this%coef%B, work%size())
288 end if
289 call field_add2(fu, work)
290
291 ! v
292 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
293 call field_col3(accumulate, chi_gl, fld_gl)
294 if (neko_bcknd_device .eq. 1) then
295 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
296 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
297 call device_invcol2(work%x_d, this%coef%B_d, work%size())
298 else
299 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
300 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
301 call invcol2(work%x, this%coef%B, work%size())
302 end if
303 call field_add2(fv, work)
304
305 ! w
306 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
307 call field_col3(accumulate, chi_gl, fld_gl)
308 if (neko_bcknd_device .eq. 1) then
309 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
310 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
311 call device_invcol2(work%x_d, this%coef%B_d, work%size())
312 else
313 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
314 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
315 call invcol2(work%x, this%coef%B, work%size())
316 end if
317 call field_add2(fw, work)
318
319 call this%scratch_GL%relinquish_field(temp_indices_gl)
320
321 else
322 ! multiple and add the RHS
323 call field_addcol3(fu, this%u, work)
324 call field_addcol3(fv, this%v, work)
325 call field_addcol3(fw, this%w, work)
326
327 end if
328
329 call neko_scratch_registry%relinquish_field(temp_indices)
330
331 end subroutine adjoint_lube_source_term_compute
332
Implements the adjoint_lube_source_term_t type.
Implements the design_t.
Definition design.f90:34
Some common Masking operations we may need.
Definition mask_ops.f90:34
A adjoint source term corresponding to an objective of.
A topology optimization design variable.
An abstract design type.
Definition design.f90:52