292 character*60 task, csave
294 integer n, m, iprint, nbd(n), index(n),
295 + iwhere(n), indx2(n), isave(23)
296 double precision f, factr, pgtol,
297 + x(n), l(n), u(n), g(n), z(n), r(n), d(n), t(n),
301 + ws(n, m), wy(n, m), sy(m, m), ss(m, m),
302 + wt(m, m), wn(2*m, 2*m), snd(2*m, 2*m), dsave(29)
484 logical prjctd,cnstnd,boxed,updatd,wrk
486 integer i,k,nintol,itfile,iback,nskip,
487 + head,col,iter,itail,iupdat,
488 + nseg,nfgv,info,ifun,
489 + iword,nfree,nact,ileave,nenter
490 double precision theta,fold,
ddot,dr,rr,tol,
491 + xstep,sbgnrm,ddum,dnorm,dtd,epsmch,
492 + cpu1,cpu2,cachyt,sbtime,lnscht,
time1,time2,
493 + gd,gdold,stp,stpmx,time
494 double precision one,zero
495 parameter(one=1.0d0,zero=0.0d0)
497 if (task .eq.
'START')
then 499 epsmch = epsilon(one)
550 if (iprint .ge. 1)
then 552 open (8, file =
'iterate.dat', status =
'unknown')
557 call errclb(n,m,factr,l,u,nbd,task,info,k)
558 if (task(1:5) .eq.
'ERROR')
then 559 call prn3lb(n,x,f,task,iprint,info,itfile,
560 + iter,nfgv,nintol,nskip,nact,sbgnrm,
561 + zero,nseg,word,iback,stp,xstep,k,
562 + cachyt,sbtime,lnscht)
566 call prn1lb(n,m,l,u,x,iprint,itfile,epsmch)
570 call active(n,l,u,nbd,x,iwhere,iprint,prjctd,cnstnd,boxed)
621 if (task(1:5) .eq.
'FG_LN')
goto 666
622 if (task(1:5) .eq.
'NEW_X')
goto 777
623 if (task(1:5) .eq.
'FG_ST')
goto 111
624 if (task(1:4) .eq.
'STOP')
then 625 if (task(7:9) .eq.
'CPU')
then 627 call dcopy(n,t,1,x,1)
628 call dcopy(n,r,1,g,1)
645 call projgr(n,l,u,nbd,x,g,sbgnrm)
647 if (iprint .ge. 1)
then 648 write (6,1002) iter,f,sbgnrm
649 write (itfile,1003) iter,nfgv,sbgnrm,f
651 if (sbgnrm .le. pgtol)
then 653 task =
'CONVERGENCE: NORM_OF_PROJECTED_GRADIENT_<=_PGTOL' 660 if (iprint .ge. 99)
write (6,1001) iter + 1
663 if (.not. cnstnd .and. col .gt. 0)
then 665 call dcopy(n,x,1,z,1)
678 call cauchy(n,x,l,u,nbd,g,indx2,iwhere,t,d,z,
679 + m,wy,ws,sy,wt,theta,col,head,
680 + wa(1),wa(2*m+1),wa(4*m+1),wa(6*m+1),nseg,
681 + iprint, sbgnrm, info, epsmch)
682 if (info .ne. 0)
then 684 if(iprint .ge. 1)
write (6, 1005)
692 cachyt = cachyt + cpu2 - cpu1
696 cachyt = cachyt + cpu2 - cpu1
697 nintol = nintol + nseg
702 call freev(n,nfree,index,nenter,ileave,indx2,
703 + iwhere,wrk,updatd,cnstnd,iprint,iter)
711 if (nfree .eq. 0 .or. col .eq. 0)
goto 555
727 if (wrk)
call formk(n,nfree,index,nenter,ileave,indx2,iupdat,
728 + updatd,wn,snd,m,ws,wy,sy,theta,col,head,info)
729 if (info .ne. 0)
then 732 if(iprint .ge. 1)
write (6, 1006)
740 sbtime = sbtime + cpu2 - cpu1
746 call cmprlb(n,m,x,g,ws,wy,sy,wt,z,r,wa,index,
747 + theta,col,head,nfree,cnstnd,info)
748 if (info .ne. 0)
goto 444
752 call subsm( n, m, nfree, index, l, u, nbd, z, r, xp, ws, wy,
753 + theta, x, g, col, head, iword, wa, wn, iprint, info)
755 if (info .ne. 0)
then 758 if(iprint .ge. 1)
write (6, 1005)
766 sbtime = sbtime + cpu2 - cpu1
771 sbtime = sbtime + cpu2 - cpu1
787 call lnsrlb(n,l,u,nbd,x,f,fold,gd,gdold,g,d,r,t,z,stp,dnorm,
788 + dtd,xstep,stpmx,iter,ifun,iback,nfgv,info,task,
789 + boxed,cnstnd,csave,isave(22),dsave(17))
790 if (info .ne. 0 .or. iback .ge. 20)
then 792 call dcopy(n,t,1,x,1)
793 call dcopy(n,r,1,g,1)
797 if (info .eq. 0)
then 804 task =
'ABNORMAL_TERMINATION_IN_LNSRCH' 809 if(iprint .ge. 1)
write (6, 1008)
810 if (info .eq. 0) nfgv = nfgv - 1
817 task =
'RESTART_FROM_LNSRCH' 819 lnscht = lnscht + cpu2 - cpu1
822 else if (task(1:5) .eq.
'FG_LN')
then 828 lnscht = lnscht + cpu2 - cpu1
833 call projgr(n,l,u,nbd,x,g,sbgnrm)
837 call prn2lb(n,x,f,g,iprint,itfile,iter,nfgv,nact,
838 + sbgnrm,nseg,word,iword,iback,stp,xstep)
845 if (sbgnrm .le. pgtol)
then 847 task =
'CONVERGENCE: NORM_OF_PROJECTED_GRADIENT_<=_PGTOL' 851 ddum =
max(abs(fold), abs(f), one)
852 if ((fold - f) .le. tol*ddum)
then 854 task =
'CONVERGENCE: REL_REDUCTION_OF_F_<=_FACTR*EPSMCH' 855 if (iback .ge. 10) info = -5
866 if (stp .eq. one)
then 870 dr = (gd - gdold)*stp
871 call dscal(n,stp,d,1)
875 if (dr .le. epsmch*ddum)
then 879 if (iprint .ge. 1)
write (6,1004) dr, ddum
894 call matupd(n,m,ws,wy,sy,ss,d,r,itail,
895 + iupdat,col,head,theta,rr,dr,stp,dtd)
902 call formt(m,wt,sy,ss,col,theta,info)
904 if (info .ne. 0)
then 907 if(iprint .ge. 1)
write (6, 1007)
930 call prn3lb(n,x,f,task,iprint,info,itfile,
931 + iter,nfgv,nintol,nskip,nact,sbgnrm,
932 + time,nseg,word,iback,stp,xstep,k,
933 + cachyt,sbtime,lnscht)
979 1001
format (//,
'ITERATION ',i5)
981 + (/,
'At iterate',i5,4x,
'f= ',1p,d12.5,4x,
'|proj g|= ',1p,d12.5)
982 1003
format (2(1x,i4),5x,
'-',5x,
'-',3x,
'-',5x,
'-',5x,
'-',8x,
'-',3x,
984 1004
format (
' ys=',1p,e10.3,
' -gs=',1p,e10.3,
' BFGS update SKIPPED')
986 +
' Singular triangular system detected;',/,
987 +
' refresh the lbfgs memory and restart the iteration.')
989 +
' Nonpositive definiteness in Cholesky factorization in formk;',/,
990 +
' refresh the lbfgs memory and restart the iteration.')
992 +
' Nonpositive definiteness in Cholesky factorization in formt;',/,
993 +
' refresh the lbfgs memory and restart the iteration.')
995 +
' Bad direction in the line search;',/,
996 +
' refresh the lbfgs memory and restart the iteration.')
subroutine timer(ttime)
Definition: timer.f:7
#define max(a, b)
Definition: cascade.c:32
subroutine dcopy(N, DX, INCX, DY, INCY)
Definition: dgmres.f:1381
subroutine formk(n, nsub, ind, nenter, ileave, indx2, iupdat, updatd, wn, wn1, m, ws, wy, sy, theta, col, head, info)
Definition: lbfgsb.f:1851
subroutine freev(n, nfree, index, nenter, ileave, indx2, iwhere, wrk, updatd, cnstnd, iprint, iter)
Definition: lbfgsb.f:2243
subroutine matupd(n, m, ws, wy, sy, ss, d, r, itail, iupdat, col, head, theta, rr, dr, stp, dtd)
Definition: lbfgsb.f:2586
subroutine active(n, l, u, nbd, x, iwhere, iprint, prjctd, cnstnd, boxed)
Definition: lbfgsb.f:1006
subroutine errclb(n, m, factr, l, u, nbd, task, info, k)
Definition: lbfgsb.f:1789
subroutine dscal(n, da, dx, incx)
Definition: dgesv.f:1746
subroutine formt(m, wt, sy, ss, col, theta, info)
Definition: lbfgsb.f:2174
subroutine cmprlb(n, m, x, g, ws, wy, sy, wt, z, r, wa, index, theta, col, head, nfree, cnstnd, info)
Definition: lbfgsb.f:1722
subroutine lnsrlb(n, l, u, nbd, x, f, fold, gd, gdold, g, d, r, t, z, stp, dnorm, dtd, xstep, stpmx, iter, ifun, iback, nfgv, info, task, boxed, cnstnd, csave, isave, dsave)
Definition: lbfgsb.f:2457
subroutine cauchy(n, x, l, u, nbd, g, iorder, iwhere, t, d, xcp, m, wy, ws, sy, wt, theta, col, head, p, c, wbp, v, nseg, iprint, sbgnrm, info, epsmch)
Definition: lbfgsb.f:1225
static double * time1
Definition: mafillsmmain.c:41
subroutine prn1lb(n, m, l, u, x, iprint, itfile, epsmch)
Definition: lbfgsb.f:2672
subroutine prn3lb(n, x, f, task, iprint, info, itfile, iter, nfgv, nintol, nskip, nact, sbgnrm, time, nseg, word, iback, stp, xstep, k, cachyt, sbtime, lnscht)
Definition: lbfgsb.f:2814
subroutine subsm(n, m, nsub, ind, l, u, nbd, x, d, xp, ws, wy, theta, xx, gg, col, head, iword, wv, wn, iprint, info)
Definition: lbfgsb.f:2994
subroutine prn2lb(n, x, f, g, iprint, itfile, iter, nfgv, nact, sbgnrm, nseg, word, iword, iback, stp, xstep)
Definition: lbfgsb.f:2744
double precision function ddot(N, DX, INCX, DY, INCY)
Definition: dgmres.f:2469
subroutine projgr(n, l, u, nbd, x, g, sbgnrm)
Definition: lbfgsb.f:2943