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
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
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
59 type(field_t),
pointer :: u => null()
61 type(field_t),
pointer :: v => null()
63 type(field_t),
pointer :: w => null()
66 type(field_t),
pointer :: adjoint_u => null()
68 type(field_t),
pointer :: adjoint_v => null()
70 type(field_t),
pointer :: adjoint_w => null()
74 type(space_t),
pointer :: xh_gll
76 type(coef_t),
pointer :: c_xh_gll
80 type(space_t),
pointer :: xh_gl
82 type(coef_t),
pointer :: c_xh_gl
84 type(interpolator_t),
pointer :: gll_to_gl
88 type(scratch_registry_t),
pointer :: scratch_gl
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
115 subroutine augmented_lagrangian_init_json_sim(this, json, design, simulation)
117 type(json_file),
intent(inout) :: json
118 class(
design_t),
intent(in) :: design
121 character(len=:),
allocatable :: name
122 character(len=:),
allocatable :: mask_name
123 real(kind=rp) :: weight
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.)
131 call this%init_from_attributes(
design, simulation, weight, name, &
133 end subroutine augmented_lagrangian_init_json_sim
143 subroutine augmented_lagrangian_init_attributes(this, design, simulation, &
144 weight, name, mask_name, dealias)
146 class(
design_t),
intent(in) :: design
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
153 call this%init_base(name,
design%size(), weight, mask_name)
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
163 this%dealias = dealias
165 this%c_Xh_GLL => simulation%neko_case%fluid%c_Xh
166 this%Xh_GLL => this%c_Xh_GLL%Xh
169 this%c_Xh_GL => simulation%adjoint_case%fluid_adj%c_Xh_GL
170 this%Xh_GL => this%c_Xh_GL%Xh
173 this%GLL_to_GL => simulation%adjoint_case%fluid_adj%GLL_to_GL
175 this%scratch_GL => simulation%adjoint_case%fluid_adj%scratch_GL
177 end subroutine augmented_lagrangian_init_attributes
180 subroutine augmented_lagrangian_free(this)
182 call this%free_base()
184 if (
associated(this%u))
nullify(this%u)
185 if (
associated(this%v))
nullify(this%v)
186 if (
associated(this%w))
nullify(this%w)
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)
192 nullify(this%scratch_GL)
195 end subroutine augmented_lagrangian_free
200 subroutine augmented_lagrangian_update_value(this, design)
202 class(
design_t),
intent(in) :: design
204 end subroutine augmented_lagrangian_update_value
210 subroutine augmented_lagrangian_update_sensitivity(this, design)
212 class(
design_t),
intent(in) :: design
213 type(field_t),
pointer :: work
214 integer :: temp_indices(1)
216 type(field_t),
pointer :: accumulate, fld_GL, adjoint_fld_GL
217 integer :: temp_indices_GL(3)
219 call neko_scratch_registry%request_field(work, temp_indices(1))
221 if (this%dealias)
then
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))
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, &
232 call field_col3(accumulate, fld_gl, adjoint_fld_gl)
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, &
237 call field_addcol3(accumulate, fld_gl, adjoint_fld_gl)
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, &
242 call field_addcol3(accumulate, fld_gl, adjoint_fld_gl)
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())
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())
255 call this%scratch_GL%relinquish_field(temp_indices_gl)
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)
262 call field_cmult(work, -1.0_rp)
264 if (neko_bcknd_device .eq. 1)
then
265 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
267 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
270 call neko_scratch_registry%relinquish_field(temp_indices)
272 end subroutine augmented_lagrangian_update_sensitivity
Implements the augmented_lagrangian_objective_t type.
Implements the objective_t type.
Implements the steady_problem_t type.
An objective function implementing our augmented lagrangian sensitivity contribution.
The abstract objective type.