36 use num_types,
only: rp
37 use field_list,
only: field_list_t
38 use field,
only: field_t
39 use scratch_registry,
only: neko_scratch_registry
40 use json_module,
only: json_file
41 use time_state,
only: time_state_t
42 use source_term,
only: source_term_t
43 use interpolation,
only: interpolator_t
44 use space,
only: space_t, gl
45 use coefs,
only: coef_t
46 use field_math,
only: field_subcol3, field_sub2, field_col3
47 use operators,
only: grad, dudxyz
48 use utils,
only: neko_error
49 use scratch_registry,
only: scratch_registry_t, neko_scratch_registry
50 use neko_config,
only: neko_bcknd_device
51 use math,
only: col2, invcol2
52 use device_math,
only: device_col2, device_invcol2
65 type,
public,
extends(source_term_t) :: &
68 type(field_t),
pointer :: s_adj => null()
70 type(field_t),
pointer :: s => null()
73 type(space_t),
pointer :: xh_gll
75 type(space_t),
pointer :: xh_gl
77 type(coef_t),
pointer :: c_xh_gl
79 type(interpolator_t),
pointer :: gll_to_gl
83 type(scratch_registry_t),
pointer :: scratch_gl
87 procedure, pass(this) :: init => &
88 adjoint_scalar_convection_source_term_init_from_json
90 procedure, pass(this) :: init_from_components => &
91 adjoint_scalar_convection_source_term_init_from_components
93 procedure, pass(this) :: free => adjoint_scalar_convection_source_term_free
95 procedure, pass(this) :: compute_ => &
96 adjoint_scalar_convection_source_term_compute
106 subroutine adjoint_scalar_convection_source_term_init_from_json(this, &
107 json, fields, coef, variable_name)
109 type(json_file),
intent(inout) :: json
110 type(field_list_t),
intent(in),
target :: fields
111 type(coef_t),
intent(in),
target :: coef
112 character(len=*),
intent(in) :: variable_name
122 end subroutine adjoint_scalar_convection_source_term_init_from_json
134 subroutine adjoint_scalar_convection_source_term_init_from_components(this,&
135 f_x, f_y, f_z, s, s_adj, coef, c_Xh_GL, GLL_to_GL, dealias, scratch_GL)
137 type(field_t),
pointer,
intent(in) :: f_x, f_y, f_z
138 type(field_t),
intent(in),
target :: s, s_adj
139 type(coef_t),
intent(in),
target :: coef
140 type(coef_t),
intent(in),
target :: c_Xh_GL
141 type(interpolator_t),
intent(in),
target :: GLL_to_GL
142 logical,
intent(in) :: dealias
143 type(scratch_registry_t),
intent(in),
target :: scratch_GL
145 type(field_list_t) :: fields
146 real(kind=rp) :: start_time
147 real(kind=rp) :: end_time
152 end_time = 100000000.0_rp
159 call fields%assign(1, f_x)
160 call fields%assign(2, f_y)
161 call fields%assign(3, f_z)
163 call this%init_base(fields, coef, start_time, end_time)
170 this%dealias = dealias
171 this%c_Xh_GL => c_xh_gl
172 this%Xh_GL => this%c_Xh_GL%Xh
173 this%Xh_GLL => this%coef%Xh
174 this%GLL_to_GL => gll_to_gl
175 this%scratch_GL => scratch_gl
177 end subroutine adjoint_scalar_convection_source_term_init_from_components
180 subroutine adjoint_scalar_convection_source_term_free(this)
183 call this%free_base()
186 nullify(this%c_Xh_GL)
189 nullify(this%GLL_to_GL)
190 nullify(this%scratch_GL)
192 end subroutine adjoint_scalar_convection_source_term_free
197 subroutine adjoint_scalar_convection_source_term_compute(this, time)
199 type(time_state_t),
intent(in) :: time
200 type(field_t),
pointer :: fu, fv, fw
201 integer :: temp_indices(4)
202 type(field_t),
pointer :: dsdx, dsdy, dsdz, work
203 type(field_t),
pointer :: accumulate, fld_GL, s_GL, s_adj_GL
204 integer :: temp_indices_GL(4)
208 call neko_scratch_registry%request_field(dsdx, temp_indices(1))
209 call neko_scratch_registry%request_field(dsdy, temp_indices(2))
210 call neko_scratch_registry%request_field(dsdz, temp_indices(3))
211 call neko_scratch_registry%request_field(work, temp_indices(4))
213 fu => this%fields%get(1)
214 fv => this%fields%get(2)
215 fw => this%fields%get(3)
218 if (this%dealias)
then
219 nel = this%coef%msh%nelv
220 n_gl = nel * this%Xh_GL%lxyz
221 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1))
222 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2))
223 call this%scratch_GL%request_field(s_gl, temp_indices_gl(3))
224 call this%scratch_GL%request_field(s_adj_gl, temp_indices_gl(4))
226 call this%GLL_to_GL%map(s_gl%x, this%s%x, nel, this%Xh_GL)
227 call this%GLL_to_GL%map(s_adj_gl%x, this%s_adj%x, nel, this%Xh_GL)
230 call dudxyz(fld_gl%x, s_gl%x, this%c_Xh_GL%drdx, &
231 this%c_Xh_GL%dsdx, this%c_Xh_GL%dtdx, this%c_Xh_GL)
232 call field_col3(accumulate, s_adj_gl, fld_gl)
234 if (neko_bcknd_device .eq. 1)
then
235 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
236 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
237 call device_invcol2(work%x_d, this%coef%B_d, work%size())
239 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
240 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
241 call invcol2(work%x, this%coef%B, work%size())
243 call field_sub2(fu, work)
246 call dudxyz(fld_gl%x, s_gl%x, this%c_Xh_GL%drdy, &
247 this%c_Xh_GL%dsdy, this%c_Xh_GL%dtdy, this%c_Xh_GL)
248 call field_col3(accumulate, s_adj_gl, fld_gl)
250 if (neko_bcknd_device .eq. 1)
then
251 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
252 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
253 call device_invcol2(work%x_d, this%coef%B_d, work%size())
255 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
256 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
257 call invcol2(work%x, this%coef%B, work%size())
259 call field_sub2(fv, work)
262 call dudxyz(fld_gl%x, s_gl%x, this%c_Xh_GL%drdz, &
263 this%c_Xh_GL%dsdz, this%c_Xh_GL%dtdz, this%c_Xh_GL)
264 call field_col3(accumulate, s_adj_gl, fld_gl)
266 if (neko_bcknd_device .eq. 1)
then
267 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
268 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
269 call device_invcol2(work%x_d, this%coef%B_d, work%size())
271 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
272 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
273 call invcol2(work%x, this%coef%B, work%size())
275 call field_sub2(fw, work)
277 call this%scratch_GL%relinquish_field(temp_indices_gl)
280 call grad(dsdx%x, dsdy%x, dsdz%x, this%s%x, this%coef)
281 call field_subcol3(fu, this%s_adj, dsdx)
282 call field_subcol3(fv, this%s_adj, dsdy)
283 call field_subcol3(fw, this%s_adj, dsdz)
287 call neko_scratch_registry%relinquish_field(temp_indices)
288 end subroutine adjoint_scalar_convection_source_term_compute
Implements the adjoint_scalar_convection_source_term type.