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