38 use num_types,
only: rp
39 use field_list,
only: field_list_t
40 use field,
only: field_t
41 use json_module,
only: json_file
42 use time_state,
only: time_state_t
43 use source_term,
only: source_term_t
44 use interpolation,
only: interpolator_t
45 use space,
only: space_t, gl
46 use coefs,
only: coef_t
47 use field_math,
only: field_subcol3, field_sub2, field_col3
48 use operators,
only: grad, dudxyz
49 use utils,
only: neko_error
50 use scratch_registry,
only: neko_scratch_registry, scratch_registry_t
51 use neko_config,
only: neko_bcknd_device
52 use math,
only: col2, invcol2
53 use device_math,
only: device_col2, device_invcol2
66 type,
public,
extends(source_term_t) :: &
69 type(field_t),
pointer :: s_adj => null()
71 type(field_t),
pointer :: s => null()
74 type(space_t),
pointer :: xh_gll
76 type(space_t),
pointer :: xh_gl
78 type(coef_t),
pointer :: c_xh_gl
80 type(interpolator_t),
pointer :: gll_to_gl
84 type(scratch_registry_t),
pointer :: scratch_gl
88 procedure, pass(this) :: init => &
91 procedure, pass(this) :: init_from_components => &
92 adjoint_scalar_convection_source_term_init_from_components
94 procedure, pass(this) :: free => adjoint_scalar_convection_source_term_free
96 procedure, pass(this) :: compute_ => &
97 adjoint_scalar_convection_source_term_compute
108 json, fields, coef, variable_name)
110 type(json_file),
intent(inout) :: json
111 type(field_list_t),
intent(in),
target :: fields
112 type(coef_t),
intent(in),
target :: coef
113 character(len=*),
intent(in) :: variable_name
135 subroutine adjoint_scalar_convection_source_term_init_from_components(this,&
136 f_x, f_y, f_z, s, s_adj, coef, c_Xh_GL, GLL_to_GL, dealias, scratch_GL)
138 type(field_t),
pointer,
intent(in) :: f_x, f_y, f_z
139 type(field_t),
intent(in),
target :: s, s_adj
140 type(coef_t),
intent(in),
target :: coef
141 type(coef_t),
intent(in),
target :: c_Xh_GL
142 type(interpolator_t),
intent(in),
target :: GLL_to_GL
143 logical,
intent(in) :: dealias
144 type(scratch_registry_t),
intent(in),
target :: scratch_GL
146 type(field_list_t) :: fields
147 real(kind=rp) :: start_time
148 real(kind=rp) :: end_time
153 end_time = 100000000.0_rp
160 call fields%assign(1, f_x)
161 call fields%assign(2, f_y)
162 call fields%assign(3, f_z)
164 call this%init_base(fields, coef, start_time, end_time)
171 this%dealias = dealias
172 this%c_Xh_GL => c_xh_gl
173 this%Xh_GL => this%c_Xh_GL%Xh
174 this%Xh_GLL => this%coef%Xh
175 this%GLL_to_GL => gll_to_gl
176 this%scratch_GL => scratch_gl
178 end subroutine adjoint_scalar_convection_source_term_init_from_components
181 subroutine adjoint_scalar_convection_source_term_free(this)
184 call this%free_base()
187 nullify(this%c_Xh_GL)
190 nullify(this%GLL_to_GL)
191 nullify(this%scratch_GL)
193 end subroutine adjoint_scalar_convection_source_term_free
198 subroutine adjoint_scalar_convection_source_term_compute(this, time)
200 type(time_state_t),
intent(in) :: time
201 type(field_t),
pointer :: fu, fv, fw
202 integer :: temp_indices(4)
203 type(field_t),
pointer :: dsdx, dsdy, dsdz, work
204 type(field_t),
pointer :: accumulate, fld_GL, s_GL, s_adj_GL
205 integer :: temp_indices_GL(4)
209 call neko_scratch_registry%request_field(dsdx, temp_indices(1), .false.)
210 call neko_scratch_registry%request_field(dsdy, temp_indices(2), .false.)
211 call neko_scratch_registry%request_field(dsdz, temp_indices(3), .false.)
212 call neko_scratch_registry%request_field(work, temp_indices(4), .false.)
214 fu => this%fields%get(1)
215 fv => this%fields%get(2)
216 fw => this%fields%get(3)
219 if (this%dealias)
then
220 nel = this%coef%msh%nelv
221 n_gl = nel * this%Xh_GL%lxyz
222 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1), .false.)
223 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2), .false.)
224 call this%scratch_GL%request_field(s_gl, temp_indices_gl(3), .false.)
225 call this%scratch_GL%request_field(s_adj_gl, temp_indices_gl(4), .false.)
227 call this%GLL_to_GL%map(s_gl%x, this%s%x, nel, this%Xh_GL)
228 call this%GLL_to_GL%map(s_adj_gl%x, this%s_adj%x, nel, this%Xh_GL)
231 call dudxyz(fld_gl%x, s_gl%x, this%c_Xh_GL%drdx, &
232 this%c_Xh_GL%dsdx, this%c_Xh_GL%dtdx, this%c_Xh_GL)
233 call field_col3(accumulate, s_adj_gl, fld_gl)
235 if (neko_bcknd_device .eq. 1)
then
236 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
237 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
238 call device_invcol2(work%x_d, this%coef%B_d, work%size())
240 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
241 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
242 call invcol2(work%x, this%coef%B, work%size())
244 call field_sub2(fu, work)
247 call dudxyz(fld_gl%x, s_gl%x, this%c_Xh_GL%drdy, &
248 this%c_Xh_GL%dsdy, this%c_Xh_GL%dtdy, this%c_Xh_GL)
249 call field_col3(accumulate, s_adj_gl, fld_gl)
251 if (neko_bcknd_device .eq. 1)
then
252 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
253 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
254 call device_invcol2(work%x_d, this%coef%B_d, work%size())
256 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
257 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
258 call invcol2(work%x, this%coef%B, work%size())
260 call field_sub2(fv, work)
263 call dudxyz(fld_gl%x, s_gl%x, this%c_Xh_GL%drdz, &
264 this%c_Xh_GL%dsdz, this%c_Xh_GL%dtdz, this%c_Xh_GL)
265 call field_col3(accumulate, s_adj_gl, fld_gl)
267 if (neko_bcknd_device .eq. 1)
then
268 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
269 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
270 call device_invcol2(work%x_d, this%coef%B_d, work%size())
272 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
273 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
274 call invcol2(work%x, this%coef%B, work%size())
276 call field_sub2(fw, work)
278 call this%scratch_GL%relinquish_field(temp_indices_gl)
281 call grad(dsdx%x, dsdy%x, dsdz%x, this%s%x, this%coef)
282 call field_subcol3(fu, this%s_adj, dsdx)
283 call field_subcol3(fv, this%s_adj, dsdy)
284 call field_subcol3(fw, this%s_adj, dsdz)
288 call neko_scratch_registry%relinquish_field(temp_indices)
289 end subroutine adjoint_scalar_convection_source_term_compute
Implements the adjoint_scalar_convection_source_term type.
subroutine adjoint_scalar_convection_source_term_init_from_json(this, json, fields, coef, variable_name)
The common constructor using a JSON object.