35 use lapack_interfaces,
only:
dgesv
36 use mpi_f08,
only: mpi_in_place, mpi_max, mpi_min
37 use comm,
only: neko_comm, pe_rank, mpi_real_precision
38 use math,
only: neko_eps
39 use profiler,
only: profiler_start_region, profiler_end_region
45 module subroutine mma_update_cpu(this, iter, x, df0dx, fval, dfdx)
53 class(mma_t),
intent(inout) :: this
54 integer,
intent(in) :: iter
55 real(kind=rp),
dimension(this%n),
intent(inout) :: x
56 real(kind=rp),
dimension(this%n),
intent(in) :: df0dx
57 real(kind=rp),
dimension(this%m),
intent(in) :: fval
58 real(kind=rp),
dimension(this%m, this%n),
intent(in) :: dfdx
60 if (.not. this%is_initialized)
then
61 call neko_error(
"The MMA object is not initialized.")
64 call profiler_start_region(
"MMA update")
67 call profiler_start_region(
"MMA gensub")
68 call mma_gensub_cpu(this, iter, x, df0dx, fval, dfdx)
69 call profiler_end_region(
"MMA gensub")
72 call profiler_start_region(
"MMA subsolve")
73 if (this%subsolver .eq.
"dip")
then
74 call mma_subsolve_dip_cpu(this, x)
76 call mma_subsolve_dpip_cpu(this, x)
78 call profiler_end_region(
"MMA subsolve")
80 call profiler_end_region(
"MMA update")
82 this%is_updated = .true.
83 end subroutine mma_update_cpu
86 module subroutine mma_kkt_cpu(this, x, df0dx, fval, dfdx)
106 class(mma_t),
intent(inout) :: this
107 real(kind=rp),
dimension(this%n),
intent(in) :: x
108 real(kind=rp),
dimension(this%n),
intent(in) :: df0dx
109 real(kind=rp),
dimension(this%m),
intent(in) :: fval
110 real(kind=rp),
dimension(this%m, this%n),
intent(in) :: dfdx
112 call profiler_start_region(
"MMA KKT computation")
114 if (this%subsolver .eq.
"dip")
then
115 call mma_dip_kkt_cpu(this, x, df0dx, fval, dfdx)
117 call mma_dpip_kkt_cpu(this, x, df0dx, fval, dfdx)
120 call profiler_end_region(
"MMA KKT computation")
121 end subroutine mma_kkt_cpu
125 module subroutine mma_dpip_kkt_cpu(this, x, df0dx, fval, dfdx)
145 class(mma_t),
intent(inout) :: this
146 real(kind=rp),
dimension(this%n),
intent(in) :: x
147 real(kind=rp),
dimension(this%n),
intent(in) :: df0dx
148 real(kind=rp),
dimension(this%m),
intent(in) :: fval
149 real(kind=rp),
dimension(this%m, this%n),
intent(in) :: dfdx
151 real(kind=rp) :: rez, rezeta
152 real(kind=rp),
dimension(this%m) :: rey, relambda, remu, res
153 real(kind=rp),
dimension(this%n) :: rex, rexsi, reeta
154 real(kind=rp),
dimension(3*this%n+4*this%m+2) :: residual
156 real(kind=rp),
dimension(4*this%m+2) :: residual_small
158 real(kind=rp) :: re_sq_norm
160 rex = df0dx + matmul(transpose(dfdx), this%lambda%x) &
161 - this%xsi%x + this%eta%x
162 rey = this%c%x + this%d%x*this%y%x - this%lambda%x - this%mu%x
163 rez = this%a0 - this%zeta - dot_product(this%lambda%x, this%a%x)
165 relambda = fval - this%a%x * this%z - this%y%x + this%s%x
166 rexsi = this%xsi%x * (x - this%xmin%x)
167 reeta = this%eta%x * (this%xmax%x - x)
168 remu = this%mu%x * this%y%x
169 rezeta = this%zeta * this%z
170 res = this%lambda%x * this%s%x
172 residual = [rex, rey, rez, relambda, rexsi, reeta, remu, rezeta, res]
173 residual_small = [rey, rez, relambda, remu, rezeta, res]
175 this%residumax = maxval(abs(residual))
176 re_sq_norm = norm2(rex)**2 + norm2(rexsi)**2 + norm2(reeta)**2
178 call mpi_allreduce(mpi_in_place, this%residumax, 1, &
179 mpi_real_precision, mpi_max, neko_comm, ierr)
181 call mpi_allreduce(mpi_in_place, re_sq_norm, 1, &
182 mpi_real_precision, mpi_sum, neko_comm, ierr)
184 this%residunorm = sqrt(norm2(residual_small)**2 + re_sq_norm)
185 end subroutine mma_dpip_kkt_cpu
189 module subroutine mma_dip_kkt_cpu(this, x, df0dx, fval, dfdx)
209 class(mma_t),
intent(inout) :: this
210 real(kind=rp),
dimension(this%n),
intent(in) :: x
211 real(kind=rp),
dimension(this%n),
intent(in) :: df0dx
212 real(kind=rp),
dimension(this%m),
intent(in) :: fval
213 real(kind=rp),
dimension(this%m, this%n),
intent(in) :: dfdx
215 real(kind=rp),
dimension(this%m) :: relambda, remu
216 real(kind=rp),
dimension(2*this%m) :: residual
219 relambda = fval - this%a%x * this%z - this%y%x + this%mu%x
221 remu = this%lambda%x * this%mu%x
223 residual = abs([relambda, remu])
224 this%residumax = maxval(residual)
225 this%residunorm = norm2(residual)
227 end subroutine mma_dip_kkt_cpu
233 subroutine mma_gensub_cpu(this, iter, x, df0dx, fval, dfdx)
239 class(mma_t),
intent(inout) :: this
240 real(kind=rp),
dimension(this%n),
intent(in) :: x
241 real(kind=rp),
dimension(this%n),
intent(in) :: df0dx
242 real(kind=rp),
dimension(this%m),
intent(in) :: fval
243 real(kind=rp),
dimension(this%m, this%n),
intent(in) :: dfdx
244 integer,
intent(in) :: iter
245 integer :: i, j, ierr
246 real(kind=rp),
dimension(this%n) :: x_diff
247 real(kind=rp) :: asy_factor
249 x_diff = this%xmax%x - this%xmin%x
253 associate(low => this%low%x, upp => this%upp%x, &
254 x_1 => this%xold1%x, x_2 => this%xold2%x)
256 if (iter .lt. 3)
then
258 low = x - this%asyinit * x_diff
259 upp = x + this%asyinit * x_diff
262 if ((x(j) - x_1(j)) * (x_1(j) - x_2(j)) .lt. 0.0_rp)
then
263 asy_factor = this%asydecr
264 else if ((x(j) - x_1(j)) * (x_1(j) - x_2(j)) .gt. 0.0_rp)
then
265 asy_factor = this%asyincr
270 low(j) = x(j) - asy_factor * (x_1(j) - low(j))
271 upp(j) = x(j) + asy_factor * (upp(j) - x_1(j))
277 low = max(low, x - 10.0_rp * x_diff)
278 low = min(low, x - 0.01_rp * x_diff)
280 upp = min(upp, x + 10.0_rp * x_diff)
281 upp = max(upp, x + 0.01_rp * x_diff)
295 associate(alpha => this%alpha%x, beta => this%beta%x, &
296 xmin => this%xmin%x, xmax => this%xmax%x, &
297 low => this%low%x, upp => this%upp%x, x => x)
299 alpha = max(xmin, low + 0.1_rp*(x - low), x - 0.5_rp*x_diff)
300 beta = min(xmax, upp - 0.1_rp*(upp - x), x + 0.5_rp*x_diff)
307 associate(p0j => this%p0j%x, q0j => this%q0j%x, &
308 pij => this%pij%x, qij => this%qij%x, &
309 low => this%low%x, upp => this%upp%x)
312 1.001_rp * max(df0dx, 0.0_rp) &
313 + 0.001_rp * max(-df0dx, 0.0_rp) &
314 + 0.00001_rp / max(x_diff, 0.00001_rp) &
318 0.001_rp * max(df0dx, 0.0_rp) &
319 + 1.001_rp * max(-df0dx, 0.0_rp) &
320 + 0.00001_rp / max(x_diff, 0.00001_rp)&
326 1.001_rp * max(dfdx(i, j), 0.0_rp) &
327 + 0.001_rp * max(-dfdx(i, j), 0.0_rp) &
328 + 0.00001_rp / max(x_diff(j), 0.00001_rp) &
329 ) * (upp(j) - x(j))**2
332 0.001_rp * max(dfdx(i, j), 0.0_rp) &
333 + 1.001_rp * max(-dfdx(i, j), 0.0_rp) &
334 + 0.00001_rp / max(x_diff(j), 0.00001_rp) &
335 ) * (x(j) - low(j))**2
344 associate(bi => this%bi%x, &
345 pij => this%pij%x, qij => this%qij%x, &
346 low => this%low%x, upp => this%upp%x)
352 + pij(i, j) / (upp(j) - x(j)) &
353 + qij(i, j) / (x(j) - low(j))
357 call mpi_allreduce(mpi_in_place, bi, this%m, &
358 mpi_real_precision, mpi_sum, neko_comm, ierr)
362 end subroutine mma_gensub_cpu
366 subroutine mma_subsolve_dpip_cpu(this, designx)
377 class(mma_t),
intent(inout) :: this
378 real(kind=rp),
dimension(this%n),
intent(inout) :: designx
381 integer :: i, j, k, iter, itto, ierr
382 real(kind=rp) :: epsi, residual_max, residual_norm, &
383 z, zeta, rez, rezeta, &
385 steg, zold, zetaold, new_residual
386 real(kind=rp),
dimension(this%m) :: y, lambda, s, mu, &
387 rey, relambda, remu, res, &
389 dy, dlambda, ds, dmu, &
390 yold, lambdaold, sold, muold
391 real(kind=rp),
dimension(this%n) :: x, xsi, eta, &
393 delx, diagx, dx, dxsi, deta, &
395 real(kind=rp),
dimension(4*this%m + 2) :: residual_small
396 real(kind=rp),
dimension(3*this%n + 4*this%m + 2) :: residual
397 real(kind=rp),
dimension(2*this%n + 4*this%m + 2) :: xx, dxx
399 real(kind=rp),
dimension(this%m, this%n) :: gg
400 real(kind=rp),
dimension(this%m+1) :: bb
401 real(kind=rp),
dimension(this%m+1, this%m+1) :: aa
402 real(kind=rp),
dimension(this%m * this%m) :: aa_buffer
407 integer,
dimension(this%m+1) :: ipiv
410 real(kind=rp) :: re_sq_norm
411 real(kind=rp) :: minimal_epsilon
418 x = 0.5_rp * (this%alpha%x + this%beta%x)
424 xsi = max(1.0_rp, 1.0_rp / (x - this%alpha%x))
425 eta = max(1.0_rp, 1.0_rp / (this%beta%x - x))
426 mu = max(1.0_rp, 0.5_rp * this%c%x)
431 minimal_epsilon = max(0.9_rp * this%epsimin, 1.0e-12_rp)
432 call mpi_allreduce(mpi_in_place, minimal_epsilon, 1, &
433 mpi_real_precision, mpi_min, neko_comm, ierr)
438 do while (epsi .gt. minimal_epsilon)
445 associate(p0j => this%p0j%x, q0j => this%q0j%x, &
446 pij => this%pij%x, qij => this%qij%x, &
447 low => this%low%x, upp => this%upp%x, &
448 alpha => this%alpha%x, beta => this%beta%x, &
449 c => this%c%x, d => this%d%x, &
450 a0 => this%a0, a => this%a%x)
452 rex = (p0j + matmul(transpose(pij), lambda)) / (upp - x)**2 &
453 - (q0j + matmul(transpose(qij), lambda)) / (x - low)**2 &
456 rey = c + d * y - lambda - mu
457 rez = a0 - zeta - dot_product(lambda, a)
463 relambda(i) = relambda(i) &
464 + pij(i, j) / (upp(j) - x(j)) &
465 + qij(i, j) / (x(j) - low(j))
475 call mpi_allreduce(mpi_in_place, relambda, this%m, &
476 mpi_real_precision, mpi_sum, neko_comm, ierr)
477 relambda = relambda - this%a%x*z - y + s - this%bi%x
479 rexsi = xsi * (x - this%alpha%x) - epsi
480 reeta = eta * (this%beta%x - x) - epsi
482 rezeta = zeta * z - epsi
483 res = lambda * s - epsi
486 residual = [rex, rey, rez, relambda, rexsi, reeta, remu, rezeta, res]
487 residual_small = [rey, rez, relambda, remu, rezeta, res]
489 residual_max = maxval(abs(residual))
490 re_sq_norm = norm2(rex)**2 + norm2(rexsi)**2 + norm2(reeta)**2
492 call mpi_allreduce(mpi_in_place, residual_max, 1, &
493 mpi_real_precision, mpi_max, neko_comm, ierr)
495 call mpi_allreduce(mpi_in_place, re_sq_norm, &
496 1, mpi_real_precision, mpi_sum, neko_comm, ierr)
498 residual_norm = sqrt(norm2(residual_small)**2 + re_sq_norm)
503 do iter = 1, this%max_iter
506 if (residual_max .lt. epsi)
exit
512 + this%pij%x(i,j) * lambda(i) / (this%upp%x(j) - x(j))**2 &
513 - this%qij%x(i,j) * lambda(i) / (x(j) - this%low%x(j))**2
518 + this%p0j%x / (this%upp%x - x)**2 &
519 - this%q0j%x / (x - this%low%x)**2 &
520 - epsi / (x - this%alpha%x) &
521 + epsi / (this%beta%x - x)
523 dely = this%c%x + this%d%x * y - lambda - epsi / y
524 delz = this%a0 - dot_product(lambda, this%a%x) - epsi / z
530 dellambda(i) = dellambda(i) &
531 + this%pij%x(i, j) / (this%upp%x(j) - x(j)) &
532 + this%qij%x(i, j) / (x(j) - this%low%x(j))
536 call mpi_allreduce(mpi_in_place, dellambda, this%m, &
537 mpi_real_precision, mpi_sum, neko_comm, ierr)
539 dellambda = dellambda - this%a%x*z - y - this%bi%x + epsi / lambda
542 gg(i,:) = this%pij%x(i,:) / (this%upp%x - x)**2 &
543 - this%qij%x(i,:) / (x - this%low%x)**2
547 (this%p0j%x + matmul(transpose(this%pij%x), lambda)) &
548 / (this%upp%x - x)**3 &
549 + (this%q0j%x + matmul(transpose(this%qij%x), lambda)) &
550 / (x - this%low%x)**3
552 diagx = 2.0_rp * diagx &
553 + xsi / (x - this%alpha%x) &
554 + eta / (this%beta%x - x)
568 bb(i) = bb(i) + gg(i, j) * (delx(j) / diagx(j))
572 call mpi_allreduce(mpi_in_place, bb, this%m, &
573 mpi_real_precision, mpi_sum, neko_comm, ierr)
575 bb(1:this%m) = dellambda + dely / (this%d%x + mu / y) - bb(1:this%m)
576 bb(this%m + 1) = delz
593 aa(i, j) = aa(i, j) &
594 + gg(i, k) * (1.0_rp / diagx(k)) * gg(j, k)
599 aa_buffer = reshape(aa(1:this%m, 1:this%m), [this%m * this%m])
601 call mpi_allreduce(mpi_in_place, aa_buffer, &
602 this%m*this%m, mpi_real_precision, mpi_sum, neko_comm, ierr)
604 aa(1:this%m, 1:this%m) = reshape(aa_buffer, [this%m, this%m])
608 aa(i, i) = aa(i, i) &
610 + 1.0_rp / (this%d%x(i) + mu(i) / y(i))
613 aa(1:this%m, this%m+1) = this%a%x
614 aa(this%m+1, 1:this%m) = this%a%x
615 aa(this%m+1, this%m+1) = - zeta/z
617 call dgesv(this%m + 1, 1, aa, this%m + 1, ipiv, bb, this%m + 1, info)
619 if (info .ne. 0)
then
620 call neko_error(
"DGESV failed to solve the linear system in " // &
621 "mma_subsolve_dpip.")
625 dlambda = bb(1:this%m)
629 dx = - delx / diagx - matmul(transpose(gg), dlambda) / diagx
630 dy = (-dely + dlambda) / (this%d%x + mu / y)
632 dxsi = -xsi + (epsi - dx * xsi) / (x - this%alpha%x)
633 deta = -eta + (epsi + dx * eta) / (this%beta%x - x)
634 dmu = -mu + (epsi - mu * dy) / y
635 dzeta = -zeta + (epsi - zeta * dz) / z
636 ds = -s + (epsi - dlambda * s) / lambda
638 dxx = [dy, dz, dlambda, dxsi, deta, dmu, dzeta, ds]
639 xx = [y, z, lambda, xsi, eta, mu, zeta, s]
641 steg = 1.0_rp / maxval([ &
643 -1.01_rp * dxx / xx, &
644 -1.01_rp * dx / (x - this%alpha%x), &
645 1.01_rp * dx / (this%beta%x - x) &
659 new_residual = 2.0_rp * residual_norm
662 call mpi_allreduce(mpi_in_place, steg, 1, &
663 mpi_real_precision, mpi_min, neko_comm, ierr)
664 call mpi_allreduce(mpi_in_place, new_residual, 1, &
665 mpi_real_precision, mpi_min, neko_comm, ierr)
670 do while ((new_residual .gt. residual_norm) .and. (itto .lt. 50))
678 lambda = lambdaold + steg*dlambda
679 xsi = xsiold + steg*dxsi
680 eta = etaold + steg*deta
681 mu = muold + steg*dmu
682 zeta = zetaold + steg*dzeta
687 rex = (this%p0j%x + matmul(transpose(this%pij%x), lambda)) &
688 / (this%upp%x - x)**2 &
689 - (this%q0j%x + matmul(transpose(this%qij%x), lambda)) &
690 / (x - this%low%x)**2 &
693 rey = this%c%x + this%d%x*y - lambda - mu
694 rez = this%a0 - zeta - dot_product(lambda, this%a%x)
700 relambda(i) = relambda(i) &
701 + this%pij%x(i, j) / (this%upp%x(j) - x(j)) &
702 + this%qij%x(i, j) / (x(j) - this%low%x(j))
706 call mpi_allreduce(mpi_in_place, relambda, this%m, &
707 mpi_real_precision, mpi_sum, neko_comm, ierr)
709 relambda = relambda - this%a%x*z - y + s - this%bi%x
711 rexsi = xsi * (x - this%alpha%x) - epsi
712 reeta = eta * (this%beta%x - x) - epsi
714 rezeta = zeta * z - epsi
715 res = lambda * s - epsi
718 re_sq_norm = norm2(rex)**2 + norm2(rexsi)**2 + norm2(reeta)**2
719 call mpi_allreduce(mpi_in_place, re_sq_norm, &
720 1, mpi_real_precision, mpi_sum, neko_comm, ierr)
722 residual_small = [rey, rez, relambda, remu, rezeta, res]
723 new_residual = sqrt(norm2(residual_small)**2 + re_sq_norm)
729 residual = [rex, rey, rez, relambda, rexsi, reeta, remu, rezeta, res]
732 residual_norm = new_residual
733 residual_max = maxval(abs(residual))
734 call mpi_allreduce(mpi_in_place, residual_max, 1, &
735 mpi_real_precision, mpi_max, neko_comm, ierr)
742 this%xold2%x = this%xold1%x
743 this%xold1%x = designx
749 this%lambda%x = lambda
756 end subroutine mma_subsolve_dpip_cpu
760 subroutine mma_subsolve_dip_cpu(this, designx)
795 class(mma_t),
intent(inout) :: this
796 real(kind=rp),
dimension(this%n),
intent(inout) :: designx
799 integer :: i, j, k, iter, ierr
800 real(kind=rp) :: epsi, residual_max, z, steg
801 real(kind=rp),
dimension(this%m) :: y, lambda, mu, &
802 relambda, remu, dlambda, dmu, gradlambda
803 real(kind=rp),
dimension(this%n) :: x, pjlambda, qjlambda
809 real(kind=rp),
dimension(this%n) :: ljjxinv
810 real(kind=rp),
dimension(this%m,this%n) :: hijx
811 real(kind=rp),
dimension(this%m,this%m) :: hess
812 real(kind=rp) :: hesstrace
818 integer,
dimension(this%m+1) :: ipiv
821 real(kind=rp) :: minimal_epsilon
831 lambda = max(1.0_rp, 0.5_rp * this%c%x)
839 minimal_epsilon = max(0.9_rp * this%epsimin, 1.0e-12_rp)
840 call mpi_allreduce(mpi_in_place, minimal_epsilon, 1, &
841 mpi_real_precision, mpi_min, neko_comm, ierr)
846 do while (epsi .gt. minimal_epsilon)
853 associate(p0j => this%p0j%x, q0j => this%q0j%x, &
854 pij => this%pij%x, qij => this%qij%x, &
855 low => this%low%x, upp => this%upp%x, &
856 alpha => this%alpha%x, beta => this%beta%x, &
857 c => this%c%x, d => this%d%x, &
858 a0 => this%a0, a => this%a%x, &
867 if (abs(d(i)) < neko_eps)
then
869 y(i) = max(0.0_rp, (lambda(i) - c(i)) / (1.0e-8_rp))
872 y(i) = max(0.0_rp, (lambda(i) - c(i)) / (d(i)))
880 z = merge(0.0_rp, 1.0_rp, a0 - dot_product(lambda, a) >= 0.0_rp)
886 pjlambda = (p0j + matmul(transpose(pij), lambda))
887 qjlambda = (q0j + matmul(transpose(qij), lambda))
888 x = (sqrt(pjlambda) * low + sqrt(qjlambda) * upp) / &
889 (sqrt(pjlambda) + sqrt(qjlambda))
892 x = merge(alpha, x, x .lt. alpha)
893 x = merge(beta, x, x .gt. beta)
896 relambda = matmul(pij, 1.0_rp / (upp - x)) + &
897 matmul(qij, 1.0_rp / (x - low))
900 call mpi_allreduce(mpi_in_place, relambda, this%m, &
901 mpi_real_precision, mpi_sum, neko_comm, ierr)
902 relambda = relambda - bi - y - a * z + mu
905 remu = mu * lambda - epsi
907 residual_max = maxval(abs([relambda, remu]))
911 do iter = 1, this%max_iter
914 if (residual_max .lt. epsi)
exit
918 gradlambda = matmul(pij, 1.0_rp / (upp - x)) + &
919 matmul(qij, 1.0_rp/(x - low))
922 call mpi_allreduce(mpi_in_place, gradlambda, this%m, &
923 mpi_real_precision, mpi_sum, neko_comm, ierr)
924 gradlambda = gradlambda - bi - y - a * z
927 gradlambda = - gradlambda - epsi / lambda
933 ljjxinv= - 1.0_rp / ( (2.0_rp * pjlambda/(upp - x)**3) + &
934 (2.0_rp * qjlambda/(x - low)**3))
938 ljjxinv = merge(0.0_rp, ljjxinv, x - alpha < neko_eps)
939 ljjxinv = merge(0.0_rp, ljjxinv, beta - x < neko_eps)
942 hijx(i,:) = pij(i,:) / (upp - x)**2 &
943 - qij(i,:) / (x - low)**2
953 hess(i, j) = hess(i, j) &
954 + hijx(i, k) * (ljjxinv(k)) * hijx(j, k)
959 call mpi_allreduce(mpi_in_place, hess, &
960 this%m*this%m, mpi_real_precision, mpi_sum, neko_comm, ierr)
971 if (y(i) .gt. 0.0_rp)
then
972 if (abs(d(i)) < neko_eps)
then
975 hess(i, i) = hess(i, i) - 1.0_rp/d(i)
979 hess(i, i) = hess(i, i) - mu(i) / lambda(i)
986 hesstrace = hesstrace + hess(i, i)
989 hess(i,i) = hess(i, i) - &
990 max(-1.0e-4_rp*hesstrace/this%m, 1.0e-7_rp)
993 call dgesv(this%m , 1, hess, this%m , ipiv, &
994 gradlambda, this%m, info)
996 if (info .ne. 0)
then
997 call neko_error(
"DGESV failed to solve the linear system in " // &
1000 dlambda = gradlambda
1003 dmu = -mu + epsi / lambda - dlambda * mu / lambda
1009 steg = merge(-1.01_rp * dlambda(i) / lambda(i), steg, &
1010 steg < -1.01_rp * dlambda(i) / lambda(i))
1011 steg = merge(-1.01_rp * dmu(i) / mu(i), &
1012 steg, steg < -1.01_rp * dmu(i) / mu(i))
1015 steg = 1.0_rp / steg
1017 lambda = lambda + steg*dlambda
1027 if (abs(d(i)) < neko_eps)
then
1029 y(i) = max(0.0_rp, (lambda(i) - c(i)) / (1.0e-8_rp))
1032 y(i) = max(0.0_rp, (lambda(i) - c(i)) / (d(i)))
1040 z = merge(0.0_rp, 1.0_rp, a0 - dot_product(lambda, a) >= 0.0_rp)
1046 pjlambda = (p0j + matmul(transpose(pij), lambda))
1047 qjlambda = (q0j + matmul(transpose(qij), lambda))
1048 x = (sqrt(pjlambda) * low + sqrt(qjlambda) * upp) / &
1049 (sqrt(pjlambda) + sqrt(qjlambda))
1052 x = merge(alpha, x, x .lt. alpha)
1053 x = merge(beta, x, x .gt. beta)
1056 relambda = matmul(pij, 1.0_rp / (upp - x)) + &
1057 matmul(qij, 1.0_rp / (x - low))
1059 call mpi_allreduce(mpi_in_place, relambda, this%m, &
1060 mpi_real_precision, mpi_sum, neko_comm, ierr)
1061 relambda = relambda - bi - y - a * z + mu
1064 remu = mu * lambda - epsi
1066 residual_max = maxval(abs([relambda, remu]))
1070 epsi = 0.1_rp * epsi
1074 this%xold2%x = this%xold1%x
1075 this%xold1%x = designx
1082 this%lambda%x = lambda
1084 end subroutine mma_subsolve_dip_cpu
1086end submodule mma_cpu