6 use num_types,
only: rp
7 use utils,
only: neko_error
8 use json_module,
only: json_file
9 use json_utils,
only: json_get, json_get_or_default
13 use field,
only: field_t
14 use field_registry,
only: neko_field_registry
15 use profiler,
only: profiler_start_region, profiler_end_region
17 use vector,
only: vector_t
18 use matrix,
only: matrix_t
21 use comm,
only: neko_comm, pe_rank
22 use mpi_f08,
only: mpi_integer, mpi_sum, mpi_allreduce
24 use neko_config,
only: neko_bcknd_device
26 use,
intrinsic :: iso_fortran_env, only: stderr => error_unit
28 use math,
only: copy, cmult
29 use device_math,
only: device_copy, device_cmult
30 use field_math,
only: field_rzero
31 use vector_math,
only: vector_cmult
34 use device,
only: device_memcpy, host_to_device, device_to_host
53 real(kind=rp) :: scale
57 logical :: enable_output
62 generic :: init => init_from_json, init_from_components
63 procedure, pass(this) :: init_from_json => mma_optimizer_init_from_json
64 procedure, pass(this) :: init_from_components => &
65 mma_optimizer_init_from_components
67 procedure, pass(this) :: run => mma_optimizer_run
68 procedure, pass(this) :: validate => mma_optimizer_validate
69 procedure, pass(this) :: free => mma_optimizer_free
76 subroutine mma_optimizer_init_from_json(this, parameters, problem, design, &
79 type(json_file),
intent(inout) :: parameters
83 logical :: enable_output
84 integer :: max_iterations
85 real(kind=rp) :: tolerance
87 character(len=1024) :: optimization_header
88 character(len=1024) :: problem_header
91 type(json_file) :: solver_parameters
92 logical :: unconstrained_problem
95 call this%logger%init(
'optimization_data.csv')
98 problem_header =
problem%get_log_header()
99 optimization_header =
'iter, ' // trim(problem_header) // &
100 ', KKTmax, KKTnorm2, scaling factor'
104 if (pe_rank .eq. 0)
then
105 print *,
"Initializing mma_optimizer with steady_state_problem_t."
108 call json_get(parameters,
"optimization.solver", &
113 unconstrained_problem = (
problem%get_n_constraints() == 0)
114 if (unconstrained_problem)
then
115 call this%mma%init(x,
design%size(), 1, &
116 solver_parameters, this%scale, this%auto_scale)
118 call this%mma%init(x,
design%size(),
problem%get_n_constraints(), &
119 solver_parameters, this%scale, this%auto_scale)
122 call json_get_or_default(parameters,
"optimization.solver.max_iterations", &
124 call json_get_or_default(parameters,
"optimization.solver.tolerance", &
125 tolerance, 1.0e-3_rp)
126 call json_get_or_default(parameters,
"optimization.solver.enable_output", &
127 enable_output, .true.)
130 max_iterations, tolerance, enable_output, simulation)
132 call this%logger%set_header(trim(optimization_header) // &
133 this%mma%get_backend_and_subsolver())
135 end subroutine mma_optimizer_init_from_json
138 subroutine mma_optimizer_init_from_components(this, problem, design, &
139 max_iterations, tolerance, enable_output, simulation)
143 integer,
intent(in) :: max_iterations
144 real(kind=rp),
intent(in) :: tolerance
146 logical,
intent(in) :: enable_output
149 this%enable_output = enable_output
151 call this%init_base(max_iterations, tolerance)
153 end subroutine mma_optimizer_init_from_components
156 subroutine mma_optimizer_run(this, problem, design, simulation)
160 type(
simulation_t),
optional,
intent(inout) :: simulation
164 integer :: iter, ierr, nglobal, n
165 real(kind=rp) :: scaling_factor
167 real(kind=rp) :: objective_value
168 type(vector_t) :: all_objectives
169 type(vector_t) :: constraint_value
170 type(vector_t) :: objective_sensitivities
171 type(matrix_t) :: constraint_sensitivities
173 type(vector_t) :: log_data
174 logical :: unconstrained_problem = .false.
176 type(json_file) :: parameters
179 call mpi_allreduce(n, nglobal, 1, mpi_integer, mpi_sum, neko_comm, ierr)
181 unconstrained_problem = (
problem%get_n_constraints() == 0)
182 if (unconstrained_problem)
then
184 call dummy_con%init(parameters,
design)
185 call problem%add_constraint(dummy_con)
190 call all_objectives%init(
problem%get_n_objectives())
191 call constraint_value%init(
problem%get_n_constraints())
192 call objective_sensitivities%init(n)
193 call constraint_sensitivities%init(
problem%get_n_constraints(), n)
196 scaling_factor = 1.0_rp
197 if (pe_rank .eq. 0)
then
198 print *,
"max_iterations for the optimization loop = ", &
201 call profiler_start_region(
"Optimizer iteration")
206 call problem%get_objective_value(objective_value)
207 call problem%get_constraint_values(constraint_value)
208 select type (des =>
design)
210 call des%get_sensitivity(objective_sensitivities)
212 call problem%get_objective_sensitivities(objective_sensitivities)
214 call problem%get_constraint_sensitivities(constraint_sensitivities)
216 call profiler_end_region(
"Optimizer iteration")
218 if (this%enable_output)
then
219 call profiler_start_region(
"Optimizer logging")
221 call problem%get_all_objective_values(all_objectives)
222 call mma_logger_assemble_data(log_data, 0, objective_value, &
223 all_objectives, constraint_value, 0.0_rp, 0.0_rp, scaling_factor, &
225 unconstrained_problem)
226 call this%logger%write(log_data)
228 if (
present(simulation))
then
229 call simulation%write(0)
232 call profiler_end_region(
"Optimizer logging")
236 do iter = 1, this%max_iterations
237 if (this%mma%get_residumax() .lt. this%tolerance)
exit
241 call profiler_start_region(
"Optimizer iteration")
244 if (this%auto_scale .eqv. .true.)
then
245 scaling_factor = abs(this%scale/constraint_value%x(1))
247 scaling_factor = abs(this%scale)
250 call vector_cmult(constraint_value, scaling_factor)
252 if (neko_bcknd_device .eq. 1)
then
253 call device_cmult(constraint_sensitivities%x_d, scaling_factor, &
254 constraint_sensitivities%size())
256 call cmult(constraint_sensitivities%x, scaling_factor, &
257 constraint_sensitivities%size())
261 call profiler_start_region(
"MMA update")
262 call this%mma%update(iter, x, objective_sensitivities, &
263 constraint_value, constraint_sensitivities)
264 call profiler_end_region(
"MMA update")
266 call design%update_design(x)
271 call problem%get_objective_value(objective_value)
272 call problem%get_constraint_values(constraint_value)
273 select type (des =>
design)
275 call des%get_sensitivity(objective_sensitivities)
277 call problem%get_objective_sensitivities(objective_sensitivities)
279 call problem%get_constraint_sensitivities(constraint_sensitivities)
281 call profiler_start_region(
"MMA KKT computation")
282 call this%mma%KKT(x, objective_sensitivities, &
283 constraint_value, constraint_sensitivities)
284 call profiler_end_region(
"MMA KKT computation")
286 call profiler_end_region(
"Optimizer iteration")
288 if (this%enable_output)
then
289 call profiler_start_region(
"Optimizer logging")
291 call problem%get_all_objective_values(all_objectives)
292 call mma_logger_assemble_data(log_data, iter, objective_value, &
293 all_objectives, constraint_value, this%mma%get_residumax(), &
294 this%mma%get_residunorm(), scaling_factor, &
296 unconstrained_problem)
297 call this%logger%write(log_data)
299 if (
present(simulation))
then
300 call simulation%write(iter)
303 call profiler_end_region(
"Optimizer logging")
311 if (pe_rank .eq. 0)
then
312 print *,
"MMA Optimization completed after", iter-1,
"iterations."
318 call all_objectives%free()
319 call constraint_value%free()
320 call objective_sensitivities%free()
321 call constraint_sensitivities%free()
323 end subroutine mma_optimizer_run
326 subroutine mma_optimizer_validate(this, problem, design)
331 type(vector_t) :: constraint_values
333 call constraint_values%init(
problem%get_n_constraints())
334 call problem%get_constraint_values(constraint_values)
335 if (neko_bcknd_device .eq. 1)
then
336 call device_memcpy(constraint_values%x, constraint_values%x_d, &
337 constraint_values%size(), host_to_device, .true.)
340 if (any(constraint_values%x .gt. 0.0_rp))
then
341 call neko_error(
"MMA optimizer validation failed: " // &
342 "Constraints are not satisfied.")
346 call constraint_values%free()
348 end subroutine mma_optimizer_validate
351 subroutine mma_optimizer_free(this)
356 end subroutine mma_optimizer_free
359 subroutine mma_logger_assemble_data(log_data, iter, objective_value, &
360 all_objectives, constraint_value, residumax, residunorm, &
361 scaling_factor, n, m, unconstrained_problem)
362 type(vector_t),
intent(out) :: log_data
363 integer,
intent(in) :: iter
364 real(kind=rp),
intent(in) ::objective_value
365 type(vector_t),
intent(in) :: all_objectives
366 type(vector_t),
intent(in) :: constraint_value
367 real(kind=rp),
intent(in) :: residumax, residunorm, scaling_factor
368 logical,
intent(in) :: unconstrained_problem
369 integer,
intent(in) :: n, m
370 integer :: i_tmp1, i_tmp2
376 if (unconstrained_problem)
then
377 call log_data%init(5 + n)
379 call log_data%init(5 + n + m)
383 log_data%x(1) = real(iter, kind=rp)
386 log_data%x(2) = objective_value
390 i_tmp2 = i_tmp1 + n - 1
391 log_data%x(i_tmp1 : i_tmp2) = all_objectives%x
394 if (.not. unconstrained_problem)
then
396 i_tmp2 = i_tmp1 + m - 1
397 log_data%x(i_tmp1 : i_tmp2) = constraint_value%x
402 log_data%x(i_tmp2 + 1) = residumax
403 log_data%x(i_tmp2 + 2) = residunorm
404 log_data%x(i_tmp2 + 3) = scaling_factor
406 end subroutine mma_logger_assemble_data
407end module mma_optimizer
Implements the constraint_t type.
Implements the dummy_constraint_t type.
Some common Masking operations we may need.
Contains extensions to the neko library required to run the topology optimization code.
subroutine, public reset(neko_case)
Reset the case data structure.
Module for handling the optimization problem.
Implements the steady_problem_t type.
A topology optimization design variable.
The abstract constraint type.
Abstract optimizer class.
The abstract problem type.