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