Version 0.1
#include <cblas.h>
#include <math.h>
#include <assert.h>
#include <string.h>
#include "bcls.h"
#include "bclib.h"
#include "bccgls.h"
#include "bclsqr.h"
#include "bcsolver.h"
Include dependency graph for bcsolver.c:
Go to the source code of this file.
Functions | |
static int | bcls_newton_step (BCLS *ls, int m, int nFree, int ix[], double damp, int itnLim, double tol, double dxFree[], double x[], double c[], double r[], int *itns, double *opt) |
Compute a Newton step on the free variables. | |
static int | bcls_proj_search (BCLS *ls, int m, int n, double x[], double dx[], double f, double g[], double aBreak[], int iBreak[], int ix[], double bl[], double bu[], double ex[], double Ad[], double Ae[], int *totHits) |
Find a minimizer along the projected search direction. | |
static int | bcls_proj_backtrack (BCLS *ls, int m, int n, double x[], double dx[], double f, double g[], double bl[], double bu[], double s[], int ix[], double As[], int *nSteps) |
Simple backtracking to satisfy an Armijo condition. | |
void | bcls_solver (BCLS *ls, int m, int n, double *bNorm, double x[], double b[], double c[], double bl[], double bu[], double r[], double g[], double dx[], double dxFree[], int ix[], double aBreak[], int iBreak[], double anorm[]) |
Implementation of the BCLS projected Newton/gradient search. |
Definition in file bcsolver.c.
static int bcls_newton_step | ( | BCLS * | ls, | |
int | m, | |||
int | nFree, | |||
int | ix[], | |||
double | damp, | |||
int | itnLim, | |||
double | tol, | |||
double | dxFree[], | |||
double | x[], | |||
double | c[], | |||
double | r[], | |||
int * | itns, | |||
double * | opt | |||
) | [static] |
Compute a Newton step on the free variables.
Compute a Newton step on the variables indexed by ix. The integer array ix (of length nFree) stores the indices of free variables. Thus, the kth free variable, with 0 <= k < nFree, is actually variable ix[k].
Conceptually, suppose that x and A are partitioned as
x = [ x_B x_N ] and A = [ B N ]
with x_N denoting the free variables. The objective function is given by
f(x) = 1/2 || Ax - b ||^2 + c'x + 1/2 damp^2 || x ||^2.
With respect to the free variables, the gradient and Hessian are
g_N(x) = N'(Ax - b) + cF + damp^2 xF = -N'r + cF + damp^2 xF and H_N(x) = N'N + damp^2 I,
where
The (regularized) Newton step dx_N is then the solution of
H_N(x) dx_N = - g_N(x) (*)
or equivalently,
(N'N + damp^2 I) dx_N = N'r - cF - damp^2 xF, (**)
which is really the normal equations for the least-squares problem
|| ( N ) ( r ) || min || ( ) dx - ( ) ||. (***) || ( damp I ) ( - 1/damp cF - damp xF ) ||
Recall that I has dimension (nFree,nFree).
Inexact Newton -------------- An inexact Newton step is based on approximately satisfying (*). The optimality criteria for this step is given by
|| H_N(x) dx_N + g_N(x) || <= rTol || g_N(x) ||, (****)
where rTol < 1. This routine (bcls_newton_step) computes the RHS of (****) and passes that to the subroutines as a required optimality.
SUBROUTINES ----------- This routine uses either one of the following
[in,out] | ls | The BCLS problem context. |
[in] | m | Number of rows in A. |
[in] | nFree | Number of columns in N. |
[in] | ix | Index of free variables (cols of N). |
[in] | damp | Sqrt of regularization parameter on x. |
[in] | itnLim | Maximum number of minor iterations allowed. |
[in] | tol | Required subproblem optimality. |
[out] | dxFree | The Newton step (solution of the LS subprob). |
[in] | x | The current iterate. |
[in] | c | The linear term. |
[in,out] | r | The current residual. |
[out] | itns | No. of iterations taken by the subprob solver. |
[out] | opt | The optimality of the subproblem solution. |
Definition at line 407 of file bcsolver.c.
References bcls_aprod(), bcls_newton_step_cgls(), bcls_newton_step_lsqr(), BCLS_NEWTON_STEP_LSQR, BCLS_PRECON_INIT, BCLS_PRECON_TERM, BCLS_PRECON_U, BCLS_PROD_INIT, BCLS_PROD_TERM, bcls_usolve(), cblas_dcopy(), BCLS::dx, and BCLS::Usolve.
Referenced by bcls_solver().
Here is the call graph for this function:
static int bcls_proj_backtrack | ( | BCLS * | ls, | |
int | m, | |||
int | n, | |||
double | x[], | |||
double | dx[], | |||
double | f, | |||
double | g[], | |||
double | bl[], | |||
double | bu[], | |||
double | s[], | |||
int | ix[], | |||
double | As[], | |||
int * | nSteps | |||
) | [static] |
Simple backtracking to satisfy an Armijo condition.
Do a simple backtracking search on the quadratic function
q(x + s) = f(x) + g(x)'s + 0.5 s' H s,
where H = A'A + damp^2 I and s = P[dx] is the projection of dx onto the box defined by bl and bu. The function value f = f(x) and the gradient g = g(x) are defined at entry and left unchanged. The step dx is reduced by a contant amount at each iteration until its projection s satisfies an Armijo condition, ie, sufficient descent is achieved when
g's + 0.5 s' H s <= mu g's,
and mu is a fixed constant.
Definition at line 846 of file bcsolver.c.
References BCLS::backtrack, BCLS::backtrack_limit, bcls_aprod(), BCLS_EXIT_LFAIL, BCLS_EXIT_UNBND, BCLS_PROD_A, bcls_project_step(), BCLS::BigNum, cblas_daxpy(), cblas_ddot(), BCLS::damp, BCLS::eps, BCLS::epsx, BCLS::mu, and PRINT3.
Referenced by bcls_solver().
Here is the call graph for this function:
static int bcls_proj_search | ( | BCLS * | ls, | |
int | m, | |||
int | n, | |||
double | x[], | |||
double | dx[], | |||
double | f, | |||
double | g[], | |||
double | aBreak[], | |||
int | iBreak[], | |||
int | ix[], | |||
double | bl[], | |||
double | bu[], | |||
double | ex[], | |||
double | Ad[], | |||
double | Ae[], | |||
int * | totHits | |||
) | [static] |
Find a minimizer along the projected search direction.
Find the exact minimizer of the quadratic function
q(t) = 1/2 || A X(t) - b ||^2 + 1/2 damp^2 || X(t) ||^2,
where X(t) is the projection of (x0 + t dx) onto the box defined by bl and bu.
The arguments y and e are workspace vectors. They must be at least length n and their contents will be overwritten.
Definition at line 491 of file bcsolver.c.
References BCLS::BigNum, BCLS::damp, BCLS::eps, BCLS::epsfixed, BCLS::epsx, PRINT6, and BCLS::print_level.
Referenced by bcls_solver().
void bcls_solver | ( | BCLS * | ls, | |
int | m, | |||
int | n, | |||
double * | bNorm, | |||
double | x[], | |||
double | b[], | |||
double | c[], | |||
double | bl[], | |||
double | bu[], | |||
double | r[], | |||
double | g[], | |||
double | dx[], | |||
double | dxFree[], | |||
int | ix[], | |||
double | aBreak[], | |||
int | iBreak[], | |||
double | anorm[] | |||
) |
Implementation of the BCLS projected Newton/gradient search.
[in,out] | ls | The BCLS problem context. |
[in] | m | Number of rows in A. |
[in] | n | Number of columns in A. |
[in,out] | bNorm | Norm of the RHS. |
[in,out] | x | The current iterate. |
[in] | b | The RHS vector. |
[in] | c | The linear-term vector. |
[in] | bl | Lower bounds on x. |
[in] | bu | Upper bounds on x. |
[in,out] | r | The residual vector r = b - Ax. |
[in,out] | g | The objective gradient. |
[in,out] | dx | The step direction in the full space. |
[in,out] | dxFree | The step direction on the free variables. |
[in,out] | ix | Indices of the free variables. |
[in,out] | aBreak | A list of breakpoints for a search direction. |
[in,out] | iBreak | The indices of the breakpoints. |
[in,out] | anorm | Column norms of A. |
Definition at line 95 of file bcsolver.c.
References bcls_aprod(), bcls_callback(), bcls_dload(), bcls_dual_inf(), BCLS_EXIT_CALLBK, BCLS_EXIT_CNVGD, BCLS_EXIT_LFAIL, BCLS_EXIT_MAJOR, BCLS_EXIT_MINOR, BCLS_EXIT_UNBND, BCLS_EXIT_UNDEF, bcls_free_vars(), bcls_mid(), bcls_newton_step(), BCLS_PROD_A, BCLS_PROD_At, bcls_proj_backtrack(), bcls_proj_search(), BCLS_PROJ_SEARCH_APPROX, BCLS_PROJ_SEARCH_EXACT, bcls_scatter(), BCLS_SOLN_OPTIM, cblas_daxpy(), cblas_dcopy(), cblas_ddot(), cblas_dnrm2(), cblas_dscal(), BCLS::damp, BCLS::itnMaj, BCLS::itnMin, BCLS::optTol, and PRINT1.
Referenced by bcls_solve_prob().
Here is the call graph for this function:
Generated on Sun Mar 4 22:50:03 2007 by Doxygen 1.5.1