Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
mma.f90
Go to the documentation of this file.
1
34
35!===========================================================================!
36! Method of Moving Asymptotes !
37! This implementation is based on the following documents: !
38! 1-https://people.kth.se/~krille/mmagcmma.pdf !
39! 2-https://people.kth.se/~krille/originalmma.pdf !
40! 2-https://comsolyar.com/wp-content/uploads/2020/03/gcmma.pdf !
41! ------------------------------------------------------------------------- !
42! !
43! This module solves the following original optimization problem: !
44! !
45! Minimize f_0(x) + a_0*z + sum( c_i*y_i + 0.5*d_i*(y_i)^2 ) !
46! subject to f_i(x) - a_i*z - y_i <= 0, i = 1,...,m !
47! xmin_j <= x_j <= xmax_j, j = 1,...,n !
48! z >= 0, y_i >= 0, i = 1,...,m !
49! !
50! by first creating the following convex approximation of the original !
51! problem: !
52! !
53! Minimize sum_{j = 1,...,n} (p0j / (upp_j-x_j) + q0j / (x_j-low_j)) + !
54! a0*z + sum_i = 1,...,m(c_i*y_i + 0.5*d_i*y_i^2) !
55! subject to sum_{j = 1,...,n} (pij / (upp_j-x_j) + qij / (x_j-low_j)) + !
56! a_i*z + y_i <= b_i, i = 1,...,m !
57! xmin_j <= alpha_j <= x_j <= beta_j <= xmax_j j = 1,...,n !
58! y_i>=0 i = 1,...,m !
59! z>=0 !
60! !
61! note that based on eq(3.5) there should be r0 in the approximated problem !
62! however since it is just a constant added to a minimization problem, it !
63! is ignored. !
64! A primal-dual algorithm is then employed to solve the aproximated problem !
65! using interior point method. !
66!===========================================================================!
67
68module mma
69 ! Inclusions from Neko
70 use num_types, only: rp
71 use json_module, only: json_file
72 use json_utils, only: json_get_or_default
73 use vector, only: vector_t
74 use matrix, only: matrix_t
75 use mpi_f08, only: mpi_allreduce, mpi_integer, mpi_sum, mpi_comm_world
76 use comm, only: pe_rank
77 use utils, only: neko_error
78 use neko_config, only: neko_bcknd_device, neko_bcknd_cuda, neko_bcknd_hip, &
79 neko_bcknd_opencl
80 use device, only: device_memcpy, host_to_device, device_to_host
81 use, intrinsic :: iso_c_binding, only: c_ptr
82 use logger, only: neko_log
83
84 implicit none
85 private
86
87 type, public :: mma_t
88 private
89 integer :: n, m, max_iter
90 real(kind=rp) :: a0, f0val, asyinit, asyincr, asydecr, epsimin, &
91 residumax, residunorm
92 type(vector_t) :: xold1, xold2, low, upp, alpha, beta, a, c, d, xmax, xmin
93 logical :: is_initialized = .false.
94 logical :: is_updated = .false.
95 character(len=:), allocatable :: bcknd, subsolver
96
97 ! Internal dummy variables for MMA
98 type(vector_t) :: p0j, q0j
99 type(matrix_t) :: pij, qij
100 type(vector_t) :: bi
101
102 !---nesessary for KKT check after updating df0dx, fval, dfdx --------
103 real(kind=rp) :: z, zeta
104 type(vector_t) :: y, lambda, s, mu
105 type(vector_t) :: xsi, eta
106 contains
108 generic, public :: init => init_from_json, init_from_components
109 procedure, public, pass(this) :: init_from_json => mma_init_from_json
110 procedure, public, pass(this) :: init_from_components => &
111 mma_init_from_components
112 procedure, public, pass(this) :: free => mma_free
113 procedure, public, pass(this) :: get_n => mma_get_n
114 procedure, public, pass(this) :: get_m => mma_get_m
115 procedure, public, pass(this) :: get_residumax => mma_get_residumax
116 procedure, public, pass(this) :: get_residunorm => mma_get_residunorm
117 procedure, public, pass(this) :: get_max_iter => mma_get_max_iter
118 procedure, public, pass(this) :: get_backend_and_subsolver => &
119 mma_get_backend_and_subsolver
120
121
122 generic, public :: update => update_vector, update_cpu, update_device
123 procedure, pass(this) :: update_vector => mma_update_vector
124 procedure, pass(this) :: update_cpu => mma_update_cpu
125 procedure, pass(this) :: update_device => mma_update_device
126
127 generic, public :: kkt => kkt_vector, kkt_cpu, kkt_device
128 procedure, pass(this) :: kkt_vector => mma_kkt_vector
129 procedure, pass(this) :: kkt_cpu => mma_kkt_cpu
130 procedure, pass(this) :: kkt_device => mma_kkt_device
131
132 end type mma_t
133
134 interface
135 ! ======================================================================= !
136 ! interface for cpu backend module subroutines
137
139 module subroutine mma_update_cpu(this, iter, x, df0dx, fval, dfdx)
140 class(mma_t), intent(inout) :: this
141 integer, intent(in) :: iter
142 real(kind=rp), dimension(this%n), intent(inout) :: x
143 real(kind=rp), dimension(this%n), intent(in) :: df0dx
144 real(kind=rp), dimension(this%m), intent(in) :: fval
145 real(kind=rp), dimension(this%m, this%n), intent(in) :: dfdx
146 end subroutine mma_update_cpu
147
149 module subroutine mma_kkt_cpu(this, x, df0dx, fval, dfdx)
150 class(mma_t), intent(inout) :: this
151 real(kind=rp), dimension(this%n), intent(in) :: x
152 real(kind=rp), dimension(this%n), intent(in) :: df0dx
153 real(kind=rp), dimension(this%m), intent(in) :: fval
154 real(kind=rp), dimension(this%m, this%n), intent(in) :: dfdx
155 end subroutine mma_kkt_cpu
156
157 ! ======================================================================= !
158 ! interface for device backend module subroutines
159
161 module subroutine mma_update_device(this, iter, x, df0dx, fval, dfdx)
162 class(mma_t), intent(inout) :: this
163 integer, intent(in) :: iter
164 type(c_ptr), intent(inout) :: x
165 type(c_ptr), intent(in) :: df0dx, fval, dfdx
166 end subroutine mma_update_device
167
169 module subroutine mma_kkt_device(this, x, df0dx, fval, dfdx)
170 class(mma_t), intent(inout) :: this
171 type(c_ptr), intent(in) :: x, df0dx, fval, dfdx
172 end subroutine mma_kkt_device
173
174 end interface
175
176contains
177
179 subroutine mma_init_from_json(this, x, n, m, json, scale, auto_scale)
180 ! ----------------------------------------------------- !
181 ! Initializing the mma object and all the parameters !
182 ! required for MMA method. (a_i, c_i, d_i, ...) !
183 ! x: the design varaibles(DV), n: number of DV, !
184 ! m: number of constraints !
185 ! !
186 ! Note that residumax & residunorm of the KKT conditions!
187 ! are initialized with 10^5. This is done to avoid !
188 ! unnecessary extera computation of KKT norms for the !
189 ! initial design. !
190 ! ----------------------------------------------------- !
191 class(mma_t), intent(inout) :: this
192 integer, intent(in) :: n, m
193 type(vector_t), intent(in) :: x
194
195 type(json_file), intent(inout) :: json
196
197 ! Read the scaling info for fval and dfdx from json
198 real(kind=rp), intent(out) :: scale
199 logical, intent(out) :: auto_scale
200 ! -------------------------------------------------------------------!
201 ! Internal parameters for MMA !
202 ! Minimize f_0(x) + a_0*z + sum( c_i*y_i + 0.5*d_i*(y_i)^2 ) !
203 ! subject to f_i(x) - a_i*z - y_i <= 0, i = 1,...,m !
204 ! xmin_j <= x_j <= xmax_j, j = 1,...,n !
205 ! z >= 0, y_i >= 0, i = 1,...,m !
206 ! -------------------------------------------------------------------!
207 real(kind=rp), dimension(n) :: xmax, xmin
208 real(kind=rp), dimension(m) :: a, c, d
209 character(len=:), allocatable :: subsolver, bcknd, bcknd_default
210
211 ! For reading the values from json and then set the value for the arrays
212 real(kind=rp) :: a0 , xmax_const, xmin_const, a_const, c_const, d_const
213
214 integer :: max_iter, n_global, ierr
215 real(kind=rp) :: epsimin, asyinit, asyincr, asydecr
216
217 call mpi_allreduce(n, n_global, 1, mpi_integer, &
218 mpi_sum, mpi_comm_world, ierr)
219
220 ! Assign default values for the backend based on the NEKO_BCKND_DEVICE
221 if (neko_bcknd_device .eq. 1) then
222 bcknd_default = "device"
223 else
224 bcknd_default = "cpu"
225 end if
226
227 ! ------------------------------------------------------------------------ !
228 ! Assign defaults if nothing is parsed
229 ! based on the Cpp Code by Niels
230 call json_get_or_default(json, 'mma.epsimin', epsimin, &
231 1.0e-9_rp * sqrt(real(m + n_global, rp)))
232 call json_get_or_default(json, 'mma.max_iter', max_iter, 100)
233
234 ! Following parameters are set based on eq.3.8:--------
235 call json_get_or_default(json, 'mma.asyinit', asyinit, 0.5_rp)
236 call json_get_or_default(json, 'mma.asyincr', asyincr, 1.2_rp)
237 call json_get_or_default(json, 'mma.asydecr', asydecr, 0.7_rp)
238
239 call json_get_or_default(json, 'mma.backend', bcknd, bcknd_default)
240 call json_get_or_default(json, 'mma.subsolver', subsolver, 'dip')
241
242 call json_get_or_default(json, 'mma.xmin', xmin_const, 0.0_rp)
243 call json_get_or_default(json, 'mma.xmax', xmax_const, 1.0_rp)
244 call json_get_or_default(json, 'mma.a0', a0, 1.0_rp)
245 call json_get_or_default(json, 'mma.a', a_const, 0.0_rp)
246 call json_get_or_default(json, 'mma.c', c_const, 100.0_rp)
247 call json_get_or_default(json, 'mma.d', d_const, 0.0_rp)
248
249 call json_get_or_default(json, 'mma.scale', scale, 10.0_rp)
250 call json_get_or_default(json, 'mma.auto_scale', auto_scale, .false.)
251
252 ! Initialize the MMA object with the parsed parameters
253 a = a_const
254 c = c_const
255 d = d_const
256 xmin = xmin_const
257 xmax = xmax_const
258
259 ! ------------------------------------------------------------------------ !
260 ! Initialize the MMA object with the parameters read from json
261
262 call this%init(x, n, m, a0, a, c, d, xmin, xmax, &
263 max_iter, epsimin, asyinit, asyincr, asydecr, bcknd, subsolver)
264
265 end subroutine mma_init_from_json
266
268 subroutine mma_free(this)
269 class(mma_t), intent(inout) :: this
270 ! Deallocate the internal vectors
271 call this%xold1%free()
272 call this%xold2%free()
273 call this%alpha%free()
274 call this%beta%free()
275 call this%a%free()
276 call this%c%free()
277 call this%d%free()
278 call this%low%free()
279 call this%upp%free()
280 call this%xmax%free()
281 call this%xmin%free()
282 call this%p0j%free()
283 call this%q0j%free()
284 call this%bi%free()
285 call this%y%free()
286 call this%lambda%free()
287 call this%s%free()
288 call this%mu%free()
289 call this%xsi%free()
290 call this%eta%free()
291
292 ! Deallocate the internal dummy matrices
293 call this%pij%free()
294 call this%qij%free()
295
296 this%is_initialized = .false.
297 this%is_updated = .false.
298 end subroutine mma_free
299
301 subroutine mma_init_from_components(this, x, n, m, a0, a, c, d, xmin, xmax, &
302 max_iter, epsimin, asyinit, asyincr, asydecr, bcknd, subsolver)
303 ! ----------------------------------------------------- !
304 ! Initializing the mma object and all the parameters !
305 ! required for MMA method. (a_i, c_i, d_i, ...) !
306 ! x: the design varaibles(DV), n: number of DV, !
307 ! m: number of constraints !
308 ! !
309 ! Note that residumax & residunorm of the KKT conditions!
310 ! are initialized with 10^5. This is done to avoid !
311 ! unnecessary extera computation of KKT norms for the !
312 ! initial design. !
313 ! ----------------------------------------------------- !
314 class(mma_t), intent(inout) :: this
315 integer, intent(in) :: n, m
316 type(vector_t), intent(in) :: x
317 ! -------------------------------------------------------------------!
318 ! Internal parameters for MMA !
319 ! Minimize f_0(x) + a_0*z + sum( c_i*y_i + 0.5*d_i*(y_i)^2 ) !
320 ! subject to f_i(x) - a_i*z - y_i <= 0, i = 1,...,m !
321 ! xmin_j <= x_j <= xmax_j, j = 1,...,n !
322 ! z >= 0, y_i >= 0, i = 1,...,m !
323 ! -------------------------------------------------------------------!
324 real(kind=rp), intent(in), dimension(n) :: xmax, xmin
325 real(kind=rp), intent(in), dimension(m) :: a, c, d
326 real(kind=rp), intent(in) :: a0
327 integer, intent(in), optional :: max_iter
328 real(kind=rp), intent(in), optional :: epsimin, asyinit, asyincr, asydecr
329 character(len=:), intent(in), allocatable :: bcknd, subsolver
330 character(len=256) :: log_msg
331 integer :: i
332
333 call this%free()
334
335 this%n = n
336 this%m = m
337
338 call this%xold1%init(n)
339 call this%xold2%init(n)
340 this%xold1 = x
341 this%xold2 = x
342
343 call this%alpha%init(n)
344 call this%beta%init(n)
345
346 call this%a%init(m)
347 call this%c%init(m)
348 call this%d%init(m)
349 call this%low%init(n)
350 call this%upp%init(n)
351 call this%xmax%init(n)
352 call this%xmin%init(n)
353
354 !internal dummy variables for MMA
355 call this%p0j%init(n)
356 call this%q0j%init(n)
357 call this%pij%init(m, n)
358 call this%qij%init(m, n)
359 call this%bi%init(m)
360
361 ! Necessary for KKT check after updating df0dx, fval, dfdx
362 call this%y%init(m)
363 call this%lambda%init(m)
364 call this%s%init(m)
365 call this%mu%init(m)
366 call this%xsi%init(n)
367 call this%eta%init(n)
368
369 this%a0 = a0
370 this%a%x = a
371 this%c%x = c
372 this%d%x = d
373
374 ! Set the bounds for the design variable based on the problem
375 this%xmax%x = xmax
376 this%xmin%x = xmin
377
378 if (neko_bcknd_device .eq. 1) then
379 call device_memcpy(this%a%x, this%a%x_d, m, host_to_device, &
380 sync = .false.)
381 call device_memcpy(this%c%x, this%c%x_d, m, host_to_device, &
382 sync = .false.)
383 call device_memcpy(this%d%x, this%d%x_d, m, host_to_device, &
384 sync = .false.)
385 call device_memcpy(this%xmax%x, this%xmax%x_d, n, host_to_device, &
386 sync = .false.)
387 call device_memcpy(this%xmin%x, this%xmin%x_d, n, host_to_device, &
388 sync = .true.)
389 end if
390
391
392 ! Set KKT norms to a large number for the initial design
393 this%residumax = huge(0.0_rp)
394 this%residunorm = huge(0.0_rp)
395
396 ! ------------------------------------------------------------------------ !
397 ! Assign defaults if nothing is parsed
398
399 ! Based on the Cpp Code by Niels
400 if (.not. present(epsimin)) this%epsimin = 1.0e-9_rp * sqrt(real(m + n, rp))
401 if (.not. present(max_iter)) this%max_iter = 100
402
403 ! Following parameters are set based on eq.3.8
404 if (.not. present(asyinit)) this%asyinit = 0.5_rp
405 if (.not. present(asyincr)) this%asyincr = 1.2_rp
406 if (.not. present(asydecr)) this%asydecr = 0.7_rp
407
408 ! Assign values from inputs when present
409 if (present(max_iter)) this%max_iter = max_iter
410 if (present(epsimin)) this%epsimin = epsimin
411 if (present(asyinit)) this%asyinit = asyinit
412 if (present(asyincr)) this%asyincr = asyincr
413 if (present(asydecr)) this%asydecr = asydecr
414 this%bcknd = bcknd
415 this%subsolver = subsolver
416
417 call neko_log%section('MMA Parameters')
418
419 write(log_msg, '(A10,1X,A)') 'backend ', trim(this%bcknd)
420 call neko_log%message(log_msg)
421 write(log_msg, '(A10,1X,A)') 'subsolver ', trim(this%subsolver)
422 call neko_log%message(log_msg)
423
424 write(log_msg, '(A10,1X,I0)') 'n ', this%n
425 call neko_log%message(log_msg)
426 write(log_msg, '(A10,1X,I0)') 'm ', this%m
427 call neko_log%message(log_msg)
428 write(log_msg, '(A10,1X,I0)') 'max_iter ', this%max_iter
429 call neko_log%message(log_msg)
430
431 write(log_msg, '(A10,1X,E11.5)') 'epsimin ', this%epsimin
432 call neko_log%message(log_msg)
433
434 write(log_msg, '(A10,1X,E11.5)') 'asyinit ', this%asyinit
435 call neko_log%message(log_msg)
436 write(log_msg, '(A10,1X,E11.5)') 'asyincr ', this%asyincr
437 call neko_log%message(log_msg)
438 write(log_msg, '(A10,1X,E11.5)') 'asydecr ', this%asydecr
439 call neko_log%message(log_msg)
440 write(log_msg, '(A10,1X,E11.5)') 'a0 ', this%a0
441 call neko_log%message(log_msg)
442
443 call neko_log%message('Parameters a')
444 do i = 1, this%m
445 write(log_msg, '(3X,A,I2,A,E11.5)') 'a(', i, ') = ', this%a%x(i)
446 call neko_log%message(log_msg)
447 end do
448 call neko_log%message('Parameters c')
449 do i = 1, this%m
450 write(log_msg, '(3X,A,I2,A,E11.5)') 'c(', i, ') = ', this%c%x(i)
451 call neko_log%message(log_msg)
452 end do
453 call neko_log%message('Parameters d')
454 do i = 1, this%m
455 write(log_msg, '(3X,A,I2,A,E11.5)') 'd(', i, ') = ', this%d%x(i)
456 call neko_log%message(log_msg)
457 end do
458
459 call neko_log%end_section()
460
461 ! The object is correctly initialized
462 this%is_initialized = .true.
463 end subroutine mma_init_from_components
464
466 subroutine mma_update_vector(this, iter, x, df0dx, fval, dfdx)
467 class(mma_t), intent(inout) :: this
468 integer, intent(in) :: iter
469 type(vector_t), intent(inout) :: x
470 type(vector_t), intent(inout) :: df0dx, fval
471 type(matrix_t), intent(inout) :: dfdx
472
473 ! Select backend type
474 select case (this%bcknd)
475 case ("cpu")
476 if (neko_bcknd_device .eq. 1) then
477 call device_memcpy(x%x, x%x_d, this%n, device_to_host, &
478 sync = .false.)
479 call device_memcpy(df0dx%x, df0dx%x_d, this%n, device_to_host, &
480 sync = .false.)
481 call device_memcpy(fval%x, fval%x_d, this%m, device_to_host, &
482 sync = .false.)
483 call device_memcpy(dfdx%x, dfdx%x_d, this%m * this%n, device_to_host,&
484 sync = .true.)
485 end if
486
487 call mma_update_cpu(this, iter, x%x, df0dx%x, fval%x, dfdx%x)
488
489 if (neko_bcknd_device .eq. 1) then
490 call device_memcpy(x%x, x%x_d, this%n, host_to_device, sync = .true.)
491 end if
492
493 case ("device")
494 call mma_update_device(this, iter, x%x_d, df0dx%x_d, fval%x_d, dfdx%x_d)
495 end select
496
497 end subroutine mma_update_vector
498
500 subroutine mma_kkt_vector(this, x, df0dx, fval, dfdx)
501 class(mma_t), intent(inout) :: this
502 type(vector_t), intent(inout) :: x, df0dx, fval
503 type(matrix_t), intent(inout) :: dfdx
504
505 ! Select backend type
506 select case (this%bcknd )
507 case ("cpu")
508 if (neko_bcknd_device .eq. 1) then
509 call device_memcpy(x%x, x%x_d, this%n, device_to_host, &
510 sync = .false.)
511 call device_memcpy(df0dx%x, df0dx%x_d, this%n, device_to_host, &
512 sync = .false.)
513 call device_memcpy(fval%x, fval%x_d, this%m, device_to_host, &
514 sync = .false.)
515 call device_memcpy(dfdx%x, dfdx%x_d, this%m * this%n, device_to_host,&
516 sync = .true.)
517 end if
518
519 call mma_kkt_cpu(this, x%x, df0dx%x, fval%x, dfdx%x)
520 case ("device")
521 call mma_kkt_device(this, x%x_d, df0dx%x_d, fval%x_d, dfdx%x_d)
522 end select
523 end subroutine mma_kkt_vector
524
525 ! ========================================================================== !
526 ! Getters and setters
527
529 pure function mma_get_n(this) result(n)
530 class(mma_t), intent(in) :: this
531 integer :: n
532 n = this%n
533 end function mma_get_n
534
536 pure function mma_get_m(this) result(m)
537 class(mma_t), intent(in) :: this
538 integer :: m
539 m = this%m
540 end function mma_get_m
541
543 pure function mma_get_residumax(this) result(residumax)
544 class(mma_t), intent(in) :: this
545 real(kind=rp) :: residumax
546 residumax = this%residumax
547 end function mma_get_residumax
548
550 pure function mma_get_residunorm(this) result(residunorm)
551 class(mma_t), intent(in) :: this
552 real(kind=rp) :: residunorm
553 residunorm = this%residunorm
554 end function mma_get_residunorm
555
557 pure function mma_get_max_iter(this) result(max_iter_value)
558 class(mma_t), intent(in) :: this
559 integer :: max_iter_value
560 max_iter_value = this%max_iter
561 end function mma_get_max_iter
562
564 pure function mma_get_backend_and_subsolver(this) result(backend_subsolver)
565 class(mma_t), intent(in) :: this
566 character(len=:), allocatable :: backend_subsolver
567 character(len=:), allocatable :: backend
568
569 if (neko_bcknd_cuda .eq. 1) then
570 backend = "cuda"
571 else if (neko_bcknd_hip .eq. 1) then
572 backend = "hip"
573 else if (neko_bcknd_opencl .eq. 1) then
574 backend = "opencl"
575 else
576 backend = "cpu"
577 end if
578
579 backend_subsolver = 'backend:' // trim(backend) // ', subsolver:' // &
580 trim(this%subsolver)
581 end function mma_get_backend_and_subsolver
582end module mma