gslmm::lq_decomposition< double > Struct Template Reference
[LQ decomposition]

#include <lq_decomposition_double.hh>

Inheritance diagram for gslmm::lq_decomposition< double >:

Inheritance graph
[legend]
Collaboration diagram for gslmm::lq_decomposition< double >:

Collaboration graph
[legend]
List of all members.

Detailed Description

template<>
struct gslmm::lq_decomposition< double >

class template specialisation for LQ decompositions of matricies of type double .

Factorise a general $ N\times M$ matrix $ A$ into

\[ A = L Q \]

where $ Q$ is orthogonal ($ M \times M$) and $ L$ is lower triangular ($ N \times M$).

$ Q$ is stored as a packed set of Householder transformations in the strict upper triangular part of the input matrix.

$ R$ is stored in the diagonal and lower triangle of the input matrix.

The full matrix for $ Q$ can be obtained as the product

\[ Q = Q_k .. Q_2 Q_1 \]

where $ k = min(M,N)$ and

\[ Q_i = (I - \tau_i * v_i * v_i') \]

and where $ v_i$ is a Householder vector

\[ v_i = \left[1, m_{i+1,i}, m_{i+2,i}, ... , m_{M,i}\right] \]

This storage scheme is the same as in LAPACK.

Todo:
gsl_linalg_L_solve_T


Public Types

typedef vector< double > vector_type
typedef matrix< double > matrix_type
typedef matrix_type::value_type value_type
typedef matrix_type::element_type element_type

Public Member Functions

 lq_decomposition (const matrix_type &m)
virtual ~lq_decomposition ()
virtual bool solve (const vector_type &b, vector_type &x) const
virtual bool solve (vector_type &x) const
virtual bool solve (const vector_type &b, vector_type &x, vector_type &r) const
virtual bool l_solve (const vector_type &b, vector_type &x) const
virtual bool l_solve (vector_type &x) const
virtual bool lq_solve (const vector_type &b, vector_type &x) const
virtual bool qt_multiply (vector_type &v) const
virtual bool update (const vector_type &w, const vector_type &v)
virtual bool q_multiply (vector_type &v) const
const matrix_typel () const
const matrix_typeq () const
const matrix_typelq () const

Protected Types

typedef matrix_type::data_type data_type

Protected Member Functions

void unpack () const
 lq_decomposition (const matrix_type &m, bool nodo)

Protected Attributes

vector_type _tau
matrix_type_q
matrix_type_l


Member Typedef Documentation

typedef matrix_type::data_type gslmm::lq_decomposition< double >::data_type [protected]
 

The low level matrix type.

Reimplemented from gslmm::matrix< double >.

Reimplemented in gslmm::lq_decomposition_pivot< double >.

typedef matrix_type::element_type gslmm::lq_decomposition< double >::element_type
 

The element type.

Reimplemented from gslmm::matrix< double >.

Reimplemented in gslmm::lq_decomposition_pivot< double >.

typedef matrix<double > gslmm::lq_decomposition< double >::matrix_type
 

The matrix type.

Reimplemented in gslmm::lq_decomposition_pivot< double >.

typedef matrix_type::value_type gslmm::lq_decomposition< double >::value_type
 

The element type.

Reimplemented from gslmm::matrix< double >.

Reimplemented in gslmm::lq_decomposition_pivot< double >.

typedef vector<double > gslmm::lq_decomposition< double >::vector_type
 

The vector type.

Reimplemented from gslmm::matrix< double >.

Reimplemented in gslmm::lq_decomposition_pivot< double >.


Constructor & Destructor Documentation

gslmm::lq_decomposition< double >::lq_decomposition const matrix_type m  )  [inline]
 

Create an LQ decomposition of the input matrix m.

A reference to the input matrix is stored for later use.

gslmm::lq_decomposition< double >::~lq_decomposition  )  [inline, virtual]
 

Destructor.

gslmm::lq_decomposition< double >::lq_decomposition const matrix_type m,
bool  nodo
[inline, protected]
 

Create an LQ decomposition of the input matrix m.

A reference to the input matrix is stored for later use.


Member Function Documentation

const lq_decomposition< double >::matrix_type & gslmm::lq_decomposition< double >::l  )  const [inline]
 

Returns:
The matrix $ R$ of the LQ decomposition of $ A$

bool gslmm::lq_decomposition< double >::l_solve vector_type x  )  const [inline, virtual]
 

This function solves the triangular system $ x^T L = b^T$ for $ x$ in-place.

On input $ x$ should contain the right-hand side $ b$ and is replaced by the solution on output. This function may be useful if the product $ b' = Q^T b$ has already been computed using qt_multiply

Parameters:
x On input, the vector $ b$ On return, the vector $ x$ that solves $ L x = b$
Returns:
true on success.

Reimplemented in gslmm::lq_decomposition_pivot< double >.

bool gslmm::lq_decomposition< double >::l_solve const vector_type b,
vector_type x
const [inline, virtual]
 

This function solves the triangular system $ x^T L = b^T$ for $ x$.

It may be useful if the product $ b' = Q^T b$ has already been computed using lq_multiply

Parameters:
b The vector $ b$
x On return, the vector $ x$ that solves $ L x = b$
Returns:
true on success.

Reimplemented in gslmm::lq_decomposition_pivot< double >.

const lq_decomposition< double >::matrix_type & gslmm::lq_decomposition< double >::lq  )  const [inline]
 

Returns:
The matrix $ LQ$ of the LQ decomposition of $ A$

bool gslmm::lq_decomposition< double >::lq_solve const vector_type b,
vector_type x
const [inline, virtual]
 

This function solves the system $ R x = Q^T b$ for $ x$.

It can be used when the LQ decomposition of a matrix is wanted in unpacked form as $(Q, R)$

Parameters:
b The vector $ b$
x The vector $ x$, the solution to $ R x = Q^T b$.
Returns:
true on success.

Reimplemented in gslmm::lq_decomposition_pivot< double >.

const lq_decomposition< double >::matrix_type & gslmm::lq_decomposition< double >::q  )  const [inline]
 

Returns:
The matrix $ Q$ of the LQ decomposition of $ A$

bool gslmm::lq_decomposition< double >::q_multiply vector_type v  )  const [inline, virtual]
 

This function applies the matrix $ Q$ encoded in the decomposition to the vector $ v$, storing the result $ v Q $ in $ v$.

The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix $ Q$.

Parameters:
v On input, the vector $ v$. On output the product $ v Q$.
Returns:
true on success.

bool gslmm::lq_decomposition< double >::qt_multiply vector_type v  )  const [inline, virtual]
 

This function applies the matrix $ Q^T$ encoded in the decomposition to the vector $ v$, storing the result $ v Q^T$ in $ v$.

The matrix multiplication is carried out directly using the encoding of the Householder vectors without needing to form the full matrix $ Q^T$.

Parameters:
v On input, the vector $ v$. On output the product $ v Q^T $.
Returns:
true on success.

bool gslmm::lq_decomposition< double >::solve const vector_type b,
vector_type x,
vector_type r
const [inline, virtual]
 

This function finds the least squares solution to the overdetermined system $ x^T A = b^T$ where the matrix $ A$ has more rows than columns.

The least squares solution minimizes the Euclidean norm of the residual, $ ||Ax - b||$. The routine uses the LQ decomposition of $ A$. The solution is returned in $ x$. The residual is computed as a by-product and stored in $ r$.

Parameters:
b The input vector $ b$. Must have the same dimension as the number of rows in $ A$.
x The output solution vector $ x$ Must have the same dimension as the number of columns in $ A$.
r The initial residual vector $ r$. Must have the same dimension as the number of rows in $ A$.
Returns:
True on success.

Reimplemented in gslmm::lq_decomposition_pivot< double >.

bool gslmm::lq_decomposition< double >::solve vector_type x  )  const [inline, virtual]
 

This member function solve the system $ x^T A = b^T$ using the LQ decomposition of $ A$,

\[ x^T L = b^T Q^T \]

The input matrix $ A$ must be square.

This can be used to find the inverse of matrix $ A$ by doing a back-subsitution of the identity matrix.

        gslmm::lq_decomposition<T> LQ(A);
        gslmm::matrix<T> Ainv(A.row_size(), B.column_size());
        for (size_t i = 0; i < A.row_size(); i++) {
          gslmm::vector<T> v(A.row_size());
          v.basis(i);
          LQ.solve(v);
          for (size_t j = 0; j < A.column_size(); j++) Ainv(i,j) = x[j];
        }
Parameters:
x On input, the vector $ b$, and on output the vector $ x$

Reimplemented in gslmm::lq_decomposition_pivot< double >.

bool gslmm::lq_decomposition< double >::solve const vector_type b,
vector_type x
const [inline, virtual]
 

This member function solve the system $ x^T A = b^T$ using the LQ decomposition of $ A$,

\[ x^T L = b^T Q^T \]

The input matrix $ A$ must be square.

This can be used to find the inverse of matrix $ A$ by doing a back-subsitution of the identity matrix.

        gslmm::lq_decomposition<T> LQ(A);
        gslmm::matrix<T> Ainv(A.row_size(), B.column_size());
        for (size_t i = 0; i < A.row_size(); i++) {
          gslmm::vector<T> v(A.row_size());
          gslmm::vector<T> x(A.row_size());
          v.basis(i);
          LQ.solve(v, x);
          for (size_t j = 0; j < A.column_size(); j++) Ainv(i,j) = x[j];
        }
Parameters:
b The vector $ b$
x The vector $ x$
Returns:
true on success

Reimplemented in gslmm::lq_decomposition_pivot< double >.

void gslmm::lq_decomposition< double >::unpack  )  const [inline, protected]
 

Unpack LQ into $ Q$ and $ R$.

bool gslmm::lq_decomposition< double >::update const vector_type w,
const vector_type v
[inline, virtual]
 

This function performs a rank-1 update $ w v^T$ of the.

Update a LQ factorisation for $ A= L Q$, $ A' = A + v u^T$,

\begin{eqnarray*} L' Q' = LQ + v u^T\\ = (L + v u^T Q^T) Q\\ = (L + v w^T) Q \end{eqnarray*}

$ where w = Q u$. Algorithm from Golub and Van Loan, "Matrix Computations", Section 12.5 (Updating Matrix Factorizations, Rank-One Changes).

Reimplemented in gslmm::lq_decomposition_pivot< double >.


Member Data Documentation

matrix_type* gslmm::lq_decomposition< double >::_l [mutable, protected]
 

The R of the LQ decomposition of the input matrix.

matrix_type* gslmm::lq_decomposition< double >::_q [mutable, protected]
 

The Q of the LQ decomposition of the input matrix.

vector_type gslmm::lq_decomposition< double >::_tau [protected]
 

The $ tau$ of the LQ decomposition.


The documentation for this struct was generated from the following file:
Top of page Last update Tue May 9 10:11:31 2006
Christian Holm
Created by DoxyGen 1.4.6