Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
mma_optimizer.f90
Go to the documentation of this file.
1
34
35module mma_optimizer
36 use optimizer, only: optimizer_t
37 use mma, only: mma_t
38 use problem, only: problem_t
39 use num_types, only: rp
40 use utils, only: neko_error
41 use json_utils, only: json_get, json_get_or_default
42 use simulation_m, only: simulation_t
43 use design, only: design_t
44 use brinkman_design, only: brinkman_design_t
45 use constraint, only: constraint_t
47
48 ! External modules
49 use json_module, only: json_file
50 use vector, only: vector_t
51 use matrix, only: matrix_t
52 use math, only: abscmp
53 use profiler, only: profiler_start_region, profiler_end_region
54 use logger, only: neko_log
55 use csv_file, only: csv_file_t
56 use vector_math, only: vector_cmult
57 use matrix_math, only: matrix_cmult
58 use device, only: device_memcpy, device_to_host
59 use scratch_registry, only: neko_scratch_registry
60 implicit none
61 private
62
63 public :: mma_optimizer_t
64
65 ! Concrete type for MMA optimizer
66 type, extends(optimizer_t) :: mma_optimizer_t
67
68 type(mma_t), private :: mma
69
76 real(kind=rp), private :: scale = 1.0_rp
77 real(kind=rp), private :: scaling_factor = 1.0_rp
78 logical, private :: auto_scale = .false.
79
80 ! Set to flags to remove logging for optimal performance
81 logical, private :: unconstrained_problem = .false.
82
84 logical, private :: enable_output = .true.
85 type(csv_file_t), private :: csv_log
86 contains
87
88 ! Override the deferred methods
89 generic :: init => init_from_json, init_from_components
90 procedure, pass(this) :: init_from_json => mma_optimizer_init_from_json
91 procedure, pass(this) :: init_from_components => &
92 mma_optimizer_init_from_components
93
94 procedure, pass(this) :: initialize => mma_optimizer_initialize
95 procedure, pass(this) :: step => mma_optimizer_step
96 procedure, pass(this) :: validate => mma_optimizer_validate
97 procedure, pass(this) :: write => mma_optimizer_write
98 procedure, pass(this) :: free => mma_optimizer_free
99
100 end type mma_optimizer_t
101
102contains
103
104 ! -------------------------------------------------------------------------- !
105 ! Allocator and deallocator methods for the MMA optimizer
106
108 subroutine mma_optimizer_init_from_json(this, parameters, problem, design, &
109 simulation)
110 class(mma_optimizer_t), intent(inout) :: this
111 type(json_file), intent(inout) :: parameters
112 class(problem_t), intent(inout) :: problem
113 class(design_t), intent(in) :: design
114 type(simulation_t), optional, intent(in) :: simulation
115
116 ! Variables for settings
117 type(json_file) :: solver_parameters
118 logical :: enable_output
119 integer :: max_iterations
120 real(kind=rp) :: tolerance
121
122 ! Read the solver properties from the JSON file
123 call json_get(parameters, 'optimization.solver', solver_parameters)
124 call json_get_or_default(solver_parameters, 'max_iterations', &
125 max_iterations, 100)
126 call json_get_or_default(solver_parameters, 'tolerance', &
127 tolerance, 1.0e-3_rp)
128 call json_get_or_default(solver_parameters, 'enable_output', &
129 enable_output, .true.)
130
131 call this%init_from_components(problem, design, max_iterations, tolerance, &
132 enable_output, solver_parameters, simulation)
133
134 end subroutine mma_optimizer_init_from_json
135
137 subroutine mma_optimizer_init_from_components(this, problem, design, &
138 max_iterations, tolerance, enable_output, solver_parameters, simulation)
139 class(mma_optimizer_t), intent(inout) :: this
140 class(problem_t), intent(inout) :: problem
141 class(design_t), intent(in) :: design
142 integer, intent(in) :: max_iterations
143 real(kind=rp), intent(in) :: tolerance
144 logical, intent(in) :: enable_output
145 type(json_file), intent(inout), optional :: solver_parameters
146 type(simulation_t), intent(in), optional :: simulation
147
148 ! Local variables
149 type(vector_t), pointer :: x
150 integer :: ind
151 character(len=1024) :: header
152 class(constraint_t), allocatable :: dummy_con
153
154 call neko_log%section('Optimizer Initialization')
155
156 ! Check if the problem is unconstrained
157 this%unconstrained_problem = problem%get_n_constraints() .eq. 0
158 if (this%unconstrained_problem) then
159 call neko_log%message('Unconstrained problem detected. ' // &
160 'Adding a dummy constraint to enable MMA optimization.')
161
162 allocate(dummy_constraint_t::dummy_con)
163 select type (con => dummy_con)
164 type is (dummy_constraint_t)
165 call con%init_from_attributes(design)
166 end select
167
168 call problem%add_constraint(dummy_con)
169 if (allocated(dummy_con)) deallocate(dummy_con)
170 end if
171
172 ! Initialize mma_t, handling the dummy_constraint added for unconstrained
173 ! problems in mma_optimizer_run()
174 call neko_scratch_registry%request(x, ind, design%size(), .false.)
175
176 call design%get_values(x)
177 call this%mma%init(x, design%size(), problem%get_n_constraints(), &
178 solver_parameters, this%scale, this%auto_scale)
179
180 call neko_scratch_registry%relinquish(ind)
181
182 !set the enable_output flag
183 this%enable_output = enable_output
184 this%scaling_factor = this%scale
185
186 ! Initialize the logger
187 if (this%enable_output) then
188 call this%csv_log%init('optimization_data.csv')
189 header = 'iter, ' // trim(problem%get_log_header()) // &
190 ', KKTmax, KKTnorm2, scaling factor, ' // &
191 trim(this%mma%get_backend_and_subsolver())
192
193 call this%csv_log%set_header(trim(header))
194 end if
195
196 call this%init_base(max_iterations, tolerance)
197
198 call neko_log%end_section()
199
200 end subroutine mma_optimizer_init_from_components
201
202 ! Free resources associated with the MMA optimizer
203 subroutine mma_optimizer_free(this)
204 class(mma_optimizer_t), intent(inout) :: this
205
206 ! Free MMA-specific data
207 call this%free_base()
208 call this%mma%free()
209 end subroutine mma_optimizer_free
210
211 ! -------------------------------------------------------------------------- !
212 ! Implementation of the deferred methods for the MMA optimizer
213
215 subroutine mma_optimizer_initialize(this, problem, design, simulation)
216 class(mma_optimizer_t), intent(inout) :: this
217 class(problem_t), intent(inout) :: problem
218 class(design_t), intent(inout) :: design
219 type(simulation_t), optional, intent(inout) :: simulation
220
221 type(vector_t), pointer :: x
222 type(vector_t), pointer :: constraint_value
223 type(vector_t), pointer :: objective_sensitivities
224 type(matrix_t), pointer :: constraint_sensitivities
225 integer :: n_design, n_constraint, indices(4)
226
227 n_design = design%size()
228 n_constraint = problem%get_n_constraints()
229
230 ! Grab some local pointers
231 call neko_scratch_registry%request(x, indices(1), n_design, .false.)
232 call neko_scratch_registry%request(constraint_value, indices(2), &
233 n_constraint, .false.)
234 call neko_scratch_registry%request(objective_sensitivities, indices(3), &
235 n_design, .false.)
236 call neko_scratch_registry%request(constraint_sensitivities, indices(4), &
237 n_constraint, n_design, .false.)
238
239 ! Evaluate the problem based on the updated design
240 call problem%compute(design, simulation)
241 call problem%compute_sensitivity(design, simulation)
242
243 ! Retrieve the updated objective and constraint values and sensitivities
244 call design%get_values(x)
245 call problem%get_constraint_values(constraint_value)
246
247 select type (des => design)
248 type is (brinkman_design_t)
249 call des%get_sensitivity(objective_sensitivities)
250 class default
251 call problem%get_objective_sensitivities(objective_sensitivities)
252 end select
253
254 call problem%get_constraint_sensitivities(constraint_sensitivities)
255
256 ! Check the KKT conditions and check for convergence
257 call this%mma%KKT(x, objective_sensitivities, &
258 constraint_value, constraint_sensitivities)
259
260 call neko_scratch_registry%relinquish(indices)
261 end subroutine mma_optimizer_initialize
262
264 function mma_optimizer_step(this, iter, problem, design, simulation) &
265 result(converged)
266 class(mma_optimizer_t), intent(inout) :: this
267 integer, intent(in) :: iter
268 class(problem_t), intent(inout) :: problem
269 class(design_t), intent(inout) :: design
270 type(simulation_t), optional, intent(inout) :: simulation
271
272 type(vector_t), pointer :: x
273 type(vector_t), pointer :: constraint_value
274 type(vector_t), pointer :: objective_sensitivities
275 type(matrix_t), pointer :: constraint_sensitivities
276 integer :: n_design, n_constraint, indices(4)
277
278 logical :: converged
279
280 n_design = design%size()
281 n_constraint = problem%get_n_constraints()
282
283 ! Grab some local pointers
284 call neko_scratch_registry%request(x, indices(1), n_design, .false.)
285 call neko_scratch_registry%request(constraint_value, indices(2), &
286 n_constraint, .false.)
287 call neko_scratch_registry%request(objective_sensitivities, indices(3), &
288 n_design, .false.)
289 call neko_scratch_registry%request(constraint_sensitivities, indices(4), &
290 n_constraint, n_design, .false.)
291
292 ! Retrieve the current objective and constraint values and sensitivities
293 call design%get_values(x)
294 call problem%get_constraint_values(constraint_value)
295
296 select type (des => design)
297 type is (brinkman_design_t)
298 call des%get_sensitivity(objective_sensitivities)
299 class default
300 call problem%get_objective_sensitivities(objective_sensitivities)
301 end select
302
303 call problem%get_constraint_sensitivities(constraint_sensitivities)
304
305 ! Execute the scaling
306 if (this%auto_scale) then
307 call constraint_value%copy_from(device_to_host, sync = .true.)
308 this%scaling_factor = abs(this%scale / constraint_value%x(1))
309 end if
310
311 if (.not. abscmp(this%scaling_factor, 1.0_rp)) then
312 call vector_cmult(constraint_value, this%scaling_factor)
313 call matrix_cmult(constraint_sensitivities, this%scaling_factor)
314 end if
315
316 ! Update the design variable
317 call this%mma%update(iter, x, objective_sensitivities, &
318 constraint_value, constraint_sensitivities)
319 call design%update_design(x)
320
321 ! Evaluate the problem based on the updated design
322 call problem%compute(design, simulation)
323 if (present(simulation) .and. this%enable_output) then
324 call simulation%write_forward(iter)
325 end if
326 call problem%compute_sensitivity(design, simulation)
327 if (present(simulation) .and. this%enable_output) then
328 call simulation%write_adjoint(iter)
329 end if
330
331 ! Retrieve the updated objective and constraint values and sensitivities
332 call problem%get_constraint_values(constraint_value)
333
334 select type (des => design)
335 type is (brinkman_design_t)
336 call des%get_sensitivity(objective_sensitivities)
337 class default
338 call problem%get_objective_sensitivities(objective_sensitivities)
339 end select
340
341 call problem%get_constraint_sensitivities(constraint_sensitivities)
342
343 ! Check the KKT conditions and check for convergence
344 call this%mma%KKT(x, objective_sensitivities, &
345 constraint_value, constraint_sensitivities)
346 converged = this%mma%get_residumax() .lt. this%tolerance
347
348 ! Free local resources
349 call neko_scratch_registry%relinquish(indices)
350
351 end function mma_optimizer_step
352
354 subroutine mma_optimizer_validate(this, problem, design)
355 class(mma_optimizer_t), intent(inout) :: this
356 class(problem_t), intent(in) :: problem
357 class(design_t), intent(in) :: design
358
359 type(vector_t), pointer :: constraint_values
360 integer :: ind
361
362 call neko_scratch_registry%request(constraint_values, ind, &
363 problem%get_n_constraints(), .false.)
364
365 call problem%get_constraint_values(constraint_values)
366 call constraint_values%copy_from(device_to_host, sync = .true.)
367
368 if (any(constraint_values%x .gt. 0.0_rp)) then
369 call neko_error('MMA optimizer validation failed: ' // &
370 'Constraints are not satisfied.')
371 end if
372
373 ! Free local resources
374 call neko_scratch_registry%relinquish(ind)
375
376 end subroutine mma_optimizer_validate
377
378 ! -------------------------------------------------------------------------- !
379 ! Logging and IO methods for the MMA optimizer
380
387 subroutine mma_optimizer_write(this, iter, problem)
388 class(mma_optimizer_t), intent(inout) :: this
389 integer, intent(in) :: iter
390 class(problem_t), intent(in) :: problem
391
392 type(vector_t), pointer :: log_data
393 type(vector_t), pointer :: all_objectives
394 type(vector_t), pointer :: constraint_value
395 real(kind=rp) :: objective_value
396
397 integer :: log_size, ind(3), n, m, i_tmp1, i_tmp2
398
399 if (.not. this%enable_output) return
400 call profiler_start_region('Optimizer logging')
401
402 n = problem%get_n_objectives()
403 m = problem%get_n_constraints()
404 if (this%unconstrained_problem) then
405 log_size = 5 + n
406 else
407 log_size = 5 + n + m
408 endif
409
410 call neko_scratch_registry%request(log_data, ind(1), log_size, .false.)
411 call neko_scratch_registry%request(all_objectives, ind(2), n, .false.)
412 call neko_scratch_registry%request(constraint_value, ind(3), m, .false.)
413
414 ! Prepare data for logging
415 call problem%get_objective_value(objective_value)
416 call problem%get_all_objective_values(all_objectives)
417 call problem%get_constraint_values(constraint_value)
418
419 call all_objectives%copy_from(device_to_host, sync = .true.)
420 call constraint_value%copy_from(device_to_host, sync = .true.)
421
422 ! Assemble the log data
423 log_data%x(1) = real(iter, kind=rp)
424
425 ! total objective
426 log_data%x(2) = objective_value
427
428 ! individual objectives
429 i_tmp1 = 3
430 i_tmp2 = i_tmp1 + n - 1
431 log_data%x(i_tmp1 : i_tmp2) = all_objectives%x
432
433 ! constraints
434 if (.not. this%unconstrained_problem) then
435 i_tmp1 = i_tmp2 + 1
436 i_tmp2 = i_tmp1 + m - 1
437 log_data%x(i_tmp1 : i_tmp2) = constraint_value%x
438 end if
439
440 ! convergence stuff
441 if (iter .eq. 0) then
442 log_data%x(i_tmp2 + 1) = 0.0_rp
443 log_data%x(i_tmp2 + 2) = 0.0_rp
444 else
445 log_data%x(i_tmp2 + 1) = this%mma%get_residumax()
446 log_data%x(i_tmp2 + 2) = this%mma%get_residunorm()
447 end if
448 log_data%x(i_tmp2 + 3) = this%scaling_factor
449
450 call this%csv_log%write(log_data)
451
452 ! Free local resources
453 call neko_scratch_registry%relinquish(ind)
454
455 call profiler_end_region('Optimizer logging')
456 end subroutine mma_optimizer_write
457
458end module mma_optimizer
459
Implements the constraint_t type.
Implements the design_t.
Definition design.f90:36
Implements the dummy_constraint_t type.
Module for handling the optimization problem.
Definition problem.f90:41
Implements the steady_problem_t type.
A topology optimization design variable.
The abstract constraint type.
An abstract design type.
Definition design.f90:54
Abstract optimizer class.
Definition optimizer.f90:54
The abstract problem type.
Definition problem.f90:67