Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_scalar_convection_source_term.f90
1! Copyright (c) 2023, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
34! this is a such a dumb name
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
53 implicit none
54 private
55
56 ! I don't know how to name this term, but when you have a passive
57 ! scalar you get an extra term in the adjoint velocity equation, which comes
58 ! from the convective term in the passive scalar equation.
59 ! Maybe this should be called just `adjoint_scalar_convection`?
60 ! but it does come in a source term...
61
62 ! In any case,
63 ! it's a source term acting on the adjoint velocity equations, of the form:
64 ! \f$\nabla s s_adj\f$
65 type, public, extends(source_term_t) :: &
68 type(field_t), pointer :: s_adj => null()
70 type(field_t), pointer :: s => null()
71 ! --- for over-integration
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
81 logical :: dealias
83 type(scratch_registry_t), pointer :: scratch_gl
84
85 contains
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
98
99contains
106 subroutine adjoint_scalar_convection_source_term_init_from_json(this, &
107 json, fields, coef, variable_name)
108 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
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
113
114 ! this is a bit weird... because I don't think this should come from the
115 ! JSON.
116 ! Maybe we should think of all these source terms as only "appendable"
117 !
118 ! Because we'll never have the whole case here, so we'll never be able
119 ! init from components anyway...
120
121
122 end subroutine adjoint_scalar_convection_source_term_init_from_json
123
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)
136 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
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
144
145 type(field_list_t) :: fields
146 real(kind=rp) :: start_time
147 real(kind=rp) :: end_time
148
149 ! I wish you didn't need a start time and end time...
150 ! but I'm just going to set a super big number...
151 start_time = 0.0_rp
152 end_time = 100000000.0_rp
153
154 call this%free()
155
156 ! this is copying the fluid source term init
157 ! We package the fields for the source term to operate on in a field list.
158 call fields%init(3)
159 call fields%assign(1, f_x)
160 call fields%assign(2, f_y)
161 call fields%assign(3, f_z)
162
163 call this%init_base(fields, coef, start_time, end_time)
164
165 ! point everything in the correct places
166 this%s_adj => s_adj
167 this%s => s
168
169 ! for over integration
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
176
177 end subroutine adjoint_scalar_convection_source_term_init_from_components
178
180 subroutine adjoint_scalar_convection_source_term_free(this)
181 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
182
183 call this%free_base()
184 nullify(this%s_adj)
185 nullify(this%s)
186 nullify(this%c_Xh_GL)
187 nullify(this%Xh_GL)
188 nullify(this%Xh_GLL)
189 nullify(this%GLL_to_GL)
190 nullify(this%scratch_GL)
191
192 end subroutine adjoint_scalar_convection_source_term_free
193
197 subroutine adjoint_scalar_convection_source_term_compute(this, time)
198 class(adjoint_scalar_convection_source_term_t), intent(inout) :: this
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)
205 integer :: n_GL, nel
206
207
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))
212
213 fu => this%fields%get(1)
214 fv => this%fields%get(2)
215 fw => this%fields%get(3)
216
217 ! we need the term \f$\nabla s s_adj\f$
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))
225
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)
228
229 ! u
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)
233 ! Evaluate term on GL and preempt the GLL premultiplication
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())
238 else
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())
242 end if
243 call field_sub2(fu, work)
244
245 ! v
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)
249 ! Evaluate term on GL and preempt the GLL premultiplication
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())
254 else
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())
258 end if
259 call field_sub2(fv, work)
260
261 ! w
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)
265 ! Evaluate term on GL and preempt the GLL premultiplication
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())
270 else
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())
274 end if
275 call field_sub2(fw, work)
276
277 call this%scratch_GL%relinquish_field(temp_indices_gl)
278
279 else
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)
284 end if
285
286 ! free the scratch
287 call neko_scratch_registry%relinquish_field(temp_indices)
288 end subroutine adjoint_scalar_convection_source_term_compute
289
Implements the adjoint_scalar_convection_source_term type.