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 aproximated problem !
65! 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
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 generic :: write => write_hdf5
134 procedure, pass(this) :: write_hdf5 => mma_write_hdf5
135
136 generic :: read => read_hdf5
137 procedure, pass(this) :: read_hdf5 => mma_read_hdf5
138
139 ! Private utilities
140 procedure, pass(this) :: copy_from => mma_copy_from
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 ! ======================================================================= !
184 ! Interface for IO routines
185
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
191
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
196 end interface
197
198contains
199
201 subroutine mma_init_from_json(this, x, n, m, json, scale, auto_scale)
202 ! ----------------------------------------------------- !
203 ! Initializing the mma object and all the parameters !
204 ! required for MMA method. (a_i, c_i, d_i, ...) !
205 ! x: the design varaibles(DV), n: number of DV, !
206 ! m: number of constraints !
207 ! !
208 ! Note that residumax & residunorm of the KKT conditions!
209 ! are initialized with 10^5. This is done to avoid !
210 ! unnecessary extera computation of KKT norms for the !
211 ! initial design. !
212 ! ----------------------------------------------------- !
213 class(mma_t), intent(inout) :: this
214 integer, intent(in) :: n, m
215 type(vector_t), intent(in) :: x
216
217 type(json_file), intent(inout) :: json
218
219 ! Read the scaling info for fval and dfdx from json
220 real(kind=rp), intent(out) :: scale
221 logical, intent(out) :: auto_scale
222 ! -------------------------------------------------------------------!
223 ! Internal parameters for MMA !
224 ! Minimize f_0(x) + a_0*z + sum( c_i*y_i + 0.5*d_i*(y_i)^2 ) !
225 ! subject to f_i(x) - a_i*z - y_i <= 0, i = 1,...,m !
226 ! xmin_j <= x_j <= xmax_j, j = 1,...,n !
227 ! z >= 0, y_i >= 0, i = 1,...,m !
228 ! -------------------------------------------------------------------!
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
232
233 ! For reading the values from json and then set the value for the arrays
234 real(kind=rp) :: a0 , xmax_const, xmin_const, a_const, c_const, d_const
235
236 integer :: max_iter, n_global, ierr
237 real(kind=rp) :: epsimin, asyinit, asyincr, asydecr
238
239 call mpi_allreduce(n, n_global, 1, mpi_integer, &
240 mpi_sum, neko_comm, ierr)
241
242 ! Assign default values for the backend based on the NEKO_BCKND_DEVICE
243 if (neko_bcknd_device .eq. 1) then
244 bcknd_default = "device"
245 else
246 bcknd_default = "cpu"
247 end if
248
249 ! ------------------------------------------------------------------------ !
250 ! Assign defaults if nothing is parsed
251 ! based on the Cpp Code by Niels
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)
255
256 ! Following parameters are set based on eq.3.8:--------
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)
260
261 call json_get_or_default(json, 'mma.backend', bcknd, bcknd_default)
262 call json_get_or_default(json, 'mma.subsolver', subsolver, 'dip')
263
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)
270
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.)
273
274 ! Initialize the MMA object with the parsed parameters
275 a = a_const
276 c = c_const
277 d = d_const
278 xmin = xmin_const
279 xmax = xmax_const
280
281 ! ------------------------------------------------------------------------ !
282 ! Initialize the MMA object with the parameters read from json
283
284 call this%init(x, n, m, a0, a, c, d, xmin, xmax, &
285 max_iter, epsimin, asyinit, asyincr, asydecr, bcknd, subsolver)
286
287 end subroutine mma_init_from_json
288
290 subroutine mma_free(this)
291 class(mma_t), intent(inout) :: this
292 ! Deallocate the internal vectors
293 call this%xold1%free()
294 call this%xold2%free()
295 call this%alpha%free()
296 call this%beta%free()
297 call this%a%free()
298 call this%c%free()
299 call this%d%free()
300 call this%low%free()
301 call this%upp%free()
302 call this%xmax%free()
303 call this%xmin%free()
304 call this%p0j%free()
305 call this%q0j%free()
306 call this%bi%free()
307 call this%y%free()
308 call this%lambda%free()
309 call this%s%free()
310 call this%mu%free()
311 call this%xsi%free()
312 call this%eta%free()
313
314 ! Deallocate the internal dummy matrices
315 call this%pij%free()
316 call this%qij%free()
317
318 this%is_initialized = .false.
319 this%is_updated = .false.
320 end subroutine mma_free
321
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)
325 ! ----------------------------------------------------- !
326 ! Initializing the mma object and all the parameters !
327 ! required for MMA method. (a_i, c_i, d_i, ...) !
328 ! x: the design varaibles(DV), n: number of DV, !
329 ! m: number of constraints !
330 ! !
331 ! Note that residumax & residunorm of the KKT conditions!
332 ! are initialized with 10^5. This is done to avoid !
333 ! unnecessary extera computation of KKT norms for the !
334 ! initial design. !
335 ! ----------------------------------------------------- !
336 class(mma_t), intent(inout) :: this
337 integer, intent(in) :: n, m
338 type(vector_t), intent(in) :: x
339 ! -------------------------------------------------------------------!
340 ! Internal parameters for MMA !
341 ! Minimize f_0(x) + a_0*z + sum( c_i*y_i + 0.5*d_i*(y_i)^2 ) !
342 ! subject to f_i(x) - a_i*z - y_i <= 0, i = 1,...,m !
343 ! xmin_j <= x_j <= xmax_j, j = 1,...,n !
344 ! z >= 0, y_i >= 0, i = 1,...,m !
345 ! -------------------------------------------------------------------!
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
353 integer :: i, ierr
354
355 call this%free()
356
357 this%n = n
358 this%m = m
359
360 call this%xold1%init(n)
361 call this%xold2%init(n)
362 this%xold1 = x
363 this%xold2 = x
364
365 call this%alpha%init(n)
366 call this%beta%init(n)
367
368 call this%a%init(m)
369 call this%c%init(m)
370 call this%d%init(m)
371 call this%low%init(n)
372 call this%upp%init(n)
373 call this%xmax%init(n)
374 call this%xmin%init(n)
375
376 !internal dummy variables for MMA
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)
381 call this%bi%init(m)
382
383 ! Necessary for KKT check after updating df0dx, fval, dfdx
384 call this%y%init(m)
385 call this%lambda%init(m)
386 call this%s%init(m)
387 call this%mu%init(m)
388 call this%xsi%init(n)
389 call this%eta%init(n)
390
391 this%a0 = a0
392 this%a%x = a
393 this%c%x = c
394 this%d%x = d
395
396 ! Set the bounds for the design variable based on the problem
397 this%xmax%x = xmax
398 this%xmin%x = xmin
399
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.)
406 end if
407
408 ! Set KKT norms to a large number for the initial design
409 this%residumax = huge(0.0_rp)
410 this%residunorm = huge(0.0_rp)
411
412 ! ------------------------------------------------------------------------ !
413 ! Assign defaults if nothing is parsed
414
415 ! Based on the Cpp Code by Niels
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))
418
419 ! Following parameters are set based on eq.3.8
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
423
424 ! Assign values from inputs when present
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
430 this%bcknd = bcknd
431 this%subsolver = subsolver
432
433 ! Sync parameters across MPI
434 call mpi_allreduce(this%n, this%n_global, 1, &
435 mpi_integer, mpi_sum, neko_comm, ierr)
436
437 call neko_log%section('MMA Parameters')
438
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)
443
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)
450
451 write(log_msg, '(A10,1X,E11.5)') 'epsimin ', this%epsimin
452 call neko_log%message(log_msg)
453
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)
462
463 call neko_log%message('Parameters a')
464 do i = 1, this%m
465 write(log_msg, '(3X,A,I2,A,E11.5)') 'a(', i, ') = ', this%a%x(i)
466 call neko_log%message(log_msg)
467 end do
468 call neko_log%message('Parameters c')
469 do i = 1, this%m
470 write(log_msg, '(3X,A,I2,A,E11.5)') 'c(', i, ') = ', this%c%x(i)
471 call neko_log%message(log_msg)
472 end do
473 call neko_log%message('Parameters d')
474 do i = 1, this%m
475 write(log_msg, '(3X,A,I2,A,E11.5)') 'd(', i, ') = ', this%d%x(i)
476 call neko_log%message(log_msg)
477 end do
478
479 call neko_log%end_section()
480
481 ! The object is correctly initialized
482 this%is_initialized = .true.
483 end subroutine mma_init_from_components
484
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
492
493 ! Select backend type
494 select case (this%bcknd)
495 case ("cpu")
496 if (neko_bcknd_device .eq. 1) then
497 call device_memcpy(x%x, x%x_d, this%n, device_to_host, &
498 sync = .false.)
499 call device_memcpy(df0dx%x, df0dx%x_d, this%n, device_to_host, &
500 sync = .false.)
501 call device_memcpy(fval%x, fval%x_d, this%m, device_to_host, &
502 sync = .false.)
503 call device_memcpy(dfdx%x, dfdx%x_d, this%m * this%n, device_to_host,&
504 sync = .true.)
505 end if
506
507 call mma_update_cpu(this, iter, x%x, df0dx%x, fval%x, dfdx%x)
508
509 if (neko_bcknd_device .eq. 1) then
510 call device_memcpy(x%x, x%x_d, this%n, host_to_device, sync = .true.)
511 end if
512
513 case ("device")
514 call mma_update_device(this, iter, x%x_d, df0dx%x_d, fval%x_d, dfdx%x_d)
515 end select
516
517 end subroutine mma_update_vector
518
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
524
525 ! Select backend type
526 select case (this%bcknd )
527 case ("cpu")
528 if (neko_bcknd_device .eq. 1) then
529 call device_memcpy(x%x, x%x_d, this%n, device_to_host, &
530 sync = .false.)
531 call device_memcpy(df0dx%x, df0dx%x_d, this%n, device_to_host, &
532 sync = .false.)
533 call device_memcpy(fval%x, fval%x_d, this%m, device_to_host, &
534 sync = .false.)
535 call device_memcpy(dfdx%x, dfdx%x_d, this%m * this%n, device_to_host,&
536 sync = .true.)
537 end if
538
539 call mma_kkt_cpu(this, x%x, df0dx%x, fval%x, dfdx%x)
540 case ("device")
541 call mma_kkt_device(this, x%x_d, df0dx%x_d, fval%x_d, dfdx%x_d)
542 end select
543 end subroutine mma_kkt_vector
544
545 ! ========================================================================== !
546 ! Getters and setters
547
549 pure function mma_get_n(this) result(n)
550 class(mma_t), intent(in) :: this
551 integer :: n
552 n = this%n
553 end function mma_get_n
554
556 pure function mma_get_m(this) result(m)
557 class(mma_t), intent(in) :: this
558 integer :: m
559 m = this%m
560 end function mma_get_m
561
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
568
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
575
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
582
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
588
589 if (neko_bcknd_cuda .eq. 1) then
590 backend = "cuda"
591 else if (neko_bcknd_hip .eq. 1) then
592 backend = "hip"
593 else if (neko_bcknd_opencl .eq. 1) then
594 backend = "opencl"
595 else
596 backend = "cpu"
597 end if
598
599 backend_subsolver = 'backend:' // trim(backend) // ', subsolver:' // &
600 trim(this%subsolver)
601 end function mma_get_backend_and_subsolver
602
603 ! ========================================================================== !
604 ! Private utilities
605
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
611
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.)
616
617 call this%low%copy_from(direction, sync = .false.)
618 call this%upp%copy_from(direction, sync = .false.)
619
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.)
625
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.)
631
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)
638
639 end subroutine mma_copy_from
640end 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:202
MMA type.
Definition mma.f90:89