Object-oriented Scientific Computing Library: Version 0.910
Data Structures | Functions
o2scl_linalg Namespace Reference

The namespace for linear algebra classes and functions. More...


Detailed Description

This namespace documentation is in the file src/base/lib_settings.h

Idea for Future:
Move this documentation to the linalg directory.

Data Structures

class  linear_solver
 A generic solver for the linear system $ A x = b $ [abstract base]. More...
class  linear_solver_lu
 Generic linear solver using LU decomposition. More...
class  linear_solver_qr
 Generic linear solver using QR decomposition. More...
class  linear_solver_hh
 Generic Householder linear solver. More...
class  gsl_solver_LU
 GSL solver by LU decomposition. More...
class  gsl_solver_QR
 GSL solver by QR decomposition. More...
class  gsl_solver_HH
 GSL Householder solver. More...
class  pointer_2_mem
 Allocation object for 2 C-style arrays of equal size. More...
class  pointer_4_mem
 Allocation object for 4 C-style arrays of equal size. More...
class  pointer_5_mem
 Allocation object for 5 C-style arrays of equal size. More...
class  uvector_2_mem
 Allocation object for 2 arrays of equal size. More...
class  uvector_4_mem
 Allocation object for 4 arrays of equal size. More...
class  uvector_5_mem
 Allocation object for 5 arrays of equal size. More...

Functions

void create_givens (const double a, const double b, double &c, double &s)
 Desc.
template<class mat1_t , class mat2_t >
void apply_givens_qr (size_t M, size_t N, mat1_t &Q, mat2_t &R, size_t i, size_t j, double c, double s)
 Apply a rotation to matrices from the QR decomposition.
template<class mat1_t , class mat2_t >
void apply_givens_lq (size_t M, size_t N, mat1_t &Q, mat2_t &L, size_t i, size_t j, double c, double s)
 Apply a rotation to matrices from the LQ decomposition.
template<class vec_t >
void apply_givens_vec (vec_t &v, size_t i, size_t j, double c, double s)
 Apply a rotation to a vector, $ v \rightarrow G^{T} v $.
template<class mat_t , class vec_t , class vec2_t >
int HH_solve (size_t n, mat_t &A, const vec_t &b, vec2_t &x)
 Solve linear system after Householder decomposition.
template<class mat_t , class vec_t >
int HH_svx (size_t N, size_t M, mat_t &A, vec_t &x)
 Solve a linear system after Householder decomposition in place.
template<class vec_t >
double householder_transform (const size_t n, vec_t &v)
 Replace the vector v with a Householder vector and a coefficient tau that annihilates v[1] through v[n-1] (inclusive)
template<class mat_t >
double householder_transform_subcol (mat_t &A, const size_t ir, const size_t ic, const size_t M)
 Compute the Householder transform of a vector formed with n rows of a column of a matrix.
template<class mat_t >
double householder_transform_subrow (mat_t &A, const size_t ir, const size_t ic, const size_t N)
 Compute the Householder transform of a vector formed with the last n columns of a row of a matrix.
template<class vec_t , class mat_t >
void householder_hm (const size_t M, const size_t N, double tau, const vec_t &v, mat_t &A)
 Apply a Householder transformation $ (v,\tau) $ to matrix $ A $ of size $ M $ by $ N $.
template<class mat_t >
void householder_hm_sub (mat_t &M, const size_t ir, const size_t ic, const size_t nr, const size_t nc, const mat_t &M2, const size_t ir2, const size_t ic2, double tau)
 Apply a Householder transformation to submatrix of a larger matrix.
template<class mat1_t , class mat2_t >
void householder_hm_sub2 (const size_t M, const size_t ic, double tau, const mat1_t &mv, mat2_t &A)
 Special version of Householder transformation for QR_unpack()
template<class vec_t , class vec2_t >
void householder_hv (const size_t N, double tau, const vec_t &v, vec2_t &w)
 Apply a Householder transformation v to vector w.
template<class mat_t , class vec_t >
void householder_hv_subcol (const mat_t &A, vec_t &w, double tau, const size_t ie, const size_t N)
 Apply a Householder transformation v to vector w where v is stored as a column in a matrix A.
template<class mat_t >
void householder_hm1 (const size_t M, const size_t N, double tau, mat_t &A)
 Apply a Householder transformation $ (v,\tau) $ to a matrix being build up from the identity matrix, using the first column of A as a Householder vector.
template<class vec_t , class mat_t >
void householder_mh (const size_t M, const size_t N, double tau, const vec_t &v, mat_t &A)
 Apply the Householder transformation (v,tau) to the right-hand side of the matrix A.
template<class mat_t , class mat2_t >
void householder_mh_sub (mat_t &M, const size_t ir, const size_t ic, const size_t nr, const size_t nc, const mat2_t &M2, const size_t ir2, const size_t ic2, double tau)
 Apply the Householder transformation (v,tau) to the right-hand side of the matrix A.
template<size_t N>
int LU_decomp_array_2d (const size_t n, double A[][N], o2scl::permutation &p, int &signum)
 Specialized version of LU_decomp for C-style 2D arrays.
template<class mat_t >
int diagonal_has_zero (const size_t N, mat_t &A)
 Return 1 if at least one of the elements in the diagonal is zero.
template<class mat_t >
int LU_decomp (const size_t N, mat_t &A, o2scl::permutation &p, int &signum)
 Compute the LU decomposition of the matrix A.
template<class mat_t , class mat_row_t >
int LU_decomp_alt (const size_t N, mat_t &A, o2scl::permutation &p, int &signum)
template<class mat_t , class vec_t , class vec2_t >
int LU_solve (const size_t N, const mat_t &LU, const o2scl::permutation &p, const vec_t &b, vec2_t &x)
 Solve a linear system after LU decomposition.
template<class mat_t , class vec_t >
int LU_svx (const size_t N, const mat_t &LU, const o2scl::permutation &p, vec_t &x)
 Solve a linear system after LU decomposition in place.
template<class mat_t , class mat2_t , class vec_t , class vec2_t , class vec3_t >
int LU_refine (const size_t N, const mat_t &A, const mat2_t &LU, const o2scl::permutation &p, const vec_t &b, vec2_t &x, vec3_t &residual)
 Refine the solution of a linear system.
template<class mat_t , class mat2_t , class mat_col_t >
int LU_invert (const size_t N, const mat_t &LU, const o2scl::permutation &p, mat2_t &inverse)
 Compute the inverse of a matrix from its LU decomposition.
template<class mat_t >
double LU_det (const size_t N, const mat_t &LU, int signum)
 Compute the determinant of a matrix from its LU decomposition.
template<class mat_t >
double LU_lndet (const size_t N, const mat_t &LU)
 Compute the logarithm of the absolute value of the determinant of a matrix from its LU decomposition.
template<class mat_t >
int LU_sgndet (const size_t N, const mat_t &LU, int signum)
 Compute the sign of the determinant of a matrix from its LU decomposition.
template<class mat_t , class vec_t >
void QR_decomp (size_t M, size_t N, mat_t &A, vec_t &tau)
 Compute the QR decomposition of matrix A.
template<class mat_t , class vec_t , class vec2_t , class vec3_t >
void QR_solve (size_t N, const mat_t &QR, const vec_t &tau, const vec2_t &b, vec3_t &x)
 Solve the system A x = b using the QR factorization.
template<class mat_t , class vec_t , class vec2_t >
void QR_svx (size_t M, size_t N, const mat_t &QR, const vec_t &tau, vec2_t &x)
 Solve the system A x = b in place using the QR factorization.
template<class mat_t , class vec_t , class vec2_t >
void QR_QTvec (const size_t M, const size_t N, const mat_t &QR, const vec_t &tau, vec2_t &v)
 Form the product Q^T v from a QR factorized matrix.
template<class mat1_t , class mat2_t , class mat3_t , class vec_t >
void QR_unpack (const size_t M, const size_t N, const mat1_t &QR, const vec_t &tau, mat2_t &Q, mat3_t &R)
 Unpack the QR matrix to the individual Q and R components.
template<class mat1_t , class mat2_t , class vec1_t , class vec2_t >
void QR_update (size_t M, size_t N, mat1_t &Q, mat2_t &R, vec1_t &w, vec2_t &v)
 Update a QR factorisation for A= Q R, A' = A + u v^T,.
template<class vec_t , class vec2_t >
int chop_small_elements (size_t N, vec_t &d, vec2_t &f)
 Desc.
template<class vec_t , class vec2_t >
double trailing_eigenvalue (size_t n, const vec_t &d, const vec_t &f)
 Desc.
int create_schur (double d0, double f0, double d1, double &c, double &s)
 Desc.
template<class vec_t , class vec2_t , class mat_t , class mat2_t >
int svd2 (size_t M, size_t N, vec_t &d, vec2_t &f, mat_t &U, mat2_t &V)
 Desc.
template<class vec_t , class vec2_t , class mat_t >
int chase_out_intermediate_zero (size_t M, size_t n, vec_t &d, vec2_t &f, mat_t &U, size_t k0)
 Desc.
template<class vec_t , class vec2_t , class mat_t >
int chase_out_trailing_zero (size_t N, size_t n, vec_t &d, vec2_t &f, mat_t &V)
 Desc.
template<class vec_t , class vec2_t , class mat_t , class mat2_t >
int qrstep (size_t M, size_t N, size_t n, vec_t &d, vec2_t &f, mat_t &U, mat2_t &V)
 Desc.
template<class vec_t , class vec2_t , class vec3_t , class vec4_t , class mem_t , class mem_vec_t >
void solve_tridiag_sym (const vec_t &diag, const vec2_t &offdiag, const vec3_t &b, vec4_t &x, size_t N, mem_t &m)
 Solve a symmetric tridiagonal linear system with user-specified memory.
template<class vec_t , class vec2_t , class vec3_t , class vec4_t , class vec5_t , class mem_t , class mem_vec_t >
void solve_tridiag_nonsym (const vec_t &diag, const vec2_t &abovediag, const vec3_t &belowdiag, const vec4_t &rhs, vec5_t &x, size_t N, mem_t &m)
 Solve an asymmetric tridiagonal linear system with user-specified memory.
template<class vec_t , class vec2_t , class vec3_t , class vec4_t , class mem_t , class mem_vec_t >
void solve_cyc_tridiag_sym (const vec_t &diag, const vec2_t &offdiag, const vec3_t &b, vec4_t &x, size_t N, mem_t &m)
 Solve a symmetric cyclic tridiagonal linear system with user specified memory.
template<class vec_t , class vec2_t , class vec3_t , class vec4_t , class vec5_t , class mem_t , class mem_vec_t >
void solve_cyc_tridiag_nonsym (const vec_t &diag, const vec2_t &abovediag, const vec3_t &belowdiag, const vec4_t &rhs, vec5_t &x, size_t N, mem_t &m)
 Solve an asymmetric cyclic tridiagonal linear system with user-specified memory.
template<class vec_t , class vec2_t , class vec3_t , class vec4_t >
void solve_tridiag_sym (const vec_t &diag, const vec2_t &offdiag, const vec3_t &b, vec4_t &x, size_t N)
 Solve a symmetric tridiagonal linear system.
template<class vec_t , class vec2_t , class vec3_t , class vec4_t , class vec5_t >
void solve_tridiag_nonsym (const vec_t &diag, const vec2_t &abovediag, const vec3_t &belowdiag, const vec4_t &rhs, vec5_t &x, size_t N)
 Solve an asymmetric tridiagonal linear system.
template<class vec_t , class vec2_t , class vec3_t , class vec4_t >
void solve_cyc_tridiag_sym (const vec_t &diag, const vec2_t &offdiag, const vec3_t &b, vec4_t &x, size_t N)
 Solve a symmetric cyclic tridiagonal linear system.
template<class vec_t , class vec2_t , class vec3_t , class vec4_t , class vec5_t >
void solve_cyc_tridiag_nonsym (const vec_t &diag, const vec2_t &abovediag, const vec3_t &belowdiag, const vec4_t &rhs, vec5_t &x, size_t N)
 Solve an asymmetric cyclic tridiagonal linear system.

Function Documentation

template<class mat1_t , class mat2_t >
void o2scl_linalg::apply_givens_qr ( size_t  M,
size_t  N,
mat1_t &  Q,
mat2_t &  R,
size_t  i,
size_t  j,
double  c,
double  s 
)

This performs $ Q \rightarrow Q G $ and $ R \rightarrow G^{T} R $.

Definition at line 58 of file givens_base.h.

template<class mat1_t , class mat2_t >
void o2scl_linalg::apply_givens_lq ( size_t  M,
size_t  N,
mat1_t &  Q,
mat2_t &  L,
size_t  i,
size_t  j,
double  c,
double  s 
)

This performs $ Q \rightarrow Q G $ and $ L \rightarrow L G^{T} $.

Definition at line 83 of file givens_base.h.

template<class mat_t , class vec_t >
int o2scl_linalg::HH_svx ( size_t  N,
size_t  M,
mat_t &  A,
vec_t &  x 
)
Idea for Future:
Handle memory allocation like the tridiagonal functions.

Definition at line 68 of file hh_base.h.

template<class vec_t >
double o2scl_linalg::householder_transform ( const size_t  n,
vec_t &  v 
)

On exit, this function returns the value of $ \tau = 2/ (v^{T} v) $. If n is less than or equal to 1 then this function returns zero without calling the error handler.

Definition at line 69 of file householder_base.h.

template<class mat_t >
double o2scl_linalg::householder_transform_subcol ( mat_t &  A,
const size_t  ir,
const size_t  ic,
const size_t  M 
)

This performs a Householder transform of a vector defined by a column of a matrix A which starts at element A[ir][ic] and ends at element A[N-1][ic].

Used in QR_decomp().

Definition at line 114 of file householder_base.h.

template<class mat_t >
double o2scl_linalg::householder_transform_subrow ( mat_t &  A,
const size_t  ir,
const size_t  ic,
const size_t  N 
)

This performs a Householder transform of a vector defined by a row of a matrix A which starts at element A[ir][ic] and ends at element A[ir][N-1]

Definition at line 156 of file householder_base.h.

template<class mat_t >
void o2scl_linalg::householder_hm_sub ( mat_t &  M,
const size_t  ir,
const size_t  ic,
const size_t  nr,
const size_t  nc,
const mat_t &  M2,
const size_t  ir2,
const size_t  ic2,
double  tau 
)

This applies a householder transformation (v,tau) to submatrix of M. The submatrix has nr rows and nc columns and starts at row ir of column ic of the original matrix M. The vector v is taken from a column of M2 starting at row ir2 and column ic2.

This function is used in QR_decomp().

Definition at line 237 of file householder_base.h.

template<class mat_t , class vec_t >
void o2scl_linalg::householder_hv_subcol ( const mat_t &  A,
vec_t &  w,
double  tau,
const size_t  ie,
const size_t  N 
)

Used in QR_QTvec().

Definition at line 337 of file householder_base.h.

template<class mat_t , class mat2_t >
void o2scl_linalg::householder_mh_sub ( mat_t &  M,
const size_t  ir,
const size_t  ic,
const size_t  nr,
const size_t  nc,
const mat2_t &  M2,
const size_t  ir2,
const size_t  ic2,
double  tau 
)

This is unfinished but needed for the bidiagonalization functions.

Definition at line 446 of file householder_base.h.

template<size_t N>
int o2scl_linalg::LU_decomp_array_2d ( const size_t  n,
double  A[][N],
o2scl::permutation &  p,
int &  signum 
)
Note:
Note that N and n must be equal, by no checking is done to ensure that this is the case

Definition at line 42 of file lu.h.

template<class mat_t >
int o2scl_linalg::LU_decomp ( const size_t  N,
mat_t &  A,
o2scl::permutation &  p,
int &  signum 
)

On output the diagonal and upper triangular part of the input matrix A contain the matrix U. The lower triangular part of the input matrix (excluding the diagonal) contains L. The diagonal elements of L are unity, and are not stored.

The permutation matrix P is encoded in the permutation p. The j-th column of the matrix P is given by the k-th column of the identity matrix, where k = p_j the j-th element of the permutation vector. The sign of the permutation is given by signum. It has the value (-1)^n, where n is the number of interchanges in the permutation.

The algorithm used in the decomposition is Gaussian Elimination with partial pivoting (Golub & Van Loan, Matrix Computations, Algorithm 3.4.1).

Idea for Future:
The "swap rows j and i_pivot" section could probably be made more efficient using a "matrix_row"-like object as done in GSL. (7/16/09 - I've tried this, and it doesn't seem to improve the speed significantly.)

Definition at line 86 of file lu_base.h.

template<class mat_t , class vec_t , class vec2_t >
int o2scl_linalg::LU_solve ( const size_t  N,
const mat_t &  LU,
const o2scl::permutation &  p,
const vec_t &  b,
vec2_t &  x 
)

This function solve the square system A x = b using the LU decomposition of A into (LU, p) given by gsl_linalg_LU_decomp or gsl_linalg_complex_LU_decomp.

Definition at line 207 of file lu_base.h.

template<class mat_t , class vec_t >
int o2scl_linalg::LU_svx ( const size_t  N,
const mat_t &  LU,
const o2scl::permutation &  p,
vec_t &  x 
)

These functions solve the square system A x = b in-place using the LU decomposition of A into (LU,p). On input x should contain the right-hand side b, which is replaced by the solution on output.

Definition at line 230 of file lu_base.h.

template<class mat_t , class mat2_t , class vec_t , class vec2_t , class vec3_t >
int o2scl_linalg::LU_refine ( const size_t  N,
const mat_t &  A,
const mat2_t &  LU,
const o2scl::permutation &  p,
const vec_t &  b,
vec2_t &  x,
vec3_t &  residual 
)

These functions apply an iterative improvement to x, the solution of A x = b, using the LU decomposition of A into (LU,p). The initial residual r = A x - b is also computed and stored in residual.

Definition at line 262 of file lu_base.h.

template<class mat_t , class mat2_t , class mat_col_t >
int o2scl_linalg::LU_invert ( const size_t  N,
const mat_t &  LU,
const o2scl::permutation &  p,
mat2_t &  inverse 
)

These functions compute the inverse of a matrix A from its LU decomposition (LU,p), storing the result in the matrix inverse. The inverse is computed by solving the system A x = b for each column of the identity matrix. It is preferable to avoid direct use of the inverse whenever possible, as the linear solver functions can obtain the same result more efficiently and reliably.

Idea for Future:
Could rewrite to avoid mat_col_t, (9/16/09 - However, the function may be faster if mat_col_t is left in, so it's unclear what's best.)

Definition at line 300 of file lu_base.h.

template<class mat_t >
double o2scl_linalg::LU_det ( const size_t  N,
const mat_t &  LU,
int  signum 
)

Compute the determinant of a matrix A from its LU decomposition, LU. The determinant is computed as the product of the diagonal elements of U and the sign of the row permutation signum.

Definition at line 339 of file lu_base.h.

template<class mat_t >
double o2scl_linalg::LU_lndet ( const size_t  N,
const mat_t &  LU 
)

Compute the logarithm of the absolute value of the determinant of a matrix A, $ \ln|\det(A)| $, from its LU decomposition, LU. This function may be useful if the direct computation of the determinant would overflow or underflow.

Definition at line 361 of file lu_base.h.

template<class mat_t >
int o2scl_linalg::LU_sgndet ( const size_t  N,
const mat_t &  LU,
int  signum 
)

Compute the sign or phase factor of the determinant of a matrix A, $ \det(A)/|\det(A)| $, from its LU decomposition, LU.

Definition at line 380 of file lu_base.h.

template<class mat1_t , class mat2_t , class vec1_t , class vec2_t >
void o2scl_linalg::QR_update ( size_t  M,
size_t  N,
mat1_t &  Q,
mat2_t &  R,
vec1_t &  w,
vec2_t &  v 
)

The parameters M and N are the number of rows and columns of the matrix R.

      * Q' R' = QR + u v^T
      *       = Q (R + Q^T u v^T)
      *       = Q (R + w v^T)
      *
      * where w = Q^T u.
      *
      * Algorithm from Golub and Van Loan, "Matrix Computations", Section
      * 12.5 (Updating Matrix Factorizations, Rank-One Changes)
      

Definition at line 161 of file qr_base.h.

template<class vec_t , class vec2_t >
int o2scl_linalg::chop_small_elements ( size_t  N,
vec_t &  d,
vec2_t &  f 
)

The parameter N is the size of d.

Definition at line 55 of file svdstep_base.h.

template<class vec_t , class vec2_t >
double o2scl_linalg::trailing_eigenvalue ( size_t  n,
const vec_t &  d,
const vec_t &  f 
)

The parameter n is the size of the vector d.

Should be finished.

Definition at line 82 of file svdstep_base.h.

int o2scl_linalg::create_schur ( double  d0,
double  f0,
double  d1,
double &  c,
double &  s 
)

Should be finished.

Definition at line 117 of file svdstep_base.h.

template<class vec_t , class vec2_t , class mat_t , class mat2_t >
int o2scl_linalg::svd2 ( size_t  M,
size_t  N,
vec_t &  d,
vec2_t &  f,
mat_t &  U,
mat2_t &  V 
)

The parameter M is the number of rows in U and N is the number of rows in V.

Definition at line 169 of file svdstep_base.h.

template<class vec_t , class vec2_t , class mat_t >
int o2scl_linalg::chase_out_intermediate_zero ( size_t  M,
size_t  n,
vec_t &  d,
vec2_t &  f,
mat_t &  U,
size_t  k0 
)

Should be finished.

Definition at line 305 of file svdstep_base.h.

template<class vec_t , class vec2_t , class mat_t >
int o2scl_linalg::chase_out_trailing_zero ( size_t  N,
size_t  n,
vec_t &  d,
vec2_t &  f,
mat_t &  V 
)

Should be finished.

Definition at line 352 of file svdstep_base.h.

template<class vec_t , class vec2_t , class mat_t , class mat2_t >
int o2scl_linalg::qrstep ( size_t  M,
size_t  N,
size_t  n,
vec_t &  d,
vec2_t &  f,
mat_t &  U,
mat2_t &  V 
)

Should be finished.

The parameter M is the number of rows in U, N is the number of rows in V, and n is the length of the vector d.

Definition at line 402 of file svdstep_base.h.

template<class vec_t , class vec2_t , class vec3_t , class vec4_t , class mem_t , class mem_vec_t >
void o2scl_linalg::solve_tridiag_sym ( const vec_t &  diag,
const vec2_t &  offdiag,
const vec3_t &  b,
vec4_t &  x,
size_t  N,
mem_t &  m 
)

This function solves the system $ A x = b $ where $ A $ is a matrix of the form

      *
      *     diag[0]  offdiag[0]             0   .....
      *  offdiag[0]     diag[1]    offdiag[1]   .....
      *           0  offdiag[1]       diag[2]
      *           0           0    offdiag[2]   .....
      

given the N diagonal elements in diag, N-1 diagonal elements in offdiag, and the N elements b from the RHS.

See EngelnMullges96 .

Definition at line 208 of file tridiag_base.h.

template<class vec_t , class vec2_t , class vec3_t , class vec4_t , class vec5_t , class mem_t , class mem_vec_t >
void o2scl_linalg::solve_tridiag_nonsym ( const vec_t &  diag,
const vec2_t &  abovediag,
const vec3_t &  belowdiag,
const vec4_t &  rhs,
vec5_t &  x,
size_t  N,
mem_t &  m 
)

This function solves the system $ A x = b $ where $ A $ is a matrix of the form

      * 
      *       diag[0]  abovediag[0]             0   .....
      *  belowdiag[0]       diag[1]  abovediag[1]   .....
      *             0  belowdiag[1]       diag[2]
      *             0             0  belowdiag[2]   .....
      

This function uses plain Gauss elimination, not bothering with the zeroes.

Idea for Future:
Offer an option to avoid throwing on divide by zero?

Definition at line 273 of file tridiag_base.h.

template<class vec_t , class vec2_t , class vec3_t , class vec4_t , class mem_t , class mem_vec_t >
void o2scl_linalg::solve_cyc_tridiag_sym ( const vec_t &  diag,
const vec2_t &  offdiag,
const vec3_t &  b,
vec4_t &  x,
size_t  N,
mem_t &  m 
)

This function solves the system $ A x = b $ where $ A $ is a matrix of the form

      *
      *      diag[0]  offdiag[0]             0   .....  offdiag[N-1]
      *   offdiag[0]     diag[1]    offdiag[1]   .....
      *            0  offdiag[1]       diag[2]
      *            0           0    offdiag[2]   .....
      *          ...         ...
      * offdiag[N-1]         ...
      

See EngelnMullges96 .

Definition at line 331 of file tridiag_base.h.

template<class vec_t , class vec2_t , class vec3_t , class vec4_t , class vec5_t , class mem_t , class mem_vec_t >
void o2scl_linalg::solve_cyc_tridiag_nonsym ( const vec_t &  diag,
const vec2_t &  abovediag,
const vec3_t &  belowdiag,
const vec4_t &  rhs,
vec5_t &  x,
size_t  N,
mem_t &  m 
)

This function solves the system $ A x = b $ where $ A $ is a matrix of the form

      *
      *        diag[0]  abovediag[0]             0   .....  belowdiag[N-1]
      *   belowdiag[0]       diag[1]  abovediag[1]   .....
      *              0  belowdiag[1]       diag[2]
      *              0             0  belowdiag[2]   .....
      *            ...           ...
      * abovediag[N-1]           ...
      

This function solves the following system without the corner elements and then use Sherman-Morrison formula to compensate for them.

Idea for Future:
Offer an option to avoid throwing on divide by zero?

Definition at line 427 of file tridiag_base.h.

 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines

Documentation generated with Doxygen. Provided under the GNU Free Documentation License (see License Information).

Get Object-oriented Scientific Computing
Lib at SourceForge.net. Fast, secure and Free Open Source software
downloads.