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
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
91 type(coef_t),
pointer :: c_xh_gll
93 type(space_t),
pointer :: xh_gll
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
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
136 type(json_file),
intent(inout) :: json
137 class(
design_t),
intent(in) :: design
140 character(len=:),
allocatable :: mask_name
141 character(len=:),
allocatable :: name
142 real(kind=rp) :: weight
143 logical :: dealias_sensitivity, dealias_forcing
144 real(kind=rp) :: start_time, end_time
146 call nekotop_continuation%json_get_or_register(json,
'weight', &
147 this%weight, weight, 1.0_rp)
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 call json_get_or_default(json,
"start_time", start_time, 0.0_rp)
155 call json_get_or_default(json,
"end_time", end_time, huge(0.0_rp))
157 call this%init_from_attributes(
design, simulation, weight, name, &
158 mask_name, dealias_sensitivity, dealias_forcing, start_time, &
173 subroutine brinkman_dissipation_init_attributes(this, design, simulation, &
174 weight, name, mask_name, dealias_sensitivity, dealias_forcing, &
175 start_time, end_time)
177 class(
design_t),
intent(in) :: design
179 real(kind=rp),
intent(in) :: weight
180 character(len=*),
intent(in) :: name
181 character(len=*),
intent(in) :: mask_name
182 logical,
intent(in) :: dealias_sensitivity
183 logical,
intent(in) :: dealias_forcing
184 real(kind=rp),
intent(in) :: start_time
185 real(kind=rp),
intent(in) :: end_time
190 call this%init_base(name,
design%size(), weight, mask_name, &
191 start_time, end_time)
193 this%dealias_forcing = dealias_forcing
194 this%dealias_sensitivity = dealias_sensitivity
199 this%brinkman_amplitude => &
200 neko_registry%get_field(
"brinkman_amplitude")
204 call neko_error(
'Brinkman dissipation only works with '// &
208 this%u => neko_registry%get_field(
'u')
209 this%v => neko_registry%get_field(
'v')
210 this%w => neko_registry%get_field(
'w')
213 this%c_Xh_GLL => simulation%neko_case%fluid%c_Xh
214 this%Xh_GLL => simulation%neko_case%fluid%c_Xh%Xh
215 this%gdim = this%c_Xh_GLL%msh%gdim
218 this%c_Xh_GL => simulation%adjoint_case%fluid_adj%c_Xh_GL
219 this%Xh_GL => this%c_Xh_GL%Xh
222 this%GLL_to_GL => simulation%adjoint_case%fluid_adj%GLL_to_GL
224 this%scratch_GL => simulation%adjoint_case%fluid_adj%scratch_GL
227 if (this%has_mask)
then
230 this%volume = this%c_Xh_GLL%volume
235 associate(f_adj_x => simulation%adjoint_fluid%f_adj_x, &
236 f_adj_y => simulation%adjoint_fluid%f_adj_y, &
237 f_adj_z => simulation%adjoint_fluid%f_adj_z)
239 call brinkman_dissipation%init_from_components(f_adj_x, f_adj_y, &
240 f_adj_z,
design, this%weight, this%u, this%v, this%w, this%mask, &
241 this%has_mask, this%c_Xh_GLL, this%c_Xh_GL, this%GLL_to_GL, &
242 this%dealias_forcing, this%volume, this%scratch_GL, this%gdim, &
243 this%start_time, this%end_time)
248 select type (f => simulation%adjoint_fluid)
250 call f%source_term%add_source_term(brinkman_dissipation)
253 end subroutine brinkman_dissipation_init_attributes
256 subroutine brinkman_dissipation_free(this)
258 call this%free_base()
263 this%c_Xh_GLL => null()
264 this%brinkman_amplitude => null()
265 nullify(this%c_Xh_GL)
268 nullify(this%GLL_to_GL)
269 nullify(this%scratch_GL)
271 end subroutine brinkman_dissipation_free
276 subroutine brinkman_dissipation_update_value(this, design)
278 class(design_t),
intent(in) :: design
279 type(field_t),
pointer :: work
280 integer :: temp_indices(1)
282 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
284 call field_col3(work, this%u, this%u)
285 call field_addcol3(work, this%v, this%v)
286 if (this%gdim .eq. 3)
then
287 call field_addcol3(work, this%w, this%w)
289 call field_col2(work, this%brinkman_amplitude)
291 if (this%has_mask)
then
292 if (neko_bcknd_device .eq. 1)
then
295 call mask_exterior_const(work, this%mask, 0.0_rp)
296 this%value = device_glsc2(work%x_d, this%c_Xh_GLL%B_d,
design%size())
298 this%value = glsc2_mask(work%x, this%c_Xh_GLL%B,
design%size(), &
299 this%mask%mask%get(), this%mask%size)
302 if (neko_bcknd_device .eq. 1)
then
303 this%value = device_glsc2(work%x_d, this%c_Xh_GLL%B_d,
design%size())
305 this%value = glsc2(work%x, this%c_Xh_GLL%B,
design%size())
308 this%value = 0.5_rp * this%value / this%volume
310 call neko_scratch_registry%relinquish_field(temp_indices)
312 end subroutine brinkman_dissipation_update_value
318 subroutine brinkman_dissipation_update_sensitivity(this, design)
320 class(design_t),
intent(in) :: design
321 type(field_t),
pointer :: work
322 integer :: temp_indices(1)
324 type(field_t),
pointer :: accumulate, fld_GL
325 integer :: temp_indices_GL(2)
329 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
331 if(this%dealias_sensitivity)
then
332 nel = this%c_Xh_GLL%msh%nelv
333 n_gl = nel * this%Xh_GL%lxyz
334 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1), &
336 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2), .false.)
338 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
339 call field_col3(accumulate, fld_gl, fld_gl)
340 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
341 call field_addcol3(accumulate, fld_gl, fld_gl)
342 if (this%gdim .eq. 3)
then
343 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
344 call field_addcol3(accumulate, fld_gl, fld_gl)
347 call field_cmult(accumulate, this%weight * 0.5_rp / this%volume)
350 if (neko_bcknd_device .eq. 1)
then
351 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
352 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
353 call device_invcol2(work%x_d, this%c_Xh_GLL%B_d, work%size())
355 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
356 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
357 call invcol2(work%x, this%c_Xh_GLL%B, work%size())
360 call this%scratch_GL%relinquish_field(temp_indices_gl)
363 call field_col3(work, this%u, this%u)
364 call field_addcol3(work, this%v, this%v)
365 if (this%gdim .eq. 3)
then
366 call field_addcol3(work, this%w, this%w)
369 call field_cmult(work, this%weight * 0.5_rp / this%volume)
372 if (neko_bcknd_device .eq. 1)
then
373 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
375 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
378 call neko_scratch_registry%relinquish_field(temp_indices)
380 end subroutine brinkman_dissipation_update_sensitivity
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.
Some common Masking operations we may need.
real(kind=rp) function, public compute_masked_volume(mask, coef)
Compute the volume of the domain contained within the mask.
Implements the objective_t type.
Implements the steady_problem_t type.
A adjoint source term corresponding to an objective of.
A topology optimization design variable.
An objective function corresponding to out of plane stresses .
The abstract objective type.