GeneralBrokenLines  trunk_rev46
 All Classes Files Functions Variables Typedefs Pages
Public Member Functions | Private Member Functions | Private Attributes
BorderedBandMatrix Class Reference

(Symmetric) Bordered Band Matrix. More...

#include <BorderedBandMatrix.h>

List of all members.

Public Member Functions

 BorderedBandMatrix ()
 Create bordered band matrix.
virtual ~BorderedBandMatrix ()
void resize (unsigned int nSize, unsigned int nBorder=1, unsigned int nBand=5)
 Resize bordered band matrix.
void solveAndInvertBorderedBand (const VVector &aRightHandSide, VVector &aSolution)
 Solve linear equation system, partially calculate inverse.
void addBlockMatrix (double aWeight, const std::vector< unsigned int > *anIndex, const std::vector< double > *aVector)
 Add symmetric block matrix.
TMatrixDSym getBlockMatrix (const std::vector< unsigned int > anIndex) const
 Retrieve symmetric block matrix.
void printMatrix () const
 Print bordered band matrix.

Private Member Functions

void decomposeBand ()
 (root free) Cholesky decomposition of band part: C=LDL^T
VVector solveBand (const VVector &aRightHandSide) const
 Solve for band part.
VMatrix solveBand (const VMatrix &aRightHandSide) const
 solve band part for mixed part (border rows).
VMatrix invertBand ()
 Invert band part.
VMatrix bandOfAVAT (const VMatrix &anArray, const VSymMatrix &aSymArray) const
 Calculate band part of: 'anArray * aSymArray * anArray.T'.

Private Attributes

unsigned int numSize
 Matrix size.
unsigned int numBorder
 Border size.
unsigned int numBand
 Band width.
unsigned int numCol
 Band matrix size.
VSymMatrix theBorder
 Border part.
VMatrix theMixed
 Mixed part.
VMatrix theBand
 Band part.

Detailed Description

(Symmetric) Bordered Band Matrix.

Separate storage of border, mixed and band parts (as vector<double>).

 *  Example for matrix size=8 with border size and band width of two
 *
 *     +-                                 -+
 *     |  B11 B12 M13 M14 M15 M16 M17 M18  |
 *     |  B12 B22 M23 M24 M25 M26 M27 M28  |
 *     |  M13 M23 C33 C34 C35  0.  0.  0.  |
 *     |  M14 M24 C34 C44 C45 C46  0.  0.  |
 *     |  M15 M25 C35 C45 C55 C56 C57  0.  |
 *     |  M16 M26  0. C46 C56 C66 C67 C68  |
 *     |  M17 M27  0.  0. C57 C67 C77 C78  |
 *     |  M18 M28  0.  0.  0. C68 C78 C88  |
 *     +-                                 -+
 *
 *  Is stored as::
 *
 *     +-         -+     +-                         -+
 *     |  B11 B12  |     |  M13 M14 M15 M16 M17 M18  |
 *     |  B12 B22  |     |  M23 M24 M25 M26 M27 M28  |
 *     +-         -+     +-                         -+
 *
 *                       +-                         -+
 *                       |  C33 C44 C55 C66 C77 C88  |
 *                       |  C34 C45 C56 C67 C78  0.  |
 *                       |  C35 C46 C57 C68  0.  0.  |
 *                       +-                         -+
 *

Definition at line 53 of file BorderedBandMatrix.h.


Constructor & Destructor Documentation

BorderedBandMatrix::BorderedBandMatrix ( )

Create bordered band matrix.

Definition at line 11 of file BorderedBandMatrix.cpp.

BorderedBandMatrix::~BorderedBandMatrix ( )
virtual

Definition at line 14 of file BorderedBandMatrix.cpp.


Member Function Documentation

void BorderedBandMatrix::addBlockMatrix ( double  aWeight,
const std::vector< unsigned int > *  anIndex,
const std::vector< double > *  aVector 
)

Add symmetric block matrix.

Add (extended) block matrix defined by 'aVector * aWeight * aVector.T' to bordered band matrix: BBmatrix(anIndex(i),anIndex(j)) += aVector(i) * aWeight * aVector(j).

Parameters:
aWeight[in] Weight
anIndex[in] List of rows/colums to be used
aVector[in] Vector

Definition at line 43 of file BorderedBandMatrix.cpp.

References numBand, numBorder, theBand, theBorder, and theMixed.

Referenced by GblTrajectory::buildLinearEquationSystem().

VMatrix BorderedBandMatrix::bandOfAVAT ( const VMatrix anArray,
const VSymMatrix aSymArray 
) const
private

Calculate band part of: 'anArray * aSymArray * anArray.T'.

Returns:
Band part of product

Definition at line 291 of file BorderedBandMatrix.cpp.

References numBand, numBorder, and numCol.

Referenced by solveAndInvertBorderedBand().

void BorderedBandMatrix::decomposeBand ( )
private

(root free) Cholesky decomposition of band part: C=LDL^T

Decompose band matrix into diagonal matrix D and lower triangular band matrix L (diagonal=1). Overwrite band matrix with D and off-diagonal part of L.

Exceptions:
2: matrix is singular.
3: matrix is not positive definite.

Definition at line 174 of file BorderedBandMatrix.cpp.

References numBand, numCol, and theBand.

Referenced by solveAndInvertBorderedBand().

TMatrixDSym BorderedBandMatrix::getBlockMatrix ( const std::vector< unsigned int >  anIndex) const

Retrieve symmetric block matrix.

Get (compressed) block from bordered band matrix: aMatrix(i,j) = BBmatrix(anIndex(i),anIndex(j)).

Parameters:
anIndex[in] List of rows/colums to be used

Definition at line 72 of file BorderedBandMatrix.cpp.

References numBorder, theBand, theBorder, and theMixed.

Referenced by GblTrajectory::getResAndErr(), and GblTrajectory::getResults().

VMatrix BorderedBandMatrix::invertBand ( )
private

Invert band part.

Returns:
Inverted band

Definition at line 267 of file BorderedBandMatrix.cpp.

References numBand, numCol, and theBand.

Referenced by solveAndInvertBorderedBand().

void BorderedBandMatrix::printMatrix ( ) const

Print bordered band matrix.

Definition at line 155 of file BorderedBandMatrix.cpp.

References VMatrix::print(), VSymMatrix::print(), theBand, theBorder, and theMixed.

void BorderedBandMatrix::resize ( unsigned int  nSize,
unsigned int  nBorder = 1,
unsigned int  nBand = 5 
)

Resize bordered band matrix.

Parameters:
nSize[in] Size of matrix
nBorder[in] Size of border (=1 for q/p + additional local parameters)
nBand[in] Band width (usually = 5, for simplified jacobians = 4)

Definition at line 23 of file BorderedBandMatrix.cpp.

References numBand, numBorder, numCol, numSize, VMatrix::resize(), VSymMatrix::resize(), theBand, theBorder, and theMixed.

Referenced by GblTrajectory::buildLinearEquationSystem().

void BorderedBandMatrix::solveAndInvertBorderedBand ( const VVector aRightHandSide,
VVector aSolution 
)

Solve linear equation system, partially calculate inverse.

Solve linear equation A*x=b system with bordered band matrix A, calculate bordered band part of inverse of A. Use decomposition in border and band part for block matrix algebra:

| A  Ct |   | x1 |   | b1 |        , A  is the border part
|       | * |    | = |    |        , Ct is the mixed part
| C  D  |   | x2 |   | b2 |        , D  is the band part

Explicit inversion of D is avoided by using solution X of D*X=C (X=D^-1*C, obtained from Cholesky decomposition and forward/backward substitution)

| x1 |   | E*b1 - E*Xt*b2 |        , E^-1 = A-Ct*D^-1*C = A-Ct*X
|    | = |                |
| x2 |   |  x   - X*x1    |        , x is solution of D*x=b2 (x=D^-1*b2)

Inverse matrix is:

|  E   -E*Xt          |
|                     |            , only band part of (D^-1 + X*E*Xt)
| -X*E  D^-1 + X*E*Xt |              is calculated
Parameters:
[in]aRightHandSideRight hand side (vector) 'b' of A*x=b
[out]aSolutionSolution (vector) x of A*x=b

Definition at line 122 of file BorderedBandMatrix.cpp.

References bandOfAVAT(), decomposeBand(), VVector::getVec(), VSymMatrix::invert(), invertBand(), numBorder, numCol, VVector::putVec(), solveBand(), theBand, theBorder, theMixed, and VMatrix::transpose().

Referenced by GblTrajectory::fit().

VVector BorderedBandMatrix::solveBand ( const VVector aRightHandSide) const
private

Solve for band part.

Solve C*x=b for band part using decomposition C=LDL^T and forward (L*z=b) and backward substitution (L^T*x=D^-1*z).

Parameters:
[in]aRightHandSideRight hand side (vector) 'b' of C*x=b
Returns:
Solution (vector) 'x' of C*x=b

Definition at line 209 of file BorderedBandMatrix.cpp.

References VMatrix::getNumCols(), VMatrix::getNumRows(), and theBand.

Referenced by solveAndInvertBorderedBand().

VMatrix BorderedBandMatrix::solveBand ( const VMatrix aRightHandSide) const
private

solve band part for mixed part (border rows).

Solve C*X=B for mixed part using decomposition C=LDL^T and forward and backward substitution.

Parameters:
[in]aRightHandSideRight hand side (matrix) 'B' of C*X=B
Returns:
Solution (matrix) 'X' of C*X=B

Definition at line 238 of file BorderedBandMatrix.cpp.

References VMatrix::getNumCols(), VMatrix::getNumRows(), numBorder, and theBand.


Member Data Documentation

unsigned int BorderedBandMatrix::numBand
private

Band width.

Definition at line 70 of file BorderedBandMatrix.h.

Referenced by addBlockMatrix(), bandOfAVAT(), decomposeBand(), invertBand(), and resize().

unsigned int BorderedBandMatrix::numBorder
private
unsigned int BorderedBandMatrix::numCol
private

Band matrix size.

Definition at line 71 of file BorderedBandMatrix.h.

Referenced by bandOfAVAT(), decomposeBand(), invertBand(), resize(), and solveAndInvertBorderedBand().

unsigned int BorderedBandMatrix::numSize
private

Matrix size.

Definition at line 68 of file BorderedBandMatrix.h.

Referenced by resize().

VMatrix BorderedBandMatrix::theBand
private
VSymMatrix BorderedBandMatrix::theBorder
private

Border part.

Definition at line 72 of file BorderedBandMatrix.h.

Referenced by addBlockMatrix(), getBlockMatrix(), printMatrix(), resize(), and solveAndInvertBorderedBand().

VMatrix BorderedBandMatrix::theMixed
private

Mixed part.

Definition at line 73 of file BorderedBandMatrix.h.

Referenced by addBlockMatrix(), getBlockMatrix(), printMatrix(), resize(), and solveAndInvertBorderedBand().


The documentation for this class was generated from the following files: