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
52 use field_math,
only: field_addcol3, field_copy, field_cmult, field_add2, &
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
68 type,
public,
extends(source_term_t) :: &
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()
82 class(point_zone_t),
pointer :: mask => null()
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
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
121 class(source_term_t),
allocatable,
intent(inout) :: obj
131 subroutine adjoint_brinkman_dissipation_source_term_init_from_json(this, &
132 json, fields, coef, variable_name)
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
146 end subroutine adjoint_brinkman_dissipation_source_term_init_from_json
163 subroutine adjoint_brinkman_dissipation_source_term_init_from_components( &
164 this, f_x, f_y, f_z, design, K, &
167 coef, c_Xh_GL, GLL_to_GL, dealias, volume, scratch_GL, gdim)
169 type(field_t),
pointer,
intent(in) :: f_x, f_y, f_z
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
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
189 end_time = 100000000.0_rp
196 call fields%assign(1, f_x)
197 call fields%assign(2, f_y)
198 call fields%assign(3, f_z)
200 call this%init_base(fields, coef, start_time, end_time)
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
218 this%chi => neko_registry%get_field(
"brinkman_amplitude")
220 call neko_error(
'Unknown design type')
224 this%if_mask = if_mask
225 if (this%if_mask)
then
229 end subroutine adjoint_brinkman_dissipation_source_term_init_from_components
232 subroutine adjoint_brinkman_dissipation_source_term_free(this)
235 call this%free_base()
239 nullify(this%c_Xh_GL)
242 nullify(this%GLL_to_GL)
244 nullify(this%scratch_GL)
246 end subroutine adjoint_brinkman_dissipation_source_term_free
251 subroutine adjoint_brinkman_dissipation_source_term_compute(this, time)
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)
261 fu => this%fields%get_by_index(1)
262 fv => this%fields%get_by_index(2)
263 fw => this%fields%get_by_index(3)
269 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
270 call field_copy(work, this%chi)
273 call field_cmult(work, this%K / this%volume)
276 if (this%if_mask)
then
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), &
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.)
288 call this%GLL_to_GL%map(chi_gl%x, work%x, nel, this%Xh_GL)
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())
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())
302 call field_add2(fu, work)
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())
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())
316 call field_add2(fv, work)
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())
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())
331 call field_add2(fw, work)
334 call this%scratch_GL%relinquish_field(temp_indices_gl)
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)
346 call neko_scratch_registry%relinquish_field(temp_indices)
348 end subroutine adjoint_brinkman_dissipation_source_term_compute
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.
Some common Masking operations we may need.
A adjoint source term corresponding to an objective of.
A topology optimization design variable.