ExaTron.dbreakpt
— MethodSubroutine dbreakpt
This subroutine computes the number of break-points, and the minimal and maximal break-points of the projection of x + alpha*w on the n-dimensional interval [xl,xu].
ExaTron.dcauchy
— MethodSubroutine dcauchy
This subroutine computes a Cauchy step that satisfies a trust region constraint and a sufficient decrease condition.
Ths Cauchy step is computed for the quadratic
q(s) = 0.5*s'*A*s + g'*s
where A is a symmetric matrix in compressed row storage, and g is a vector. Given a parameter alpha, the Cauchy step is
s[alpha] = P[x - alpha*g] - x,
with P the projection onto the n-dimensional interval [xl,xu]. The Cauchy step satisfies the trust region constraint and the sufficient decrease condition
|| s || <= delta, q(s) <= mu_0*(g'*s),
where mu_0 is a constant in (0,1).
MINPACK-2 Project. March 1999. Argonne National Laboratory. Chih-Jen Lin and Jorge J. More'.
ExaTron.dgpstep
— MethodSubroutine dgpstep
This subroutine computes the gradient projection step
s = P[x + alpha*w] - x,
where P is the projection on the n-dimensional interval [xl,xu].
ExaTron.dicf
— MethodSubroutine dicf
Given a sparse symmetric matrix A in compressed row storage, this subroutine computes an incomplete Cholesky factorization.
Implementation of dicf is based on the Jones-Plassmann code. Arrays indf and list define the data structure. At the beginning of the computation of the j-th column,
For k < j, indf[k] is the index of A for the first
nonzero l[i,k] in the k-th column with i >= j.
For k < j, list[i] is a pointer to a linked list of column
indices k with i = L.rowval[indf[k]].
For the computation of the j-th column, the array indr records the row indices. Hence, if nlj is the number of nonzeros in the j-th column, then indr[1],...,indr[nlj] are the row indices. Also, for i > j, indf[i] marks the row indices in the j-th column so that indf[i] = 1 if l[i,j] is not zero.
MINPACK-2 Project. May 1998. Argonne National Laboratory. Chih-Jen Lin and Jorge J. More'.
ExaTron.dicfs
— MethodSubroutine dicfs
Given a symmetric matrix A in compreessed column storage, this subroutine computes an incomplete Cholesky factor of A + alpha*D, where alpha is a shift and D is the diagonal matrix with entries set to the l2 norms of the columns of A.
MINPACK-2 Project. October 1998. Argonne National Laboratory. Chih-Jen Lin and Jorge J. More'.
ExaTron.dmid
— MethodSubroutine dmid
This subroutine computes the projection of x on the n-dimensional interval [xl,xu].
MINPACK-2 Project. March 1999. Argonne National Laboratory. Chih-Jen Lin and Jorge J. More'.
ExaTron.dprsrch
— MethodSubroutine dprsrch
This subroutine uses a projected search to compute a step that satisfies a sufficient decrease condition for the quadratic
q(s) = 0.5*s'*A*s + g'*s,
where A is a symmetric matrix in compressed column storage, and g is a vector. Given the parameter alpha, the step is
s[alpha] = P[x + alpha*w] - x,
where w is the search direction and P the projection onto the n-dimensional interval [xl,xu]. The final step s = s[alpha] satisfies the sufficient decrease condition
q(s) <= mu_0*(g'*s),
where mu_0 is a constant in (0,1).
The search direction w must be a descent direction for the quadratic q at x such that the quadratic is decreasing in the ray x + alpha*w for 0 <= alpha <= 1.
MINPACK-2 Project. March 1999. Argonne National Laboratory. Chih-Jen Lin and Jorge J. More'.
ExaTron.dsel2
— MethodSubroutine dsel2
Given an array x, this subroutine permutes the elements of the array keys so that
abs(x(keys(i))) <= abs(x(keys(k))), 1 <= i <= k, abs(x(keys(k))) <= abs(x(keys(i))), k <= i <= n.
In other words, the smallest k elements of x in absolute value are x(keys(i)), i = 1,...,k, and x(keys(k)) is the kth smallest element.
MINPACK-2 Project. March 1998. Argonne National Laboratory. William D. Kastak, Chih-Jen Lin, and Jorge J. More'.
Revised October 1999. Length of x was incorrectly set to n.
ExaTron.dspcg
— MethodSubroutine dspcg
This subroutine generates a sequence of approximate minimizers for the subproblem
min { q(x) : xl <= x <= xu }.
The quadratic is defined by
q(x[0]+s) = 0.5*s'*A*s + g'*s,
where x[0] is a base point provided by the user, A is a symmetric matrix in compressed column storage, and g is a vector.
At each stage we have an approximate minimizer x[k], and generate a direction p[k] by using a preconditioned conjugate gradient method on the subproblem
min { q(x[k]+p) : || L'*p || <= delta, s(fixed) = 0 },
where fixed is the set of variables fixed at x[k], delta is the trust region bound, and L is an incomplete Cholesky factorization of the submatrix
B = A(free:free),
where free is the set of free variables at x[k]. Given p[k], the next minimizer x[k+1] is generated by a projected search.
The starting point for this subroutine is x[1] = x[0] + s, where x[0] is a base point and s is the Cauchy step.
The subroutine converges when the step s satisfies
|| (g + A*s)[free] || <= rtol*|| g[free] ||
In this case the final x is an approximate minimizer in the face defined by the free variables.
The subroutine terminates when the trust region bound does not allow further progress, that is, || L'*p[k] || = delta. In this case the final x satisfies q(x) < q(x[k]).
MINPACK-2 Project. March 1999. Argonne National Laboratory. Chih-Jen Lin and Jorge J. More'.
March 2000
Clarified documentation of nv variable. Eliminated the nnz = max(nnz,1) statement.
ExaTron.dssyax
— MethodSubroutine dssyax
This subroutine computes the matrix-vector product y = A*x, where A is a symmetric matrix with the strict lower triangular part in compressed column storage.
ExaTron.dtron
— MethodSubroutine dtron
This subroutine implements a trust region Newton method for the solution of large bound-constrained optimization problems
min { f(x) : xl <= x <= xu }
where the Hessian matrix is sparse. The user must evaluate the function, gradient, and the Hessian matrix.
ExaTron.dtrpcg
— MethodSubroutine dtrpcg
Given a sparse symmetric matrix A in compressed column storage, this subroutine uses a preconditioned conjugate gradient method to find an approximate minimizer of the trust region subproblem
min { q(s) : || L'*s || <= delta }.
where q is the quadratic
q(s) = 0.5s'As + g's,
A is a symmetric matrix in compressed column storage, L is a lower triangular matrix in compressed column storage, and g is a vector.
This subroutine generates the conjugate gradient iterates for the equivalent problem
min { Q(w) : || w || <= delta },
where Q is the quadratic defined by
Q(w) = q(s), w = L'*s.
Termination occurs if the conjugate gradient iterates leave the trust regoin, a negative curvature direction is generated, or one of the following two convergence tests is satisfied.
Convergence in the original variables:
|| grad q(s) || <= tol
Convergence in the scaled variables:
|| grad Q(w) || <= stol
Note that if w = L's, then Lgrad Q(w) = grad q(s).
MINPACK-2 Project. March 1999. Argonne National Laboratory. Chih-Jen Lin and Jorge J. More'.
August 1999
Corrected documentation for l, ldiag, lcolptr, and lrowind.
February 2001
We now set iters = 0 in the special case g = 0.
ExaTron.dtrqsol
— MethodSubroutine dtrqsol
This subroutine computes the largest (non-negative) solution of the quadratic trust region equation
||x + sigma*p|| = delta.
The code is only guaranteed to produce a non-negative solution if ||x|| <= delta, and p != 0. If the trust region equation has no solution, sigma = 0.
MINPACK-2 Project. March 1999. Argonne National Laboratory. Chih-Jen Lin and Jorge J. More'.
ExaTron.ihsort
— MethodSubroutine ihsort
Given an integer array keys of length n, this subroutine uses a heap sort to sort the keys in increasing order.
This subroutine is a minor modification of code written by Mark Jones and Paul Plassmann.
MINPACK-2 Project. March 1998. Argonne National Laboratory. Chih-Jen Lin and Jorge J. More'.
ExaTron.insort
— MethodSubroutine insort
Given an integer array keys of length n, this subroutine uses an insertion sort to sort the keys in increasing order.
MINPACK-2 Project. March 1998. Argonne National Laboratory. Chih-Jen Lin and Jorge J. More'.
ExaTron.tron_daxpy
— MethodSubroutine daxpy
This subroutine computes constant times a vector plus a vector. It uses unrolled loops for increments equal to one. Jack Dongarra, LINPACK, 3/11/78.
ExaTron.tron_ddot
— MethodSubroutine ddot
This subroutine forms the dot product of two vectors. It uses unrolled loops for increments equal to one. Jack Dongarra, LINPACK, 3/11/78.
ExaTron.tron_dnrm2
— MethodDNRM2 returns the euclidean norm of a vector via the function name, so that
DNRM2 := sqrt( x'*x )
– This version written on 25-October-1982. Modified on 14-October-1993 to inline the call to DLASSQ. Sven Hammarling, Nag Ltd.