29#include "../base/base_uses.f90"
35 CHARACTER(len=*),
PARAMETER,
PRIVATE :: moduleN =
'cp_lbfgs'
186 SUBROUTINE setulb(n, m, x, lower_bound, upper_bound, nbd, f, g, factr, pgtol, wa, iwa, &
187 task, iprint, csave, lsave, isave, dsave, trust_radius, spgr, iwunit)
189 INTEGER,
INTENT(in) :: n, m
190 REAL(kind=
dp),
INTENT(inout) :: x(n)
191 REAL(kind=
dp) :: lower_bound(n), upper_bound(n)
193 REAL(kind=
dp) :: f, g(n)
194 REAL(kind=
dp),
INTENT(in) :: factr, pgtol
195 REAL(kind=
dp) :: wa(2*m*n + 5*n + 11*m*m + 8*m)
197 CHARACTER(LEN=60) :: task
199 CHARACTER(LEN=60) :: csave
202 REAL(kind=
dp) :: dsave(29)
203 REAL(kind=
dp),
INTENT(in) :: trust_radius
204 TYPE(
spgr_type),
OPTIONAL,
POINTER :: spgr
205 INTEGER,
OPTIONAL :: iwunit
207 INTEGER :: i, ld, lr, lsnd, lss, lsy, lt, lwa, lwn, &
208 lws, lwt, lwy, lxp, lz, wunit
227 IF (
PRESENT(iwunit) .AND. iwunit > 0) wunit = iwunit
229 IF (task ==
'START')
THEN
237 isave(5) = isave(4) + isave(1)
239 isave(6) = isave(5) + isave(1)
241 isave(7) = isave(6) + isave(2)
243 isave(8) = isave(7) + isave(2)
245 isave(9) = isave(8) + isave(2)
247 isave(10) = isave(9) + isave(3)
249 isave(11) = isave(10) + isave(3)
251 isave(12) = isave(11) + n
253 isave(13) = isave(12) + n
255 isave(14) = isave(13) + n
257 isave(15) = isave(14) + n
259 isave(16) = isave(15) + n
279 IF (trust_radius >= 0)
THEN
281 lower_bound(i) = x(i) - trust_radius
282 upper_bound(i) = x(i) + trust_radius
288 CALL mainlb(n, m, x, lower_bound, upper_bound, nbd, f, g, factr, pgtol, &
289 wa(lws), wa(lwy), wa(lsy), wa(lss), wa(lwt), &
290 wa(lwn), wa(lsnd), wa(lz), wa(lr), wa(ld), wa(lt), wa(lxp), &
292 iwa(1), iwa(n + 1), iwa(2*n + 1), task, iprint, &
293 csave, lsave, isave(22), dsave, spgr, wunit)
399 SUBROUTINE mainlb(n, m, x, lower_bound, upper_bound, nbd, f, g, factr, pgtol, ws, wy, &
400 sy, ss, wt, wn, snd, z, r, d, t, xp, wa, &
401 index, iwhere, indx2, task, &
402 iprint, csave, lsave, isave, dsave, spgr, iwunit)
403 INTEGER,
INTENT(in) :: n, m
404 REAL(kind=
dp),
INTENT(inout) :: x(n)
405 REAL(kind=
dp),
INTENT(in) :: lower_bound(n), upper_bound(n)
407 REAL(kind=
dp) :: f, g(n), factr, pgtol, ws(n, m), wy(n, m), sy(m, m), ss(m, m), wt(m, m), &
408 wn(2*m, 2*m), snd(2*m, 2*m), z(n), r(n), d(n), t(n), xp(n), wa(8*m)
409 INTEGER :: index(n), iwhere(n), indx2(n)
410 CHARACTER(LEN=60) :: task
412 CHARACTER(LEN=60) :: csave
415 REAL(kind=
dp) :: dsave(29)
416 TYPE(
spgr_type),
OPTIONAL,
POINTER :: spgr
417 INTEGER,
OPTIONAL :: iwunit
419 REAL(kind=
dp),
PARAMETER :: one = 1.0_dp, zero = 0.0_dp
421 CHARACTER(LEN=3) :: word
422 INTEGER :: col, head, i, iback, ifun, ileave, info, &
423 itail, iter, itfile, iupdat, iword, k, &
424 nact, nenter, nfgv, nfree, nintol, &
426 LOGICAL :: boxed, constrained, first, &
427 keep_space_group, updatd, wrk, &
429 REAL(kind=
dp) :: cachyt, cpu1, cpu2, ddot, ddum, dnorm, dr, dtd, epsmch, fold, g_inf_norm, &
430 gd, gdold, lnscht, rr, sbtime, step_max, stp, theta, time, time1, time2, tol, xstep
453 IF (
PRESENT(iwunit) .AND. iwunit > 0) wunit = iwunit
455 keep_space_group = .false.
456 IF (
PRESENT(spgr))
THEN
457 IF (
ASSOCIATED(spgr)) keep_space_group = spgr%keep_space_group
460 IF (task ==
'START')
THEN
462 epsmch = epsilon(one)
513 IF (iprint >= 1)
THEN
515 CALL open_file(file_name=
'iterate.dat', unit_number=itfile, file_action=
'WRITE', file_status=
'UNKNOWN')
520 CALL errclb(n, m, factr, lower_bound, upper_bound, nbd, task, info, k)
521 IF (task(1:5) ==
'ERROR')
THEN
522 CALL prn3lb(n, x, f, task, iprint, info, itfile, &
523 iter, nfgv, nintol, nskip, nact, g_inf_norm, &
524 zero, nseg, word, iback, stp, xstep, k, &
525 cachyt, sbtime, lnscht, wunit)
529 CALL prn1lb(n, m, lower_bound, upper_bound, x, iprint, itfile, epsmch, wunit)
533 CALL active(n, lower_bound, upper_bound, nbd, x, iwhere, iprint, x_projected, constrained, boxed, wunit)
535 IF (keep_space_group)
THEN
542 CALL save_local(lsave, isave, dsave, x_projected, constrained, boxed, updatd, nintol, itfile, iback, nskip, head, col, itail, &
543 iter, iupdat, nseg, nfgv, info, ifun, iword, nfree, nact, ileave, nenter, theta, fold, tol, dnorm, epsmch, &
544 cpu1, cachyt, sbtime, lnscht, time1, gd, step_max, g_inf_norm, stp, gdold, dtd)
548 IF (keep_space_group)
THEN
555 x_projected = lsave(1)
556 constrained = lsave(2)
591 g_inf_norm = dsave(13)
599 IF (task(1:4) ==
'STOP')
THEN
600 IF (task(7:9) ==
'CPU')
THEN
602 CALL dcopy(n, t, 1, x, 1)
603 CALL dcopy(n, r, 1, g, 1)
608 CALL prn3lb(n, x, f, task, iprint, info, itfile, &
609 iter, nfgv, nintol, nskip, nact, g_inf_norm, &
610 time, nseg, word, iback, stp, xstep, k, &
611 cachyt, sbtime, lnscht, wunit)
612 CALL save_local(lsave, isave, dsave, x_projected, constrained, boxed, updatd, nintol, itfile, iback, nskip, head, col, itail, &
613 iter, iupdat, nseg, nfgv, info, ifun, iword, nfree, nact, ileave, nenter, theta, fold, tol, dnorm, epsmch, &
614 cpu1, cachyt, sbtime, lnscht, time1, gd, step_max, g_inf_norm, stp, gdold, dtd)
619 IF (.NOT. (task(1:5) ==
'FG_LN' .OR. task(1:5) ==
'NEW_X'))
THEN
626 CALL projgr(n, lower_bound, upper_bound, nbd, x, g, g_inf_norm)
628 IF (iprint >= 1)
THEN
629 WRITE (wunit, 1002) iter, f, g_inf_norm
630 WRITE (itfile, 1003) iter, nfgv, g_inf_norm, f
632 IF (g_inf_norm <= pgtol)
THEN
634 task =
'CONVERGENCE: NORM_OF_PROJECTED_GRADIENT_<=_PGTOL'
637 CALL prn3lb(n, x, f, task, iprint, info, itfile, &
638 iter, nfgv, nintol, nskip, nact, g_inf_norm, &
639 time, nseg, word, iback, stp, xstep, k, &
640 cachyt, sbtime, lnscht, wunit)
641 CALL save_local(lsave, isave, dsave, x_projected, constrained, boxed, updatd, nintol, itfile, iback, nskip, head, col, itail, &
642 iter, iupdat, nseg, nfgv, info, ifun, iword, nfree, nact, ileave, nenter, theta, fold, tol, dnorm, epsmch, &
643 cpu1, cachyt, sbtime, lnscht, time1, gd, step_max, g_inf_norm, stp, gdold, dtd)
650 IF (.NOT. first .OR. .NOT. (task(1:5) ==
'FG_LN' .OR. task(1:5) ==
'NEW_X'))
THEN
651 IF (iprint >= 99)
WRITE (wunit, 1001) iter + 1
654 IF (.NOT. constrained .AND. col > 0)
THEN
656 CALL dcopy(n, x, 1, z, 1)
668 CALL cauchy(n, x, lower_bound, upper_bound, nbd, g, indx2, iwhere, t, d, z, &
669 m, wy, ws, sy, wt, theta, col, head, &
670 wa(1), wa(2*m + 1), wa(4*m + 1), wa(6*m + 1), nseg, &
671 iprint, g_inf_norm, info, epsmch, wunit)
673 IF (keep_space_group)
THEN
678 IF (iprint >= 1)
WRITE (wunit, 1005)
686 cachyt = cachyt + cpu2 - cpu1
691 cachyt = cachyt + cpu2 - cpu1
692 nintol = nintol + nseg
697 CALL freev(n, nfree, index, nenter, ileave, indx2, &
698 iwhere, wrk, updatd, constrained, iprint, iter, wunit)
706 IF (.NOT. (nfree == 0 .OR. col == 0))
THEN
722 IF (wrk)
CALL formk(n, nfree, index, nenter, ileave, indx2, iupdat, &
723 updatd, wn, snd, m, ws, wy, sy, theta, col, head, info)
727 IF (iprint >= 1)
WRITE (wunit, 1006)
735 sbtime = sbtime + cpu2 - cpu1
742 CALL cmprlb(n, m, x, g, ws, wy, sy, wt, z, r, wa, index, &
743 theta, col, head, nfree, constrained, info)
745 IF (keep_space_group)
THEN
752 CALL subsm(n, m, nfree, index, lower_bound, upper_bound, nbd, z, r, xp, ws, wy, &
753 theta, x, g, col, head, iword, wa, wn, iprint, info, wunit)
755 IF (keep_space_group)
THEN
763 IF (iprint >= 1)
WRITE (wunit, 1005)
771 sbtime = sbtime + cpu2 - cpu1
777 sbtime = sbtime + cpu2 - cpu1
788 IF (keep_space_group)
THEN
797 IF (.NOT. first .OR. .NOT. (task(1:5) ==
'NEW_X'))
THEN
799 IF (keep_space_group)
THEN
806 CALL lnsrlb(n, lower_bound, upper_bound, nbd, x, f, fold, gd, gdold, g, d, r, t, z, stp, dnorm, &
807 dtd, xstep, step_max, iter, ifun, iback, nfgv, info, task, &
808 boxed, constrained, csave, isave(22), dsave(17), wunit)
810 IF (keep_space_group)
THEN
814 IF (info /= 0 .OR. iback >= 20)
THEN
816 CALL dcopy(n, t, 1, x, 1)
817 CALL dcopy(n, r, 1, g, 1)
828 task =
'ABNORMAL_TERMINATION_IN_LNSRCH'
832 CALL prn3lb(n, x, f, task, iprint, info, itfile, &
833 iter, nfgv, nintol, nskip, nact, g_inf_norm, &
834 time, nseg, word, iback, stp, xstep, k, &
835 cachyt, sbtime, lnscht, wunit)
836 CALL save_local(lsave, isave, dsave, x_projected, constrained, boxed, updatd, nintol, itfile, iback, nskip, head, col, itail, &
837 iter, iupdat, nseg, nfgv, info, ifun, iword, nfree, nact, ileave, nenter, theta, fold, tol, dnorm, epsmch, &
838 cpu1, cachyt, sbtime, lnscht, time1, gd, step_max, g_inf_norm, stp, gdold, dtd)
842 IF (iprint >= 1)
WRITE (wunit, 1008)
843 IF (info == 0) nfgv = nfgv - 1
850 task =
'RESTART_FROM_LNSRCH'
852 lnscht = lnscht + cpu2 - cpu1
856 ELSE IF (task(1:5) ==
'FG_LN')
THEN
858 CALL save_local(lsave, isave, dsave, x_projected, constrained, boxed, updatd, nintol, itfile, iback, nskip, head, col, itail, &
859 iter, iupdat, nseg, nfgv, info, ifun, iword, nfree, nact, ileave, nenter, theta, fold, tol, dnorm, epsmch, &
860 cpu1, cachyt, sbtime, lnscht, time1, gd, step_max, g_inf_norm, stp, gdold, dtd)
865 lnscht = lnscht + cpu2 - cpu1
870 CALL projgr(n, lower_bound, upper_bound, nbd, x, g, g_inf_norm)
874 CALL prn2lb(n, x, f, g, iprint, itfile, iter, nfgv, nact, &
875 g_inf_norm, nseg, word, iword, iback, stp, xstep, wunit)
876 CALL save_local(lsave, isave, dsave, x_projected, constrained, boxed, updatd, nintol, itfile, iback, nskip, head, col, itail, &
877 iter, iupdat, nseg, nfgv, info, ifun, iword, nfree, nact, ileave, nenter, theta, fold, tol, dnorm, epsmch, &
878 cpu1, cachyt, sbtime, lnscht, time1, gd, step_max, g_inf_norm, stp, gdold, dtd)
885 IF (g_inf_norm <= pgtol)
THEN
887 task =
'CONVERGENCE: NORM_OF_PROJECTED_GRADIENT_<=_PGTOL'
890 CALL prn3lb(n, x, f, task, iprint, info, itfile, &
891 iter, nfgv, nintol, nskip, nact, g_inf_norm, &
892 time, nseg, word, iback, stp, xstep, k, &
893 cachyt, sbtime, lnscht, wunit)
894 CALL save_local(lsave, isave, dsave, x_projected, constrained, boxed, updatd, nintol, itfile, iback, nskip, head, col, itail, &
895 iter, iupdat, nseg, nfgv, info, ifun, iword, nfree, nact, ileave, nenter, theta, fold, tol, dnorm, epsmch, &
896 cpu1, cachyt, sbtime, lnscht, time1, gd, step_max, g_inf_norm, stp, gdold, dtd)
900 ddum = max(abs(fold), abs(f), one)
901 IF ((fold - f) <= tol*ddum)
THEN
903 task =
'CONVERGENCE: REL_REDUCTION_OF_F_<=_FACTR*EPSMCH'
904 IF (iback >= 10) info = -5
908 CALL prn3lb(n, x, f, task, iprint, info, itfile, &
909 iter, nfgv, nintol, nskip, nact, g_inf_norm, &
910 time, nseg, word, iback, stp, xstep, k, &
911 cachyt, sbtime, lnscht, wunit)
912 CALL save_local(lsave, isave, dsave, x_projected, constrained, boxed, updatd, nintol, itfile, iback, nskip, head, col, itail, &
913 iter, iupdat, nseg, nfgv, info, ifun, iword, nfree, nact, ileave, nenter, theta, fold, tol, dnorm, epsmch, &
914 cpu1, cachyt, sbtime, lnscht, time1, gd, step_max, g_inf_norm, stp, gdold, dtd)
919 IF (keep_space_group)
THEN
926 rr = ddot(n, r, 1, r, 1)
931 dr = (gd - gdold)*stp
932 CALL dscal(n, stp, d, 1)
936 IF (dr <= epsmch*ddum)
THEN
940 IF (iprint >= 1)
WRITE (wunit, 1004) dr, ddum
956 CALL matupd(n, m, ws, wy, sy, ss, d, r, itail, &
957 iupdat, col, head, theta, rr, dr, stp, dtd)
964 CALL formt(m, wt, sy, ss, col, theta, info)
969 IF (iprint >= 1)
WRITE (wunit, 1007)
9861001
FORMAT(//,
' L-BFGS| ITERATION ', i5)
988 (/,
' L-BFGS| At iterate', i5, 4x,
'f= ', 1p, d12.5, 4x,
'|proj g|= ', 1p, d12.5)
9891003
FORMAT(2(1x, i4), 5x,
'-', 5x,
'-', 3x,
'-', 5x,
'-', 5x,
'-', 8x,
'-', 3x, &
9911004
FORMAT(
' L-BFGS| ys=', 1p, e10.3,
' -gs=', 1p, e10.3,
' BFGS update SKIPPED')
993 ' L-BFGS| Singular triangular system detected;', /, &
994 ' L-BFGS| refresh the lbfgs memory and restart the iteration.')
996 ' L-BFGS| Nonpositive definiteness in Cholesky factorization in formk;', /, &
997 ' L-BFGS| refresh the lbfgs memory and restart the iteration.')
999 ' L-BFGS| Nonpositive definiteness in Cholesky factorization in formt;', /, &
1000 ' L-BFGS| refresh the lbfgs memory and restart the iteration.')
1002 ' L-BFGS| Bad direction in the line search;', /, &
1003 ' L-BFGS| refresh the lbfgs memory and restart the iteration.')
1007 END SUBROUTINE mainlb
1033 SUBROUTINE active(n, lower_bound, upper_bound, nbd, x, iwhere, iprint, &
1034 x_projected, constrained, boxed, iwunit)
1036 INTEGER,
INTENT(in) :: n
1037 REAL(kind=
dp),
INTENT(in) :: lower_bound(n), upper_bound(n)
1039 REAL(kind=
dp) :: x(n)
1040 INTEGER,
INTENT(out) :: iwhere(n)
1042 LOGICAL :: x_projected, constrained, boxed
1043 INTEGER,
OPTIONAL :: iwunit
1045 REAL(kind=
dp),
PARAMETER :: zero = 0.0_dp
1047 INTEGER :: i, nbdd, wunit
1050 IF (
PRESENT(iwunit) .AND. iwunit > 0) wunit = iwunit
1056 x_projected = .false.
1057 constrained = .false.
1063 IF (nbd(i) > 0)
THEN
1064 IF (nbd(i) <= 2 .AND. x(i) <= lower_bound(i))
THEN
1065 IF (x(i) < lower_bound(i))
THEN
1066 x_projected = .true.
1067 x(i) = lower_bound(i)
1070 ELSE IF (nbd(i) >= 2 .AND. x(i) >= upper_bound(i))
THEN
1071 IF (x(i) > upper_bound(i))
THEN
1072 x_projected = .true.
1073 x(i) = upper_bound(i)
1083 IF (nbd(i) /= 2) boxed = .false.
1084 IF (nbd(i) == 0)
THEN
1090 constrained = .true.
1091 IF (nbd(i) == 2 .AND. upper_bound(i) - lower_bound(i) <= zero)
THEN
1100 IF (iprint >= 0)
THEN
1101 IF (x_projected)
WRITE (wunit, 2001)
1102 IF (.NOT. constrained)
WRITE (wunit, 3001)
1105 IF (iprint > 0)
WRITE (wunit, 1001) nbdd
11071001
FORMAT(/,
' L-BFGS| At X0 ', i9,
' variables are exactly at the bounds')
11082001
FORMAT(
' L-BFGS| The initial X is infeasible. Restart with its projection.')
11093001
FORMAT(
' L-BFGS| This problem is unconstrained.')
1113 END SUBROUTINE active
1137 SUBROUTINE bmv(m, sy, wt, col, v, p, info)
1140 REAL(kind=
dp) :: sy(m, m), wt(m, m)
1142 REAL(kind=
dp),
INTENT(in) :: v(2*col)
1143 REAL(kind=
dp),
INTENT(out) :: p(2*col)
1144 INTEGER,
INTENT(out) :: info
1147 REAL(kind=
dp) :: sum
1149 IF (col == 0)
RETURN
1155 p(col + 1) = v(col + 1)
1160 sum = sum + sy(i, k)*v(k)/sy(k, k)
1165 CALL dtrsl(wt, m, col, p(col + 1), 11, info)
1166 IF (info /= 0)
RETURN
1170 p(i) = v(i)/sqrt(sy(i, i))
1177 CALL dtrsl(wt, m, col, p(col + 1), 01, info)
1178 IF (info /= 0)
RETURN
1183 p(i) = -p(i)/sqrt(sy(i, i))
1188 sum = sum + sy(k, i)*p(col + k)/sy(i, i)
1278 SUBROUTINE cauchy(n, x, lower_bound, upper_bound, nbd, g, iorder, iwhere, t, d, xcp, &
1279 m, wy, ws, sy, wt, theta, col, head, p, c, wbp, &
1280 v, nseg, iprint, g_inf_norm, info, epsmch, iwunit)
1281 INTEGER,
INTENT(in) :: n
1282 REAL(kind=
dp),
INTENT(in) :: x(n), lower_bound(n), upper_bound(n)
1283 INTEGER,
INTENT(in) :: nbd(n)
1284 REAL(kind=
dp),
INTENT(in) :: g(n)
1285 INTEGER :: iorder(n)
1286 INTEGER,
INTENT(inout) :: iwhere(n)
1287 REAL(kind=
dp) :: t(n), d(n), xcp(n)
1288 INTEGER,
INTENT(in) :: m
1289 REAL(kind=
dp),
INTENT(in) :: sy(m, m), wt(m, m), theta
1290 INTEGER,
INTENT(in) :: col
1291 REAL(kind=
dp),
INTENT(in) :: ws(n, col), wy(n, col)
1292 INTEGER,
INTENT(in) :: head
1293 REAL(kind=
dp) :: p(2*m), c(2*m), wbp(2*m), v(2*m)
1294 INTEGER :: nseg, iprint
1295 REAL(kind=
dp),
INTENT(in) :: g_inf_norm
1296 INTEGER,
INTENT(inout) :: info
1297 REAL(kind=
dp) :: epsmch
1298 INTEGER,
OPTIONAL :: iwunit
1300 REAL(kind=
dp),
PARAMETER :: one = 1.0_dp, zero = 0.0_dp
1302 INTEGER :: col2, i, ibkmin, ibp, iter, j, nbreak, &
1303 nfree, nleft, pointr, wunit
1304 LOGICAL :: bnded, xlower, xupper
1305 REAL(kind=
dp) :: bkmin, ddot, dibp, dibp2, dt, dtm, f1, &
1306 f2, f2_org, neggi, tj, tj0, tl, tsum, &
1307 tu, wmc, wmp, wmw, zibp
1310 IF (
PRESENT(iwunit) .AND. iwunit > 0) wunit = iwunit
1331 IF (g_inf_norm <= zero)
THEN
1332 IF (iprint >= 0)
WRITE (wunit, 7010)
1333 CALL dcopy(n, x, 1, xcp, 1)
1343 IF (iprint >= 99)
WRITE (wunit, 3010)
1357 IF (iwhere(i) /= 3 .AND. iwhere(i) /= -1)
THEN
1360 IF (nbd(i) <= 2) tl = x(i) - lower_bound(i)
1361 IF (nbd(i) >= 2) tu = upper_bound(i) - x(i)
1365 xlower = nbd(i) <= 2 .AND. tl <= zero
1366 xupper = nbd(i) >= 2 .AND. tu <= zero
1371 IF (neggi <= zero) iwhere(i) = 1
1372 ELSE IF (xupper)
THEN
1373 IF (neggi >= zero) iwhere(i) = 2
1375 IF (abs(neggi) <= zero) iwhere(i) = -3
1379 IF (iwhere(i) /= 0 .AND. iwhere(i) /= -1)
THEN
1383 f1 = f1 - neggi*neggi
1386 p(j) = p(j) + wy(i, pointr)*neggi
1387 p(col + j) = p(col + j) + ws(i, pointr)*neggi
1388 pointr = mod(pointr, m) + 1
1390 IF (nbd(i) <= 2 .AND. nbd(i) /= 0 &
1391 & .AND. neggi < zero)
THEN
1395 t(nbreak) = tl/(-neggi)
1396 IF (nbreak == 1 .OR. t(nbreak) < bkmin)
THEN
1400 ELSE IF (nbd(i) >= 2 .AND. neggi > zero)
THEN
1404 t(nbreak) = tu/neggi
1405 IF (nbreak == 1 .OR. t(nbreak) < bkmin)
THEN
1413 IF (abs(neggi) > zero) bnded = .false.
1422 IF (theta /= one)
THEN
1424 CALL dscal(col, theta, p(col + 1), 1)
1429 CALL dcopy(n, x, 1, xcp, 1)
1431 IF (nbreak == 0 .AND. nfree == n + 1)
THEN
1433 IF (iprint > 100)
WRITE (wunit, 1010) (xcp(i), i=1, n)
1448 CALL bmv(m, sy, wt, col, p, v, info)
1449 IF (info /= 0)
RETURN
1450 f2 = f2 - ddot(col2, v, 1, p, 1)
1456 WRITE (wunit, 1011) nbreak
1465 IF (nleft == 0)
THEN
1466 IF (iprint >= 99)
THEN
1468 WRITE (wunit, 4010) nseg, f1, f2
1469 WRITE (wunit, 6010) dtm
1471 IF (dtm <= zero) dtm = zero
1477 CALL daxpy(n, tsum, d, 1, xcp, 1)
1480 DO WHILE (nleft > 0)
1491 ibp = iorder(ibkmin)
1496 IF (ibkmin /= nbreak)
THEN
1497 t(ibkmin) = t(nbreak)
1498 iorder(ibkmin) = iorder(nbreak)
1503 CALL hpsolb(nleft, t, iorder, iter - 2)
1510 IF (dt /= zero .AND. iprint >= 100)
THEN
1511 WRITE (wunit, 4011) nseg, f1, f2
1512 WRITE (wunit, 5010) dt
1513 WRITE (wunit, 6010) dtm
1519 IF (iprint >= 99)
THEN
1521 WRITE (wunit, 4010) nseg, f1, f2
1522 WRITE (wunit, 6010) dtm
1524 IF (dtm <= zero) dtm = zero
1530 CALL daxpy(n, tsum, d, 1, xcp, 1)
1542 IF (dibp > zero)
THEN
1543 zibp = upper_bound(ibp) - x(ibp)
1544 xcp(ibp) = upper_bound(ibp)
1547 zibp = lower_bound(ibp) - x(ibp)
1548 xcp(ibp) = lower_bound(ibp)
1551 IF (iprint >= 100)
WRITE (wunit, 8010) ibp
1552 IF (nleft == 0 .AND. nbreak == n)
THEN
1567 f1 = f1 + dt*f2 + dibp2 - theta*dibp*zibp
1568 f2 = f2 - theta*dibp2
1572 CALL daxpy(col2, dt, p, 1, c, 1)
1578 wbp(j) = wy(ibp, pointr)
1579 wbp(col + j) = theta*ws(ibp, pointr)
1580 pointr = mod(pointr, m) + 1
1584 CALL bmv(m, sy, wt, col, wbp, v, info)
1585 IF (info /= 0)
RETURN
1586 wmc = ddot(col2, c, 1, v, 1)
1587 wmp = ddot(col2, p, 1, v, 1)
1588 wmw = ddot(col2, wbp, 1, v, 1)
1591 CALL daxpy(col2, -dibp, wbp, 1, p, 1)
1595 f2 = f2 + 2.0_dp*dibp*wmp - dibp2*wmw
1598 f2 = max(epsmch*f2_org, f2)
1611 IF (iprint >= 99)
THEN
1613 WRITE (wunit, 4010) nseg, f1, f2
1614 WRITE (wunit, 6010) dtm
1616 IF (dtm <= zero) dtm = zero
1622 CALL daxpy(n, tsum, d, 1, xcp, 1)
1630 IF (col > 0)
CALL daxpy(col2, dtm, p, 1, c, 1)
1631 IF (iprint > 100)
WRITE (wunit, 1010) (xcp(i), i=1, n)
1632 IF (iprint >= 99)
WRITE (wunit, 2010)
16341010
FORMAT(
' L-BFGS| Cauchy X = ', /, (4x, 1p, 6(1x, d11.4)))
16351011
FORMAT(/,
' L-BFGS| There are ', i12,
' breakpoints ')
16362010
FORMAT(/,
' L-BFGS| ---------------- exit CAUCHY-----------------')
16373010
FORMAT(/,
' L-BFGS| ---------------- enter CAUCHY ---------------')
16384010
FORMAT(
' L-BFGS| Piece ', i3,
' --f1, f2 at start point ', 1p, 2(1x, d11.4))
16394011
FORMAT(/,
' L-BFGS| Piece ', i3,
' --f1, f2 at start point ', &
16414012
FORMAT(/,
' L-BFGS| GCP found in this segment')
16425010
FORMAT(
' L-BFGS| Distance to the next break point = ', 1p, d11.4)
16436010
FORMAT(
' L-BFGS| Distance to the stationary point = ', 1p, d11.4)
16447010
FORMAT(
' L-BFGS| Subgnorm = 0. GCP = X.')
16458010
FORMAT(
' L-BFGS| Variable ', i12,
' is fixed.')
1649 END SUBROUTINE cauchy
1679 SUBROUTINE cmprlb(n, m, x, g, ws, wy, sy, wt, z, r, wa, index, &
1680 theta, col, head, nfree, constrained, info)
1682 INTEGER,
INTENT(in) :: n, m
1683 REAL(kind=
dp),
INTENT(in) :: x(n), g(n), ws(n, m), wy(n, m), &
1684 sy(m, m), wt(m, m), z(n)
1685 REAL(kind=
dp),
INTENT(out) :: r(n), wa(4*m)
1686 INTEGER,
INTENT(in) :: index(n)
1687 REAL(kind=
dp),
INTENT(in) :: theta
1688 INTEGER,
INTENT(in) :: col, head, nfree
1689 LOGICAL,
INTENT(in) :: constrained
1692 INTEGER :: i, j, k, pointr
1693 REAL(kind=
dp) :: a1, a2
1695 IF (.NOT. constrained .AND. col > 0)
THEN
1702 r(i) = -theta*(z(k) - x(k)) - g(k)
1704 CALL bmv(m, sy, wt, col, wa(2*m + 1), wa(1), info)
1712 a2 = theta*wa(col + j)
1715 r(i) = r(i) + wy(k, pointr)*a1 + ws(k, pointr)*a2
1717 pointr = mod(pointr, m) + 1
1723 END SUBROUTINE cmprlb
1743 SUBROUTINE errclb(n, m, factr, lower_bound, upper_bound, nbd, task, info, k)
1745 INTEGER,
INTENT(in) :: n, m
1746 REAL(kind=
dp),
INTENT(in) :: factr, lower_bound(n), upper_bound(n)
1748 CHARACTER(LEN=60) :: task
1751 REAL(kind=
dp),
PARAMETER :: zero = 0.0_dp
1757 IF (n <= 0) task =
'ERROR: N <= 0'
1758 IF (m <= 0) task =
'ERROR: M <= 0'
1759 IF (factr < zero) task =
'ERROR: FACTR < 0'
1764 IF (nbd(i) < 0 .OR. nbd(i) > 3)
THEN
1766 task =
'ERROR: INVALID NBD'
1770 IF (nbd(i) == 2)
THEN
1771 IF (lower_bound(i) > upper_bound(i))
THEN
1773 task =
'ERROR: NO FEASIBLE SOLUTION'
1782 END SUBROUTINE errclb
1831 SUBROUTINE formk(n, nsub, ind, nenter, ileave, indx2, iupdat, &
1832 updatd, wn, wn1, m, ws, wy, sy, theta, col, &
1835 INTEGER,
INTENT(in) :: n, nsub, ind(n), nenter, ileave, &
1838 INTEGER,
INTENT(in) :: m
1839 REAL(kind=
dp) :: wn1(2*m, 2*m)
1840 REAL(kind=
dp),
INTENT(out) :: wn(2*m, 2*m)
1841 REAL(kind=
dp),
INTENT(in) :: ws(n, m), wy(n, m), sy(m, m), theta
1842 INTEGER,
INTENT(in) :: col, head
1843 INTEGER,
INTENT(out) :: info
1845 REAL(kind=
dp),
PARAMETER :: zero = 0.0_dp
1847 INTEGER :: col2, dbegin, dend, i, ipntr, is, is1, &
1848 iy, jpntr, js, js1, jy, k, k1, m2, &
1850 REAL(kind=
dp) :: ddot, temp1, temp2, temp3, temp4
1873 IF (iupdat > m)
THEN
1877 CALL dcopy(m - jy, wn1(jy + 1, jy + 1), 1, wn1(jy, jy), 1)
1878 CALL dcopy(m - jy, wn1(js + 1, js + 1), 1, wn1(js, js), 1)
1879 CALL dcopy(m - 1, wn1(m + 2, jy + 1), 1, wn1(m + 1, jy), 1)
1890 ipntr = head + col - 1
1891 IF (ipntr > m) ipntr = ipntr - m
1901 temp1 = temp1 + wy(k1, ipntr)*wy(k1, jpntr)
1906 temp2 = temp2 + ws(k1, ipntr)*ws(k1, jpntr)
1907 temp3 = temp3 + ws(k1, ipntr)*wy(k1, jpntr)
1912 jpntr = mod(jpntr, m) + 1
1917 jpntr = head + col - 1
1918 IF (jpntr > m) jpntr = jpntr - m
1926 temp3 = temp3 + ws(k1, ipntr)*wy(k1, jpntr)
1928 ipntr = mod(ipntr, m) + 1
1950 temp1 = temp1 + wy(k1, ipntr)*wy(k1, jpntr)
1951 temp2 = temp2 + ws(k1, ipntr)*ws(k1, jpntr)
1955 temp3 = temp3 + wy(k1, ipntr)*wy(k1, jpntr)
1956 temp4 = temp4 + ws(k1, ipntr)*ws(k1, jpntr)
1958 wn1(iy, jy) = wn1(iy, jy) + temp1 - temp3
1959 wn1(is, js) = wn1(is, js) - temp2 + temp4
1960 jpntr = mod(jpntr, m) + 1
1962 ipntr = mod(ipntr, m) + 1
1967 DO is = m + 1, m + upcl
1974 temp1 = temp1 + ws(k1, ipntr)*wy(k1, jpntr)
1978 temp3 = temp3 + ws(k1, ipntr)*wy(k1, jpntr)
1980 IF (is <= jy + m)
THEN
1981 wn1(is, jy) = wn1(is, jy) + temp1 - temp3
1983 wn1(is, jy) = wn1(is, jy) - temp1 + temp3
1985 jpntr = mod(jpntr, m) + 1
1987 ipntr = mod(ipntr, m) + 1
2000 wn(jy, iy) = wn1(iy, jy)/theta
2001 wn(js, is) = wn1(is1, js1)*theta
2004 wn(jy, is) = -wn1(is1, jy)
2007 wn(jy, is) = wn1(is1, jy)
2009 wn(iy, iy) = wn(iy, iy) + sy(iy, iy)
2017 CALL dpofa(wn, m2, col, info)
2024 DO js = col + 1, col2
2025 CALL dtrsl(wn, m2, col, wn(1, js), 11, info)
2031 DO is = col + 1, col2
2033 wn(is, js) = wn(is, js) + ddot(col, wn(1, is), 1, wn(1, js), 1)
2039 CALL dpofa(wn(col + 1, col + 1), m2, col, info)
2047 END SUBROUTINE formk
2068 SUBROUTINE formt(m, wt, sy, ss, col, theta, info)
2071 REAL(kind=
dp) :: wt(m, m), sy(m, m), ss(m, m)
2073 REAL(kind=
dp) :: theta
2076 REAL(kind=
dp),
PARAMETER :: zero = 0.0_dp
2078 INTEGER :: i, j, k, k1
2079 REAL(kind=
dp) :: ddum
2085 wt(1, j) = theta*ss(1, j)
2092 ddum = ddum + sy(i, k)*sy(j, k)/sy(k, k)
2094 wt(i, j) = ddum + theta*ss(i, j)
2101 CALL dpofa(wt, m, col, info)
2108 END SUBROUTINE formt
2143 SUBROUTINE freev(n, nfree, index, nenter, ileave, indx2, &
2144 iwhere, wrk, updatd, constrained, iprint, iter, iwunit)
2147 INTEGER,
INTENT(inout) :: index(n)
2148 INTEGER :: nenter, ileave
2149 INTEGER,
INTENT(out) :: indx2(n)
2150 INTEGER :: iwhere(n)
2151 LOGICAL :: wrk, updatd, constrained
2152 INTEGER :: iprint, iter
2153 INTEGER,
OPTIONAL :: iwunit
2155 INTEGER :: i, iact, k, wunit
2158 IF (
PRESENT(iwunit) .AND. iwunit > 0) wunit = iwunit
2162 IF (iter > 0 .AND. constrained)
THEN
2167 IF (iwhere(k) > 0)
THEN
2170 IF (iprint >= 100)
WRITE (wunit, 1030) k
2175 IF (iwhere(k) <= 0)
THEN
2178 IF (iprint >= 100)
WRITE (wunit, 2030) k
2181 IF (iprint >= 99)
WRITE (wunit, 3030) n + 1 - ileave, nenter
2183 wrk = (ileave < n + 1) .OR. (nenter > 0) .OR. updatd
2190 IF (iwhere(i) <= 0)
THEN
2198 IF (iprint >= 99)
WRITE (wunit, 4030) nfree, iter + 1
22001030
FORMAT(
' L-BFGS| Variable ', i12,
' leaves the set of free variables')
22012030
FORMAT(
' L-BFGS| Variable ', i12,
' enters the set of free variables')
22023030
FORMAT(
' L-BFGS| ', i12,
' variables leave; ', i12,
' variables enter')
22034030
FORMAT(
' L-BFGS| ', i12,
' variables are free at GCP ', i12)
2207 END SUBROUTINE freev
2229 SUBROUTINE hpsolb(n, t, iorder, iheap)
2230 INTEGER,
INTENT(in) :: n
2231 REAL(kind=
dp),
INTENT(inout) :: t(n)
2232 INTEGER,
INTENT(inout) :: iorder(n)
2233 INTEGER,
INTENT(in) :: iheap
2235 INTEGER :: i, indxin, indxou, j, k
2236 REAL(kind=
dp) :: ddum, out
2244 IF (iheap == 0)
THEN
2256 IF (ddum < t(j))
THEN
2258 iorder(i) = iorder(j)
2282 DO WHILE (j <= n - 1)
2283 IF (t(j + 1) < t(j)) j = j + 1
2284 IF (t(j) < ddum)
THEN
2286 iorder(i) = iorder(j)
2304 END SUBROUTINE hpsolb
2349 SUBROUTINE lnsrlb(n, lower_bound, upper_bound, nbd, x, f, fold, gd, gdold, g, d, r, t, &
2350 z, stp, dnorm, dtd, xstep, step_max, iter, ifun, &
2351 iback, nfgv, info, task, boxed, constrained, csave, &
2352 isave, dsave, iwunit)
2354 INTEGER,
INTENT(in) :: n
2355 REAL(kind=
dp),
INTENT(in) :: lower_bound(n), upper_bound(n)
2357 REAL(kind=
dp) :: x(n), f, fold, gd, gdold, g(n), d(n), &
2358 r(n), t(n), z(n), stp, dnorm, dtd, &
2360 INTEGER :: iter, ifun, iback, nfgv, info
2361 CHARACTER(LEN=60) :: task
2362 LOGICAL :: boxed, constrained
2363 CHARACTER(LEN=60) :: csave
2365 REAL(kind=
dp) :: dsave(13)
2366 INTEGER,
OPTIONAL :: iwunit
2368 REAL(kind=
dp),
PARAMETER :: big = 1.0e10_dp, ftol = 1.0e-3_dp, &
2369 gtol = 0.9_dp, one = 1.0_dp, &
2370 xtol = 0.1_dp, zero = 0.0_dp
2373 REAL(kind=
dp) :: a1, a2, ddot
2376 IF (
PRESENT(iwunit) .AND. iwunit > 0) wunit = iwunit
2378 IF (.NOT. (task(1:5) ==
'FG_LN'))
THEN
2380 dtd = ddot(n, d, 1, d, 1)
2386 IF (constrained)
THEN
2392 IF (nbd(i) /= 0)
THEN
2393 IF (a1 < zero .AND. nbd(i) <= 2)
THEN
2394 a2 = lower_bound(i) - x(i)
2395 IF (a2 >= zero)
THEN
2397 ELSE IF (a1*step_max < a2)
THEN
2400 ELSE IF (a1 > zero .AND. nbd(i) >= 2)
THEN
2401 a2 = upper_bound(i) - x(i)
2402 IF (a2 <= zero)
THEN
2404 ELSE IF (a1*step_max > a2)
THEN
2413 IF (iter == 0 .AND. .NOT. boxed)
THEN
2414 stp = min(one/dnorm, step_max)
2419 CALL dcopy(n, x, 1, t, 1)
2420 CALL dcopy(n, g, 1, r, 1)
2426 gd = ddot(n, g, 1, d, 1)
2429 IF (gd >= zero)
THEN
2432 WRITE (wunit, 1020) gd
2438 CALL dcsrch(f, gd, stp, ftol, gtol, xtol, zero, step_max, csave, isave, dsave)
2441 IF (csave(1:4) /=
'CONV' .AND. csave(1:4) /=
'WARN')
THEN
2446 IF (stp == one)
THEN
2447 CALL dcopy(n, z, 1, x, 1)
2450 x(i) = stp*d(i) + t(i)
24571020
FORMAT(
' L-BFGS| ascent direction in projection gd = ', d12.5)
2461 END SUBROUTINE lnsrlb
2489 SUBROUTINE matupd(n, m, ws, wy, sy, ss, d, r, itail, &
2490 iupdat, col, head, theta, rr, dr, stp, dtd)
2493 REAL(kind=
dp) :: ws(n, m), wy(n, m), sy(m, m), ss(m, m), &
2495 INTEGER :: itail, iupdat, col, head
2496 REAL(kind=
dp) :: theta, rr, dr, stp, dtd
2498 REAL(kind=
dp),
PARAMETER :: one = 1.0_dp
2500 INTEGER :: j, pointr
2501 REAL(kind=
dp) :: ddot
2506 IF (iupdat <= m)
THEN
2508 itail = mod(head + iupdat - 2, m) + 1
2510 itail = mod(itail, m) + 1
2511 head = mod(head, m) + 1
2516 CALL dcopy(n, d, 1, ws(1, itail), 1)
2517 CALL dcopy(n, r, 1, wy(1, itail), 1)
2527 IF (iupdat > m)
THEN
2530 CALL dcopy(j, ss(2, j + 1), 1, ss(1, j), 1)
2531 CALL dcopy(col - j, sy(j + 1, j + 1), 1, sy(j, j), 1)
2538 sy(col, j) = ddot(n, d, 1, wy(1, pointr), 1)
2539 ss(j, col) = ddot(n, ws(1, pointr), 1, d, 1)
2540 pointr = mod(pointr, m) + 1
2542 IF (stp == one)
THEN
2545 ss(col, col) = stp*stp*dtd
2551 END SUBROUTINE matupd
2575 SUBROUTINE prn1lb(n, m, lower_bound, upper_bound, x, iprint, itfile, epsmch, iwunit)
2577 INTEGER,
INTENT(in) :: n, m
2578 REAL(kind=
dp),
INTENT(in) :: lower_bound(n), upper_bound(n), x(n)
2579 INTEGER :: iprint, itfile
2580 REAL(kind=
dp) :: epsmch
2581 INTEGER,
OPTIONAL :: iwunit
2586 IF (
PRESENT(iwunit) .AND. iwunit > 0) wunit = iwunit
2588 IF (iprint >= 0)
THEN
2589 WRITE (wunit, 7001) epsmch
2590 WRITE (wunit, 7002) n, m
2591 IF (iprint >= 1)
THEN
2592 WRITE (itfile, 2001) epsmch
2593 WRITE (itfile, 7003) n, m
2594 WRITE (itfile, 9001)
2595 IF (iprint > 100)
THEN
2596 WRITE (wunit, 1004)
' L-BFGS| L =', (lower_bound(i), i=1, n)
2597 WRITE (wunit, 1004)
' L-BFGS| X0 =', (x(i), i=1, n)
2598 WRITE (wunit, 1004)
' L-BFGS| U =', (upper_bound(i), i=1, n)
26031004
FORMAT(/, a13, 1p, /, (4x, 1p, 6(1x, d11.4)))
26042001
FORMAT(
'RUNNING THE L-BFGS-B CODE', /, /, &
2605 'it = iteration number', /, &
2606 'nf = number of function evaluations', /, &
2607 'nseg = number of segments explored during the Cauchy search', /, &
2608 'nact = number of active bounds at the generalized Cauchy point' &
2610 'sub = manner in which the subspace minimization terminated:' &
2611 , /,
' con = converged, bnd = a bound was reached', /, &
2612 'itls = number of iterations performed in the line search', /, &
2613 'stepl = step length used', /, &
2614 'tstep = norm of the displacement (total step)', /, &
2615 'projg = norm of the projected gradient', /, &
2616 'f = function value', /, /, &
2618 'Machine precision =', 1p, d10.3)
26197001
FORMAT(/,
' L-BFGS| RUNNING THE L-BFGS-B CODE', /, &
2620 ' L-BFGS| Machine precision =', 1p, d10.3)
26217002
FORMAT(/,
' L-BFGS| N = ', i12,
' M = ', i12)
26227003
FORMAT(
' N = ', i12,
' M = ', i12)
26239001
FORMAT(/, 3x,
'it', 3x,
'nf', 2x,
'nseg', 2x,
'nact', 2x,
'sub', 2x,
'itls', &
2624 2x,
'stepl', 4x,
'tstep', 5x,
'projg', 8x,
'f')
2628 END SUBROUTINE prn1lb
2657 SUBROUTINE prn2lb(n, x, f, g, iprint, itfile, iter, nfgv, nact, &
2658 g_inf_norm, nseg, word, iword, iback, stp, xstep, iwunit)
2660 INTEGER,
INTENT(in) :: n
2661 REAL(kind=
dp),
INTENT(in) :: x(n), f, g(n)
2662 INTEGER,
INTENT(in) :: iprint, itfile, iter, nfgv, nact
2663 REAL(kind=
dp),
INTENT(in) :: g_inf_norm
2664 INTEGER,
INTENT(in) :: nseg
2665 CHARACTER(LEN=3) :: word
2666 INTEGER :: iword, iback
2667 REAL(kind=
dp) :: stp, xstep
2668 INTEGER,
OPTIONAL :: iwunit
2670 INTEGER :: i, imod, wunit
2673 IF (
PRESENT(iwunit) .AND. iwunit > 0) wunit = iwunit
2677 IF (iword == 0)
THEN
2680 ELSE IF (iword == 1)
THEN
2683 ELSE IF (iword == 5)
THEN
2689 IF (iprint >= 99)
THEN
2690 WRITE (wunit, 2002) iback, xstep
2691 WRITE (wunit, 2001) iter, f, g_inf_norm
2692 IF (iprint > 100)
THEN
2693 WRITE (wunit, 1004)
' L-BFGS| X =', (x(i), i=1, n)
2694 WRITE (wunit, 1004)
' L-BFGS| G =', (g(i), i=1, n)
2696 ELSE IF (iprint > 0)
THEN
2697 imod = mod(iter, iprint)
2698 IF (imod == 0)
WRITE (wunit, 2001) iter, f, g_inf_norm
2700 IF (iprint >= 1)
WRITE (itfile, 3001) &
2701 iter, nfgv, nseg, nact, word, iback, stp, xstep, g_inf_norm, f
27031004
FORMAT(/, a12, 1p, /, (4x, 1p, 6(1x, d11.4)))
2705 (/,
' L-BFGS| At iterate', i5, 4x,
'f= ', 1p, d12.5, 4x,
'|proj g|= ', 1p, d12.5)
27062002
FORMAT(/,
' L-BFGS| LINE SEARCH ', i12,
' times; norm of step = ', 1p, d24.15)
27073001
FORMAT(2(1x, i4), 2(1x, i5), 2x, a3, 1x, i4, 1p, 2(2x, d7.1), 1p, 2(1x, d10.3))
2711 END SUBROUTINE prn2lb
2749 SUBROUTINE prn3lb(n, x, f, task, iprint, info, itfile, &
2750 iter, nfgv, nintol, nskip, nact, g_inf_norm, &
2751 time, nseg, word, iback, stp, xstep, k, &
2752 cachyt, sbtime, lnscht, iwunit)
2754 INTEGER,
INTENT(in) :: n
2755 REAL(kind=
dp),
INTENT(in) :: x(n), f
2756 CHARACTER(LEN=60),
INTENT(in) :: task
2757 INTEGER,
INTENT(in) :: iprint, info, itfile, iter, nfgv, &
2759 REAL(kind=
dp),
INTENT(in) :: g_inf_norm, time
2760 INTEGER,
INTENT(in) :: nseg
2761 CHARACTER(LEN=3) :: word
2763 REAL(kind=
dp) :: stp, xstep
2765 REAL(kind=
dp) :: cachyt, sbtime, lnscht
2766 INTEGER,
OPTIONAL :: iwunit
2771 IF (
PRESENT(iwunit) .AND. iwunit > 0) wunit = iwunit
2773 IF (iprint >= 0 .AND. .NOT. (task(1:5) ==
'ERROR'))
THEN
2776 WRITE (wunit, 3005) n, iter, nfgv, nintol, nskip, nact, g_inf_norm, f
2777 IF (iprint >= 100)
THEN
2778 WRITE (wunit, 1004)
' L-BFGS| X =', (x(i), i=1, n)
2780 IF (iprint >= 1)
WRITE (wunit, 3006) f
2782 IF (iprint >= 0)
THEN
2785 WRITE (wunit, 3009) task
2787 IF (info == -1)
WRITE (wunit, 9011)
2788 IF (info == -2)
WRITE (wunit, 9012)
2789 IF (info == -3)
WRITE (wunit, 9013)
2790 IF (info == -4)
WRITE (wunit, 9014)
2791 IF (info == -5)
WRITE (wunit, 9015)
2792 IF (info == -6)
WRITE (wunit, 9016) k
2793 IF (info == -7)
WRITE (wunit, 9017) k, k
2794 IF (info == -8)
WRITE (wunit, 9018)
2795 IF (info == -9)
WRITE (wunit, 9019)
2797 IF (iprint >= 1)
WRITE (wunit, 3007) cachyt, sbtime, lnscht
2798 WRITE (wunit, 3008) time
2801 IF (iprint >= 1)
THEN
2802 IF (info == -4 .OR. info == -9)
THEN
2803 WRITE (itfile, 3002) &
2804 iter, nfgv, nseg, nact, word, iback, stp, xstep
2806 WRITE (itfile, 4009) task
2808 IF (info == -1)
WRITE (itfile, 9011)
2809 IF (info == -2)
WRITE (itfile, 9012)
2810 IF (info == -3)
WRITE (itfile, 9013)
2811 IF (info == -4)
WRITE (itfile, 9014)
2812 IF (info == -5)
WRITE (itfile, 9015)
2813 IF (info == -8)
WRITE (itfile, 9018)
2814 IF (info == -9)
WRITE (itfile, 9019)
2816 WRITE (itfile, 3008) time
28201004
FORMAT(/, a12, 1p, /, (4x, 1p, 6(1x, d11.4)))
28213001
FORMAT(/,
' L-BFGS| ---------------- Information ----------------')
28223002
FORMAT(2(1x, i4), 2(1x, i5), 2x, a3, 1x, i4, 1p, 2(2x, d7.1), 6x,
'-', 10x,
'-')
2824 ' L-BFGS| * * *', /, /, &
2825 ' L-BFGS| Tit = total number of iterations', /, &
2826 ' L-BFGS| Tnf = total number of function evaluations', /, &
2827 ' L-BFGS| Tnint = total number of segments explored during', &
2828 ' L-BFGS| Cauchy searches', /, &
2829 ' L-BFGS| Skip = number of BFGS updates skipped', /, &
2830 ' L-BFGS| Nact = number of active bounds at final generalized', &
2831 ' L-BFGS| Cauchy point', /, &
2832 ' L-BFGS| Projg = norm of the final projected gradient', /, &
2833 ' L-BFGS| F = final function value', /, /, &
28353004
FORMAT(/,
' L-BFGS| ', 3x,
'N', 4x,
'Tit', 5x,
'Tnf', 2x,
'Tnint', 2x, &
2836 'Skip', 2x,
'Nact', 5x,
'Projg', 8x,
'F')
28373005
FORMAT(
' L-BFGS| ', i5, 2(1x, i6), (1x, i6), (2x, i4), (1x, i5), 1p, 2(2x, d10.3))
28383006
FORMAT(
' L-BFGS| F =', d12.5)
2840 ' L-BFGS| Cauchy time', 1p, e10.3,
' seconds.', / &
2841 ' L-BFGS| Subspace minimization time', 1p, e10.3,
' seconds.', / &
2842 ' L-BFGS| Line search time', 1p, e10.3,
' seconds.')
28433008
FORMAT(/,
' Total User time', 1p, e10.3,
' seconds.',/)
28443009
FORMAT(/,
' L-BFGS| ', a60)
2847 ' Matrix in 1st Cholesky factorization in formk is not Pos. Def.')
2849 ' Matrix in 2st Cholesky factorization in formk is not Pos. Def.')
2851 ' Matrix in the Cholesky factorization in formt is not Pos. Def.')
2853 ' Derivative >= 0, backtracking line search impossible.', /, &
2854 ' Previous x, f and g restored.', /, &
2855 ' Possible causes: 1 error in function or gradient evaluation;', /, &
2856 ' 2 rounding errors dominate computation.')
2858 ' Warning: more than 10 function and gradient', /, &
2859 ' evaluations in the last line search. Termination', /, &
2860 ' may possibly be caused by a bad search direction.')
28619016
FORMAT(
' Input nbd(', i12,
') is invalid.')
28629017
FORMAT(
' l(', i12,
') > u(', i12,
'). No feasible solution.')
28639018
FORMAT(/,
' The triangular system is singular.')
2865 ' Line search cannot locate an adequate point after 20 function', /, &
2866 ' and gradient evaluations. Previous x, f and g restored.', /, &
2867 ' Possible causes: 1 error in function or gradient evaluation;', /, &
2868 ' 2 rounding error dominate computation.')
2872 END SUBROUTINE prn3lb
2890 SUBROUTINE projgr(n, lower_bound, upper_bound, nbd, x, g, g_inf_norm)
2892 INTEGER,
INTENT(in) :: n
2893 REAL(kind=
dp),
INTENT(in) :: lower_bound(n), upper_bound(n)
2894 INTEGER,
INTENT(in) :: nbd(n)
2895 REAL(kind=
dp),
INTENT(in) :: x(n), g(n)
2896 REAL(kind=
dp) :: g_inf_norm
2898 REAL(kind=
dp),
PARAMETER :: zero = 0.0_dp
2906 IF (nbd(i) /= 0)
THEN
2908 IF (nbd(i) >= 2) gi = max((x(i) - upper_bound(i)), gi)
2910 IF (nbd(i) <= 2) gi = min((x(i) - lower_bound(i)), gi)
2913 g_inf_norm = max(g_inf_norm, abs(gi))
2918 END SUBROUTINE projgr
3027 SUBROUTINE subsm(n, m, nsub, ind, lower_bound, upper_bound, nbd, x, d, xp, ws, wy, &
3029 col, head, iword, wv, wn, iprint, info, iwunit)
3030 INTEGER,
INTENT(in) :: n, m, nsub, ind(nsub)
3031 REAL(kind=
dp),
INTENT(in) :: lower_bound(n), upper_bound(n)
3032 INTEGER,
INTENT(in) :: nbd(n)
3033 REAL(kind=
dp),
INTENT(inout) :: x(n), d(n)
3034 REAL(kind=
dp) :: xp(n)
3035 REAL(kind=
dp),
INTENT(in) :: ws(n, m), wy(n, m), theta, xx(n), gg(n)
3036 INTEGER,
INTENT(in) :: col, head
3037 INTEGER,
INTENT(out) :: iword
3038 REAL(kind=
dp) :: wv(2*m)
3039 REAL(kind=
dp),
INTENT(in) :: wn(2*m, 2*m)
3041 INTEGER,
INTENT(out) :: info
3042 INTEGER,
OPTIONAL :: iwunit
3044 REAL(kind=
dp),
PARAMETER :: one = 1.0_dp, zero = 0.0_dp
3046 INTEGER :: col2, i, ibd, j, js, jy, k, m2, pointr, &
3048 REAL(kind=
dp) :: alpha, dd_p, dk, temp1, temp2, xk
3062 IF (
PRESENT(iwunit) .AND. iwunit > 0) wunit = iwunit
3064 IF (nsub <= 0)
RETURN
3065 IF (iprint >= 99)
WRITE (wunit, 4001)
3075 temp1 = temp1 + wy(k, pointr)*d(j)
3076 temp2 = temp2 + ws(k, pointr)*d(j)
3079 wv(col + i) = theta*temp2
3080 pointr = mod(pointr, m) + 1
3087 CALL dtrsl(wn, m2, col2, wv, 11, info)
3088 IF (info /= 0)
RETURN
3092 CALL dtrsl(wn, m2, col2, wv, 01, info)
3093 IF (info /= 0)
RETURN
3102 d(i) = d(i) + wy(k, pointr)*wv(jy)/theta &
3103 & + ws(k, pointr)*wv(js)
3105 pointr = mod(pointr, m) + 1
3108 CALL dscal(nsub, one/theta, d, 1)
3115 CALL dcopy(n, x, 1, xp, 1)
3121 IF (nbd(k) /= 0)
THEN
3124 IF (nbd(k) == 1)
THEN
3125 x(k) = max(lower_bound(k), xk + dk)
3126 IF (x(k) == lower_bound(k)) iword = 1
3130 IF (nbd(k) == 2)
THEN
3131 xk = max(lower_bound(k), xk + dk)
3132 x(k) = min(upper_bound(k), xk)
3133 IF (x(k) == lower_bound(k) .OR. x(k) == upper_bound(k)) iword = 1
3137 IF (nbd(k) == 3)
THEN
3138 x(k) = min(upper_bound(k), xk + dk)
3139 IF (x(k) == upper_bound(k)) iword = 1
3150 IF (.NOT. (iword == 0))
THEN
3156 dd_p = dd_p + (x(i) - xx(i))*gg(i)
3158 IF (dd_p > zero)
THEN
3159 CALL dcopy(n, xp, 1, x, 1)
3160 IF (iprint > 0)
WRITE (wunit, 4002)
3161 IF (iprint > 0)
WRITE (wunit, 4003)
3168 IF (nbd(k) /= 0)
THEN
3169 IF (dk < zero .AND. nbd(k) <= 2)
THEN
3170 temp2 = lower_bound(k) - x(k)
3171 IF (temp2 >= zero)
THEN
3173 ELSE IF (dk*alpha < temp2)
THEN
3176 ELSE IF (dk > zero .AND. nbd(k) >= 2)
THEN
3177 temp2 = upper_bound(k) - x(k)
3178 IF (temp2 <= zero)
THEN
3180 ELSE IF (dk*alpha > temp2)
THEN
3184 IF (temp1 < alpha)
THEN
3191 IF (alpha < one)
THEN
3195 x(k) = upper_bound(k)
3197 ELSE IF (dk < zero)
THEN
3198 x(k) = lower_bound(k)
3204 x(k) = x(k) + alpha*d(i)
3209 IF (iprint >= 99)
WRITE (wunit, 4004)
32114001
FORMAT(/,
' L-BFGS| ---------------- enter SUBSM ----------------',/)
32124002
FORMAT(
' L-BFGS| Positive dir derivative in projection ')
32134003
FORMAT(
' L-BFGS| Using the backtracking step ')
32144004
FORMAT(/,
' L-BFGS| ---------------- exit SUBSM -----------------',/)
3218 END SUBROUTINE subsm
3306 SUBROUTINE dcsrch(f, g, stp, ftol, gtol, xtol, stpmin, stpmax, &
3308 REAL(kind=
dp) :: f, g
3309 REAL(kind=
dp),
INTENT(inout) :: stp
3310 REAL(kind=
dp) :: ftol, gtol, xtol, stpmin, stpmax
3311 CHARACTER(LEN=*) :: task
3313 REAL(kind=
dp) :: dsave(13)
3315 REAL(kind=
dp),
PARAMETER :: p5 = 0.5_dp, p66 = 0.66_dp, &
3316 xtrapl = 1.1_dp, xtrapu = 4.0_dp, &
3321 REAL(kind=
dp) :: finit, fm, ftest, fx, fxm, fy, fym, &
3322 ginit, gm, gtest, gx, gxm, gy, gym, &
3323 stmax, stmin, stx, sty, width, width1
3340 IF (task(1:5) ==
'START')
THEN
3344 IF (stp < stpmin) task =
'ERROR: STP < STPMIN'
3345 IF (stp > stpmax) task =
'ERROR: STP > STPMAX'
3346 IF (g >= zero) task =
'ERROR: INITIAL G >= ZERO'
3347 IF (ftol < zero) task =
'ERROR: FTOL < ZERO'
3348 IF (gtol < zero) task =
'ERROR: GTOL < ZERO'
3349 IF (xtol < zero) task =
'ERROR: XTOL < ZERO'
3350 IF (stpmin < zero) task =
'ERROR: STPMIN < ZERO'
3351 IF (stpmax < stpmin) task =
'ERROR: STPMAX < STPMIN'
3355 IF (task(1:5) ==
'ERROR')
RETURN
3364 width = stpmax - stpmin
3381 stmax = stp + xtrapu*stp
3388 IF (isave(1) == 1)
THEN
3411 ftest = finit + stp*gtest
3412 IF (stage == 1 .AND. f <= ftest .AND. g >= zero) &
3417 IF (brackt .AND. (stp <= stmin .OR. stp >= stmax)) &
3418 task =
'WARNING: ROUNDING ERRORS PREVENT PROGRESS'
3419 IF (brackt .AND. stmax - stmin <= xtol*stmax) &
3420 task =
'WARNING: XTOL TEST SATISFIED'
3421 IF (stp == stpmax .AND. f <= ftest .AND. g <= gtest) &
3422 task =
'WARNING: STP = STPMAX'
3423 IF (stp == stpmin .AND. (f > ftest .OR. g >= gtest)) &
3424 task =
'WARNING: STP = STPMIN'
3428 IF (f <= ftest .AND. abs(g) <= gtol*(-ginit)) &
3429 task =
'CONVERGENCE'
3433 IF (.NOT. (task(1:4) ==
'WARN' .OR. task(1:4) ==
'CONV'))
THEN
3439 IF (stage == 1 .AND. f <= fx .AND. f > ftest)
THEN
3444 fxm = fx - stx*gtest
3445 fym = fy - sty*gtest
3452 CALL dcstep(stx, fxm, gxm, sty, fym, gym, stp, fm, gm, &
3453 brackt, stmin, stmax)
3457 fx = fxm + stx*gtest
3458 fy = fym + sty*gtest
3466 CALL dcstep(stx, fx, gx, sty, fy, gy, stp, f, g, &
3467 brackt, stmin, stmax)
3474 IF (abs(sty - stx) >= p66*width1) stp = stx + p5*(sty - stx)
3476 width = abs(sty - stx)
3482 stmin = min(stx, sty)
3483 stmax = max(stx, sty)
3485 stmin = stp + xtrapl*(stp - stx)
3486 stmax = stp + xtrapu*(stp - stx)
3491 stp = max(stp, stpmin)
3492 stp = min(stp, stpmax)
3497 IF (brackt .AND. (stp <= stmin .OR. stp >= stmax) &
3498 .OR. (brackt .AND. stmax - stmin <= xtol*stmax)) stp = stx
3530 END SUBROUTINE dcsrch
3575 SUBROUTINE dcstep(stx, fx, dx, sty, fy, dy, stp, fp, dp_loc, brackt, &
3577 REAL(kind=
dp),
INTENT(inout) :: stx, fx, dx, sty, fy, dy, stp
3578 REAL(kind=
dp),
INTENT(in) :: fp, dp_loc
3579 LOGICAL,
INTENT(inout) :: brackt
3580 REAL(kind=
dp),
INTENT(in) :: stpmin, stpmax
3582 REAL(kind=
dp),
PARAMETER :: p66 = 0.66_dp, three = 3.0_dp, &
3583 two = 2.0_dp, zero = 0.0_dp
3585 REAL(kind=
dp) ::
gamma, p, q, r, s, sgnd, stpc, stpf, &
3599 sgnd = dp_loc*sign(1.0_dp, dx)
3607 theta = three*(fx - fp)/(stp - stx) + dx + dp_loc
3608 s = max(abs(theta), abs(dx), abs(dp_loc))
3609 gamma = s*sqrt((theta/s)**2 - (dx/s)*(dp_loc/s))
3611 p = (
gamma - dx) + theta
3614 stpc = stx + r*(stp - stx)
3615 stpq = stx + ((dx/((fx - fp)/(stp - stx) + dx))/two)* &
3617 IF (abs(stpc - stx) < abs(stpq - stx))
THEN
3620 stpf = stpc + (stpq - stpc)/two
3629 ELSE IF (sgnd < zero)
THEN
3630 theta = three*(fx - fp)/(stp - stx) + dx + dp_loc
3631 s = max(abs(theta), abs(dx), abs(dp_loc))
3632 gamma = s*sqrt((theta/s)**2 - (dx/s)*(dp_loc/s))
3634 p = (
gamma - dp_loc) + theta
3637 stpc = stp + r*(stx - stp)
3638 stpq = stp + (dp_loc/(dp_loc - dx))*(stx - stp)
3639 IF (abs(stpc - stp) > abs(stpq - stp))
THEN
3649 ELSE IF (abs(dp_loc) < abs(dx))
THEN
3656 theta = three*(fx - fp)/(stp - stx) + dx + dp_loc
3657 s = max(abs(theta), abs(dx), abs(dp_loc))
3662 gamma = s*sqrt(max(zero, (theta/s)**2 - (dx/s)*(dp_loc/s)))
3664 p = (
gamma - dp_loc) + theta
3667 IF (r < zero .AND.
gamma /= zero)
THEN
3668 stpc = stp + r*(stx - stp)
3669 ELSE IF (stp > stx)
THEN
3674 stpq = stp + (dp_loc/(dp_loc - dx))*(stx - stp)
3682 IF (abs(stpc - stp) < abs(stpq - stp))
THEN
3688 stpf = min(stp + p66*(sty - stp), stpf)
3690 stpf = max(stp + p66*(sty - stp), stpf)
3698 IF (abs(stpc - stp) > abs(stpq - stp))
THEN
3703 stpf = min(stpmax, stpf)
3704 stpf = max(stpmin, stpf)
3714 theta = three*(fp - fy)/(sty - stp) + dy + dp_loc
3715 s = max(abs(theta), abs(dy), abs(dp_loc))
3716 gamma = s*sqrt((theta/s)**2 - (dy/s)*(dp_loc/s))
3718 p = (
gamma - dp_loc) + theta
3721 stpc = stp + r*(sty - stp)
3723 ELSE IF (stp > stx)
THEN
3737 IF (sgnd < zero)
THEN
3752 END SUBROUTINE dcstep
3776 SUBROUTINE dpofa(a, lda, n, info)
3777 INTEGER,
INTENT(in) :: lda
3778 REAL(kind=
dp) :: a(lda, *)
3779 INTEGER,
INTENT(in) :: n
3782 INTEGER :: j, jm1, k
3783 REAL(kind=
dp) :: ddot, s, t
3797 IF (.NOT. (jm1 < 1))
THEN
3799 t = a(k, j) - ddot(k - 1, a(1, k), 1, a(1, j), 1)
3807 IF (s <= 0.0_dp)
EXIT
3812 END SUBROUTINE dpofa
3844 SUBROUTINE dtrsl(t, ldt, n, b, job, info)
3845 INTEGER,
INTENT(in) :: ldt
3846 REAL(kind=
dp),
INTENT(in) :: t(ldt, *)
3847 INTEGER,
INTENT(in) :: n
3848 REAL(kind=
dp),
INTENT(inout) :: b(*)
3849 INTEGER,
INTENT(in) :: job
3850 INTEGER,
INTENT(out) :: info
3852 INTEGER :: case, j, jj
3853 REAL(kind=
dp) :: ddot, temp
3865 IF (t(info, info) == 0.0_dp)
RETURN
3872 IF (mod(job, 10) /= 0)
CASE = 2
3873 IF (mod(job, 100)/10 /= 0)
CASE =
CASE + 2
3884 CALL daxpy(n - j + 1, temp, t(j, j - 1), 1, b(j), 1)
3897 CALL daxpy(j, temp, t(1, j + 1), 1, b(1), 1)
3909 b(j) = b(j) - ddot(jj - 1, t(j + 1, j), 1, b(j + 1), 1)
3918 IF (.NOT. (n < 2))
THEN
3920 b(j) = b(j) - ddot(j - 1, t(1, j), 1, b(1), 1)
3925 cpabort(
"unexpected case")
3929 END SUBROUTINE dtrsl
3939 SUBROUTINE timer(ttime)
3940 REAL(kind=
dp) :: ttime
3961 END SUBROUTINE timer
4048 SUBROUTINE save_local(lsave,isave,dsave,x_projected,constrained,boxed,updatd,nintol,itfile,iback,nskip,head,col,itail,&
4049 iter, iupdat, nseg, nfgv, info, ifun, iword, nfree, nact, ileave, nenter, theta, fold, tol, dnorm, epsmch, cpu1, &
4050 cachyt, sbtime, lnscht, time1, gd, step_max, g_inf_norm, stp, gdold, dtd)
4051 LOGICAL,
INTENT(out) :: lsave(4)
4052 INTEGER,
INTENT(out) :: isave(23)
4053 REAL(kind=
dp),
INTENT(out) :: dsave(29)
4054 LOGICAL,
INTENT(in) :: x_projected, constrained, boxed, updatd
4055 INTEGER,
INTENT(in) :: nintol, itfile, iback, nskip, head, col, &
4056 itail, iter, iupdat, nseg, nfgv, info, &
4057 ifun, iword, nfree, nact, ileave, &
4059 REAL(kind=
dp),
INTENT(in) :: theta, fold, tol, dnorm, epsmch, cpu1, &
4060 cachyt, sbtime, lnscht, time1, gd, &
4061 step_max, g_inf_norm, stp, gdold, dtd
4063 lsave(1) = x_projected
4064 lsave(2) = constrained
4098 dsave(12) = step_max
4099 dsave(13) = g_inf_norm
4104 END SUBROUTINE save_local
collects all references to literature in CP2K as new algorithms / method are included from literature...
integer, save, public byrd1995
Utility routines to open and close files. Tracking of preconnections.
subroutine, public open_file(file_name, file_status, file_form, file_action, file_position, file_pad, unit_number, debug, skip_get_unit_number, file_access)
Opens the requested file using a free unit number.
LBFGS-B routine (version 3.0, April 25, 2011)
subroutine, public setulb(n, m, x, lower_bound, upper_bound, nbd, f, g, factr, pgtol, wa, iwa, task, iprint, csave, lsave, isave, dsave, trust_radius, spgr, iwunit)
This subroutine partitions the working arrays wa and iwa, and then uses the limited memory BFGS metho...
Calculation of the incomplete Gamma function F_n(t) for multi-center integrals over Cartesian Gaussia...
Defines the basic variable types.
integer, parameter, public dp
Machine interface based on Fortran 2003 and POSIX.
integer, parameter, public default_output_unit
real(kind=dp) function, public m_walltime()
returns time from a real-time clock, protected against rolling early/easily
Space Group Symmetry Type Module (version 1.0, Ferbruary 12, 2021)
Space Group Symmetry Module (version 1.0, January 16, 2020)
subroutine, public spgr_apply_rotations_coord(spgr, coord)
routine applies the rotation matrices to the coordinates.
subroutine, public spgr_apply_rotations_force(spgr, force)
routine applies the rotation matrices to the forces.