Version 0.1
#include <cblas.h>
#include <string.h>
#include <assert.h>
#include <float.h>
#include "bcls.h"
#include "bclib.h"
#include "bccgls.h"
Include dependency graph for bccgls.c:
Go to the source code of this file.
Functions | |
static void | aprod_free (BCLS *ls, const int mode, const int m, const int nFree, double dxFree[], double y[]) |
Mat-vec routine called by CGLS. | |
static int | bcls_cgls (BCLS *ls, int m, int n, int kmax, double tol, double damp, double c[], double x[], double r[], double s[], double p[], int *itns, double *opt, FILE *nout) |
CG for least squares. | |
int | bcls_newton_step_cgls (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 using CGLS. |
Definition in file bccgls.c.
static void aprod_free | ( | BCLS * | ls, | |
const int | mode, | |||
const int | m, | |||
const int | nFree, | |||
double | dxFree[], | |||
double | y[] | |||
) | [static] |
Mat-vec routine called by CGLS.
CGLS call this routine, which in turn calls bcls_aprod:
This routine is declared "static" so that it won't be confused with the user's own Aprod routine.
[in,out] | ls | BCLS solver context. |
[in] | mode | Specifies if the product with A or A' is required. |
[in] | m | Number of rows in A. |
[in] | nFree | Length of dxFree (and ix which is recovered via ls). |
[in,out] | dxFree | RHS or LHS (depending on the mode). |
[in,out] | y | RHS or LHS (depending on the mode). |
Definition at line 67 of file bccgls.c.
References bcls_aprod(), bcls_gather(), BCLS_PROD_A, BCLS_PROD_At, bcls_scatter(), BCLS::dx, and BCLS::ix.
Referenced by bcls_cgls().
Here is the call graph for this function:
static int bcls_cgls | ( | BCLS * | ls, | |
int | m, | |||
int | n, | |||
int | kmax, | |||
double | tol, | |||
double | damp, | |||
double | c[], | |||
double | x[], | |||
double | r[], | |||
double | s[], | |||
double | p[], | |||
int * | itns, | |||
double * | opt, | |||
FILE * | nout | |||
) | [static] |
CG for least squares.
This routine implements CG for solving the symmetric linear system
( N'N + damp I ) x = N'r + c
where N is an m-by-n submatrix of the general matrix A. This implementation is adapted from Hestenes and Stiefel's CG for least-squares method to accomodate the damping term and the additional vector on the RHS. (See, eg, Bjork, Numerical Methods for Least-Squares Problems, 1996, p. 289.)
[in,out] | ls | BCLS solver context. |
[in] | m | No. of rows in A. |
[in] | n | No. of cols in A. |
[in] | kmax | Maximum no. of CG iterations. |
[in] | tol | Required optimality tolerance:
|
[in] | damp | Regularization term on x. |
[in] | c | Additional vector (set NULL if not present). |
[out] | x | Solution vector. |
[in,out] | r | On entry, r must be defined. On exit, r is the residual. |
[in,out] | s | Work vector. |
[in,out] | p | Work vector. |
[out] | itns | No. of CG iterations. |
[out] | opt | Relative residual (see tol, above). |
[in] | nout | File for log output (set NULL if not needed). |
Definition at line 129 of file bccgls.c.
References aprod_free(), bcls_dload(), BCLS_PROD_A, BCLS_PROD_At, cblas_daxpy(), cblas_dcopy(), cblas_dnrm2(), and cblas_dscal().
Referenced by bcls_newton_step_cgls().
Here is the call graph for this function:
int bcls_newton_step_cgls | ( | 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 using CGLS.
See the description of bcls_newton_step.
[in,out] | ls | BCLS solver context |
[in] | m | No. of rows in A. |
[in] | nFree | No. of free columns, i.e., no. of cols in N. |
[in] | ix | Index of vars for which a Newton step is computed. |
[in] | damp | Regularization term on dxFree. |
[in] | itnLim | Maximum no. of CG iterations. |
[in] | tol | Required optimality tolerance:
|
[in] | dxFree | Newton direction in the free vars (len = nFree). |
[in] | x | Current iterate. |
[in] | c | Linear term. (Set to NULL if not present.) |
[in,out] | r | On entry, r is current residual: r = b - Ax. Will be used as workspace. |
[in] | itns | No. of CG iterations. |
[in] | opt | Relative residual (see tol, above). |
Definition at line 258 of file bccgls.c.
References bcls_cgls(), bcls_gather(), cblas_dscal(), BCLS::minor_file, BCLS::wrk_u, BCLS::wrk_v, and BCLS::wrk_w.
Referenced by bcls_newton_step().
Here is the call graph for this function:
Generated on Sun Mar 4 22:50:03 2007 by Doxygen 1.5.1