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
56
57 ! I don't know how to name this term, but when you have a passive
58 ! scalar you get an extra term in the adjoint velocity equation, which comes
59 ! from the convective term in the passive scalar equation.
60 ! Maybe this should be called just `adjoint_scalar_convection`?
61 ! but it does come in a source term...
62
63 ! In any case,
64 ! it's a source term acting on the adjoint velocity equations, of the form:
65 ! \f$\nabla s s_adj\f$
66 type, public, extends(source_term_t) :: &
69 type(field_t), pointer :: s_adj => null()
71 type(field_t), pointer :: s => null()
72 ! --- for over-integration
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
82 logical :: dealias
84 type(scratch_registry_t), pointer :: scratch_gl
85
86 contains
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
99
100contains
108 json, fields, coef, variable_name)
109 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
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
114
115 ! this is a bit weird... because I don't think this should come from the
116 ! JSON.
117 ! Maybe we should think of all these source terms as only "appendable"
118 !
119 ! Because we'll never have the whole case here, so we'll never be able
120 ! init from components anyway...
121
122
124
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)
137 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
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
145
146 type(field_list_t) :: fields
147 real(kind=rp) :: start_time
148 real(kind=rp) :: end_time
149
150 ! I wish you didn't need a start time and end time...
151 ! but I'm just going to set a super big number...
152 start_time = 0.0_rp
153 end_time = 100000000.0_rp
154
155 call this%free()
156
157 ! this is copying the fluid source term init
158 ! We package the fields for the source term to operate on in a field list.
159 call fields%init(3)
160 call fields%assign(1, f_x)
161 call fields%assign(2, f_y)
162 call fields%assign(3, f_z)
163
164 call this%init_base(fields, coef, start_time, end_time)
165
166 ! point everything in the correct places
167 this%s_adj => s_adj
168 this%s => s
169
170 ! for over integration
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
177
178 end subroutine adjoint_scalar_convection_source_term_init_from_components
179
181 subroutine adjoint_scalar_convection_source_term_free(this)
182 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
183
184 call this%free_base()
185 nullify(this%s_adj)
186 nullify(this%s)
187 nullify(this%c_Xh_GL)
188 nullify(this%Xh_GL)
189 nullify(this%Xh_GLL)
190 nullify(this%GLL_to_GL)
191 nullify(this%scratch_GL)
192
193 end subroutine adjoint_scalar_convection_source_term_free
194
198 subroutine adjoint_scalar_convection_source_term_compute(this, time)
199 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
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)
206 integer :: n_GL, nel
207
208
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.)
213
214 fu => this%fields%get(1)
215 fv => this%fields%get(2)
216 fw => this%fields%get(3)
217
218 ! we need the term \f$\nabla s s_adj\f$
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.)
226
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)
229
230 ! u
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)
234 ! Evaluate term on GL and preempt the GLL premultiplication
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())
239 else
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())
243 end if
244 call field_sub2(fu, work)
245
246 ! v
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)
250 ! Evaluate term on GL and preempt the GLL premultiplication
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())
255 else
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())
259 end if
260 call field_sub2(fv, work)
261
262 ! w
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)
266 ! Evaluate term on GL and preempt the GLL premultiplication
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())
271 else
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())
275 end if
276 call field_sub2(fw, work)
277
278 call this%scratch_GL%relinquish_field(temp_indices_gl)
279
280 else
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)
285 end if
286
287 ! free the scratch
288 call neko_scratch_registry%relinquish_field(temp_indices)
289 end subroutine adjoint_scalar_convection_source_term_compute
290
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.