Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_scalar_convection_source_term.f90
Go to the documentation of this file.
1
34!
36! this is a such a dumb name
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
54 implicit none
55 private
57
58 ! I don't know how to name this term, but when you have a passive
59 ! scalar you get an extra term in the adjoint velocity equation, which comes
60 ! from the convective term in the passive scalar equation.
61 ! Maybe this should be called just `adjoint_scalar_convection`?
62 ! but it does come in a source term...
63
64 ! In any case,
65 ! it's a source term acting on the adjoint velocity equations, of the form:
66 ! \f$\nabla s s_adj\f$
67 type, public, extends(source_term_t) :: &
70 type(field_t), pointer :: s_adj => null()
72 type(field_t), pointer :: s => null()
73 ! --- for over-integration
75 type(space_t), pointer :: xh_gll
77 type(space_t), pointer :: xh_gl
79 type(coef_t), pointer :: c_xh_gl
81 type(interpolator_t), pointer :: gll_to_gl
83 logical :: dealias
85 type(scratch_registry_t), pointer :: scratch_gl
86
87 contains
89 procedure, pass(this) :: init => &
90 adjoint_scalar_convection_source_term_init_from_json
92 procedure, pass(this) :: init_from_components => &
93 adjoint_scalar_convection_source_term_init_from_components
95 procedure, pass(this) :: free => adjoint_scalar_convection_source_term_free
97 procedure, pass(this) :: compute_ => &
98 adjoint_scalar_convection_source_term_compute
100
101contains
102
105 class(source_term_t), allocatable, intent(inout) :: obj
108
115 subroutine adjoint_scalar_convection_source_term_init_from_json(this, &
116 json, fields, coef, variable_name)
117 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
118 type(json_file), intent(inout) :: json
119 type(field_list_t), intent(in), target :: fields
120 type(coef_t), intent(in), target :: coef
121 character(len=*), intent(in) :: variable_name
122
123 ! this is a bit weird... because I don't think this should come from the
124 ! JSON.
125 ! Maybe we should think of all these source terms as only "appendable"
126 !
127 ! Because we'll never have the whole case here, so we'll never be able
128 ! init from components anyway...
129
130
131 end subroutine adjoint_scalar_convection_source_term_init_from_json
132
143 subroutine adjoint_scalar_convection_source_term_init_from_components(this,&
144 f_x, f_y, f_z, s, s_adj, coef, c_Xh_GL, GLL_to_GL, dealias, scratch_GL)
145 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
146 type(field_t), pointer, intent(in) :: f_x, f_y, f_z
147 type(field_t), intent(in), target :: s, s_adj
148 type(coef_t), intent(in), target :: coef
149 type(coef_t), intent(in), target :: c_xh_gl
150 type(interpolator_t), intent(in), target :: gll_to_gl
151 logical, intent(in) :: dealias
152 type(scratch_registry_t), intent(in), target :: scratch_gl
153
154 type(field_list_t) :: fields
155 real(kind=rp) :: start_time
156 real(kind=rp) :: end_time
157
158 ! I wish you didn't need a start time and end time...
159 ! but I'm just going to set a super big number...
160 start_time = 0.0_rp
161 end_time = 100000000.0_rp
162
163 call this%free()
164
165 ! this is copying the fluid source term init
166 ! We package the fields for the source term to operate on in a field list.
167 call fields%init(3)
168 call fields%assign(1, f_x)
169 call fields%assign(2, f_y)
170 call fields%assign(3, f_z)
171
172 call this%init_base(fields, coef, start_time, end_time)
173
174 ! point everything in the correct places
175 this%s_adj => s_adj
176 this%s => s
177
178 ! for over integration
179 this%dealias = dealias
180 this%c_Xh_GL => c_xh_gl
181 this%Xh_GL => this%c_Xh_GL%Xh
182 this%Xh_GLL => this%coef%Xh
183 this%GLL_to_GL => gll_to_gl
184 this%scratch_GL => scratch_gl
185
186 end subroutine adjoint_scalar_convection_source_term_init_from_components
187
189 subroutine adjoint_scalar_convection_source_term_free(this)
190 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
191
192 call this%free_base()
193 nullify(this%s_adj)
194 nullify(this%s)
195 nullify(this%c_Xh_GL)
196 nullify(this%Xh_GL)
197 nullify(this%Xh_GLL)
198 nullify(this%GLL_to_GL)
199 nullify(this%scratch_GL)
200
201 end subroutine adjoint_scalar_convection_source_term_free
202
206 subroutine adjoint_scalar_convection_source_term_compute(this, time)
207 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
208 type(time_state_t), intent(in) :: time
209 type(field_t), pointer :: fu, fv, fw
210 integer :: temp_indices(4)
211 type(field_t), pointer :: dsdx, dsdy, dsdz, work
212 type(field_t), pointer :: accumulate, fld_gl, s_gl, s_adj_gl
213 integer :: temp_indices_gl(4)
214 integer :: n_gl, nel
215
216
217 call neko_scratch_registry%request_field(dsdx, temp_indices(1), .false.)
218 call neko_scratch_registry%request_field(dsdy, temp_indices(2), .false.)
219 call neko_scratch_registry%request_field(dsdz, temp_indices(3), .false.)
220 call neko_scratch_registry%request_field(work, temp_indices(4), .false.)
221
222 fu => this%fields%get(1)
223 fv => this%fields%get(2)
224 fw => this%fields%get(3)
225
226 ! we need the term \f$\nabla s s_adj\f$
227 if (this%dealias) then
228 nel = this%coef%msh%nelv
229 n_gl = nel * this%Xh_GL%lxyz
230 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1), .false.)
231 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2), .false.)
232 call this%scratch_GL%request_field(s_gl, temp_indices_gl(3), .false.)
233 call this%scratch_GL%request_field(s_adj_gl, temp_indices_gl(4), .false.)
234
235 call this%GLL_to_GL%map(s_gl%x, this%s%x, nel, this%Xh_GL)
236 call this%GLL_to_GL%map(s_adj_gl%x, this%s_adj%x, nel, this%Xh_GL)
237
238 ! u
239 call dudxyz(fld_gl%x, s_gl%x, this%c_Xh_GL%drdx, &
240 this%c_Xh_GL%dsdx, this%c_Xh_GL%dtdx, this%c_Xh_GL)
241 call field_col3(accumulate, s_adj_gl, fld_gl)
242 ! Evaluate term on GL and preempt the GLL premultiplication
243 if (neko_bcknd_device .eq. 1) then
244 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
245 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
246 call device_invcol2(work%x_d, this%coef%B_d, work%size())
247 else
248 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
249 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
250 call invcol2(work%x, this%coef%B, work%size())
251 end if
252 call field_sub2(fu, work)
253
254 ! v
255 call dudxyz(fld_gl%x, s_gl%x, this%c_Xh_GL%drdy, &
256 this%c_Xh_GL%dsdy, this%c_Xh_GL%dtdy, this%c_Xh_GL)
257 call field_col3(accumulate, s_adj_gl, fld_gl)
258 ! Evaluate term on GL and preempt the GLL premultiplication
259 if (neko_bcknd_device .eq. 1) then
260 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
261 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
262 call device_invcol2(work%x_d, this%coef%B_d, work%size())
263 else
264 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
265 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
266 call invcol2(work%x, this%coef%B, work%size())
267 end if
268 call field_sub2(fv, work)
269
270 ! w
271 call dudxyz(fld_gl%x, s_gl%x, this%c_Xh_GL%drdz, &
272 this%c_Xh_GL%dsdz, this%c_Xh_GL%dtdz, this%c_Xh_GL)
273 call field_col3(accumulate, s_adj_gl, fld_gl)
274 ! Evaluate term on GL and preempt the GLL premultiplication
275 if (neko_bcknd_device .eq. 1) then
276 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
277 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
278 call device_invcol2(work%x_d, this%coef%B_d, work%size())
279 else
280 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
281 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
282 call invcol2(work%x, this%coef%B, work%size())
283 end if
284 call field_sub2(fw, work)
285
286 call this%scratch_GL%relinquish_field(temp_indices_gl)
287
288 else
289 call grad(dsdx%x, dsdy%x, dsdz%x, this%s%x, this%coef)
290 call field_subcol3(fu, this%s_adj, dsdx)
291 call field_subcol3(fv, this%s_adj, dsdy)
292 call field_subcol3(fw, this%s_adj, dsdz)
293 end if
294
295 ! free the scratch
296 call neko_scratch_registry%relinquish_field(temp_indices)
297 end subroutine adjoint_scalar_convection_source_term_compute
298
Implements the adjoint_scalar_convection_source_term type.
subroutine, public adjoint_scalar_convection_source_term_allocate(obj)
Allocator for the adjoint scalar convection source term.