Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
augmented_lagrangian_objective.f90
1! Copyright (c) 2024-25, 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!
35 use num_types, only: rp
36 use field, only: field_t
37 use field_math, only: field_col3, field_addcol3, field_cmult
38 use scratch_registry, only: scratch_registry_t, neko_scratch_registry
39 use objective, only: objective_t
40 use simulation_m, only: simulation_t
41 use neko_config, only: neko_bcknd_device
42 use math, only: copy, col2, invcol2
43 use device_math, only: device_copy, device_col2, device_invcol2
44 use design, only: design_t
45 use json_module, only: json_file
46 use json_utils, only: json_get_or_default
47 use interpolation, only: interpolator_t
48 use space, only: space_t, gl
49 use coefs, only: coef_t
50 implicit none
51 private
52
56 private
57
59 type(field_t), pointer :: u => null()
61 type(field_t), pointer :: v => null()
63 type(field_t), pointer :: w => null()
64
66 type(field_t), pointer :: adjoint_u => null()
68 type(field_t), pointer :: adjoint_v => null()
70 type(field_t), pointer :: adjoint_w => null()
71
72 ! ---- everything GLL ----
74 type(space_t), pointer :: xh_gll
76 type(coef_t), pointer :: c_xh_gll
77
78 ! ---- everything for GL ----
80 type(space_t), pointer :: xh_gl
82 type(coef_t), pointer :: c_xh_gl
84 type(interpolator_t), pointer :: gll_to_gl
86 logical :: dealias
88 type(scratch_registry_t), pointer :: scratch_gl
89
90 contains
92 procedure, public, pass(this) :: init_json_sim => &
93 augmented_lagrangian_init_json_sim
95 procedure, public, pass(this) :: init_from_attributes => &
96 augmented_lagrangian_init_attributes
98 procedure, public, pass(this) :: free => augmented_lagrangian_free
100 procedure, public, pass(this) :: update_value => &
101 augmented_lagrangian_update_value
103 procedure, public, pass(this) :: update_sensitivity => &
104 augmented_lagrangian_update_sensitivity
105
107
108contains
109
115 subroutine augmented_lagrangian_init_json_sim(this, json, design, simulation)
116 class(augmented_lagrangian_objective_t), intent(inout) :: this
117 type(json_file), intent(inout) :: json
118 class(design_t), intent(in) :: design
119 type(simulation_t), target, intent(inout) :: simulation
120
121 character(len=:), allocatable :: name
122 character(len=:), allocatable :: mask_name
123 real(kind=rp) :: weight
124 logical :: dealias
125
126 call json_get_or_default(json, "weight", weight, 1.0_rp)
127 call json_get_or_default(json, "mask_name", mask_name, "")
128 call json_get_or_default(json, "name", name, "Augmented Lagrangian")
129 call json_get_or_default(json, "dealias", dealias, .true.)
130
131 call this%init_from_attributes(design, simulation, weight, name, &
132 mask_name, dealias)
133 end subroutine augmented_lagrangian_init_json_sim
134
143 subroutine augmented_lagrangian_init_attributes(this, design, simulation, &
144 weight, name, mask_name, dealias)
145 class(augmented_lagrangian_objective_t), intent(inout) :: this
146 class(design_t), intent(in) :: design
147 type(simulation_t), target, intent(inout) :: simulation
148 real(kind=rp), intent(in) :: weight
149 character(len=*), intent(in) :: name
150 character(len=*), intent(in) :: mask_name
151 logical, intent(in) :: dealias
152
153 call this%init_base(name, design%size(), weight, mask_name)
154
155 ! Save the simulation and design
156 this%u => simulation%neko_case%fluid%u
157 this%v => simulation%neko_case%fluid%v
158 this%w => simulation%neko_case%fluid%w
159 this%adjoint_u => simulation%adjoint_case%fluid_adj%u_adj
160 this%adjoint_v => simulation%adjoint_case%fluid_adj%v_adj
161 this%adjoint_w => simulation%adjoint_case%fluid_adj%w_adj
162
163 this%dealias = dealias
164 ! GLL
165 this%c_Xh_GLL => simulation%neko_case%fluid%c_Xh
166 this%Xh_GLL => this%c_Xh_GLL%Xh
167
168 ! GL
169 this%c_Xh_GL => simulation%adjoint_case%fluid_adj%c_Xh_GL
170 this%Xh_GL => this%c_Xh_GL%Xh
171
172 ! GLL to GL
173 this%GLL_to_GL => simulation%adjoint_case%fluid_adj%GLL_to_GL
174
175 this%scratch_GL => simulation%adjoint_case%fluid_adj%scratch_GL
176
177 end subroutine augmented_lagrangian_init_attributes
178
180 subroutine augmented_lagrangian_free(this)
181 class(augmented_lagrangian_objective_t), intent(inout) :: this
182 call this%free_base()
183
184 if (associated(this%u)) nullify(this%u)
185 if (associated(this%v)) nullify(this%v)
186 if (associated(this%w)) nullify(this%w)
187
188 if (associated(this%adjoint_u)) nullify(this%adjoint_u)
189 if (associated(this%adjoint_v)) nullify(this%adjoint_v)
190 if (associated(this%adjoint_w)) nullify(this%adjoint_w)
191
192 nullify(this%scratch_GL)
193
194
195 end subroutine augmented_lagrangian_free
196
200 subroutine augmented_lagrangian_update_value(this, design)
201 class(augmented_lagrangian_objective_t), intent(inout) :: this
202 class(design_t), intent(in) :: design
203
204 end subroutine augmented_lagrangian_update_value
205
210 subroutine augmented_lagrangian_update_sensitivity(this, design)
211 class(augmented_lagrangian_objective_t), intent(inout) :: this
212 class(design_t), intent(in) :: design
213 type(field_t), pointer :: work
214 integer :: temp_indices(1)
215 integer :: n_GL, nel
216 type(field_t), pointer :: accumulate, fld_GL, adjoint_fld_GL
217 integer :: temp_indices_GL(3)
218
219 call neko_scratch_registry%request_field(work, temp_indices(1))
220
221 if (this%dealias) then
222
223 nel = this%c_Xh_GLL%msh%nelv
224 n_gl = nel * this%Xh_GL%lxyz
225 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1))
226 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2))
227 call this%scratch_GL%request_field(adjoint_fld_gl, temp_indices_gl(3))
228
229 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
230 call this%GLL_to_GL%map(adjoint_fld_gl%x, this%adjoint_u%x, nel, &
231 this%Xh_GL)
232 call field_col3(accumulate, fld_gl, adjoint_fld_gl)
233
234 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
235 call this%GLL_to_GL%map(adjoint_fld_gl%x, this%adjoint_v%x, nel, &
236 this%Xh_GL)
237 call field_addcol3(accumulate, fld_gl, adjoint_fld_gl)
238
239 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
240 call this%GLL_to_GL%map(adjoint_fld_gl%x, this%adjoint_w%x, nel, &
241 this%Xh_GL)
242 call field_addcol3(accumulate, fld_gl, adjoint_fld_gl)
243
244 ! Evaluate term on GL and preempt the GLL premultiplication
245 if (neko_bcknd_device .eq. 1) then
246 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
247 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
248 call device_invcol2(work%x_d, this%c_Xh_GLL%B_d, work%size())
249 else
250 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
251 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
252 call invcol2(work%x, this%c_Xh_GLL%B, work%size())
253 end if
254
255 call this%scratch_GL%relinquish_field(temp_indices_gl)
256 else
257 call field_col3(work, this%u, this%adjoint_u)
258 call field_addcol3(work, this%v, this%adjoint_v)
259 call field_addcol3(work, this%w, this%adjoint_w)
260 end if
261 ! but negative
262 call field_cmult(work, -1.0_rp)
263
264 if (neko_bcknd_device .eq. 1) then
265 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
266 else
267 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
268 end if
269
270 call neko_scratch_registry%relinquish_field(temp_indices)
271
272 end subroutine augmented_lagrangian_update_sensitivity
273
Implements the augmented_lagrangian_objective_t type.
Implements the design_t.
Definition design.f90:34
Implements the objective_t type.
Definition objective.f90:35
Implements the steady_problem_t type.
An objective function implementing our augmented lagrangian sensitivity contribution.
An abstract design type.
Definition design.f90:52
The abstract objective type.
Definition objective.f90:51