Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_minimum_dissipation_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!
35!
36! If the objective function $\int |\nabla u|^2$,
37! the corresponding adjoint forcing is $ \int \nabla v \cdot \nabla u $ in weak
38! form.
40 use num_types, only: rp
41 use field_list, only: field_list_t
42 use json_module, only: json_file
43 use json_utils, only: json_get, json_get_or_default
44 use source_term, only: source_term_t
45 use coefs, only: coef_t
46 use neko_config, only: neko_bcknd_device
47 use utils, only: neko_error
48 use field, only: field_t
49 use field_math, only: field_subcol3, field_add2, field_add2s2, field_rzero
50 use json_module, only: json_file
51 use time_state, only: time_state_t
53 use simcomp_executor, only: neko_simcomps
54 use user_source_term, only: user_source_term_t
55 use num_types, only: rp
56 use field, only: field_t
57 use field_registry, only: neko_field_registry
58 use math, only: rzero, copy, chsign, cfill, invcol2
59 use device_math, only: device_copy, device_cmult, device_cfill, device_invcol2
60 use neko_config, only: neko_bcknd_device
61 use scratch_registry, only: neko_scratch_registry
63 use point_zone, only: point_zone_t
64 use ax_product, only : ax_t, ax_helm_factory
65
66 implicit none
67 private
68
70 ! $\int \nabla v \cdot \nabla u $
71 type, public, extends(source_term_t) :: &
74 type(field_t), pointer :: u => null()
76 type(field_t), pointer :: v => null()
78 type(field_t), pointer :: w => null()
80 real(kind=rp) :: obj_scale
82 class(point_zone_t), pointer :: mask => null()
84 logical :: if_mask
86 class(ax_t), allocatable :: ax
88 real(kind=rp) :: volume
89
90 contains
92 procedure, pass(this) :: init => &
93 adjoint_minimum_dissipation_source_term_init_from_json
95 procedure, pass(this) :: init_from_components => &
96 adjoint_minimum_dissipation_source_term_init_from_components
98 procedure, pass(this) :: free => &
99 adjoint_minimum_dissipation_source_term_free
101 procedure, pass(this) :: compute_ => &
102 adjoint_minimum_dissipation_source_term_compute
104
105contains
112 subroutine adjoint_minimum_dissipation_source_term_init_from_json(this, &
113 json, fields, coef, variable_name)
114 class(adjoint_minimum_dissipation_source_term_t), intent(inout) :: this
115 type(json_file), intent(inout) :: json
116 type(field_list_t), intent(in), target :: fields
117 type(coef_t), intent(in), target :: coef
118 character(len=*), intent(in) :: variable_name
119
120 ! we shouldn't be initializing this from JSON
121 ! maybe throw an error?
122
123
124 end subroutine adjoint_minimum_dissipation_source_term_init_from_json
125
135 subroutine adjoint_minimum_dissipation_source_term_init_from_components(this,&
136 f_x, f_y, f_z, &
137 u, v, w, obj_scale, &
138 mask, if_mask, &
139 coef, volume)
140 class(adjoint_minimum_dissipation_source_term_t), intent(inout) :: this
141 type(field_t), pointer, intent(in) :: f_x, f_y, f_z
142 type(field_list_t) :: fields
143 type(coef_t) :: coef
144 real(kind=rp) :: start_time
145 real(kind=rp) :: end_time
146 real(kind=rp) :: obj_scale
147 type(field_t), intent(in), target :: u, v, w
148 class(point_zone_t), intent(in), target :: mask
149 real(kind=rp), intent(in) :: volume
150 logical :: if_mask
151
152 ! I wish you didn't need a start time and end time...
153 ! but I'm just going to set a super big number...
154 start_time = 0.0_rp
155 end_time = 100000000.0_rp
156
157 call this%free()
158
159 ! this is copying the fluid source term init
160 ! We package the fields for the source term to operate on in a field list.
161 call fields%init(3)
162 call fields%assign(1, f_x)
163 call fields%assign(2, f_y)
164 call fields%assign(3, f_z)
165
166 call this%init_base(fields, coef, start_time, end_time)
167
168 ! point everything in the correct places
169 this%u => u
170 this%v => v
171 this%w => w
172
173 this%obj_scale = obj_scale
174 this%volume = volume
175
176 this%if_mask = if_mask
177 if (this%if_mask) then
178 this%mask => mask
179 end if
180
181 ! Initialize the ax_helm object
182 call ax_helm_factory(this%Ax, full_formulation = .false.)
183
184 end subroutine adjoint_minimum_dissipation_source_term_init_from_components
185
187 subroutine adjoint_minimum_dissipation_source_term_free(this)
188 class(adjoint_minimum_dissipation_source_term_t), intent(inout) :: this
189
190 call this%free_base()
191 nullify(this%u)
192 nullify(this%v)
193 nullify(this%w)
194 nullify(this%mask)
195 if (allocated(this%Ax)) then
196 deallocate(this%Ax)
197 end if
198
199 end subroutine adjoint_minimum_dissipation_source_term_free
200
204 subroutine adjoint_minimum_dissipation_source_term_compute(this, time)
205 class(adjoint_minimum_dissipation_source_term_t), intent(inout) :: this
206 type(time_state_t), intent(in) :: time
207 type(field_t), pointer :: u, v, w
208 type(field_t), pointer :: fu, fv, fw
209 type(field_t), pointer :: work
210 integer :: temp_indices(1)
211 integer n
212
213 fu => this%fields%get_by_index(1)
214 fv => this%fields%get_by_index(2)
215 fw => this%fields%get_by_index(3)
216
217 n = fu%size()
218
219 u => this%u
220 v => this%v
221 w => this%w
222
223 call neko_scratch_registry%request_field(work, temp_indices(1))
224
225 associate(coef => this%coef)
226
227 ! Note that axhelm computes h1 * lap u + h2 * u
228 ! we set h1 = 1 and h2 = 0 to compute the weak laplacian.
229 if (neko_bcknd_device .eq. 1) then
230 call device_cfill(coef%h1_d, 1.0_rp, n)
231 call device_cfill(coef%h2_d, 0.0_rp, n)
232 else
233 call cfill(coef%h1, 1.0_rp, n)
234 call cfill(coef%h2, 0.0_rp, n)
235 end if
236 coef%ifh2 = .false.
237
238 ! ------------------------------------------------------------------------
239 ! u
240
241 call this%Ax%compute(work%x, u%x, coef, coef%msh, coef%xh)
242
243 ! pre-divide out the mass matrix to counteract it's multiplication
244 if (neko_bcknd_device .eq. 1) then
245 call device_invcol2(work%x_d, coef%B_d, work%size())
246 else
247 call invcol2(work%x, coef%B, work%size())
248 end if
249
250 ! mask
251 if (this%if_mask) then
252 call mask_exterior_const(work, this%mask, 0.0_rp)
253 end if
254
255 ! add to RHS
256 call field_add2s2(fu, work, this%obj_scale / this%volume)
257
258 ! ------------------------------------------------------------------------
259 ! v
260
261 call this%Ax%compute(work%x, v%x, coef, coef%msh, coef%xh)
262
263 ! pre-divide out the mass matrix to counteract it's multiplication
264 if (neko_bcknd_device .eq. 1) then
265 call device_invcol2(work%x_d, coef%B_d, work%size())
266 else
267 call invcol2(work%x, coef%B, work%size())
268 end if
269
270 ! mask
271 if (this%if_mask) then
272 call mask_exterior_const(work, this%mask, 0.0_rp)
273 end if
274
275 ! add to RHS
276 call field_add2s2(fv, work, this%obj_scale / this%volume)
277
278 ! ------------------------------------------------------------------------
279 ! w
280
281 call this%Ax%compute(work%x, w%x, coef, coef%msh, coef%xh)
282
283 ! pre-divide out the mass matrix to counteract it's multiplication
284 if (neko_bcknd_device .eq. 1) then
285 call device_invcol2(work%x_d, coef%B_d, work%size())
286 else
287 call invcol2(work%x, coef%B, work%size())
288 end if
289
290 ! mask
291 if (this%if_mask) then
292 call mask_exterior_const(work, this%mask, 0.0_rp)
293 end if
294
295 ! add to RHS
296 call field_add2s2(fw, work, this%obj_scale / this%volume)
297
298 call neko_scratch_registry%relinquish_field(temp_indices)
299
300 end associate
301
302 end subroutine adjoint_minimum_dissipation_source_term_compute
303
Implements the adjoint_minimum_dissipation_source_term_t type.
Some common Masking operations we may need.
Definition mask_ops.f90:34
Implements the steady_simcomp_t type.
The steady_simcomp_t type is a simulation component that terminates a simulation when the normed diff...