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