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
71
73 ! $K \int_\Omega \frac{1}{2}\chi|\mathbf{u}|^2$.
74 type, public, extends(source_term_t) :: adjoint_lube_source_term_t
75
77 type(field_t), pointer :: u => null()
79 type(field_t), pointer :: v => null()
81 type(field_t), pointer :: w => null()
83 type(field_t), pointer :: chi => null()
85 real(kind=rp) :: k
87 class(point_zone_t), pointer :: mask => null()
89 logical :: if_mask
91 logical :: dealias
93 type(space_t), pointer :: xh_gll
95 type(space_t), pointer :: xh_gl
97 type(coef_t), pointer :: c_xh_gl
99 type(interpolator_t), pointer :: gll_to_gl
101 type(scratch_registry_t), pointer :: scratch_gl
103 real(kind=rp) :: volume
104
105 contains
107 procedure, pass(this) :: init => adjoint_lube_source_term_init_from_json
109 procedure, pass(this) :: init_from_components => &
110 adjoint_lube_source_term_init_from_components
112 procedure, pass(this) :: free => adjoint_lube_source_term_free
114 procedure, pass(this) :: compute_ => adjoint_lube_source_term_compute
116
117contains
124 subroutine adjoint_lube_source_term_init_from_json(this, json, fields, coef, &
125 variable_name)
126 class(adjoint_lube_source_term_t), intent(inout) :: this
127 type(json_file), intent(inout) :: json
128 type(field_list_t), intent(in), target :: fields
129 type(coef_t), intent(in), target :: coef
130 character(len=*), intent(in) :: variable_name
131 ! real(kind=rp), allocatable :: values(:)
132 ! real(kind=rp) :: start_time, end_time
133
134
135 ! we shouldn't be initializing this from JSON
136 ! maybe throw an error?
137
138
140
155 subroutine adjoint_lube_source_term_init_from_components(this, &
156 f_x, f_y, f_z, design, K, &
157 u, v, w, &
158 mask, if_mask, &
159 coef, c_Xh_GL, GLL_to_GL, dealias, volume, scratch_GL)
160 class(adjoint_lube_source_term_t), intent(inout) :: this
161 type(field_t), pointer, intent(in) :: f_x, f_y, f_z
162 class(design_t), intent(in), target :: design
163 real(kind=rp), intent(in) :: k
164 type(field_t), intent(in), target :: u, v, w
165 class(point_zone_t), intent(in), target :: mask
166 logical :: if_mask
167 type(coef_t), intent(in) :: coef
168 type(coef_t), intent(in), target :: c_Xh_GL
169 type(interpolator_t), intent(in), target :: GLL_to_GL
170 logical, intent(in) :: dealias
171 real(kind=rp), intent(in) :: volume
172 type(scratch_registry_t), intent(in), target :: scratch_gl
173 real(kind=rp) :: start_time
174 real(kind=rp) :: end_time
175 type(field_list_t) :: fields
176
177 ! I wish you didn't need a start time and end time...
178 ! but I'm just going to set a super big number...
179 start_time = 0.0_rp
180 end_time = 100000000.0_rp
181
182 call this%free()
183
184 ! this is copying the fluid source term init
185 ! We package the fields for the source term to operate on in a field list.
186 call fields%init(3)
187 call fields%assign(1, f_x)
188 call fields%assign(2, f_y)
189 call fields%assign(3, f_z)
190
191 call this%init_base(fields, coef, start_time, end_time)
192
193 ! point everything in the correct places
194 this%c_Xh_GL => c_xh_gl
195 this%Xh_GL => this%c_Xh_GL%Xh
196 this%Xh_GLL => this%coef%Xh
197 this%GLL_to_GL => gll_to_gl
198 this%dealias = dealias
199 this%scratch_GL => scratch_gl
200 this%u => u
201 this%v => v
202 this%w => w
203
204 this%volume = volume
205
206 select type (design)
207 type is (brinkman_design_t)
208 this%chi => neko_registry%get_field("brinkman_amplitude")
209 class default
210 call neko_error('Unknown design type')
211 end select
212
213 this%K = k
214 this%if_mask = if_mask
215 if (this%if_mask) then
216 this%mask => mask
217 end if
218
219 end subroutine adjoint_lube_source_term_init_from_components
220
222 subroutine adjoint_lube_source_term_free(this)
223 class(adjoint_lube_source_term_t), intent(inout) :: this
224
225 call this%free_base()
226 nullify(this%u)
227 nullify(this%v)
228 nullify(this%w)
229 nullify(this%c_Xh_GL)
230 nullify(this%Xh_GL)
231 nullify(this%Xh_GLL)
232 nullify(this%GLL_to_GL)
233 nullify(this%mask)
234 nullify(this%scratch_GL)
235
236 end subroutine adjoint_lube_source_term_free
237
241 subroutine adjoint_lube_source_term_compute(this, time)
242 class(adjoint_lube_source_term_t), intent(inout) :: this
243 type(time_state_t), intent(in) :: time
244 type(field_t), pointer :: fu, fv, fw
245 type(field_t), pointer :: work
246 integer :: temp_indices(1)
247 type(field_t), pointer :: accumulate, fld_GL, chi_GL
248 integer :: temp_indices_GL(3)
249 integer :: n_GL, nel
250
251 fu => this%fields%get_by_index(1)
252 fv => this%fields%get_by_index(2)
253 fw => this%fields%get_by_index(3)
254
255 ! BE SO CAREFUL HERE
256 ! It make look the same as the Brinkman term, but it's assumed that
257 ! this source term acts on the adjoint, and the u,v,w here come from
258 ! the primal
259 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
260 call field_copy(work, this%chi)
261
262 ! scale by K and volume
263 call field_cmult(work, this%K / this%volume)
264
265 ! mask
266 if (this%if_mask) then
267 call mask_exterior_const(work, this%mask, 0.0_rp)
268 end if
269
270 if (this%dealias) then
271 nel = this%coef%msh%nelv
272 n_gl = nel * this%Xh_GL%lxyz
273 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1), &
274 .false.)
275 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2), .false.)
276 call this%scratch_GL%request_field(chi_gl, temp_indices_gl(3), .false.)
277
278 call this%GLL_to_GL%map(chi_gl%x, work%x, nel, this%Xh_GL)
279
280 ! u
281 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
282 call field_col3(accumulate, chi_gl, fld_gl)
283 if (neko_bcknd_device .eq. 1) then
284 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
285 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
286 call device_invcol2(work%x_d, this%coef%B_d, work%size())
287 else
288 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
289 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
290 call invcol2(work%x, this%coef%B, work%size())
291 end if
292 call field_add2(fu, work)
293
294 ! v
295 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
296 call field_col3(accumulate, chi_gl, fld_gl)
297 if (neko_bcknd_device .eq. 1) then
298 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
299 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
300 call device_invcol2(work%x_d, this%coef%B_d, work%size())
301 else
302 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
303 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
304 call invcol2(work%x, this%coef%B, work%size())
305 end if
306 call field_add2(fv, work)
307
308 ! w
309 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
310 call field_col3(accumulate, chi_gl, fld_gl)
311 if (neko_bcknd_device .eq. 1) then
312 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
313 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
314 call device_invcol2(work%x_d, this%coef%B_d, work%size())
315 else
316 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
317 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
318 call invcol2(work%x, this%coef%B, work%size())
319 end if
320 call field_add2(fw, work)
321
322 call this%scratch_GL%relinquish_field(temp_indices_gl)
323
324 else
325 ! multiple and add the RHS
326 call field_addcol3(fu, this%u, work)
327 call field_addcol3(fv, this%v, work)
328 call field_addcol3(fw, this%w, work)
329
330 end if
331
332 call neko_scratch_registry%relinquish_field(temp_indices)
333
334 end subroutine adjoint_lube_source_term_compute
335
Implements the adjoint_lube_source_term_t type.
subroutine adjoint_lube_source_term_init_from_json(this, json, fields, coef, variable_name)
The common constructor using a JSON object.
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