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) :: xmin_eff, xmax_eff
247 real(kind=rp),
dimension(this%n) :: x_diff
248 real(kind=rp) :: asy_factor
250 xmin_eff = this%xmin%x
251 xmax_eff = this%xmax%x
252 if (this%move_limit .gt. 0.0_rp)
then
253 xmin_eff = max(xmin_eff, x - this%move_limit)
254 xmax_eff = min(xmax_eff, x + this%move_limit)
257 x_diff = max(xmax_eff - xmin_eff, 1.0e-5_rp)
261 associate(low => this%low%x, upp => this%upp%x, &
262 x_1 => this%xold1%x, x_2 => this%xold2%x)
264 if (iter .lt. 3)
then
266 low = x - this%asyinit * x_diff
267 upp = x + this%asyinit * x_diff
270 if ((x(j) - x_1(j)) * (x_1(j) - x_2(j)) .lt. 0.0_rp)
then
271 asy_factor = this%asydecr
272 else if ((x(j) - x_1(j)) * (x_1(j) - x_2(j)) .gt. 0.0_rp)
then
273 asy_factor = this%asyincr
278 low(j) = x(j) - asy_factor * (x_1(j) - low(j))
279 upp(j) = x(j) + asy_factor * (upp(j) - x_1(j))
285 low = max(low, x - 10.0_rp * x_diff)
286 low = min(low, x - 0.01_rp * x_diff)
288 upp = min(upp, x + 10.0_rp * x_diff)
289 upp = max(upp, x + 0.01_rp * x_diff)
303 associate(alpha => this%alpha%x, beta => this%beta%x, &
304 xmin => xmin_eff, xmax => xmax_eff, &
305 low => this%low%x, upp => this%upp%x, x => x)
307 alpha = max(xmin, low + 0.1_rp*(x - low), x - 0.5_rp*x_diff)
308 beta = min(xmax, upp - 0.1_rp*(upp - x), x + 0.5_rp*x_diff)
315 associate(p0j => this%p0j%x, q0j => this%q0j%x, &
316 pij => this%pij%x, qij => this%qij%x, &
317 low => this%low%x, upp => this%upp%x)
320 1.001_rp * max(df0dx, 0.0_rp) &
321 + 0.001_rp * max(-df0dx, 0.0_rp) &
322 + 0.00001_rp / max(x_diff, 0.00001_rp) &
326 0.001_rp * max(df0dx, 0.0_rp) &
327 + 1.001_rp * max(-df0dx, 0.0_rp) &
328 + 0.00001_rp / max(x_diff, 0.00001_rp)&
334 1.001_rp * max(dfdx(i, j), 0.0_rp) &
335 + 0.001_rp * max(-dfdx(i, j), 0.0_rp) &
336 + 0.00001_rp / max(x_diff(j), 0.00001_rp) &
337 ) * (upp(j) - x(j))**2
340 0.001_rp * max(dfdx(i, j), 0.0_rp) &
341 + 1.001_rp * max(-dfdx(i, j), 0.0_rp) &
342 + 0.00001_rp / max(x_diff(j), 0.00001_rp) &
343 ) * (x(j) - low(j))**2
352 associate(bi => this%bi%x, &
353 pij => this%pij%x, qij => this%qij%x, &
354 low => this%low%x, upp => this%upp%x)
360 + pij(i, j) / (upp(j) - x(j)) &
361 + qij(i, j) / (x(j) - low(j))
365 call mpi_allreduce(mpi_in_place, bi, this%m, &
366 mpi_real_precision, mpi_sum, neko_comm, ierr)
370 end subroutine mma_gensub_cpu
374 subroutine mma_subsolve_dpip_cpu(this, designx)
385 class(mma_t),
intent(inout) :: this
386 real(kind=rp),
dimension(this%n),
intent(inout) :: designx
389 integer :: i, j, k, iter, itto, ierr
390 real(kind=rp) :: epsi, residual_max, residual_norm, &
391 z, zeta, rez, rezeta, &
393 steg, zold, zetaold, new_residual
394 real(kind=rp),
dimension(this%m) :: y, lambda, s, mu, &
395 rey, relambda, remu, res, &
397 dy, dlambda, ds, dmu, &
398 yold, lambdaold, sold, muold
399 real(kind=rp),
dimension(this%n) :: x, xsi, eta, &
401 delx, diagx, dx, dxsi, deta, &
403 real(kind=rp),
dimension(4*this%m + 2) :: residual_small
404 real(kind=rp),
dimension(3*this%n + 4*this%m + 2) :: residual
405 real(kind=rp),
dimension(2*this%n + 4*this%m + 2) :: xx, dxx
407 real(kind=rp),
dimension(this%m, this%n) :: gg
408 real(kind=rp),
dimension(this%m+1) :: bb
409 real(kind=rp),
dimension(this%m+1, this%m+1) :: aa
410 real(kind=rp),
dimension(this%m * this%m) :: aa_buffer
415 integer,
dimension(this%m+1) :: ipiv
418 real(kind=rp) :: re_sq_norm
419 real(kind=rp) :: minimal_epsilon
426 x = 0.5_rp * (this%alpha%x + this%beta%x)
432 xsi = max(1.0_rp, 1.0_rp / (x - this%alpha%x))
433 eta = max(1.0_rp, 1.0_rp / (this%beta%x - x))
434 mu = max(1.0_rp, 0.5_rp * this%c%x)
439 minimal_epsilon = max(0.9_rp * this%epsimin, 1.0e-12_rp)
440 call mpi_allreduce(mpi_in_place, minimal_epsilon, 1, &
441 mpi_real_precision, mpi_min, neko_comm, ierr)
446 do while (epsi .gt. minimal_epsilon)
453 associate(p0j => this%p0j%x, q0j => this%q0j%x, &
454 pij => this%pij%x, qij => this%qij%x, &
455 low => this%low%x, upp => this%upp%x, &
456 alpha => this%alpha%x, beta => this%beta%x, &
457 c => this%c%x, d => this%d%x, &
458 a0 => this%a0, a => this%a%x)
460 rex = (p0j + matmul(transpose(pij), lambda)) / (upp - x)**2 &
461 - (q0j + matmul(transpose(qij), lambda)) / (x - low)**2 &
464 rey = c + d * y - lambda - mu
465 rez = a0 - zeta - dot_product(lambda, a)
471 relambda(i) = relambda(i) &
472 + pij(i, j) / (upp(j) - x(j)) &
473 + qij(i, j) / (x(j) - low(j))
483 call mpi_allreduce(mpi_in_place, relambda, this%m, &
484 mpi_real_precision, mpi_sum, neko_comm, ierr)
485 relambda = relambda - this%a%x*z - y + s - this%bi%x
487 rexsi = xsi * (x - this%alpha%x) - epsi
488 reeta = eta * (this%beta%x - x) - epsi
490 rezeta = zeta * z - epsi
491 res = lambda * s - epsi
494 residual = [rex, rey, rez, relambda, rexsi, reeta, remu, rezeta, res]
495 residual_small = [rey, rez, relambda, remu, rezeta, res]
497 residual_max = maxval(abs(residual))
498 re_sq_norm = norm2(rex)**2 + norm2(rexsi)**2 + norm2(reeta)**2
500 call mpi_allreduce(mpi_in_place, residual_max, 1, &
501 mpi_real_precision, mpi_max, neko_comm, ierr)
503 call mpi_allreduce(mpi_in_place, re_sq_norm, &
504 1, mpi_real_precision, mpi_sum, neko_comm, ierr)
506 residual_norm = sqrt(norm2(residual_small)**2 + re_sq_norm)
511 do iter = 1, this%max_iter
514 if (residual_max .lt. epsi)
exit
520 + this%pij%x(i,j) * lambda(i) / (this%upp%x(j) - x(j))**2 &
521 - this%qij%x(i,j) * lambda(i) / (x(j) - this%low%x(j))**2
526 + this%p0j%x / (this%upp%x - x)**2 &
527 - this%q0j%x / (x - this%low%x)**2 &
528 - epsi / (x - this%alpha%x) &
529 + epsi / (this%beta%x - x)
531 dely = this%c%x + this%d%x * y - lambda - epsi / y
532 delz = this%a0 - dot_product(lambda, this%a%x) - epsi / z
538 dellambda(i) = dellambda(i) &
539 + this%pij%x(i, j) / (this%upp%x(j) - x(j)) &
540 + this%qij%x(i, j) / (x(j) - this%low%x(j))
544 call mpi_allreduce(mpi_in_place, dellambda, this%m, &
545 mpi_real_precision, mpi_sum, neko_comm, ierr)
547 dellambda = dellambda - this%a%x*z - y - this%bi%x + epsi / lambda
550 gg(i,:) = this%pij%x(i,:) / (this%upp%x - x)**2 &
551 - this%qij%x(i,:) / (x - this%low%x)**2
555 (this%p0j%x + matmul(transpose(this%pij%x), lambda)) &
556 / (this%upp%x - x)**3 &
557 + (this%q0j%x + matmul(transpose(this%qij%x), lambda)) &
558 / (x - this%low%x)**3
560 diagx = 2.0_rp * diagx &
561 + xsi / (x - this%alpha%x) &
562 + eta / (this%beta%x - x)
576 bb(i) = bb(i) + gg(i, j) * (delx(j) / diagx(j))
580 call mpi_allreduce(mpi_in_place, bb, this%m, &
581 mpi_real_precision, mpi_sum, neko_comm, ierr)
583 bb(1:this%m) = dellambda + dely / (this%d%x + mu / y) - bb(1:this%m)
584 bb(this%m + 1) = delz
601 aa(i, j) = aa(i, j) &
602 + gg(i, k) * (1.0_rp / diagx(k)) * gg(j, k)
607 aa_buffer = reshape(aa(1:this%m, 1:this%m), [this%m * this%m])
609 call mpi_allreduce(mpi_in_place, aa_buffer, &
610 this%m*this%m, mpi_real_precision, mpi_sum, neko_comm, ierr)
612 aa(1:this%m, 1:this%m) = reshape(aa_buffer, [this%m, this%m])
616 aa(i, i) = aa(i, i) &
618 + 1.0_rp / (this%d%x(i) + mu(i) / y(i))
621 aa(1:this%m, this%m+1) = this%a%x
622 aa(this%m+1, 1:this%m) = this%a%x
623 aa(this%m+1, this%m+1) = - zeta/z
625 call dgesv(this%m + 1, 1, aa, this%m + 1, ipiv, bb, this%m + 1, info)
627 if (info .ne. 0)
then
628 call neko_error(
"DGESV failed to solve the linear system in " // &
629 "mma_subsolve_dpip.")
633 dlambda = bb(1:this%m)
637 dx = - delx / diagx - matmul(transpose(gg), dlambda) / diagx
638 dy = (-dely + dlambda) / (this%d%x + mu / y)
640 dxsi = -xsi + (epsi - dx * xsi) / (x - this%alpha%x)
641 deta = -eta + (epsi + dx * eta) / (this%beta%x - x)
642 dmu = -mu + (epsi - mu * dy) / y
643 dzeta = -zeta + (epsi - zeta * dz) / z
644 ds = -s + (epsi - dlambda * s) / lambda
646 dxx = [dy, dz, dlambda, dxsi, deta, dmu, dzeta, ds]
647 xx = [y, z, lambda, xsi, eta, mu, zeta, s]
649 steg = 1.0_rp / maxval([ &
651 -1.01_rp * dxx / xx, &
652 -1.01_rp * dx / (x - this%alpha%x), &
653 1.01_rp * dx / (this%beta%x - x) &
667 new_residual = 2.0_rp * residual_norm
670 call mpi_allreduce(mpi_in_place, steg, 1, &
671 mpi_real_precision, mpi_min, neko_comm, ierr)
672 call mpi_allreduce(mpi_in_place, new_residual, 1, &
673 mpi_real_precision, mpi_min, neko_comm, ierr)
678 do while ((new_residual .gt. residual_norm) .and. (itto .lt. 50))
686 lambda = lambdaold + steg*dlambda
687 xsi = xsiold + steg*dxsi
688 eta = etaold + steg*deta
689 mu = muold + steg*dmu
690 zeta = zetaold + steg*dzeta
695 rex = (this%p0j%x + matmul(transpose(this%pij%x), lambda)) &
696 / (this%upp%x - x)**2 &
697 - (this%q0j%x + matmul(transpose(this%qij%x), lambda)) &
698 / (x - this%low%x)**2 &
701 rey = this%c%x + this%d%x*y - lambda - mu
702 rez = this%a0 - zeta - dot_product(lambda, this%a%x)
708 relambda(i) = relambda(i) &
709 + this%pij%x(i, j) / (this%upp%x(j) - x(j)) &
710 + this%qij%x(i, j) / (x(j) - this%low%x(j))
714 call mpi_allreduce(mpi_in_place, relambda, this%m, &
715 mpi_real_precision, mpi_sum, neko_comm, ierr)
717 relambda = relambda - this%a%x*z - y + s - this%bi%x
719 rexsi = xsi * (x - this%alpha%x) - epsi
720 reeta = eta * (this%beta%x - x) - epsi
722 rezeta = zeta * z - epsi
723 res = lambda * s - epsi
726 re_sq_norm = norm2(rex)**2 + norm2(rexsi)**2 + norm2(reeta)**2
727 call mpi_allreduce(mpi_in_place, re_sq_norm, &
728 1, mpi_real_precision, mpi_sum, neko_comm, ierr)
730 residual_small = [rey, rez, relambda, remu, rezeta, res]
731 new_residual = sqrt(norm2(residual_small)**2 + re_sq_norm)
732 call mpi_allreduce(mpi_in_place, new_residual, 1, &
733 mpi_real_precision, mpi_min, neko_comm, ierr)
739 residual = [rex, rey, rez, relambda, rexsi, reeta, remu, rezeta, res]
742 residual_norm = new_residual
743 residual_max = maxval(abs(residual))
744 call mpi_allreduce(mpi_in_place, residual_max, 1, &
745 mpi_real_precision, mpi_max, neko_comm, ierr)
752 this%xold2%x = this%xold1%x
753 this%xold1%x = designx
759 this%lambda%x = lambda
766 end subroutine mma_subsolve_dpip_cpu
770 subroutine mma_subsolve_dip_cpu(this, designx)
808 class(mma_t),
intent(inout) :: this
809 real(kind=rp),
dimension(this%n),
intent(inout) :: designx
812 integer :: i, j, k, iter, ierr
813 real(kind=rp) :: epsi, residual_max, z, steg
814 real(kind=rp),
dimension(this%m) :: y, lambda, mu, &
815 relambda, remu, dlambda, dmu, gradlambda
816 real(kind=rp),
dimension(this%n) :: x, pjlambda, qjlambda
822 real(kind=rp),
dimension(this%n) :: ljjxinv
823 real(kind=rp),
dimension(this%m,this%n) :: hijx
824 real(kind=rp),
dimension(this%m,this%m) :: hess
825 real(kind=rp) :: hesstrace
831 integer,
dimension(this%m+1) :: ipiv
834 real(kind=rp) :: minimal_epsilon
844 lambda = max(1.0_rp, 0.5_rp * this%c%x)
852 minimal_epsilon = max(0.9_rp * this%epsimin, 1.0e-12_rp)
853 call mpi_allreduce(mpi_in_place, minimal_epsilon, 1, &
854 mpi_real_precision, mpi_min, neko_comm, ierr)
859 do while (epsi .gt. minimal_epsilon)
866 associate(p0j => this%p0j%x, q0j => this%q0j%x, &
867 pij => this%pij%x, qij => this%qij%x, &
868 low => this%low%x, upp => this%upp%x, &
869 alpha => this%alpha%x, beta => this%beta%x, &
870 c => this%c%x, d => this%d%x, &
871 a0 => this%a0, a => this%a%x, &
879 y = max(0.0_rp, lambda - c)
884 z = max(0.0_rp,dot_product(lambda, a) - a0)
890 pjlambda = (p0j + matmul(transpose(pij), lambda))
891 qjlambda = (q0j + matmul(transpose(qij), lambda))
892 x = (sqrt(pjlambda) * low + sqrt(qjlambda) * upp) / &
893 (sqrt(pjlambda) + sqrt(qjlambda))
896 x = merge(alpha, x, x .lt. alpha)
897 x = merge(beta, x, x .gt. beta)
900 relambda = matmul(pij, 1.0_rp / (upp - x)) + &
901 matmul(qij, 1.0_rp / (x - low))
904 call mpi_allreduce(mpi_in_place, relambda, this%m, &
905 mpi_real_precision, mpi_sum, neko_comm, ierr)
906 relambda = relambda - bi - y - a * z + mu
909 remu = mu * lambda - epsi
911 residual_max = maxval(abs([relambda, remu]))
912 call mpi_allreduce(mpi_in_place, residual_max, 1, &
913 mpi_real_precision, mpi_max, neko_comm, ierr)
917 do iter = 1, this%max_iter
920 if (residual_max .lt. epsi)
exit
924 gradlambda = matmul(pij, 1.0_rp / (upp - x)) + &
925 matmul(qij, 1.0_rp/(x - low))
928 call mpi_allreduce(mpi_in_place, gradlambda, this%m, &
929 mpi_real_precision, mpi_sum, neko_comm, ierr)
930 gradlambda = gradlambda - bi - y - a * z
933 gradlambda = - gradlambda - epsi / lambda
939 ljjxinv= - 1.0_rp / ( (2.0_rp * pjlambda/(upp - x)**3) + &
940 (2.0_rp * qjlambda/(x - low)**3))
944 ljjxinv = merge(0.0_rp, ljjxinv, x - alpha < neko_eps)
945 ljjxinv = merge(0.0_rp, ljjxinv, beta - x < neko_eps)
948 hijx(i,:) = pij(i,:) / (upp - x)**2 &
949 - qij(i,:) / (x - low)**2
959 hess(i, j) = hess(i, j) &
960 + hijx(i, k) * (ljjxinv(k)) * hijx(j, k)
965 call mpi_allreduce(mpi_in_place, hess, &
966 this%m * this%m, mpi_real_precision, mpi_sum, neko_comm, ierr)
971 if (dot_product(lambda, a) .gt. 0.0_rp)
then
974 hess(i, j) = hess(i, j) - a(i) * a(j)
982 if (y(i) .gt. 0.0_rp)
then
983 hess(i, i) = hess(i, i) - 1.0_rp
986 hess(i, i) = hess(i, i) - mu(i) / lambda(i)
993 hesstrace = hesstrace + hess(i, i)
996 hess(i,i) = hess(i, i) - &
997 max(-1.0e-4_rp*hesstrace/this%m, 1.0e-7_rp)
1000 call dgesv(this%m , 1, hess, this%m , ipiv, &
1001 gradlambda, this%m, info)
1003 if (info .ne. 0)
then
1004 call neko_error(
"DGESV failed to solve the linear system in " // &
1005 "mma_subsolve_dip.")
1007 dlambda = gradlambda
1010 dmu = -mu + epsi / lambda - dlambda * mu / lambda
1016 steg = merge(-1.01_rp * dlambda(i) / lambda(i), steg, &
1017 steg < -1.01_rp * dlambda(i) / lambda(i))
1018 steg = merge(-1.01_rp * dmu(i) / mu(i), &
1019 steg, steg < -1.01_rp * dmu(i) / mu(i))
1022 steg = 1.0_rp / steg
1024 lambda = lambda + steg*dlambda
1033 y = max(0.0_rp, lambda - c)
1039 z = max(0.0_rp,dot_product(lambda, a) - a0)
1045 pjlambda = (p0j + matmul(transpose(pij), lambda))
1046 qjlambda = (q0j + matmul(transpose(qij), lambda))
1047 x = (sqrt(pjlambda) * low + sqrt(qjlambda) * upp) / &
1048 (sqrt(pjlambda) + sqrt(qjlambda))
1051 x = merge(alpha, x, x .lt. alpha)
1052 x = merge(beta, x, x .gt. beta)
1055 relambda = matmul(pij, 1.0_rp / (upp - x)) + &
1056 matmul(qij, 1.0_rp / (x - low))
1058 call mpi_allreduce(mpi_in_place, relambda, this%m, &
1059 mpi_real_precision, mpi_sum, neko_comm, ierr)
1060 relambda = relambda - bi - y - a * z + mu
1063 remu = mu * lambda - epsi
1065 residual_max = maxval(abs([relambda, remu]))
1066 call mpi_allreduce(mpi_in_place, residual_max, 1, &
1067 mpi_real_precision, mpi_max, neko_comm, ierr)
1071 epsi = 0.1_rp * epsi
1075 this%xold2%x = this%xold1%x
1076 this%xold1%x = designx
1083 this%lambda%x = lambda
1085 end subroutine mma_subsolve_dip_cpu
1087end submodule mma_cpu