Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_minimum_dissipation_source_term.f90
Go to the documentation of this file.
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
50 use user_intf, only: user_t, simulation_component_user_settings
51 use json_module, only: json_file
53 use simcomp_executor, only: neko_simcomps
54 use fluid_user_source_term, only: fluid_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
59 use device_math, only: device_copy, device_cmult
60 use neko_config, only: neko_bcknd_device
61 use operators, only: curl
62 use scratch_registry, only: neko_scratch_registry
64 use point_zone, only: point_zone_t
65 implicit none
66 private
67
69 ! $\int \nabla v \cdot \nabla u $
70 type, public, extends(source_term_t) :: &
73 type(field_t), pointer :: u,v,w
75 real(kind=rp) :: obj_scale
77 class(point_zone_t), pointer :: mask
79 logical :: if_mask
80
81 contains
83 procedure, pass(this) :: init => &
86 procedure, pass(this) :: init_from_components => &
89 procedure, pass(this) :: free => &
92 procedure, pass(this) :: compute_ => &
95
96contains
102 json, fields, coef)
103 class(adjoint_minimum_dissipation_source_term_t), intent(inout) :: this
104 type(json_file), intent(inout) :: json
105 type(field_list_t), intent(in), target :: fields
106 type(coef_t), intent(in), target :: coef
107
108 ! we shouldn't be initializing this from JSON
109 ! maybe throw an error?
110
111
113
118 ! u,v,w reffer to the primal, not the adjoint
120 f_x, f_y, f_z, &
121 u, v, w, obj_scale, &
122 mask, if_mask, &
123 coef)
124 class(adjoint_minimum_dissipation_source_term_t), intent(inout) :: this
125 type(field_t), pointer, intent(in) :: f_x, f_y, f_z
126 type(field_list_t) :: fields
127 type(coef_t) :: coef
128 real(kind=rp) :: start_time
129 real(kind=rp) :: end_time
130 real(kind=rp) :: obj_scale
131 type(field_t), intent(in), target :: u, v, w
132 class(point_zone_t), intent(in), target :: mask
133 logical :: if_mask
134
135 ! I wish you didn't need a start time and end time...
136 ! but I'm just going to set a super big number...
137 start_time = 0.0_rp
138 end_time = 100000000.0_rp
139
140 call this%free()
141
142 ! this is copying the fluid source term init
143 ! We package the fields for the source term to operate on in a field list.
144 call fields%init(3)
145 call fields%assign(1, f_x)
146 call fields%assign(2, f_y)
147 call fields%assign(3, f_z)
148
149 call this%init_base(fields, coef, start_time, end_time)
150
151 ! point everything in the correct places
152 this%u => u
153 this%v => v
154 this%w => w
155
156 this%obj_scale = obj_scale
157
158 this%if_mask = if_mask
159 if (this%if_mask) then
160 this%mask => mask
161 end if
162
164
167 class(adjoint_minimum_dissipation_source_term_t), intent(inout) :: this
168
169 call this%free_base()
171
176 class(adjoint_minimum_dissipation_source_term_t), intent(inout) :: this
177 real(kind=rp), intent(in) :: t
178 integer, intent(in) :: tstep
179 type(field_t), pointer :: u, v, w
180 type(field_t), pointer :: fu, fv, fw
181 !type(field_t), pointer :: dudx, dudy, dudz
182 !type(field_t), pointer :: dvdx, dvdy, dvdz
183 !type(field_t), pointer :: dwdx, dwdy, dwdz
184 type(field_t), pointer :: wo1, wo2, wo3, wo4, wo5, wo6
185 type(field_t), pointer :: t1 , t2
186 integer :: temp_indices(8)
187 integer n
188
189
190 fu => this%fields%get_by_index(1)
191 fv => this%fields%get_by_index(2)
192 fw => this%fields%get_by_index(3)
193
194 n = fu%size()
195
196 ! fuck I'm not sure about this... I need a pen and paper
197 ! also there should be a way to pre-process this forcing term...
198 ! instead of recalculating it every time
199 u => this%u
200 v => this%v
201 w => this%w
202 call neko_scratch_registry%request_field(wo1, temp_indices(1))
203 call neko_scratch_registry%request_field(wo2, temp_indices(2))
204 call neko_scratch_registry%request_field(wo3, temp_indices(3))
205 call neko_scratch_registry%request_field(wo4, temp_indices(4))
206 call neko_scratch_registry%request_field(wo5, temp_indices(5))
207 call neko_scratch_registry%request_field(wo6, temp_indices(6))
208 call neko_scratch_registry%request_field(t1 , temp_indices(7))
209 call neko_scratch_registry%request_field(t2 , temp_indices(8))
210
211 ! TODO
212 ! ok we're computing gradients at every timestep... which is stupid...
213 ! BUT
214 ! if this was unsteady we would have to do this.
215
216 ! TODO
217 ! this is cheating a little bit...
218 ! in strong form, \nabla u . \nabla v => - v . \nabla^2 u + bdry
219 !
220 ! we can do this properly later in weak form, ideally using ax_helm or so
221 !
222 ! for now we'll work in strong form and ignore the bdry
223 ! and suffer the double derivative :/
224 !
225 ! in fact, we'll go even quicker and use
226 ! \nabla ^2 u = grad (div (u)) - curl ( curl (u ))
227 ! and assume divergence free u
228
229 ! => - \nabla ^2 u = curl ( curl (u ))
230
231 call curl(wo1, wo2, wo3, u, v, w, t1, t2, this%coef)
232 call curl(wo4, wo5, wo6, wo1, wo2, wo3, t1, t2, this%coef)
233
234 ! mask
235 if (this%if_mask) then
236 call mask_exterior_const(wo4, this%mask, 0.0_rp)
237 call mask_exterior_const(wo5, this%mask, 0.0_rp)
238 call mask_exterior_const(wo6, this%mask, 0.0_rp)
239 end if
240
241 call field_add2s2(fu, wo4, this%obj_scale)
242 call field_add2s2(fv, wo5, this%obj_scale)
243 call field_add2s2(fw, wo6, this%obj_scale)
244
245 ! don't worry... we'll write this MUCH cleaner in the final version
246
247 call neko_scratch_registry%relinquish_field(temp_indices)
248
250
Implements the adjoint_minimum_dissipation_source_term_t type.
subroutine adjoint_minimum_dissipation_source_term_init_from_components(this, f_x, f_y, f_z, u, v, w, obj_scale, mask, if_mask, coef)
The constructor from type components.
subroutine adjoint_minimum_dissipation_source_term_init_from_json(this, json, fields, coef)
The common constructor using a JSON object.
subroutine adjoint_minimum_dissipation_source_term_compute(this, t, tstep)
Computes the source term and adds the result to fields.
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...