CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
gbl::BorderedBandMatrix Class Reference

(Symmetric) Bordered Band Matrix. More...

#include <BorderedBandMatrix.h>

Public Member Functions

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

Private Member Functions

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

Private Attributes

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

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 56 of file BorderedBandMatrix.h.

Constructor & Destructor Documentation

gbl::BorderedBandMatrix::BorderedBandMatrix ( )

Create bordered band matrix.

Definition at line 14 of file BorderedBandMatrix.cc.

14  : numSize(0), numBorder(0), numBand(0), numCol(0) {
15 }
unsigned int numCol
Band matrix size.
unsigned int numBorder
Border size.
unsigned int numSize
Matrix size.
unsigned int numBand
Band width.
gbl::BorderedBandMatrix::~BorderedBandMatrix ( )
virtual

Definition at line 17 of file BorderedBandMatrix.cc.

17  {
18 }

Member Function Documentation

void gbl::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 46 of file BorderedBandMatrix.cc.

References i, j, bookConverter::max, numBand, numBorder, theBand, theBorder, and theMixed.

Referenced by gbl::GblTrajectory::buildLinearEquationSystem().

48  {
49  int nBorder = numBorder;
50  for (unsigned int i = 0; i < anIndex->size(); ++i) {
51  int iIndex = (*anIndex)[i] - 1; // anIndex has to be sorted
52  for (unsigned int j = 0; j <= i; ++j) {
53  int jIndex = (*anIndex)[j] - 1;
54  if (iIndex < nBorder) {
55  theBorder(iIndex, jIndex) += (*aVector)[i] * aWeight
56  * (*aVector)[j];
57  } else if (jIndex < nBorder) {
58  theMixed(jIndex, iIndex - nBorder) += (*aVector)[i] * aWeight
59  * (*aVector)[j];
60  } else {
61  unsigned int nBand = iIndex - jIndex;
62  theBand(nBand, jIndex - nBorder) += (*aVector)[i] * aWeight
63  * (*aVector)[j];
64  numBand = std::max(numBand, nBand); // update band width
65  }
66  }
67  }
68 }
int i
Definition: DBlmapReader.cc:9
VMatrix theMixed
Mixed part.
int j
Definition: DBlmapReader.cc:9
unsigned int numBorder
Border size.
VSymMatrix theBorder
Border part.
VMatrix theBand
Band part.
unsigned int numBand
Band width.
VMatrix gbl::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 294 of file BorderedBandMatrix.cc.

References i, j, relval_steps::k, prof2calltree::l, bookConverter::max, numBand, numBorder, and numCol.

Referenced by solveAndInvertBorderedBand().

295  {
296  int nBand = numBand;
297  int nCol = numCol;
298  int nBorder = numBorder;
299  double sum;
300  VMatrix aBand((nBand + 1), nCol);
301  for (int i = 0; i < nCol; ++i) {
302  for (int j = std::max(0, i - nBand); j <= i; ++j) {
303  sum = 0.;
304  for (int l = 0; l < nBorder; ++l) { // diagonal
305  sum += anArray(i, l) * aSymArray(l, l) * anArray(j, l);
306  for (int k = 0; k < l; ++k) { // off diagonal
307  sum += anArray(i, l) * aSymArray(l, k) * anArray(j, k)
308  + anArray(i, k) * aSymArray(l, k) * anArray(j, l);
309  }
310  }
311  aBand(i - j, j) = sum;
312  }
313  }
314  return aBand;
315 }
int i
Definition: DBlmapReader.cc:9
unsigned int numCol
Band matrix size.
int j
Definition: DBlmapReader.cc:9
unsigned int numBorder
Border size.
unsigned int numBand
Band width.
void gbl::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 177 of file BorderedBandMatrix.cc.

References i, j, relval_steps::k, min(), numBand, numCol, and theBand.

Referenced by solveAndInvertBorderedBand().

177  {
178 
179  int nRow = numBand + 1;
180  int nCol = numCol;
181  VVector auxVec(nCol);
182  for (int i = 0; i < nCol; ++i) {
183  auxVec(i) = theBand(0, i) * 16.0; // save diagonal elements
184  }
185  for (int i = 0; i < nCol; ++i) {
186  if ((theBand(0, i) + auxVec(i)) != theBand(0, i)) {
187  theBand(0, i) = 1.0 / theBand(0, i);
188  if (theBand(0, i) < 0.) {
189  throw 3; // not positive definite
190  }
191  } else {
192  theBand(0, i) = 0.0;
193  throw 2; // singular
194  }
195  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
196  double rxw = theBand(j, i) * theBand(0, i);
197  for (int k = 0; k < std::min(nRow, nCol - i) - j; ++k) {
198  theBand(k, i + j) -= theBand(k + j, i) * rxw;
199  }
200  theBand(j, i) = rxw;
201  }
202  }
203 }
int i
Definition: DBlmapReader.cc:9
unsigned int numCol
Band matrix size.
int j
Definition: DBlmapReader.cc:9
T min(T a, T b)
Definition: MathUtil.h:58
VMatrix theBand
Band part.
unsigned int numBand
Band width.
TMatrixDSym gbl::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 75 of file BorderedBandMatrix.cc.

References i, j, numBorder, theBand, theBorder, and theMixed.

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

76  {
77 
78  TMatrixDSym aMatrix(anIndex.size());
79  int nBorder = numBorder;
80  for (unsigned int i = 0; i < anIndex.size(); ++i) {
81  int iIndex = anIndex[i] - 1; // anIndex has to be sorted
82  for (unsigned int j = 0; j <= i; ++j) {
83  int jIndex = anIndex[j] - 1;
84  if (iIndex < nBorder) {
85  aMatrix(i, j) = theBorder(iIndex, jIndex); // border part of inverse
86  } else if (jIndex < nBorder) {
87  aMatrix(i, j) = -theMixed(jIndex, iIndex - nBorder); // mixed part of inverse
88  } else {
89  unsigned int nBand = iIndex - jIndex;
90  aMatrix(i, j) = theBand(nBand, jIndex - nBorder); // band part of inverse
91  }
92  aMatrix(j, i) = aMatrix(i, j);
93  }
94  }
95  return aMatrix;
96 }
int i
Definition: DBlmapReader.cc:9
VMatrix theMixed
Mixed part.
int j
Definition: DBlmapReader.cc:9
unsigned int numBorder
Border size.
VSymMatrix theBorder
Border part.
VMatrix theBand
Band part.
VMatrix gbl::BorderedBandMatrix::invertBand ( )
private

Invert band part.

Returns
Inverted band

Definition at line 270 of file BorderedBandMatrix.cc.

References funct::abs(), i, j, relval_steps::k, bookConverter::max, min(), numBand, numCol, and theBand.

Referenced by solveAndInvertBorderedBand().

270  {
271 
272  int nRow = numBand + 1;
273  int nCol = numCol;
274  VMatrix inverseBand(nRow, nCol);
275 
276  for (int i = nCol - 1; i >= 0; i--) {
277  double rxw = theBand(0, i);
278  for (int j = i; j >= std::max(0, i - nRow + 1); j--) {
279  for (int k = j + 1; k < std::min(nCol, j + nRow); ++k) {
280  rxw -= inverseBand(abs(i - k), std::min(i, k))
281  * theBand(k - j, j);
282  }
283  inverseBand(i - j, j) = rxw;
284  rxw = 0.;
285  }
286  }
287  return inverseBand;
288 }
int i
Definition: DBlmapReader.cc:9
unsigned int numCol
Band matrix size.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int j
Definition: DBlmapReader.cc:9
T min(T a, T b)
Definition: MathUtil.h:58
VMatrix theBand
Band part.
unsigned int numBand
Band width.
void gbl::BorderedBandMatrix::printMatrix ( ) const

Print bordered band matrix.

Definition at line 158 of file BorderedBandMatrix.cc.

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

Referenced by gbl::GblTrajectory::printTrajectory().

158  {
159  std::cout << "Border part " << std::endl;
160  theBorder.print();
161  std::cout << "Mixed part " << std::endl;
162  theMixed.print();
163  std::cout << "Band part " << std::endl;
164  theBand.print();
165 }
void print() const
Print matrix.
Definition: VMatrix.cc:167
VMatrix theMixed
Mixed part.
VSymMatrix theBorder
Border part.
void print() const
Print matrix.
Definition: VMatrix.cc:70
VMatrix theBand
Band part.
tuple cout
Definition: gather_cfg.py:121
void gbl::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 26 of file BorderedBandMatrix.cc.

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

Referenced by Vispa.Gui.TextDialog.TextDialog::__init__(), Vispa.Plugins.ConfigEditor.ToolDialog.ToolDialog::__init__(), Vispa.Main.MainWindow.MainWindow::_loadIni(), gbl::GblTrajectory::buildLinearEquationSystem(), and Vispa.Gui.PortConnection.PointToPointConnection::updateConnection().

27  {
28  numSize = nSize;
29  numBorder = nBorder;
30  numCol = nSize - nBorder;
31  numBand = 0;
34  theBand.resize((nBand + 1), numCol);
35 }
void resize(const unsigned int nRows, const unsigned int nCols)
Resize Matrix.
Definition: VMatrix.cc:33
unsigned int numCol
Band matrix size.
void resize(const unsigned int nRows)
Resize symmetric matrix.
Definition: VMatrix.cc:153
VMatrix theMixed
Mixed part.
unsigned int numBorder
Border size.
VSymMatrix theBorder
Border part.
VMatrix theBand
Band part.
unsigned int numSize
Matrix size.
unsigned int numBand
Band width.
void gbl::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 125 of file BorderedBandMatrix.cc.

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

Referenced by gbl::GblTrajectory::fit().

126  {
127 
128  // decompose band
129  decomposeBand();
130  // invert band
131  VMatrix inverseBand = invertBand();
132  if (numBorder > 0) { // need to use block matrix decomposition to solve
133  // solve for mixed part
134  const VMatrix auxMat = solveBand(theMixed); // = Xt
135  const VMatrix auxMatT = auxMat.transpose(); // = X
136  // solve for border part
137  const VVector auxVec = aRightHandSide.getVec(numBorder)
138  - auxMat * aRightHandSide.getVec(numCol, numBorder); // = b1 - Xt*b2
139  VSymMatrix inverseBorder = theBorder - theMixed * auxMatT;
140  inverseBorder.invert(); // = E
141  const VVector borderSolution = inverseBorder * auxVec; // = x1
142  // solve for band part
143  const VVector bandSolution = solveBand(
144  aRightHandSide.getVec(numCol, numBorder)); // = x
145  aSolution.putVec(borderSolution);
146  aSolution.putVec(bandSolution - auxMatT * borderSolution, numBorder); // = x2
147  // parts of inverse
148  theBorder = inverseBorder; // E
149  theMixed = inverseBorder * auxMat; // E*Xt (-mixed part of inverse) !!!
150  theBand = inverseBand + bandOfAVAT(auxMatT, inverseBorder); // band(D^-1 + X*E*Xt)
151  } else {
152  aSolution.putVec(solveBand(aRightHandSide));
153  theBand = inverseBand;
154  }
155 }
void putVec(const VVector &aVector, unsigned int start=0)
Put part of vector.
Definition: VMatrix.cc:262
unsigned int numCol
Band matrix size.
VMatrix invertBand()
Invert band part.
VMatrix theMixed
Mixed part.
void decomposeBand()
(root free) Cholesky decomposition of band part: C=LDL^T
unsigned int numBorder
Border size.
VSymMatrix theBorder
Border part.
VMatrix theBand
Band part.
VMatrix bandOfAVAT(const VMatrix &anArray, const VSymMatrix &aSymArray) const
Calculate band part of: &#39;anArray * aSymArray * anArray.T&#39;.
VVector solveBand(const VVector &aRightHandSide) const
Solve for band part.
VVector gbl::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 212 of file BorderedBandMatrix.cc.

References gbl::VMatrix::getNumCols(), gbl::VMatrix::getNumRows(), i, j, min(), and theBand.

Referenced by solveAndInvertBorderedBand().

212  {
213 
214  int nRow = theBand.getNumRows();
215  int nCol = theBand.getNumCols();
216  VVector aSolution(aRightHandSide);
217  for (int i = 0; i < nCol; ++i) // forward substitution
218  {
219  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
220  aSolution(j + i) -= theBand(j, i) * aSolution(i);
221  }
222  }
223  for (int i = nCol - 1; i >= 0; i--) // backward substitution
224  {
225  double rxw = theBand(0, i) * aSolution(i);
226  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
227  rxw -= theBand(j, i) * aSolution(j + i);
228  }
229  aSolution(i) = rxw;
230  }
231  return aSolution;
232 }
int i
Definition: DBlmapReader.cc:9
unsigned int getNumRows() const
Get number of rows.
Definition: VMatrix.cc:57
int j
Definition: DBlmapReader.cc:9
T min(T a, T b)
Definition: MathUtil.h:58
unsigned int getNumCols() const
Get number of columns.
Definition: VMatrix.cc:65
VMatrix theBand
Band part.
VMatrix gbl::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 241 of file BorderedBandMatrix.cc.

References gbl::VMatrix::getNumCols(), gbl::VMatrix::getNumRows(), i, j, min(), numBorder, and theBand.

241  {
242 
243  int nRow = theBand.getNumRows();
244  int nCol = theBand.getNumCols();
245  VMatrix aSolution(aRightHandSide);
246  for (unsigned int iBorder = 0; iBorder < numBorder; iBorder++) {
247  for (int i = 0; i < nCol; ++i) // forward substitution
248  {
249  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
250  aSolution(iBorder, j + i) -= theBand(j, i)
251  * aSolution(iBorder, i);
252  }
253  }
254  for (int i = nCol - 1; i >= 0; i--) // backward substitution
255  {
256  double rxw = theBand(0, i) * aSolution(iBorder, i);
257  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
258  rxw -= theBand(j, i) * aSolution(iBorder, j + i);
259  }
260  aSolution(iBorder, i) = rxw;
261  }
262  }
263  return aSolution;
264 }
int i
Definition: DBlmapReader.cc:9
unsigned int getNumRows() const
Get number of rows.
Definition: VMatrix.cc:57
int j
Definition: DBlmapReader.cc:9
unsigned int numBorder
Border size.
T min(T a, T b)
Definition: MathUtil.h:58
unsigned int getNumCols() const
Get number of columns.
Definition: VMatrix.cc:65
VMatrix theBand
Band part.

Member Data Documentation

unsigned int gbl::BorderedBandMatrix::numBand
private

Band width.

Definition at line 73 of file BorderedBandMatrix.h.

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

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

Band matrix size.

Definition at line 74 of file BorderedBandMatrix.h.

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

unsigned int gbl::BorderedBandMatrix::numSize
private

Matrix size.

Definition at line 71 of file BorderedBandMatrix.h.

Referenced by resize().

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

Border part.

Definition at line 75 of file BorderedBandMatrix.h.

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

VMatrix gbl::BorderedBandMatrix::theMixed
private

Mixed part.

Definition at line 76 of file BorderedBandMatrix.h.

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