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
78 use neko_config,
only: neko_bcknd_device, neko_bcknd_cuda, neko_bcknd_hip, &
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
91 integer :: n, m, n_global, max_iter
92 real(kind=rp) :: a0, asyinit, asyincr, asydecr, epsimin, &
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
100 type(vector_t) :: p0j, q0j
101 type(matrix_t) :: pij, qij
105 real(kind=rp) :: z, zeta
106 type(vector_t) :: y, lambda, s, mu
107 type(vector_t) :: xsi, eta
110 generic,
public :: init => init_from_json, init_from_components
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
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
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
133 generic ::
write => write_hdf5
134 procedure, pass(this) :: write_hdf5 => mma_write_hdf5
136 generic ::
read => read_hdf5
137 procedure, pass(this) :: read_hdf5 => mma_read_hdf5
140 procedure, pass(this) :: copy_from => mma_copy_from
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
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
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
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
186 module subroutine mma_write_hdf5(this, filename, overwrite)
187 class(mma_t),
intent(inout) :: this
188 character(len=*),
intent(in) :: filename
189 logical,
intent(in),
optional :: overwrite
190 end subroutine mma_write_hdf5
192 module subroutine mma_read_hdf5(this, filename)
193 class(mma_t),
intent(inout) :: this
194 character(len=*),
intent(in) :: filename
195 end subroutine mma_read_hdf5
201 subroutine mma_init_from_json(this, x, n, m, json, scale, auto_scale)
213 class(mma_t),
intent(inout) :: this
214 integer,
intent(in) :: n, m
215 type(vector_t),
intent(in) :: x
217 type(json_file),
intent(inout) :: json
220 real(kind=rp),
intent(out) :: scale
221 logical,
intent(out) :: auto_scale
229 real(kind=rp),
dimension(n) :: xmax, xmin
230 real(kind=rp),
dimension(m) :: a, c, d
231 character(len=:),
allocatable :: subsolver, bcknd, bcknd_default
234 real(kind=rp) :: a0 , xmax_const, xmin_const, a_const, c_const, d_const
236 integer :: max_iter, n_global, ierr
237 real(kind=rp) :: epsimin, asyinit, asyincr, asydecr
239 call mpi_allreduce(n, n_global, 1, mpi_integer, &
240 mpi_sum, neko_comm, ierr)
243 if (neko_bcknd_device .eq. 1)
then
244 bcknd_default =
"device"
246 bcknd_default =
"cpu"
252 call json_get_or_default(json,
'mma.epsimin', epsimin, &
253 1.0e-9_rp * sqrt(real(m + n_global, rp)))
254 call json_get_or_default(json,
'mma.max_iter', max_iter, 100)
257 call json_get_or_default(json,
'mma.asyinit', asyinit, 0.5_rp)
258 call json_get_or_default(json,
'mma.asyincr', asyincr, 1.2_rp)
259 call json_get_or_default(json,
'mma.asydecr', asydecr, 0.7_rp)
261 call json_get_or_default(json,
'mma.backend', bcknd, bcknd_default)
262 call json_get_or_default(json,
'mma.subsolver', subsolver,
'dip')
264 call json_get_or_default(json,
'mma.xmin', xmin_const, 0.0_rp)
265 call json_get_or_default(json,
'mma.xmax', xmax_const, 1.0_rp)
266 call json_get_or_default(json,
'mma.a0', a0, 1.0_rp)
267 call json_get_or_default(json,
'mma.a', a_const, 0.0_rp)
268 call json_get_or_default(json,
'mma.c', c_const, 100.0_rp)
269 call json_get_or_default(json,
'mma.d', d_const, 0.0_rp)
271 call json_get_or_default(json,
'mma.scale', scale, 10.0_rp)
272 call json_get_or_default(json,
'mma.auto_scale', auto_scale, .false.)
284 call this%init(x, n, m, a0, a, c, d, xmin, xmax, &
285 max_iter, epsimin, asyinit, asyincr, asydecr, bcknd, subsolver)
287 end subroutine mma_init_from_json
290 subroutine mma_free(this)
291 class(mma_t),
intent(inout) :: this
293 call this%xold1%free()
294 call this%xold2%free()
295 call this%alpha%free()
296 call this%beta%free()
302 call this%xmax%free()
303 call this%xmin%free()
308 call this%lambda%free()
318 this%is_initialized = .false.
319 this%is_updated = .false.
320 end subroutine mma_free
323 subroutine mma_init_from_components(this, x, n, m, a0, a, c, d, xmin, xmax, &
324 max_iter, epsimin, asyinit, asyincr, asydecr, bcknd, subsolver)
336 class(mma_t),
intent(inout) :: this
337 integer,
intent(in) :: n, m
338 type(vector_t),
intent(in) :: x
346 real(kind=rp),
intent(in),
dimension(n) :: xmax, xmin
347 real(kind=rp),
intent(in),
dimension(m) :: a, c, d
348 real(kind=rp),
intent(in) :: a0
349 integer,
intent(in),
optional :: max_iter
350 real(kind=rp),
intent(in),
optional :: epsimin, asyinit, asyincr, asydecr
351 character(len=:),
intent(in),
allocatable :: bcknd, subsolver
352 character(len=256) :: log_msg
360 call this%xold1%init(n)
361 call this%xold2%init(n)
365 call this%alpha%init(n)
366 call this%beta%init(n)
371 call this%low%init(n)
372 call this%upp%init(n)
373 call this%xmax%init(n)
374 call this%xmin%init(n)
377 call this%p0j%init(n)
378 call this%q0j%init(n)
379 call this%pij%init(m, n)
380 call this%qij%init(m, n)
385 call this%lambda%init(m)
388 call this%xsi%init(n)
389 call this%eta%init(n)
400 if (neko_bcknd_device .eq. 1)
then
401 call this%a%copy_from(host_to_device, sync = .false.)
402 call this%c%copy_from(host_to_device, sync = .false.)
403 call this%d%copy_from(host_to_device, sync = .false.)
404 call this%xmax%copy_from(host_to_device, sync = .false.)
405 call this%xmin%copy_from(host_to_device, sync = .true.)
409 this%residumax = huge(0.0_rp)
410 this%residunorm = huge(0.0_rp)
416 if (.not.
present(max_iter)) this%max_iter = 100
417 if (.not.
present(epsimin)) this%epsimin = 1.0e-9_rp * sqrt(real(m + n, rp))
420 if (.not.
present(asyinit)) this%asyinit = 0.5_rp
421 if (.not.
present(asyincr)) this%asyincr = 1.2_rp
422 if (.not.
present(asydecr)) this%asydecr = 0.7_rp
425 if (
present(max_iter)) this%max_iter = max_iter
426 if (
present(epsimin)) this%epsimin = epsimin
427 if (
present(asyinit)) this%asyinit = asyinit
428 if (
present(asyincr)) this%asyincr = asyincr
429 if (
present(asydecr)) this%asydecr = asydecr
431 this%subsolver = subsolver
434 call mpi_allreduce(this%n, this%n_global, 1, &
435 mpi_integer, mpi_sum, neko_comm, ierr)
437 call neko_log%section(
'MMA Parameters')
439 write(log_msg,
'(A10,1X,A)')
'backend ', trim(this%bcknd)
440 call neko_log%message(log_msg)
441 write(log_msg,
'(A10,1X,A)')
'subsolver ', trim(this%subsolver)
442 call neko_log%message(log_msg)
444 write(log_msg,
'(A10,1X,I0)')
'n ', this%n
445 call neko_log%message(log_msg)
446 write(log_msg,
'(A10,1X,I0)')
'm ', this%m
447 call neko_log%message(log_msg)
448 write(log_msg,
'(A10,1X,I0)')
'max_iter ', this%max_iter
449 call neko_log%message(log_msg)
451 write(log_msg,
'(A10,1X,E11.5)')
'epsimin ', this%epsimin
452 call neko_log%message(log_msg)
454 write(log_msg,
'(A10,1X,E11.5)')
'asyinit ', this%asyinit
455 call neko_log%message(log_msg)
456 write(log_msg,
'(A10,1X,E11.5)')
'asyincr ', this%asyincr
457 call neko_log%message(log_msg)
458 write(log_msg,
'(A10,1X,E11.5)')
'asydecr ', this%asydecr
459 call neko_log%message(log_msg)
460 write(log_msg,
'(A10,1X,E11.5)')
'a0 ', this%a0
461 call neko_log%message(log_msg)
463 call neko_log%message(
'Parameters a')
465 write(log_msg,
'(3X,A,I2,A,E11.5)')
'a(', i,
') = ', this%a%x(i)
466 call neko_log%message(log_msg)
468 call neko_log%message(
'Parameters c')
470 write(log_msg,
'(3X,A,I2,A,E11.5)')
'c(', i,
') = ', this%c%x(i)
471 call neko_log%message(log_msg)
473 call neko_log%message(
'Parameters d')
475 write(log_msg,
'(3X,A,I2,A,E11.5)')
'd(', i,
') = ', this%d%x(i)
476 call neko_log%message(log_msg)
479 call neko_log%end_section()
482 this%is_initialized = .true.
483 end subroutine mma_init_from_components
486 subroutine mma_update_vector(this, iter, x, df0dx, fval, dfdx)
487 class(mma_t),
intent(inout) :: this
488 integer,
intent(in) :: iter
489 type(vector_t),
intent(inout) :: x
490 type(vector_t),
intent(inout) :: df0dx, fval
491 type(matrix_t),
intent(inout) :: dfdx
494 select case (this%bcknd)
496 if (neko_bcknd_device .eq. 1)
then
497 call device_memcpy(x%x, x%x_d, this%n, device_to_host, &
499 call device_memcpy(df0dx%x, df0dx%x_d, this%n, device_to_host, &
501 call device_memcpy(fval%x, fval%x_d, this%m, device_to_host, &
503 call device_memcpy(dfdx%x, dfdx%x_d, this%m * this%n, device_to_host,&
507 call mma_update_cpu(this, iter, x%x, df0dx%x, fval%x, dfdx%x)
509 if (neko_bcknd_device .eq. 1)
then
510 call device_memcpy(x%x, x%x_d, this%n, host_to_device, sync = .true.)
514 call mma_update_device(this, iter, x%x_d, df0dx%x_d, fval%x_d, dfdx%x_d)
517 end subroutine mma_update_vector
520 subroutine mma_kkt_vector(this, x, df0dx, fval, dfdx)
521 class(mma_t),
intent(inout) :: this
522 type(vector_t),
intent(inout) :: x, df0dx, fval
523 type(matrix_t),
intent(inout) :: dfdx
526 select case (this%bcknd )
528 if (neko_bcknd_device .eq. 1)
then
529 call device_memcpy(x%x, x%x_d, this%n, device_to_host, &
531 call device_memcpy(df0dx%x, df0dx%x_d, this%n, device_to_host, &
533 call device_memcpy(fval%x, fval%x_d, this%m, device_to_host, &
535 call device_memcpy(dfdx%x, dfdx%x_d, this%m * this%n, device_to_host,&
539 call mma_kkt_cpu(this, x%x, df0dx%x, fval%x, dfdx%x)
541 call mma_kkt_device(this, x%x_d, df0dx%x_d, fval%x_d, dfdx%x_d)
543 end subroutine mma_kkt_vector
549 pure function mma_get_n(this)
result(n)
550 class(mma_t),
intent(in) :: this
553 end function mma_get_n
556 pure function mma_get_m(this)
result(m)
557 class(mma_t),
intent(in) :: this
560 end function mma_get_m
563 pure function mma_get_residumax(this)
result(residumax)
564 class(mma_t),
intent(in) :: this
565 real(kind=rp) :: residumax
566 residumax = this%residumax
567 end function mma_get_residumax
570 pure function mma_get_residunorm(this)
result(residunorm)
571 class(mma_t),
intent(in) :: this
572 real(kind=rp) :: residunorm
573 residunorm = this%residunorm
574 end function mma_get_residunorm
577 pure function mma_get_max_iter(this)
result(max_iter_value)
578 class(mma_t),
intent(in) :: this
579 integer :: max_iter_value
580 max_iter_value = this%max_iter
581 end function mma_get_max_iter
584 pure function mma_get_backend_and_subsolver(this)
result(backend_subsolver)
585 class(mma_t),
intent(in) :: this
586 character(len=:),
allocatable :: backend_subsolver
587 character(len=:),
allocatable :: backend
589 if (neko_bcknd_cuda .eq. 1)
then
591 else if (neko_bcknd_hip .eq. 1)
then
593 else if (neko_bcknd_opencl .eq. 1)
then
599 backend_subsolver =
'backend:' // trim(backend) //
', subsolver:' // &
601 end function mma_get_backend_and_subsolver
607 subroutine mma_copy_from(this, direction, sync)
608 class(mma_t),
intent(inout) :: this
609 integer,
intent(in) :: direction
610 logical,
intent(in) :: sync
612 call this%xold1%copy_from(direction, sync = .false.)
613 call this%xold2%copy_from(direction, sync = .false.)
614 call this%xmax%copy_from(direction, sync = .false.)
615 call this%xmin%copy_from(direction, sync = .false.)
617 call this%low%copy_from(direction, sync = .false.)
618 call this%upp%copy_from(direction, sync = .false.)
620 call this%a%copy_from(direction, sync = .false.)
621 call this%c%copy_from(direction, sync = .false.)
622 call this%d%copy_from(direction, sync = .false.)
623 call this%y%copy_from(direction, sync = .false.)
624 call this%s%copy_from(direction, sync = .false.)
626 call this%p0j%copy_from(direction, sync = .false.)
627 call this%q0j%copy_from(direction, sync = .false.)
628 call this%pij%copy_from(direction, sync = .false.)
629 call this%qij%copy_from(direction, sync = .false.)
630 call this%bi%copy_from(direction, sync = .false.)
632 call this%alpha%copy_from(direction, sync = .false.)
633 call this%beta%copy_from(direction, sync = .false.)
634 call this%lambda%copy_from(direction, sync = .false.)
635 call this%mu%copy_from(direction, sync = .false.)
636 call this%xsi%copy_from(direction, sync = .false.)
637 call this%eta%copy_from(direction, sync = sync)
639 end subroutine mma_copy_from
subroutine mma_init_from_json(this, x, n, m, json, scale, auto_scale)
Device KKT check for convergence.