48 use num_types,
only: rp
49 use field_list,
only: field_list_t
50 use json_module,
only: json_file
51 use source_term,
only: source_term_t
52 use coefs,
only: coef_t
53 use interpolation,
only: interpolator_t
54 use space,
only: space_t, gl
55 use field,
only: field_t
56 use time_state,
only: time_state_t
59 use field_math,
only: field_addcol3, field_copy, field_cmult, field_add2, &
62 use point_zone,
only: point_zone_t
63 use utils,
only: neko_error
64 use registry,
only : neko_registry
65 use scratch_registry,
only: neko_scratch_registry, scratch_registry_t
66 use neko_config,
only: neko_bcknd_device
67 use math,
only: col2, invcol2
68 use device_math,
only: device_col2, device_invcol2
78 type(field_t),
pointer :: u => null()
80 type(field_t),
pointer :: v => null()
82 type(field_t),
pointer :: w => null()
84 type(field_t),
pointer :: chi => null()
88 class(point_zone_t),
pointer :: mask => null()
94 type(space_t),
pointer :: xh_gll
96 type(space_t),
pointer :: xh_gl
98 type(coef_t),
pointer :: c_xh_gl
100 type(interpolator_t),
pointer :: gll_to_gl
102 type(scratch_registry_t),
pointer :: scratch_gl
104 real(kind=rp) :: volume
108 procedure, pass(this) :: init => adjoint_lube_source_term_init_from_json
110 procedure, pass(this) :: init_from_components => &
111 adjoint_lube_source_term_init_from_components
113 procedure, pass(this) :: free => adjoint_lube_source_term_free
115 procedure, pass(this) :: compute_ => adjoint_lube_source_term_compute
122 class(source_term_t),
allocatable,
intent(inout) :: obj
132 subroutine adjoint_lube_source_term_init_from_json(this, json, fields, coef, &
135 type(json_file),
intent(inout) :: json
136 type(field_list_t),
intent(in),
target :: fields
137 type(coef_t),
intent(in),
target :: coef
138 character(len=*),
intent(in) :: variable_name
147 end subroutine adjoint_lube_source_term_init_from_json
163 subroutine adjoint_lube_source_term_init_from_components(this, &
164 f_x, f_y, f_z, design, K, &
167 coef, c_Xh_GL, GLL_to_GL, dealias, volume, scratch_GL)
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 real(kind=rp) :: start_time
182 real(kind=rp) :: end_time
183 type(field_list_t) :: fields
188 end_time = 100000000.0_rp
195 call fields%assign(1, f_x)
196 call fields%assign(2, f_y)
197 call fields%assign(3, f_z)
199 call this%init_base(fields, coef, start_time, end_time)
202 this%c_Xh_GL => c_xh_gl
203 this%Xh_GL => this%c_Xh_GL%Xh
204 this%Xh_GLL => this%coef%Xh
205 this%GLL_to_GL => gll_to_gl
206 this%dealias = dealias
207 this%scratch_GL => scratch_gl
216 this%chi => neko_registry%get_field(
"brinkman_amplitude")
218 call neko_error(
'Unknown design type')
222 this%if_mask = if_mask
223 if (this%if_mask)
then
227 end subroutine adjoint_lube_source_term_init_from_components
230 subroutine adjoint_lube_source_term_free(this)
233 call this%free_base()
237 nullify(this%c_Xh_GL)
240 nullify(this%GLL_to_GL)
242 nullify(this%scratch_GL)
244 end subroutine adjoint_lube_source_term_free
249 subroutine adjoint_lube_source_term_compute(this, time)
251 type(time_state_t),
intent(in) :: time
252 type(field_t),
pointer :: fu, fv, fw
253 type(field_t),
pointer :: work
254 integer :: temp_indices(1)
255 type(field_t),
pointer :: accumulate, fld_gl, chi_gl
256 integer :: temp_indices_gl(3)
259 fu => this%fields%get_by_index(1)
260 fv => this%fields%get_by_index(2)
261 fw => this%fields%get_by_index(3)
267 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
268 call field_copy(work, this%chi)
271 call field_cmult(work, this%K / this%volume)
274 if (this%if_mask)
then
278 if (this%dealias)
then
279 nel = this%coef%msh%nelv
280 n_gl = nel * this%Xh_GL%lxyz
281 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1), &
283 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2), .false.)
284 call this%scratch_GL%request_field(chi_gl, temp_indices_gl(3), .false.)
286 call this%GLL_to_GL%map(chi_gl%x, work%x, nel, this%Xh_GL)
289 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
290 call field_col3(accumulate, chi_gl, fld_gl)
291 if (neko_bcknd_device .eq. 1)
then
292 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
293 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
294 call device_invcol2(work%x_d, this%coef%B_d, work%size())
296 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
297 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
298 call invcol2(work%x, this%coef%B, work%size())
300 call field_add2(fu, work)
303 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
304 call field_col3(accumulate, chi_gl, fld_gl)
305 if (neko_bcknd_device .eq. 1)
then
306 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
307 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
308 call device_invcol2(work%x_d, this%coef%B_d, work%size())
310 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
311 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
312 call invcol2(work%x, this%coef%B, work%size())
314 call field_add2(fv, work)
317 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
318 call field_col3(accumulate, chi_gl, fld_gl)
319 if (neko_bcknd_device .eq. 1)
then
320 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
321 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
322 call device_invcol2(work%x_d, this%coef%B_d, work%size())
324 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
325 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
326 call invcol2(work%x, this%coef%B, work%size())
328 call field_add2(fw, work)
330 call this%scratch_GL%relinquish_field(temp_indices_gl)
334 call field_addcol3(fu, this%u, work)
335 call field_addcol3(fv, this%v, work)
336 call field_addcol3(fw, this%w, work)
340 call neko_scratch_registry%relinquish_field(temp_indices)
342 end subroutine adjoint_lube_source_term_compute
Implements the adjoint_lube_source_term_t type.
subroutine, public adjoint_lube_source_term_allocate(obj)
Allocator for the adjoint lube source term.
Some common Masking operations we may need.
A adjoint source term corresponding to an objective of.
A topology optimization design variable.