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 approximated !
65! problem using interior point method. !
66!===========================================================================!
67
69module mma
70 ! Inclusions from Neko
71 use num_types, only: rp, dp, sp
72 use json_module, only: json_file
73 use json_utils, only: json_get_or_default
74 use vector, only: vector_t
75 use matrix, only: matrix_t
76 use comm, only: pe_rank, neko_comm, pe_size, mpi_real_precision
77 use utils, only: neko_error, filename_suffix
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 use mpi_f08, only: mpi_sum, mpi_allreduce, mpi_integer
84
85 implicit none
86 private
87
89 type, public :: mma_t
90 private
91 integer :: n, m, n_global, max_iter
92 real(kind=rp) :: a0, asyinit, asyincr, asydecr, epsimin, &
93 residumax, residunorm
94 type(vector_t) :: xold1, xold2, low, upp, alpha, beta, a, c, d, xmax, xmin
95 logical :: is_initialized = .false.
96 logical :: is_updated = .false.
97 character(len=:), allocatable :: subsolver, bcknd
98
99 ! Internal dummy variables for MMA
100 type(vector_t) :: p0j, q0j
101 type(matrix_t) :: pij, qij
102 type(vector_t) :: bi
103
104 !---nesessary for KKT check after updating df0dx, fval, dfdx --------
105 real(kind=rp) :: z, zeta
106 type(vector_t) :: y, lambda, s, mu
107 type(vector_t) :: xsi, eta
108 contains
110 generic, public :: init => init_from_json, init_from_components
111 procedure, public, pass(this) :: init_from_json => mma_init_from_json
112 procedure, public, pass(this) :: init_from_components => &
113 mma_init_from_components
114 procedure, public, pass(this) :: free => mma_free
115 procedure, public, pass(this) :: get_n => mma_get_n
116 procedure, public, pass(this) :: get_m => mma_get_m
117 procedure, public, pass(this) :: get_residumax => mma_get_residumax
118 procedure, public, pass(this) :: get_residunorm => mma_get_residunorm
119 procedure, public, pass(this) :: get_max_iter => mma_get_max_iter
120 procedure, public, pass(this) :: get_backend_and_subsolver => &
121 mma_get_backend_and_subsolver
122
123 generic, public :: update => update_vector, update_cpu, update_device
124 procedure, pass(this) :: update_vector => mma_update_vector
125 procedure, pass(this) :: update_cpu => mma_update_cpu
126 procedure, pass(this) :: update_device => mma_update_device
127
128 generic, public :: kkt => kkt_vector, kkt_cpu, kkt_device
129 procedure, pass(this) :: kkt_vector => mma_kkt_vector
130 procedure, pass(this) :: kkt_cpu => mma_kkt_cpu
131 procedure, pass(this) :: kkt_device => mma_kkt_device
132
133 procedure, pass(this) :: save_checkpoint => mma_save_checkpoint
134 procedure, pass(this) :: load_checkpoint => mma_load_checkpoint
135
136 ! Private utilities
137 procedure, pass(this) :: copy_from => mma_copy_from
138
139 end type mma_t
140
141 ! ========================================================================== !
142 ! interface for cpu backend module subroutines
143
144 interface
145 ! CPU update function, runs one iteration of MMA
146 module subroutine mma_update_cpu(this, iter, x, df0dx, fval, dfdx)
147 class(mma_t), intent(inout) :: this
148 integer, intent(in) :: iter
149 real(kind=rp), dimension(this%n), intent(inout) :: x
150 real(kind=rp), dimension(this%n), intent(in) :: df0dx
151 real(kind=rp), dimension(this%m), intent(in) :: fval
152 real(kind=rp), dimension(this%m, this%n), intent(in) :: dfdx
153 end subroutine mma_update_cpu
154
155 ! CPU KKT check for convergence
156 module subroutine mma_kkt_cpu(this, x, df0dx, fval, dfdx)
157 class(mma_t), intent(inout) :: this
158 real(kind=rp), dimension(this%n), intent(in) :: x
159 real(kind=rp), dimension(this%n), intent(in) :: df0dx
160 real(kind=rp), dimension(this%m), intent(in) :: fval
161 real(kind=rp), dimension(this%m, this%n), intent(in) :: dfdx
162 end subroutine mma_kkt_cpu
163
164 ! ========================================================================== !
165 ! interface for device backend module subroutines
166
167 ! Device update function, runs one iteration of MMA
168 module subroutine mma_update_device(this, iter, x, df0dx, fval, dfdx)
169 class(mma_t), intent(inout) :: this
170 integer, intent(in) :: iter
171 type(c_ptr), intent(inout) :: x
172 type(c_ptr), intent(in) :: df0dx, fval, dfdx
173 end subroutine mma_update_device
174
176 module subroutine mma_kkt_device(this, x, df0dx, fval, dfdx)
177 class(mma_t), intent(inout) :: this
178 type(c_ptr), intent(in) :: x, df0dx, fval, dfdx
179 end subroutine mma_kkt_device
180
181 end interface
182
183 ! ========================================================================== !
184 ! Interface for IO routines
185
186 interface
187 module subroutine mma_save_checkpoint_hdf5(object, filename, overwrite)
188 class(mma_t), intent(inout) :: object
189 character(len=*), intent(in) :: filename
190 logical, intent(in), optional :: overwrite
191 end subroutine mma_save_checkpoint_hdf5
192
193 module subroutine mma_load_checkpoint_hdf5(object, filename)
194 class(mma_t), intent(inout) :: object
195 character(len=*), intent(in) :: filename
196 end subroutine mma_load_checkpoint_hdf5
197 end interface
198
199contains
200
201 ! ========================================================================== !
202 ! Initializers and destructors
203
205 subroutine mma_init_from_json(this, x, n, m, json, scale, auto_scale)
206 ! ----------------------------------------------------- !
207 ! Initializing the mma object and all the parameters !
208 ! required for MMA method. (a_i, c_i, d_i, ...) !
209 ! x: the design varaibles(DV), n: number of DV, !
210 ! m: number of constraints !
211 ! !
212 ! Note that residumax & residunorm of the KKT conditions!
213 ! are initialized with 10^5. This is done to avoid !
214 ! unnecessary extera computation of KKT norms for the !
215 ! initial design. !
216 ! ----------------------------------------------------- !
217 class(mma_t), intent(inout) :: this
218 integer, intent(in) :: n, m
219 type(vector_t), intent(in) :: x
220
221 type(json_file), intent(inout) :: json
222
223 ! Read the scaling info for fval and dfdx from json
224 real(kind=rp), intent(out) :: scale
225 logical, intent(out) :: auto_scale
226 ! -------------------------------------------------------------------!
227 ! Internal parameters for MMA !
228 ! Minimize f_0(x) + a_0*z + sum( c_i*y_i + 0.5*d_i*(y_i)^2 ) !
229 ! subject to f_i(x) - a_i*z - y_i <= 0, i = 1,...,m !
230 ! xmin_j <= x_j <= xmax_j, j = 1,...,n !
231 ! z >= 0, y_i >= 0, i = 1,...,m !
232 ! -------------------------------------------------------------------!
233 real(kind=rp), dimension(n) :: xmax, xmin
234 real(kind=rp), dimension(m) :: a, c, d
235 character(len=:), allocatable :: subsolver, bcknd, bcknd_default
236
237 ! For reading the values from json and then set the value for the arrays
238 real(kind=rp) :: a0 , xmax_const, xmin_const, a_const, c_const, d_const
239
240 integer :: max_iter, n_global, ierr
241 real(kind=rp) :: epsimin, asyinit, asyincr, asydecr
242
243 call mpi_allreduce(n, n_global, 1, mpi_integer, &
244 mpi_sum, neko_comm, ierr)
245
246 ! Assign default values for the backend based on the NEKO_BCKND_DEVICE
247 if (neko_bcknd_device .eq. 1) then
248 bcknd_default = "device"
249 else
250 bcknd_default = "cpu"
251 end if
252
253 ! ------------------------------------------------------------------------ !
254 ! Assign defaults if nothing is parsed
255 ! based on the Cpp Code by Niels
256 call json_get_or_default(json, 'mma.epsimin', epsimin, &
257 1.0e-9_rp * sqrt(real(m + n_global, rp)))
258 call json_get_or_default(json, 'mma.max_iter', max_iter, 100)
259
260 ! Following parameters are set based on eq.3.8:--------
261 call json_get_or_default(json, 'mma.asyinit', asyinit, 0.5_rp)
262 call json_get_or_default(json, 'mma.asyincr', asyincr, 1.2_rp)
263 call json_get_or_default(json, 'mma.asydecr', asydecr, 0.7_rp)
264
265 call json_get_or_default(json, 'mma.backend', bcknd, bcknd_default)
266 call json_get_or_default(json, 'mma.subsolver', subsolver, 'dip')
267
268 call json_get_or_default(json, 'mma.xmin', xmin_const, 0.0_rp)
269 call json_get_or_default(json, 'mma.xmax', xmax_const, 1.0_rp)
270 call json_get_or_default(json, 'mma.a0', a0, 1.0_rp)
271 call json_get_or_default(json, 'mma.a', a_const, 0.0_rp)
272 call json_get_or_default(json, 'mma.c', c_const, 100.0_rp)
273 call json_get_or_default(json, 'mma.d', d_const, 0.0_rp)
274
275 call json_get_or_default(json, 'mma.scale', scale, 10.0_rp)
276 call json_get_or_default(json, 'mma.auto_scale', auto_scale, .false.)
277
278 ! Initialize the MMA object with the parsed parameters
279 a = a_const
280 c = c_const
281 d = d_const
282 xmin = xmin_const
283 xmax = xmax_const
284
285 ! ------------------------------------------------------------------------ !
286 ! Initialize the MMA object with the parameters read from json
287
288 call this%init(x, n, m, a0, a, c, d, xmin, xmax, &
289 max_iter, epsimin, asyinit, asyincr, asydecr, bcknd, subsolver)
290
291 end subroutine mma_init_from_json
292
294 subroutine mma_free(this)
295 class(mma_t), intent(inout) :: this
296 ! Deallocate the internal vectors
297 call this%xold1%free()
298 call this%xold2%free()
299 call this%alpha%free()
300 call this%beta%free()
301 call this%a%free()
302 call this%c%free()
303 call this%d%free()
304 call this%low%free()
305 call this%upp%free()
306 call this%xmax%free()
307 call this%xmin%free()
308 call this%p0j%free()
309 call this%q0j%free()
310 call this%bi%free()
311 call this%y%free()
312 call this%lambda%free()
313 call this%s%free()
314 call this%mu%free()
315 call this%xsi%free()
316 call this%eta%free()
317
318 ! Deallocate the internal dummy matrices
319 call this%pij%free()
320 call this%qij%free()
321
322 this%is_initialized = .false.
323 this%is_updated = .false.
324 end subroutine mma_free
325
327 subroutine mma_init_from_components(this, x, n, m, a0, a, c, d, xmin, xmax, &
328 max_iter, epsimin, asyinit, asyincr, asydecr, bcknd, subsolver)
329 ! ----------------------------------------------------- !
330 ! Initializing the mma object and all the parameters !
331 ! required for MMA method. (a_i, c_i, d_i, ...) !
332 ! x: the design varaibles(DV), n: number of DV, !
333 ! m: number of constraints !
334 ! !
335 ! Note that residumax & residunorm of the KKT conditions!
336 ! are initialized with 10^5. This is done to avoid !
337 ! unnecessary extera computation of KKT norms for the !
338 ! initial design. !
339 ! ----------------------------------------------------- !
340 class(mma_t), intent(inout) :: this
341 integer, intent(in) :: n, m
342 type(vector_t), intent(in) :: x
343 ! -------------------------------------------------------------------!
344 ! Internal parameters for MMA !
345 ! Minimize f_0(x) + a_0*z + sum( c_i*y_i + 0.5*d_i*(y_i)^2 ) !
346 ! subject to f_i(x) - a_i*z - y_i <= 0, i = 1,...,m !
347 ! xmin_j <= x_j <= xmax_j, j = 1,...,n !
348 ! z >= 0, y_i >= 0, i = 1,...,m !
349 ! -------------------------------------------------------------------!
350 real(kind=rp), intent(in), dimension(n) :: xmax, xmin
351 real(kind=rp), intent(in), dimension(m) :: a, c, d
352 real(kind=rp), intent(in) :: a0
353 integer, intent(in), optional :: max_iter
354 real(kind=rp), intent(in), optional :: epsimin, asyinit, asyincr, asydecr
355 character(len=*), intent(in), optional :: bcknd, subsolver
356 character(len=256) :: log_msg
357 integer :: i, ierr
358
359 call this%free()
360
361 this%n = n
362 this%m = m
363
364 call this%xold1%init(n)
365 call this%xold2%init(n)
366 this%xold1 = x
367 this%xold2 = x
368
369 call this%alpha%init(n)
370 call this%beta%init(n)
371
372 call this%a%init(m)
373 call this%c%init(m)
374 call this%d%init(m)
375 call this%low%init(n)
376 call this%upp%init(n)
377 call this%xmax%init(n)
378 call this%xmin%init(n)
379
380 !internal dummy variables for MMA
381 call this%p0j%init(n)
382 call this%q0j%init(n)
383 call this%pij%init(m, n)
384 call this%qij%init(m, n)
385 call this%bi%init(m)
386
387 ! Necessary for KKT check after updating df0dx, fval, dfdx
388 call this%y%init(m)
389 call this%lambda%init(m)
390 call this%s%init(m)
391 call this%mu%init(m)
392 call this%xsi%init(n)
393 call this%eta%init(n)
394
395 this%a0 = a0
396 this%a%x = a
397 this%c%x = c
398 this%d%x = d
399
400 ! Set the bounds for the design variable based on the problem
401 this%xmax%x = xmax
402 this%xmin%x = xmin
403
404 if (neko_bcknd_device .eq. 1) then
405 call this%a%copy_from(host_to_device, sync = .false.)
406 call this%c%copy_from(host_to_device, sync = .false.)
407 call this%d%copy_from(host_to_device, sync = .false.)
408 call this%xmax%copy_from(host_to_device, sync = .false.)
409 call this%xmin%copy_from(host_to_device, sync = .true.)
410 end if
411
412 ! Set KKT norms to a large number for the initial design
413 this%residumax = huge(0.0_rp)
414 this%residunorm = huge(0.0_rp)
415
416 ! Sync parameters across MPI
417 call mpi_allreduce(n, this%n_global, 1, mpi_integer, mpi_sum, neko_comm, &
418 ierr)
419
420 ! ------------------------------------------------------------------------ !
421 ! Assign defaults if nothing is parsed
422
423 ! Based on the Cpp Code by Niels
424 if (.not. present(max_iter)) this%max_iter = 100
425 if (.not. present(epsimin)) then
426 this%epsimin = 1.0e-9_rp * sqrt(real(this%m + this%n_global, rp))
427 end if
428
429 ! Following parameters are set based on eq.3.8
430 if (.not. present(asyinit)) this%asyinit = 0.5_rp
431 if (.not. present(asyincr)) this%asyincr = 1.2_rp
432 if (.not. present(asydecr)) this%asydecr = 0.7_rp
433
434 ! Set default backend based on NEKO_BCKND_DEVICE
435 if (.not. present(bcknd) .and. neko_bcknd_device .eq. 0) then
436 this%bcknd = "cpu"
437 else if (.not. present(bcknd)) then
438 this%bcknd = "device"
439 end if
440
441 ! Set default subsolver
442 if (.not. present(subsolver)) this%subsolver = "dip"
443
444 ! Assign values from inputs when present
445 if (present(max_iter)) this%max_iter = max_iter
446 if (present(epsimin)) this%epsimin = epsimin
447 if (present(asyinit)) this%asyinit = asyinit
448 if (present(asyincr)) this%asyincr = asyincr
449 if (present(asydecr)) this%asydecr = asydecr
450 if (present(bcknd)) this%bcknd = bcknd
451 if (present(subsolver)) this%subsolver = subsolver
452
453 call neko_log%section('MMA Parameters')
454
455 write(log_msg, '(A10,1X,A)') 'backend ', trim(this%bcknd)
456 call neko_log%message(log_msg)
457 write(log_msg, '(A10,1X,A)') 'subsolver ', trim(this%subsolver)
458 call neko_log%message(log_msg)
459
460 write(log_msg, '(A10,1X,I0)') 'n ', this%n_global
461 call neko_log%message(log_msg)
462 write(log_msg, '(A10,1X,I0)') 'm ', this%m
463 call neko_log%message(log_msg)
464 write(log_msg, '(A10,1X,I0)') 'max_iter ', this%max_iter
465 call neko_log%message(log_msg)
466
467 write(log_msg, '(A10,1X,E11.5)') 'epsimin ', this%epsimin
468 call neko_log%message(log_msg)
469
470 write(log_msg, '(A10,1X,E11.5)') 'asyinit ', this%asyinit
471 call neko_log%message(log_msg)
472 write(log_msg, '(A10,1X,E11.5)') 'asyincr ', this%asyincr
473 call neko_log%message(log_msg)
474 write(log_msg, '(A10,1X,E11.5)') 'asydecr ', this%asydecr
475 call neko_log%message(log_msg)
476 write(log_msg, '(A10,1X,E11.5)') 'a0 ', this%a0
477 call neko_log%message(log_msg)
478
479 call neko_log%message('Parameters a')
480 do i = 1, this%m
481 write(log_msg, '(3X,A,I2,A,E11.5)') 'a(', i, ') = ', this%a%x(i)
482 call neko_log%message(log_msg)
483 end do
484 call neko_log%message('Parameters c')
485 do i = 1, this%m
486 write(log_msg, '(3X,A,I2,A,E11.5)') 'c(', i, ') = ', this%c%x(i)
487 call neko_log%message(log_msg)
488 end do
489 call neko_log%message('Parameters d')
490 do i = 1, this%m
491 write(log_msg, '(3X,A,I2,A,E11.5)') 'd(', i, ') = ', this%d%x(i)
492 call neko_log%message(log_msg)
493 end do
494
495 call neko_log%end_section()
496
497 ! The object is correctly initialized
498 this%is_initialized = .true.
499 end subroutine mma_init_from_components
500
501 ! ========================================================================== !
502 ! Updator and KKT checker
503
505 subroutine mma_update_vector(this, iter, x, df0dx, fval, dfdx)
506 class(mma_t), intent(inout) :: this
507 integer, intent(in) :: iter
508 type(vector_t), intent(inout) :: x
509 type(vector_t), intent(inout) :: df0dx, fval
510 type(matrix_t), intent(inout) :: dfdx
511
512 ! Select backend type
513 select case (this%bcknd)
514 case ("cpu")
515 if (neko_bcknd_device .eq. 1) then
516 call x%copy_from(device_to_host, sync = .false.)
517 call df0dx%copy_from(device_to_host, sync = .false.)
518 call fval%copy_from(device_to_host, sync = .false.)
519 call dfdx%copy_from(device_to_host, sync = .true.)
520 end if
521
522 call mma_update_cpu(this, iter, x%x, df0dx%x, fval%x, dfdx%x)
523
524 if (neko_bcknd_device .eq. 1) then
525 call x%copy_from(host_to_device, sync = .true.)
526 end if
527
528 case ("device")
529 call mma_update_device(this, iter, x%x_d, df0dx%x_d, fval%x_d, dfdx%x_d)
530 end select
531
532 end subroutine mma_update_vector
533
535 subroutine mma_kkt_vector(this, x, df0dx, fval, dfdx)
536 class(mma_t), intent(inout) :: this
537 type(vector_t), intent(inout) :: x, df0dx, fval
538 type(matrix_t), intent(inout) :: dfdx
539
540 ! Select backend type
541 select case (this%bcknd )
542 case ("cpu")
543 if (neko_bcknd_device .eq. 1) then
544 call device_memcpy(x%x, x%x_d, this%n, device_to_host, &
545 sync = .false.)
546 call device_memcpy(df0dx%x, df0dx%x_d, this%n, device_to_host, &
547 sync = .false.)
548 call device_memcpy(fval%x, fval%x_d, this%m, device_to_host, &
549 sync = .false.)
550 call device_memcpy(dfdx%x, dfdx%x_d, this%m * this%n, device_to_host,&
551 sync = .true.)
552 end if
553
554 call mma_kkt_cpu(this, x%x, df0dx%x, fval%x, dfdx%x)
555 case ("device")
556 call mma_kkt_device(this, x%x_d, df0dx%x_d, fval%x_d, dfdx%x_d)
557 end select
558 end subroutine mma_kkt_vector
559
560 ! ========================================================================== !
561 ! IO Functions
562
567 subroutine mma_save_checkpoint(this, filename, overwrite)
568 class(mma_t), intent(inout) :: this
569 character(len=*), intent(in) :: filename
570 logical, intent(in), optional :: overwrite
571 character(len=12) :: file_ext
572
573 ! Get the file extension
574 call filename_suffix(filename, file_ext)
575
576 select case (trim(file_ext))
577 case ('h5', 'hdf5', 'hf5')
578 call mma_save_checkpoint_hdf5(this, filename, overwrite)
579 case default
580 call neko_error('mma_save_checkpoint: Unsupported file format: ' // &
581 trim(file_ext))
582 end select
583
584 end subroutine mma_save_checkpoint
585
589 subroutine mma_load_checkpoint(this, filename)
590 class(mma_t), intent(inout) :: this
591 character(len=*), intent(in) :: filename
592 character(len=12) :: file_ext
593
594 ! Get the file extension
595 call filename_suffix(filename, file_ext)
596
597 select case (trim(file_ext))
598 case ('h5', 'hdf5', 'hf5')
599 call mma_load_checkpoint_hdf5(this, filename)
600 case default
601 call neko_error('mma_load_checkpoint: Unsupported file format: ' // &
602 trim(file_ext))
603 end select
604 end subroutine mma_load_checkpoint
605
606 ! ========================================================================== !
607 ! Getters and setters
608
610 pure function mma_get_n(this) result(n)
611 class(mma_t), intent(in) :: this
612 integer :: n
613 n = this%n
614 end function mma_get_n
615
617 pure function mma_get_m(this) result(m)
618 class(mma_t), intent(in) :: this
619 integer :: m
620 m = this%m
621 end function mma_get_m
622
624 pure function mma_get_residumax(this) result(residumax)
625 class(mma_t), intent(in) :: this
626 real(kind=rp) :: residumax
627 residumax = this%residumax
628 end function mma_get_residumax
629
631 pure function mma_get_residunorm(this) result(residunorm)
632 class(mma_t), intent(in) :: this
633 real(kind=rp) :: residunorm
634 residunorm = this%residunorm
635 end function mma_get_residunorm
636
638 pure function mma_get_max_iter(this) result(max_iter_value)
639 class(mma_t), intent(in) :: this
640 integer :: max_iter_value
641 max_iter_value = this%max_iter
642 end function mma_get_max_iter
643
645 pure function mma_get_backend_and_subsolver(this) result(backend_subsolver)
646 class(mma_t), intent(in) :: this
647 character(len=:), allocatable :: backend_subsolver
648 character(len=:), allocatable :: backend
649
650 if (neko_bcknd_cuda .eq. 1) then
651 backend = "cuda"
652 else if (neko_bcknd_hip .eq. 1) then
653 backend = "hip"
654 else if (neko_bcknd_opencl .eq. 1) then
655 backend = "opencl"
656 else
657 backend = "cpu"
658 end if
659
660 backend_subsolver = 'backend:' // trim(backend) // ', subsolver:' // &
661 trim(this%subsolver)
662 end function mma_get_backend_and_subsolver
663
664 ! ========================================================================== !
665 ! Private utilities
666
668 subroutine mma_copy_from(this, direction, sync)
669 class(mma_t), intent(inout) :: this
670 integer, intent(in) :: direction
671 logical, intent(in) :: sync
672
673 call this%xold1%copy_from(direction, sync = .false.)
674 call this%xold2%copy_from(direction, sync = .false.)
675 call this%xmax%copy_from(direction, sync = .false.)
676 call this%xmin%copy_from(direction, sync = .false.)
677
678 call this%low%copy_from(direction, sync = .false.)
679 call this%upp%copy_from(direction, sync = .false.)
680
681 call this%a%copy_from(direction, sync = .false.)
682 call this%c%copy_from(direction, sync = .false.)
683 call this%d%copy_from(direction, sync = .false.)
684 call this%y%copy_from(direction, sync = .false.)
685 call this%s%copy_from(direction, sync = .false.)
686
687 call this%p0j%copy_from(direction, sync = .false.)
688 call this%q0j%copy_from(direction, sync = .false.)
689 call this%pij%copy_from(direction, sync = .false.)
690 call this%qij%copy_from(direction, sync = .false.)
691 call this%bi%copy_from(direction, sync = .false.)
692
693 call this%alpha%copy_from(direction, sync = .false.)
694 call this%beta%copy_from(direction, sync = .false.)
695 call this%lambda%copy_from(direction, sync = .false.)
696 call this%mu%copy_from(direction, sync = .false.)
697 call this%xsi%copy_from(direction, sync = .false.)
698 call this%eta%copy_from(direction, sync = sync)
699
700 end subroutine mma_copy_from
701
702 ! ========================================================================== !
703 ! Dummy implementations for module procedures
704
705#if !HAVE_HDF5
706 module subroutine mma_save_checkpoint_hdf5(object, filename, overwrite)
707 class(mma_t), intent(inout) :: object
708 character(len=*), intent(in) :: filename
709 logical, intent(in), optional :: overwrite
710 call neko_error('mma: HDF5 support not enabled rebuild with HAVE_HDF5')
711 end subroutine mma_save_checkpoint_hdf5
712
713 module subroutine mma_load_checkpoint_hdf5(object, filename)
714 class(mma_t), intent(inout) :: object
715 character(len=*), intent(in) :: filename
716 call neko_error('mma: HDF5 support not enabled rebuild with HAVE_HDF5')
717 end subroutine mma_load_checkpoint_hdf5
718#endif
719
720end module mma
MMA module.
Definition mma.f90:69
subroutine mma_init_from_json(this, x, n, m, json, scale, auto_scale)
Device KKT check for convergence.
Definition mma.f90:206
MMA type.
Definition mma.f90:89