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