70 use num_types,
only: rp
71 use field,
only: field_t
72 use scratch_registry,
only: scratch_registry_t, neko_scratch_registry
73 use neko_config,
only: neko_bcknd_device
75 use utils,
only: neko_error
76 use json_module,
only: json_file
77 use json_utils,
only: json_get_or_default
78 use field_registry,
only: neko_field_registry
79 use interpolation,
only: interpolator_t
80 use space,
only: space_t, gl
81 use coefs,
only: coef_t
82 use math,
only: glsc2, copy, col2, invcol2
83 use device_math,
only: device_copy, device_glsc2, device_col2, device_invcol2
84 use math_ext,
only: glsc2_mask
85 use field_math,
only: field_col3, field_addcol3, field_cmult, field_col2
95 type(field_t),
pointer :: u => null()
97 type(field_t),
pointer :: v => null()
99 type(field_t),
pointer :: w => null()
101 type(field_t),
pointer :: brinkman_amplitude => null()
103 logical :: dealias_sensitivity
105 logical :: dealias_forcing
107 real(kind=rp) :: volume
111 type(coef_t),
pointer :: c_xh_gll
113 type(space_t),
pointer :: xh_gll
117 type(space_t),
pointer :: xh_gl
119 type(coef_t),
pointer :: c_xh_gl
121 type(interpolator_t),
pointer :: gll_to_gl
123 type(scratch_registry_t),
pointer :: scratch_gl
128 procedure,
public, pass(this) :: init_json_sim => lube_term_init_json_sim
130 procedure,
public, pass(this) :: init_from_attributes => &
131 lube_term_init_attributes
133 procedure,
public, pass(this) :: free => lube_term_free
135 procedure,
public, pass(this) :: update_value => &
136 lube_term_update_value
138 procedure,
public, pass(this) :: update_sensitivity => &
139 lube_term_update_sensitivity
150 subroutine lube_term_init_json_sim(this, json, design, simulation)
152 type(json_file),
intent(inout) :: json
153 class(
design_t),
intent(in) :: design
156 character(len=:),
allocatable :: mask_name
157 character(len=:),
allocatable :: name
158 real(kind=rp) :: weight
159 logical :: dealias_sensitivity, dealias_forcing
161 call json_get_or_default(json,
"weight", weight, 1.0_rp)
162 call json_get_or_default(json,
"mask_name", mask_name,
"")
163 call json_get_or_default(json,
"name", name,
"Out of plane stresses")
164 call json_get_or_default(json,
"dealias_sensitivity", &
165 dealias_sensitivity, .true.)
166 call json_get_or_default(json,
"dealias_forcing", &
167 dealias_forcing, .true.)
169 call this%init_from_attributes(
design, simulation, weight, name, &
170 mask_name, dealias_sensitivity, dealias_forcing)
171 end subroutine lube_term_init_json_sim
182 subroutine lube_term_init_attributes(this, design, simulation, weight, &
183 name, mask_name, dealias_sensitivity, dealias_forcing)
185 class(
design_t),
intent(in) :: design
187 real(kind=rp),
intent(in) :: weight
188 character(len=*),
intent(in) :: mask_name
189 character(len=*),
intent(in) :: name
190 logical,
intent(in) :: dealias_sensitivity
191 logical,
intent(in) :: dealias_forcing
195 call this%init_base(name,
design%size(), weight, mask_name)
197 this%dealias_forcing = dealias_forcing
198 this%dealias_sensitivity = dealias_sensitivity
203 this%brinkman_amplitude => &
204 neko_field_registry%get_field(
"brinkman_amplitude")
208 call neko_error(
'Minimum dissipation only works with brinkman_design')
211 this%u => neko_field_registry%get_field(
'u')
212 this%v => neko_field_registry%get_field(
'v')
213 this%w => neko_field_registry%get_field(
'w')
216 this%c_Xh_GLL => simulation%neko_case%fluid%c_Xh
217 this%Xh_GLL => simulation%neko_case%fluid%c_Xh%Xh
220 this%c_Xh_GL => simulation%adjoint_case%fluid_adj%c_Xh_GL
221 this%Xh_GL => this%c_Xh_GL%Xh
224 this%GLL_to_GL => simulation%adjoint_case%fluid_adj%GLL_to_GL
226 this%scratch_GL => simulation%adjoint_case%fluid_adj%scratch_GL
229 if (this%has_mask)
then
230 this%volume = compute_masked_volume(this%mask, this%c_Xh_GLL)
232 this%volume = this%c_Xh_GLL%volume
237 associate(f_adj_x => simulation%adjoint_fluid%f_adj_x, &
238 f_adj_y => simulation%adjoint_fluid%f_adj_y, &
239 f_adj_z => simulation%adjoint_fluid%f_adj_z)
241 call lube_term%init_from_components(f_adj_x, f_adj_y, f_adj_z,
design, &
242 this%weight, this%u, this%v, this%w, this%mask, &
243 this%has_mask, this%c_Xh_GLL, this%c_Xh_GL, this%GLL_to_GL, &
244 this%dealias_forcing, this%volume, &
250 select type (f => simulation%adjoint_fluid)
252 call f%source_term%add_source_term(lube_term)
255 end subroutine lube_term_init_attributes
258 subroutine lube_term_free(this)
260 call this%free_base()
265 this%c_Xh_GLL => null()
266 this%brinkman_amplitude => null()
267 nullify(this%c_Xh_GL)
270 nullify(this%GLL_to_GL)
271 nullify(this%scratch_GL)
273 end subroutine lube_term_free
278 subroutine lube_term_update_value(this, design)
280 class(design_t),
intent(in) :: design
281 type(field_t),
pointer :: work
282 integer :: temp_indices(1)
284 call neko_scratch_registry%request_field(work, temp_indices(1))
286 call field_col3(work, this%u, this%u)
287 call field_addcol3(work, this%v, this%v)
288 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 lube_term_update_value
318 subroutine lube_term_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))
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))
335 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2))
337 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
338 call field_col3(accumulate, fld_gl, fld_gl)
339 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
340 call field_addcol3(accumulate, fld_gl, fld_gl)
341 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
342 call field_addcol3(accumulate, fld_gl, fld_gl)
344 call field_cmult(accumulate, this%weight * 0.5_rp / this%volume)
347 if (neko_bcknd_device .eq. 1)
then
348 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
349 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
350 call device_invcol2(work%x_d, this%c_Xh_GLL%B_d, work%size())
352 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
353 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
354 call invcol2(work%x, this%c_Xh_GLL%B, work%size())
357 call this%scratch_GL%relinquish_field(temp_indices_gl)
360 call field_col3(work, this%u, this%u)
361 call field_addcol3(work, this%v, this%v)
362 call field_addcol3(work, this%w, this%w)
364 call field_cmult(work, this%weight * 0.5_rp / this%volume)
367 if (neko_bcknd_device .eq. 1)
then
368 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
370 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
373 call neko_scratch_registry%relinquish_field(temp_indices)
375 end subroutine lube_term_update_sensitivity
Adjoint Pn/Pn formulation.
Implements the adjoint_lube_source_term_t type.
Implements the lube_term_objective_t type.
Some common Masking operations we may need.
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.