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
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
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
61 type(field_t),
pointer :: u => null()
63 type(field_t),
pointer :: v => null()
65 type(field_t),
pointer :: w => null()
68 type(field_t),
pointer :: adjoint_u => null()
70 type(field_t),
pointer :: adjoint_v => null()
72 type(field_t),
pointer :: adjoint_w => null()
76 type(space_t),
pointer :: xh_gll
78 type(coef_t),
pointer :: c_xh_gll
82 type(space_t),
pointer :: xh_gl
84 type(coef_t),
pointer :: c_xh_gl
86 type(interpolator_t),
pointer :: gll_to_gl
90 type(scratch_registry_t),
pointer :: scratch_gl
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
108 procedure,
public, pass(this) :: get_log_size => &
109 augmented_lagrangian_get_log_size
111 procedure,
public, pass(this) :: get_log_headers => &
112 augmented_lagrangian_get_log_headers
114 procedure,
public, pass(this) :: get_log_values => &
115 augmented_lagrangian_get_log_values
128 type(json_file),
intent(inout) :: json
129 class(
design_t),
intent(in) :: design
132 character(len=:),
allocatable :: name
133 character(len=:),
allocatable :: mask_name
134 real(kind=rp) :: weight
137 call json_get_or_default(json,
"weight", weight, 1.0_rp)
138 call json_get_or_default(json,
"mask_name", mask_name,
"")
139 call json_get_or_default(json,
"name", name,
"Augmented Lagrangian")
140 call json_get_or_default(json,
"dealias", dealias, .true.)
142 call this%init_from_attributes(
design, simulation, weight, name, &
154 subroutine augmented_lagrangian_init_attributes(this, design, simulation, &
155 weight, name, mask_name, dealias)
157 class(
design_t),
intent(in) :: design
159 real(kind=rp),
intent(in) :: weight
160 character(len=*),
intent(in) :: name
161 character(len=*),
intent(in) :: mask_name
162 logical,
intent(in) :: dealias
164 call this%init_base(name,
design%size(), weight, mask_name)
167 this%u => simulation%neko_case%fluid%u
168 this%v => simulation%neko_case%fluid%v
169 this%w => simulation%neko_case%fluid%w
170 this%adjoint_u => simulation%adjoint_case%fluid_adj%u_adj
171 this%adjoint_v => simulation%adjoint_case%fluid_adj%v_adj
172 this%adjoint_w => simulation%adjoint_case%fluid_adj%w_adj
174 this%dealias = dealias
176 this%c_Xh_GLL => simulation%neko_case%fluid%c_Xh
177 this%Xh_GLL => this%c_Xh_GLL%Xh
180 this%c_Xh_GL => simulation%adjoint_case%fluid_adj%c_Xh_GL
181 this%Xh_GL => this%c_Xh_GL%Xh
184 this%GLL_to_GL => simulation%adjoint_case%fluid_adj%GLL_to_GL
186 this%scratch_GL => simulation%adjoint_case%fluid_adj%scratch_GL
188 end subroutine augmented_lagrangian_init_attributes
191 subroutine augmented_lagrangian_free(this)
193 call this%free_base()
195 if (
associated(this%u))
nullify(this%u)
196 if (
associated(this%v))
nullify(this%v)
197 if (
associated(this%w))
nullify(this%w)
199 if (
associated(this%adjoint_u))
nullify(this%adjoint_u)
200 if (
associated(this%adjoint_v))
nullify(this%adjoint_v)
201 if (
associated(this%adjoint_w))
nullify(this%adjoint_w)
203 nullify(this%scratch_GL)
206 end subroutine augmented_lagrangian_free
211 subroutine augmented_lagrangian_update_value(this, design)
213 class(
design_t),
intent(in) :: design
215 end subroutine augmented_lagrangian_update_value
221 subroutine augmented_lagrangian_update_sensitivity(this, design)
223 class(
design_t),
intent(in) :: design
224 type(field_t),
pointer :: work
225 integer :: temp_indices(1)
227 type(field_t),
pointer :: accumulate, fld_GL, adjoint_fld_GL
228 integer :: temp_indices_GL(3)
230 call neko_scratch_registry%request_field(work, temp_indices(1), .false.)
232 if (this%dealias)
then
234 nel = this%c_Xh_GLL%msh%nelv
235 n_gl = nel * this%Xh_GL%lxyz
236 call this%scratch_GL%request_field(accumulate, temp_indices_gl(1), &
238 call this%scratch_GL%request_field(fld_gl, temp_indices_gl(2), &
240 call this%scratch_GL%request_field(adjoint_fld_gl, temp_indices_gl(3), &
243 call this%GLL_to_GL%map(fld_gl%x, this%u%x, nel, this%Xh_GL)
244 call this%GLL_to_GL%map(adjoint_fld_gl%x, this%adjoint_u%x, nel, &
246 call field_col3(accumulate, fld_gl, adjoint_fld_gl)
248 call this%GLL_to_GL%map(fld_gl%x, this%v%x, nel, this%Xh_GL)
249 call this%GLL_to_GL%map(adjoint_fld_gl%x, this%adjoint_v%x, nel, &
251 call field_addcol3(accumulate, fld_gl, adjoint_fld_gl)
253 call this%GLL_to_GL%map(fld_gl%x, this%w%x, nel, this%Xh_GL)
254 call this%GLL_to_GL%map(adjoint_fld_gl%x, this%adjoint_w%x, nel, &
256 call field_addcol3(accumulate, fld_gl, adjoint_fld_gl)
259 if (neko_bcknd_device .eq. 1)
then
260 call device_col2(accumulate%x_d, this%c_Xh_GL%B_d, n_gl)
261 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
262 call device_invcol2(work%x_d, this%c_Xh_GLL%B_d, work%size())
264 call col2(accumulate%x, this%c_Xh_GL%B, n_gl)
265 call this%GLL_to_GL%map(work%x, accumulate%x, nel, this%Xh_GLL)
266 call invcol2(work%x, this%c_Xh_GLL%B, work%size())
269 call this%scratch_GL%relinquish_field(temp_indices_gl)
271 call field_col3(work, this%u, this%adjoint_u)
272 call field_addcol3(work, this%v, this%adjoint_v)
273 call field_addcol3(work, this%w, this%adjoint_w)
276 call field_cmult(work, -1.0_rp)
278 if (neko_bcknd_device .eq. 1)
then
279 call device_copy(this%sensitivity%x_d, work%x_d, this%sensitivity%size())
281 call copy(this%sensitivity%x, work%x, this%sensitivity%size())
284 call neko_scratch_registry%relinquish_field(temp_indices)
286 end subroutine augmented_lagrangian_update_sensitivity
291 function augmented_lagrangian_get_log_size(this)
result(n)
296 end function augmented_lagrangian_get_log_size
301 subroutine augmented_lagrangian_get_log_headers(this, headers)
303 character(len=*),
intent(out) :: headers(:)
305 if (
size(headers) .eq. 0)
return
307 end subroutine augmented_lagrangian_get_log_headers
312 subroutine augmented_lagrangian_get_log_values(this, values)
314 real(kind=rp),
intent(out) :: values(:)
316 if (
size(values) .eq. 0)
return
318 end subroutine augmented_lagrangian_get_log_values
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 objective_t type.
Implements the steady_problem_t type.
An objective function implementing our augmented lagrangian sensitivity contribution.
The abstract objective type.