Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
brinkman_dissipation_objective.f90
Go to the documentation of this file.
1
34!
36!
37! F = \int |\nabla u|^2 + K \int \chi \u^2
38!
39! |viscous dissipation | |"Brinkman dissipation"|
40!
42 use objective, only: objective_t
43 use design, only: design_t
44 use brinkman_design, only: brinkman_design_t
45 use simulation_m, only: simulation_t
49 use num_types, only: rp
50 use field, only: field_t
51 use scratch_registry, only: neko_scratch_registry, scratch_registry_t
52 use neko_config, only: neko_bcknd_device
54 use utils, only: neko_error
55 use json_module, only: json_file
56 use json_utils, only: json_get_or_default
57 use registry, only: neko_registry
58 use interpolation, only: interpolator_t
59 use space, only: space_t, gl
60 use coefs, only: coef_t
61 use math, only: glsc2, copy, col2, invcol2
62 use device_math, only: device_copy, device_glsc2, device_col2, device_invcol2
63 use math_ext, only: glsc2_mask
64 use field_math, only: field_col3, field_addcol3, field_cmult, field_col2
65 use continuation_scheduler, only: nekotop_continuation
66 implicit none
67 private
68
72 private
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 :: brinkman_amplitude => null()
83 logical :: dealias_sensitivity
85 logical :: dealias_forcing
87 real(kind=rp) :: volume
88
89 ! ---- everything GLL ----
91 type(coef_t), pointer :: c_xh_gll
93 type(space_t), pointer :: xh_gll
94
95 ! ---- everything for GL ----
97 type(space_t), pointer :: xh_gl
99 type(coef_t), pointer :: c_xh_gl
101 type(interpolator_t), pointer :: gll_to_gl
103 type(scratch_registry_t), pointer :: scratch_gl
105 integer :: gdim
106
107
108 contains
109
111 procedure, public, pass(this) :: init_json_sim => &
114 procedure, public, pass(this) :: init_from_attributes => &
115 brinkman_dissipation_init_attributes
117 procedure, public, pass(this) :: free => brinkman_dissipation_free
119 procedure, public, pass(this) :: update_value => &
120 brinkman_dissipation_update_value
122 procedure, public, pass(this) :: update_sensitivity => &
123 brinkman_dissipation_update_sensitivity
124
126
127contains
128
134 subroutine brinkman_dissipation_init_json_sim(this, json, design, simulation)
135 class(brinkman_dissipation_objective_t), intent(inout) :: this
136 type(json_file), intent(inout) :: json
137 class(design_t), intent(in) :: design
138 type(simulation_t), target, intent(inout) :: simulation
139
140 character(len=:), allocatable :: mask_name
141 character(len=:), allocatable :: name
142 real(kind=rp) :: weight
143 logical :: dealias_sensitivity, dealias_forcing
144
145 weight = 1.0_rp
146 call nekotop_continuation%json_get_or_register(json, 'weight', &
147 this%weight, weight)
148 call json_get_or_default(json, "mask_name", mask_name, "")
149 call json_get_or_default(json, "name", name, "Out of plane stresses")
150 call json_get_or_default(json, "dealias_sensitivity", &
151 dealias_sensitivity, .true.)
152 call json_get_or_default(json, "dealias_forcing", &
153 dealias_forcing, .true.)
154
155 call this%init_from_attributes(design, simulation, weight, name, &
156 mask_name, dealias_sensitivity, dealias_forcing)
158
168 subroutine brinkman_dissipation_init_attributes(this, design, simulation, &
169 weight, name, mask_name, dealias_sensitivity, dealias_forcing)
170 class(brinkman_dissipation_objective_t), intent(inout) :: this
171 class(design_t), intent(in) :: design
172 type(simulation_t), target, intent(inout) :: simulation
173 real(kind=rp), intent(in) :: weight
174 character(len=*), intent(in) :: mask_name
175 character(len=*), intent(in) :: name
176 logical, intent(in) :: dealias_sensitivity
177 logical, intent(in) :: dealias_forcing
178 type(adjoint_brinkman_dissipation_source_term_t) :: brinkman_dissipation
179
180 ! Call the base initializer
181 call this%init_base(name, design%size(), weight, mask_name)
182
183 this%dealias_forcing = dealias_forcing
184 this%dealias_sensitivity = dealias_sensitivity
185
186 ! Grab the Brinkman amplitude field.
187 select type (design)
188 type is (brinkman_design_t)
189 this%brinkman_amplitude => &
190 neko_registry%get_field("brinkman_amplitude")
191
192
193 class default
194 call neko_error('Brinkman dissipation only works with '// &
195 'brinkman_design')
196 end select
197
198 this%u => neko_registry%get_field('u')
199 this%v => neko_registry%get_field('v')
200 this%w => neko_registry%get_field('w')
201
202 ! GLL
203 this%c_Xh_GLL => simulation%neko_case%fluid%c_Xh
204 this%Xh_GLL => simulation%neko_case%fluid%c_Xh%Xh
205 this%gdim = this%c_Xh_GLL%msh%gdim
206
207 ! GL
208 this%c_Xh_GL => simulation%adjoint_case%fluid_adj%c_Xh_GL
209 this%Xh_GL => this%c_Xh_GL%Xh
210
211 ! GLL to GL
212 this%GLL_to_GL => simulation%adjoint_case%fluid_adj%GLL_to_GL
213
214 this%scratch_GL => simulation%adjoint_case%fluid_adj%scratch_GL
215
216 ! compute the volume of the objective domain
217 if (this%has_mask) then
218 this%volume = compute_masked_volume(this%mask, this%c_Xh_GLL)
219 else
220 this%volume = this%c_Xh_GLL%volume
221 end if
222
223 ! if we have the Brinkman dissipation we initialize and append it too
224
225 associate(f_adj_x => simulation%adjoint_fluid%f_adj_x, &
226 f_adj_y => simulation%adjoint_fluid%f_adj_y, &
227 f_adj_z => simulation%adjoint_fluid%f_adj_z)
228
229 call brinkman_dissipation%init_from_components(f_adj_x, f_adj_y, &
230 f_adj_z, design, this%weight, this%u, this%v, this%w, this%mask, &
231 this%has_mask, this%c_Xh_GLL, this%c_Xh_GL, this%GLL_to_GL, &
232 this%dealias_forcing, this%volume, &
233 this%scratch_GL, this%gdim)
234
235 end associate
236
237 ! append adjoint forcing term based on objective function
238 select type (f => simulation%adjoint_fluid)
239 type is (adjoint_fluid_pnpn_t)
240 call f%source_term%add_source_term(brinkman_dissipation)
241 end select
242
243 end subroutine brinkman_dissipation_init_attributes
244
246 subroutine brinkman_dissipation_free(this)
247 class(brinkman_dissipation_objective_t), intent(inout) :: this
248 call this%free_base()
249
250 this%u => null()
251 this%v => null()
252 this%w => null()
253 this%c_Xh_GLL => null()
254 this%brinkman_amplitude => null()
255 nullify(this%c_Xh_GL)
256 nullify(this%Xh_GL)
257 nullify(this%Xh_GLL)
258 nullify(this%GLL_to_GL)
259 nullify(this%scratch_GL)
260
261 end subroutine brinkman_dissipation_free
262
266 subroutine brinkman_dissipation_update_value(this, design)
267 class(brinkman_dissipation_objective_t), intent(inout) :: this
268 class(design_t), intent(in) :: design
269 type(field_t), pointer :: work
270 integer :: temp_indices(1)
271
272 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
273
274 call field_col3(work, this%u, this%u)
275 call field_addcol3(work, this%v, this%v)
276 if (this%gdim .eq. 3) then
277 call field_addcol3(work, this%w, this%w)
278 end if
279 call field_col2(work, this%brinkman_amplitude)
280
281 if (this%has_mask) then
282 if (neko_bcknd_device .eq. 1) then
283 ! note, this could be done more elagantly by writing
284 ! device_glsc2_mask
285 call mask_exterior_const(work, this%mask, 0.0_rp)
286 this%value = device_glsc2(work%x_d, this%c_Xh_GLL%B_d, design%size())
287 else
288 this%value = glsc2_mask(work%x, this%c_Xh_GLL%B, design%size(), &
289 this%mask%mask%get(), this%mask%size)
290 end if
291 else
292 if (neko_bcknd_device .eq. 1) then
293 this%value = device_glsc2(work%x_d, this%c_Xh_GLL%B_d, design%size())
294 else
295 this%value = glsc2(work%x, this%c_Xh_GLL%B, design%size())
296 end if
297 end if
298 this%value = 0.5_rp * this%value / this%volume
299
300 call neko_scratch_registry%relinquish_field(temp_indices)
301
302 end subroutine brinkman_dissipation_update_value
303
308 subroutine brinkman_dissipation_update_sensitivity(this, design)
309 class(brinkman_dissipation_objective_t), intent(inout) :: this
310 class(design_t), intent(in) :: design
311 type(field_t), pointer :: work
312 integer :: temp_indices(1)
313 integer :: n_GL, nel
314 type(field_t), pointer :: accumulate, fld_GL
315 integer :: temp_indices_GL(2)
316
317 ! The Brinkman dissipation adds an extra term in the sensitivity.
318
319 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
320
321 if(this%dealias_sensitivity) then
322 nel = this%c_Xh_GLL%msh%nelv
323 n_gl = nel * this%Xh_GL%lxyz
324 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1), &
325 .false.)
326 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2), .false.)
327
328 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
329 call field_col3(accumulate, fld_gl, fld_gl)
330 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
331 call field_addcol3(accumulate, fld_gl, fld_gl)
332 if (this%gdim .eq. 3) then
333 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
334 call field_addcol3(accumulate, fld_gl, fld_gl)
335 end if
336 ! scale
337 call field_cmult(accumulate, this%weight * 0.5_rp / this%volume)
338
339 ! Evaluate term on GL and preempt the GLL premultiplication
340 if (neko_bcknd_device .eq. 1) then
341 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
342 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
343 call device_invcol2(work%x_d, this%c_Xh_GLL%B_d, work%size())
344 else
345 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
346 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
347 call invcol2(work%x, this%c_Xh_GLL%B, work%size())
348 end if
349
350 call this%scratch_GL%relinquish_field(temp_indices_gl)
351
352 else
353 call field_col3(work, this%u, this%u)
354 call field_addcol3(work, this%v, this%v)
355 if (this%gdim .eq. 3) then
356 call field_addcol3(work, this%w, this%w)
357 end if
358 ! scale
359 call field_cmult(work, this%weight * 0.5_rp / this%volume)
360 end if
361
362 if (neko_bcknd_device .eq. 1) then
363 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
364 else
365 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
366 end if
367
368 call neko_scratch_registry%relinquish_field(temp_indices)
369
370 end subroutine brinkman_dissipation_update_sensitivity
371
Implements the adjoint_brinkman_dissipation_source_term_t type.
Adjoint Pn/Pn formulation.
Implements the brinkman_dissipation_objective_t type.
subroutine brinkman_dissipation_init_json_sim(this, json, design, simulation)
The common constructor using a JSON object.
Continuation scheduler for the optimization loop.
Implements the design_t.
Definition design.f90:36
Some common Masking operations we may need.
Definition mask_ops.f90:36
real(kind=rp) function, public compute_masked_volume(mask, coef)
Compute the volume of the domain contained within the mask.
Definition mask_ops.f90:182
Implements the objective_t type.
Definition objective.f90:36
Implements the steady_problem_t type.
A topology optimization design variable.
An objective function corresponding to out of plane stresses .
An abstract design type.
Definition design.f90:53
The abstract objective type.
Definition objective.f90:52