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