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
45 use device, only: device_memcpy, host_to_device, device_to_host
46 use, intrinsic :: iso_c_binding, only: c_ptr
47
48 implicit none
49 private
50
51 type, public :: mma_t
52 private
53 integer :: n, m, max_iter
54 real(kind=rp) :: a0, f0val, asyinit, asyincr, asydecr, epsimin, &
55 residumax, residunorm
56 type(vector_t) :: xold1, xold2, low, upp, alpha, beta, a, c, d, xmax, xmin
57 logical :: is_initialized = .false.
58 logical :: is_updated = .false.
59 character(len=:), allocatable :: bcknd, subsolver
60
61 ! Internal dummy variables for MMA
62 type(vector_t) :: p0j, q0j
63 type(matrix_t) :: pij, qij
64 type(vector_t) :: bi
65
66 !---nesessary for KKT check after updating df0dx, fval, dfdx --------
67 real(kind=rp) :: z, zeta
68 type(vector_t) :: y, lambda, s, mu
69 type(vector_t) :: xsi, eta
70 contains
72 generic, public :: init => init_from_json, init_from_components
73 procedure, public, pass(this) :: init_from_json => mma_init_from_json
74 procedure, public, pass(this) :: init_from_components => &
75 mma_init_from_components
76 procedure, public, pass(this) :: free => mma_free
77 procedure, public, pass(this) :: get_n => mma_get_n
78 procedure, public, pass(this) :: get_m => mma_get_m
79 procedure, public, pass(this) :: get_residumax => mma_get_residumax
80 procedure, public, pass(this) :: get_residunorm => mma_get_residunorm
81 procedure, public, pass(this) :: get_max_iter => mma_get_max_iter
82
83 generic, public :: update => update_vector, update_cpu, update_device
84 procedure, pass(this) :: update_vector => mma_update_vector
85 procedure, pass(this) :: update_cpu => mma_update_cpu
86 procedure, pass(this) :: update_device => mma_update_device
87
88 generic, public :: kkt => kkt_vector, kkt_cpu, kkt_device
89 procedure, pass(this) :: kkt_vector => mma_kkt_vector
90 procedure, pass(this) :: kkt_cpu => mma_kkt_cpu
91 procedure, pass(this) :: kkt_device => mma_kkt_device
92
93 end type mma_t
94
95 interface
96 ! ======================================================================= !
97 ! interface for cpu backend module subroutines
98
100 module subroutine mma_update_cpu(this, iter, x, df0dx, fval, dfdx)
101 class(mma_t), intent(inout) :: this
102 integer, intent(in) :: iter
103 real(kind=rp), dimension(this%n), intent(inout) :: x
104 real(kind=rp), dimension(this%n), intent(in) :: df0dx
105 real(kind=rp), dimension(this%m), intent(in) :: fval
106 real(kind=rp), dimension(this%m, this%n), intent(in) :: dfdx
107 end subroutine mma_update_cpu
108
110 module subroutine mma_kkt_cpu(this, x, df0dx, fval, dfdx)
111 class(mma_t), intent(inout) :: this
112 real(kind=rp), dimension(this%n), intent(in) :: x
113 real(kind=rp), dimension(this%n), intent(in) :: df0dx
114 real(kind=rp), dimension(this%m), intent(in) :: fval
115 real(kind=rp), dimension(this%m, this%n), intent(in) :: dfdx
116 end subroutine mma_kkt_cpu
117
118 ! ======================================================================= !
119 ! interface for device backend module subroutines
120
122 module subroutine mma_update_device(this, iter, x, df0dx, fval, dfdx)
123 class(mma_t), intent(inout) :: this
124 integer, intent(in) :: iter
125 type(c_ptr), intent(inout) :: x
126 type(c_ptr), intent(in) :: df0dx, fval, dfdx
127 end subroutine mma_update_device
128
130 module subroutine mma_kkt_device(this, x, df0dx, fval, dfdx)
131 class(mma_t), intent(inout) :: this
132 type(c_ptr), intent(in) :: x, df0dx, fval, dfdx
133 end subroutine mma_kkt_device
134
135 end interface
136
137contains
138
140 subroutine mma_init_from_json(this, x, n, m, json, scale, auto_scale)
141 ! ----------------------------------------------------- !
142 ! Initializing the mma object and all the parameters !
143 ! required for MMA method. (a_i, c_i, d_i, ...) !
144 ! x: the design varaibles(DV), n: number of DV, !
145 ! m: number of constraints !
146 ! !
147 ! Note that residumax & residunorm of the KKT conditions!
148 ! are initialized with 10^5. This is done to avoid !
149 ! unnecessary extera computation of KKT norms for the !
150 ! initial design. !
151 ! ----------------------------------------------------- !
152 class(mma_t), intent(inout) :: this
153 integer, intent(in) :: n, m
154 real(kind=rp), intent(in), dimension(n) :: x
155
156 type(json_file), intent(inout) :: json
157
158 ! Read the scaling info for fval and dfdx from json
159 real(kind=rp), intent(out) :: scale
160 logical, intent(out) :: auto_scale
161 ! -------------------------------------------------------------------!
162 ! Internal parameters for MMA !
163 ! Minimize f_0(x) + a_0*z + sum( c_i*y_i + 0.5*d_i*(y_i)^2 ) !
164 ! subject to f_i(x) - a_i*z - y_i <= 0, i = 1,...,m !
165 ! xmin_j <= x_j <= xmax_j, j = 1,...,n !
166 ! z >= 0, y_i >= 0, i = 1,...,m !
167 ! -------------------------------------------------------------------!
168 real(kind=rp), dimension(n) :: xmax, xmin
169 real(kind=rp), dimension(m) :: a, c, d
170 character(len=:), allocatable :: subsolver, bcknd, bcknd_default
171
172 ! For reading the values from json and then set the value for the arrays
173 real(kind=rp) :: a0 , xmax_const, xmin_const, a_const, c_const, d_const
174
175 integer :: max_iter, n_global, ierr
176 real(kind=rp) :: epsimin, asyinit, asyincr, asydecr
177
178 call mpi_allreduce(n, n_global, 1, mpi_integer, &
179 mpi_sum, mpi_comm_world, ierr)
180
181 ! Assign default values for the backend based on the NEKO_BCKND_DEVICE
182 if (neko_bcknd_cuda .eq. 1) then
183 bcknd_default = "device"
184 else
185 bcknd_default = "cpu"
186 end if
187
188 ! ------------------------------------------------------------------------ !
189 ! Assign defaults if nothing is parsed
190 ! based on the Cpp Code by Niels
191 call json_get_or_default(json, 'mma.epsimin', epsimin, &
192 1.0e-9_rp * sqrt(real(m + n_global, rp)))
193 call json_get_or_default(json, 'mma.max_iter', max_iter, 100)
194
195 ! Following parameters are set based on eq.3.8:--------
196 call json_get_or_default(json, 'mma.asyinit', asyinit, 0.5_rp)
197 call json_get_or_default(json, 'mma.asyincr', asyincr, 1.2_rp)
198 call json_get_or_default(json, 'mma.asydecr', asydecr, 0.7_rp)
199
200 call json_get_or_default(json, 'mma.backend', bcknd, bcknd_default)
201 call json_get_or_default(json, 'mma.subsolver', subsolver, 'dip')
202
203 call json_get_or_default(json, 'mma.xmin', xmin_const, 0.0_rp)
204 call json_get_or_default(json, 'mma.xmax', xmax_const, 1.0_rp)
205 call json_get_or_default(json, 'mma.a0', a0, 1.0_rp)
206 call json_get_or_default(json, 'mma.a', a_const, 0.0_rp)
207 call json_get_or_default(json, 'mma.c', c_const, 100.0_rp)
208 call json_get_or_default(json, 'mma.d', d_const, 0.0_rp)
209
210 call json_get_or_default(json, 'mma.scale', scale, 10.0_rp)
211 call json_get_or_default(json, 'mma.auto_scale', auto_scale, .false.)
212
213 call json_get_or_default(json, 'mma.epsimin', epsimin, 1.0e-9_rp)
214 ! Initialize the MMA object with the parsed parameters
215 a = a_const
216 c = c_const
217 d = d_const
218 xmin = xmin_const
219 xmax = xmax_const
220 ! initializing the mma concrete type (mma_cpu_t or mma_device_t)
221 if (pe_rank .eq. 0) then
222 print *, "Initializing MMA backend to >>> ", bcknd
223 end if
224
225 ! ------------------------------------------------------------------------ !
226 ! Initialize the MMA object with the parameters read from json
227 ! call this%init(x, n, m, a0, a, c, d, xmin, xmax, &
228 ! max_iter, epsimin, asyinit, asyincr, asydecr, bcknd)
229 call this%init(x, n, m, a0, a, c, d, xmin, xmax, &
230 max_iter, epsimin, asyinit, asyincr, asydecr, bcknd, subsolver)
231
232 end subroutine mma_init_from_json
233
235 subroutine mma_free(this)
236 class(mma_t), intent(inout) :: this
237 ! Deallocate the internal vectors
238 call this%xold1%free()
239 call this%xold2%free()
240 call this%alpha%free()
241 call this%beta%free()
242 call this%a%free()
243 call this%c%free()
244 call this%d%free()
245 call this%low%free()
246 call this%upp%free()
247 call this%xmax%free()
248 call this%xmin%free()
249 call this%p0j%free()
250 call this%q0j%free()
251 call this%bi%free()
252 call this%y%free()
253 call this%lambda%free()
254 call this%s%free()
255 call this%mu%free()
256 call this%xsi%free()
257 call this%eta%free()
258
259 ! Deallocate the internal dummy matrices
260 call this%pij%free()
261 call this%qij%free()
262
263 this%is_initialized = .false.
264 this%is_updated = .false.
265 end subroutine mma_free
266
268 subroutine mma_init_from_components(this, x, n, m, a0, a, c, d, xmin, xmax, &
269 max_iter, epsimin, asyinit, asyincr, asydecr, bcknd, subsolver)
270 ! ----------------------------------------------------- !
271 ! Initializing the mma object and all the parameters !
272 ! required for MMA method. (a_i, c_i, d_i, ...) !
273 ! x: the design varaibles(DV), n: number of DV, !
274 ! m: number of constraints !
275 ! !
276 ! Note that residumax & residunorm of the KKT conditions!
277 ! are initialized with 10^5. This is done to avoid !
278 ! unnecessary extera computation of KKT norms for the !
279 ! initial design. !
280 ! ----------------------------------------------------- !
281 class(mma_t), intent(inout) :: this
282 integer, intent(in) :: n, m
283 real(kind=rp), intent(in), dimension(n) :: x
284 ! -------------------------------------------------------------------!
285 ! Internal parameters for MMA !
286 ! Minimize f_0(x) + a_0*z + sum( c_i*y_i + 0.5*d_i*(y_i)^2 ) !
287 ! subject to f_i(x) - a_i*z - y_i <= 0, i = 1,...,m !
288 ! xmin_j <= x_j <= xmax_j, j = 1,...,n !
289 ! z >= 0, y_i >= 0, i = 1,...,m !
290 ! -------------------------------------------------------------------!
291 real(kind=rp), intent(in), dimension(n) :: xmax, xmin
292 real(kind=rp), intent(in), dimension(m) :: a, c, d
293 real(kind=rp), intent(in) :: a0
294 integer, intent(in), optional :: max_iter
295 real(kind=rp), intent(in), optional :: epsimin, asyinit, asyincr, asydecr
296 character(len=:), intent(in), allocatable :: bcknd, subsolver
297
298 call this%free()
299
300 this%n = n
301 this%m = m
302
303 call this%xold1%init(n)
304 call this%xold2%init(n)
305 this%xold1%x = x
306 this%xold2%x = x
307
308 call this%alpha%init(n)
309 call this%beta%init(n)
310
311 call this%a%init(m)
312 call this%c%init(m)
313 call this%d%init(m)
314 call this%low%init(n)
315 call this%upp%init(n)
316 call this%xmax%init(n)
317 call this%xmin%init(n)
318
319 !internal dummy variables for MMA
320 call this%p0j%init(n)
321 call this%q0j%init(n)
322 call this%pij%init(m, n)
323 call this%qij%init(m, n)
324 call this%bi%init(m)
325
326 !---nesessary for KKT check after updating df0dx, fval, dfdx --------
327 call this%y%init(m)
328 call this%lambda%init(m)
329 call this%s%init(m)
330 call this%mu%init(m)
331 call this%xsi%init(n)
332 call this%eta%init(n)
333
334 this%a0 = a0
335 this%a%x = a
336 this%c%x = c
337 this%d%x = d
338
339 !setting the bounds for the design variable based on the problem
340 this%xmax%x = xmax
341 this%xmin%x = xmin
342
343 this%low%x(:) = minval(x)
344 this%upp%x(:) = maxval(x)
345
346 !setting KKT norms to a large number for the initial design
347 this%residumax = huge(0.0_rp)
348 this%residunorm = huge(0.0_rp)
349
350 ! Select backend type
351 select case (bcknd)
352 case ("cpu")
353 if (pe_rank == 0) then
354 print *, "MMA initialized with CPU backend!"
355 end if
356 case ("device")
357 ! upload all init values to device pointers
358 call device_memcpy(this%xold1%x, this%xold1%x_d, this%n, &
359 host_to_device, sync = .false.)
360 call device_memcpy(this%xold1%x, this%xold2%x_d, this%n, &
361 host_to_device, sync = .false.)
362 call device_memcpy(this%a%x, this%a%x_d, this%m, host_to_device, &
363 sync = .false.)
364 call device_memcpy(this%c%x, this%c%x_d, this%m, host_to_device, &
365 sync = .false.)
366 call device_memcpy(this%d%x, this%d%x_d, this%m, host_to_device, &
367 sync = .false.)
368 call device_memcpy(this%xmax%x, this%xmax%x_d, this%n, host_to_device, &
369 sync = .false.)
370 call device_memcpy(this%xmin%x, this%xmin%x_d, this%n, host_to_device, &
371 sync = .false.)
372 if (pe_rank == 0) then
373 print *, "MMA initialized with CUDA backend!"
374 end if
375 case default
376 call neko_error('Unknown backend in mma_init_attributes')
377 end select
378
379 ! ------------------------------------------------------------------------ !
380 ! Assign defaults if nothing is parsed
381
382 ! based on the Cpp Code by Niels
383 if (.not. present(epsimin)) this%epsimin = 1.0e-9_rp * sqrt(real(m + n, rp))
384 if (.not. present(max_iter)) this%max_iter = 100
385
386 ! Following parameters are set based on eq.3.8:--------
387 if (.not. present(asyinit)) this%asyinit = 0.5_rp
388 if (.not. present(asyincr)) this%asyincr = 1.2_rp
389 if (.not. present(asydecr)) this%asydecr = 0.7_rp
390
391 ! Assign values from inputs when present
392 if (present(max_iter)) this%max_iter = max_iter
393 if (present(epsimin)) this%epsimin = epsimin
394 if (present(asyinit)) this%asyinit = asyinit
395 if (present(asyincr)) this%asyincr = asyincr
396 if (present(asydecr)) this%asydecr = asydecr
397 this%bcknd = bcknd
398 this%subsolver = subsolver
399 if (pe_rank == 0) then
400 if (this%subsolver .eq. "dip") then
401 print *, "Using dual solver for MMA subsolve."
402 elseif (this%subsolver .eq. "dpip") then
403 print *, "Using dual-primal solver for MMA subsolve."
404 else
405 call neko_error('Unknown subsolver for MMA, mma_init_from_components')
406 end if
407 end if
408
409 if (pe_rank .eq. 0) then
410 print *, "MMA is initialized with a0 = ", a0, ", a = ", a, ", c = ", c, &
411 ", d = ", d, "epsimin = ", this%epsimin
412 end if
413 !the object is correctly initialized
414 this%is_initialized = .true.
415 end subroutine mma_init_from_components
416
418 subroutine mma_update_vector(this, iter, x, df0dx, fval, dfdx)
419 class(mma_t), intent(inout) :: this
420 integer, intent(in) :: iter
421 type(vector_t), intent(inout) :: x
422 type(vector_t), intent(inout) :: df0dx, fval
423 type(matrix_t), intent(inout) :: dfdx
424
425 ! Select backend type
426 select case (this%bcknd)
427 case ("cpu")
428 if (neko_bcknd_device .eq. 1) then
429 call device_memcpy(x%x, x%x_d, this%n, device_to_host, &
430 sync = .false.)
431 call device_memcpy(df0dx%x, df0dx%x_d, this%n, device_to_host, &
432 sync = .false.)
433 call device_memcpy(fval%x, fval%x_d, this%m, device_to_host, &
434 sync = .false.)
435 call device_memcpy(dfdx%x, dfdx%x_d, this%m * this%n, device_to_host,&
436 sync = .true.)
437 end if
438
439 call mma_update_cpu(this, iter, x%x, df0dx%x, fval%x, dfdx%x)
440
441 if (neko_bcknd_device .eq. 1) then
442 call device_memcpy(x%x, x%x_d, this%n, host_to_device, sync = .true.)
443 end if
444
445 case ("device")
446 call mma_update_device(this, iter, x%x_d, df0dx%x_d, fval%x_d, dfdx%x_d)
447 call device_memcpy(x%x, x%x_d, this%n, device_to_host, sync = .true.)
448 end select
449
450 end subroutine mma_update_vector
451
453 subroutine mma_kkt_vector(this, x, df0dx, fval, dfdx)
454 class(mma_t), intent(inout) :: this
455 type(vector_t), intent(inout) :: x, df0dx, fval
456 type(matrix_t), intent(inout) :: dfdx
457
458 ! Select backend type
459 select case (this%bcknd )
460 case ("cpu")
461 if (neko_bcknd_device .eq. 1) then
462 call device_memcpy(x%x, x%x_d, this%n, device_to_host, &
463 sync = .false.)
464 call device_memcpy(df0dx%x, df0dx%x_d, this%n, device_to_host, &
465 sync = .false.)
466 call device_memcpy(fval%x, fval%x_d, this%m, device_to_host, &
467 sync = .false.)
468 call device_memcpy(dfdx%x, dfdx%x_d, this%m * this%n, device_to_host,&
469 sync = .true.)
470 end if
471
472 call mma_kkt_cpu(this, x%x, df0dx%x, fval%x, dfdx%x)
473 case ("device")
474 call mma_kkt_device(this, x%x_d, df0dx%x_d, fval%x_d, dfdx%x_d)
475 end select
476 end subroutine mma_kkt_vector
477
478 ! ========================================================================== !
479 ! Getters and setters
480
482 pure function mma_get_n(this) result(n)
483 class(mma_t), intent(in) :: this
484 integer :: n
485 n = this%n
486 end function mma_get_n
487
489 pure function mma_get_m(this) result(m)
490 class(mma_t), intent(in) :: this
491 integer :: m
492 m = this%m
493 end function mma_get_m
494
496 pure function mma_get_residumax(this) result(residumax)
497 class(mma_t), intent(in) :: this
498 real(kind=rp) :: residumax
499 residumax = this%residumax
500 end function mma_get_residumax
501
503 pure function mma_get_residunorm(this) result(residunorm)
504 class(mma_t), intent(in) :: this
505 real(kind=rp) :: residunorm
506 residunorm = this%residunorm
507 end function mma_get_residunorm
508
510 pure function mma_get_max_iter(this) result(max_iter_value)
511 class(mma_t), intent(in) :: this
512 integer :: max_iter_value
513 max_iter_value = this%max_iter
514 end function mma_get_max_iter
515
516end module mma
517