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