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
163 subroutine adjoint_brinkman_dissipation_source_term_init_from_components( &
164 this, f_x, f_y, f_z, design, K, &
165 u, v, w, &
166 mask, if_mask, &
167 coef, c_Xh_GL, GLL_to_GL, dealias, volume, scratch_GL, gdim)
168 class(adjoint_brinkman_dissipation_source_term_t), intent(inout) :: this
169 type(field_t), pointer, intent(in) :: f_x, f_y, f_z
170 class(design_t), intent(in), target :: design
171 real(kind=rp), intent(in) :: k
172 type(field_t), intent(in), target :: u, v, w
173 class(point_zone_t), intent(in), target :: mask
174 logical :: if_mask
175 type(coef_t), intent(in) :: coef
176 type(coef_t), intent(in), target :: c_xh_gl
177 type(interpolator_t), intent(in), target :: gll_to_gl
178 logical, intent(in) :: dealias
179 real(kind=rp), intent(in) :: volume
180 type(scratch_registry_t), intent(in), target :: scratch_gl
181 integer, intent(in) :: gdim
182 real(kind=rp) :: start_time
183 real(kind=rp) :: end_time
184 type(field_list_t) :: fields
185
186 ! I wish you didn't need a start time and end time...
187 ! but I'm just going to set a super big number...
188 start_time = 0.0_rp
189 end_time = 100000000.0_rp
190
191 call this%free()
192
193 ! this is copying the fluid source term init
194 ! We package the fields for the source term to operate on in a field list.
195 call fields%init(3)
196 call fields%assign(1, f_x)
197 call fields%assign(2, f_y)
198 call fields%assign(3, f_z)
199
200 call this%init_base(fields, coef, start_time, end_time)
201
202 ! point everything in the correct places
203 this%c_Xh_GL => c_xh_gl
204 this%Xh_GL => this%c_Xh_GL%Xh
205 this%Xh_GLL => this%coef%Xh
206 this%GLL_to_GL => gll_to_gl
207 this%dealias = dealias
208 this%scratch_GL => scratch_gl
209 this%u => u
210 this%v => v
211 this%w => w
212 this%gdim = gdim
213
214 this%volume = volume
215
216 select type (design)
217 type is (brinkman_design_t)
218 this%chi => neko_registry%get_field("brinkman_amplitude")
219 class default
220 call neko_error('Unknown design type')
221 end select
222
223 this%K = k
224 this%if_mask = if_mask
225 if (this%if_mask) then
226 this%mask => mask
227 end if
228
229 end subroutine adjoint_brinkman_dissipation_source_term_init_from_components
230
232 subroutine adjoint_brinkman_dissipation_source_term_free(this)
233 class(adjoint_brinkman_dissipation_source_term_t), intent(inout) :: this
234
235 call this%free_base()
236 nullify(this%u)
237 nullify(this%v)
238 nullify(this%w)
239 nullify(this%c_Xh_GL)
240 nullify(this%Xh_GL)
241 nullify(this%Xh_GLL)
242 nullify(this%GLL_to_GL)
243 nullify(this%mask)
244 nullify(this%scratch_GL)
245
246 end subroutine adjoint_brinkman_dissipation_source_term_free
247
251 subroutine adjoint_brinkman_dissipation_source_term_compute(this, time)
252 class(adjoint_brinkman_dissipation_source_term_t), intent(inout) :: this
253 type(time_state_t), intent(in) :: time
254 type(field_t), pointer :: fu, fv, fw
255 type(field_t), pointer :: work
256 integer :: temp_indices(1)
257 type(field_t), pointer :: accumulate, fld_gl, chi_gl
258 integer :: temp_indices_gl(3)
259 integer :: n_gl, nel
260
261 fu => this%fields%get_by_index(1)
262 fv => this%fields%get_by_index(2)
263 fw => this%fields%get_by_index(3)
264
265 ! BE SO CAREFUL HERE
266 ! It make look the same as the Brinkman term, but it's assumed that
267 ! this source term acts on the adjoint, and the u,v,w here come from
268 ! the primal
269 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
270 call field_copy(work, this%chi)
271
272 ! scale by K and volume
273 call field_cmult(work, this%K / this%volume)
274
275 ! mask
276 if (this%if_mask) then
277 call mask_exterior_const(work, this%mask, 0.0_rp)
278 end if
279
280 if (this%dealias) then
281 nel = this%coef%msh%nelv
282 n_gl = nel * this%Xh_GL%lxyz
283 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1), &
284 .false.)
285 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2), .false.)
286 call this%scratch_GL%request_field(chi_gl, temp_indices_gl(3), .false.)
287
288 call this%GLL_to_GL%map(chi_gl%x, work%x, nel, this%Xh_GL)
289
290 ! u
291 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
292 call field_col3(accumulate, chi_gl, fld_gl)
293 if (neko_bcknd_device .eq. 1) then
294 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
295 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
296 call device_invcol2(work%x_d, this%coef%B_d, work%size())
297 else
298 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
299 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
300 call invcol2(work%x, this%coef%B, work%size())
301 end if
302 call field_add2(fu, work)
303
304 ! v
305 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
306 call field_col3(accumulate, chi_gl, fld_gl)
307 if (neko_bcknd_device .eq. 1) then
308 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
309 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
310 call device_invcol2(work%x_d, this%coef%B_d, work%size())
311 else
312 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
313 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
314 call invcol2(work%x, this%coef%B, work%size())
315 end if
316 call field_add2(fv, work)
317
318 ! w
319 if (this%gdim .eq. 3) then
320 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
321 call field_col3(accumulate, chi_gl, fld_gl)
322 if (neko_bcknd_device .eq. 1) then
323 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
324 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
325 call device_invcol2(work%x_d, this%coef%B_d, work%size())
326 else
327 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
328 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
329 call invcol2(work%x, this%coef%B, work%size())
330 end if
331 call field_add2(fw, work)
332 end if
333
334 call this%scratch_GL%relinquish_field(temp_indices_gl)
335
336 else
337 ! multiple and add the RHS
338 call field_addcol3(fu, this%u, work)
339 call field_addcol3(fv, this%v, work)
340 if (this%gdim .eq. 3) then
341 call field_addcol3(fw, this%w, work)
342 end if
343
344 end if
345
346 call neko_scratch_registry%relinquish_field(temp_indices)
347
348 end subroutine adjoint_brinkman_dissipation_source_term_compute
349
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