CMS 3D CMS Logo

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...
 
void addBlockMatrix (double aWeight, unsigned int nSimple, unsigned int *anIndex, double *aVector)
 Add symmetric block matrix. More...
 
 BorderedBandMatrix ()
 Create bordered band matrix. More...
 
Eigen::MatrixXd getBlockMatrix (const std::vector< unsigned int > anIndex) const
 Retrieve symmetric block matrix. More...
 
Eigen::MatrixXd getBlockMatrix (unsigned int aSize, 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 76 of file BorderedBandMatrix.h.

Constructor & Destructor Documentation

gbl::BorderedBandMatrix::BorderedBandMatrix ( )

Create bordered band matrix.

Definition at line 37 of file BorderedBandMatrix.cc.

37  :
38  numSize(0), numBorder(0), numBand(0), numCol(0) {
39  }
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 41 of file BorderedBandMatrix.cc.

41  {
42  }

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 70 of file BorderedBandMatrix.cc.

References mps_fire::i, hpstanc_transforms::max, numBand, numBorder, theBand, theBorder, and theMixed.

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

72  {
73  int nBorder = numBorder;
74  for (unsigned int i = 0; i < anIndex->size(); ++i) {
75  int iIndex = (*anIndex)[i] - 1; // anIndex has to be sorted
76  for (unsigned int j = 0; j <= i; ++j) {
77  int jIndex = (*anIndex)[j] - 1;
78  if (iIndex < nBorder) {
79  theBorder(iIndex, jIndex) += (*aVector)[i] * aWeight
80  * (*aVector)[j];
81  } else if (jIndex < nBorder) {
82  theMixed(jIndex, iIndex - nBorder) += (*aVector)[i] * aWeight
83  * (*aVector)[j];
84  } else {
85  unsigned int nBand = iIndex - jIndex;
86  theBand(nBand, jIndex - nBorder) += (*aVector)[i] * aWeight
87  * (*aVector)[j];
88  numBand = std::max(numBand, nBand); // update band width
89  }
90  }
91  }
92  }
VMatrix theMixed
Mixed part.
unsigned int numBorder
Border size.
VSymMatrix theBorder
Border part.
VMatrix theBand
Band part.
unsigned int numBand
Band width.
void gbl::BorderedBandMatrix::addBlockMatrix ( double  aWeight,
unsigned int  aSize,
unsigned int *  anIndex,
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
aSize[in] Size of block matrix
anIndex[in] List of rows/colums to be used
aVector[in] Vector

Definition at line 104 of file BorderedBandMatrix.cc.

References mps_fire::i, hpstanc_transforms::max, numBand, numBorder, theBand, theBorder, and theMixed.

105  {
106  int nBorder = numBorder;
107  for (unsigned int i = 0; i < aSize; ++i) {
108  int iIndex = anIndex[i] - 1; // anIndex has to be sorted
109  for (unsigned int j = 0; j <= i; ++j) {
110  int jIndex = anIndex[j] - 1;
111  if (iIndex < nBorder) {
112  theBorder(iIndex, jIndex) += aVector[i] * aWeight * aVector[j];
113  } else if (jIndex < nBorder) {
114  theMixed(jIndex, iIndex - nBorder) += aVector[i] * aWeight
115  * aVector[j];
116  } else {
117  unsigned int nBand = iIndex - jIndex;
118  theBand(nBand, jIndex - nBorder) += aVector[i] * aWeight
119  * aVector[j];
120  numBand = std::max(numBand, nBand); // update band width
121  }
122  }
123  }
124  }
VMatrix theMixed
Mixed part.
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 379 of file BorderedBandMatrix.cc.

References mps_fire::i, gen::k, checklumidiff::l, hpstanc_transforms::max, numBand, numBorder, and numCol.

Referenced by solveAndInvertBorderedBand().

380  {
381  int nBand = numBand;
382  int nCol = numCol;
383  int nBorder = numBorder;
384  double sum;
385  VMatrix aBand((nBand + 1), nCol);
386  for (int i = 0; i < nCol; ++i) {
387  for (int j = std::max(0, i - nBand); j <= i; ++j) {
388  sum = 0.;
389  for (int l = 0; l < nBorder; ++l) { // diagonal
390  sum += anArray(i, l) * aSymArray(l, l) * anArray(j, l);
391  for (int k = 0; k < l; ++k) { // off diagonal
392  sum += anArray(i, l) * aSymArray(l, k) * anArray(j, k)
393  + anArray(i, k) * aSymArray(l, k) * anArray(j, l);
394  }
395  }
396  aBand(i - j, j) = sum;
397  }
398  }
399  return aBand;
400  }
unsigned int numCol
Band matrix size.
unsigned int numBorder
Border size.
int k[5][pyjets_maxn]
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 262 of file BorderedBandMatrix.cc.

References mps_fire::i, gen::k, min(), numBand, numCol, and theBand.

Referenced by solveAndInvertBorderedBand().

262  {
263 
264  int nRow = numBand + 1;
265  int nCol = numCol;
266  VVector auxVec(nCol);
267  for (int i = 0; i < nCol; ++i) {
268  auxVec(i) = theBand(0, i) * 16.0; // save diagonal elements
269  }
270  for (int i = 0; i < nCol; ++i) {
271  if ((theBand(0, i) + auxVec(i)) != theBand(0, i)) {
272  theBand(0, i) = 1.0 / theBand(0, i);
273  if (theBand(0, i) < 0.) {
274  throw 3; // not positive definite
275  }
276  } else {
277  theBand(0, i) = 0.0;
278  throw 2; // singular
279  }
280  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
281  double rxw = theBand(j, i) * theBand(0, i);
282  for (int k = 0; k < std::min(nRow, nCol - i) - j; ++k) {
283  theBand(k, i + j) -= theBand(k + j, i) * rxw;
284  }
285  theBand(j, i) = rxw;
286  }
287  }
288  }
unsigned int numCol
Band matrix size.
T min(T a, T b)
Definition: MathUtil.h:58
int k[5][pyjets_maxn]
VMatrix theBand
Band part.
unsigned int numBand
Band width.
MatrixXd 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 131 of file BorderedBandMatrix.cc.

References mps_fire::i, numBorder, theBand, theBorder, and theMixed.

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

132  {
133 
134  MatrixXd aMatrix(anIndex.size(), anIndex.size());
135  int nBorder = numBorder;
136  for (unsigned int i = 0; i < anIndex.size(); ++i) {
137  int iIndex = anIndex[i] - 1; // anIndex has to be sorted
138  for (unsigned int j = 0; j <= i; ++j) {
139  int jIndex = anIndex[j] - 1;
140  if (iIndex < nBorder) {
141  aMatrix(i, j) = theBorder(iIndex, jIndex); // border part of inverse
142  } else if (jIndex < nBorder) {
143  aMatrix(i, j) = -theMixed(jIndex, iIndex - nBorder); // mixed part of inverse
144  } else {
145  unsigned int nBand = iIndex - jIndex;
146  aMatrix(i, j) = theBand(nBand, jIndex - nBorder); // band part of inverse
147  }
148  aMatrix(j, i) = aMatrix(i, j);
149  }
150  }
151  return aMatrix;
152  }
VMatrix theMixed
Mixed part.
unsigned int numBorder
Border size.
VSymMatrix theBorder
Border part.
VMatrix theBand
Band part.
MatrixXd gbl::BorderedBandMatrix::getBlockMatrix ( unsigned int  aSize,
unsigned int *  anIndex 
) const

Retrieve symmetric block matrix.

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

Parameters
aSize[in] Matrix size
anIndex[in] Array of rows/colums to be used

Definition at line 160 of file BorderedBandMatrix.cc.

References mps_fire::i, numBorder, theBand, theBorder, and theMixed.

161  {
162 
163  MatrixXd aMatrix(aSize, aSize);
164  int nBorder = numBorder;
165  for (unsigned int i = 0; i < aSize; ++i) {
166  int iIndex = anIndex[i] - 1; // anIndex has to be sorted
167  for (unsigned int j = 0; j <= i; ++j) {
168  int jIndex = anIndex[j] - 1;
169  if (iIndex < nBorder) {
170  aMatrix(i, j) = theBorder(iIndex, jIndex); // border part of inverse
171  } else if (jIndex < nBorder) {
172  aMatrix(i, j) = -theMixed(jIndex, iIndex - nBorder); // mixed part of inverse
173  } else {
174  unsigned int nBand = iIndex - jIndex;
175  aMatrix(i, j) = theBand(nBand, jIndex - nBorder); // band part of inverse
176  }
177  aMatrix(j, i) = aMatrix(i, j);
178  }
179  }
180  return aMatrix;
181  }
VMatrix theMixed
Mixed part.
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 355 of file BorderedBandMatrix.cc.

References funct::abs(), mps_fire::i, gen::k, hpstanc_transforms::max, min(), numBand, numCol, and theBand.

Referenced by solveAndInvertBorderedBand().

355  {
356 
357  int nRow = numBand + 1;
358  int nCol = numCol;
359  VMatrix inverseBand(nRow, nCol);
360 
361  for (int i = nCol - 1; i >= 0; i--) {
362  double rxw = theBand(0, i);
363  for (int j = i; j >= std::max(0, i - nRow + 1); j--) {
364  for (int k = j + 1; k < std::min(nCol, j + nRow); ++k) {
365  rxw -= inverseBand(abs(i - k), std::min(i, k))
366  * theBand(k - j, j);
367  }
368  inverseBand(i - j, j) = rxw;
369  rxw = 0.;
370  }
371  }
372  return inverseBand;
373  }
unsigned int numCol
Band matrix size.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
T min(T a, T b)
Definition: MathUtil.h:58
int k[5][pyjets_maxn]
VMatrix theBand
Band part.
unsigned int numBand
Band width.
void gbl::BorderedBandMatrix::printMatrix ( ) const

Print bordered band matrix.

Definition at line 243 of file BorderedBandMatrix.cc.

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

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

243  {
244  std::cout << "Border part " << std::endl;
245  theBorder.print();
246  std::cout << "Mixed part " << std::endl;
247  theMixed.print();
248  std::cout << "Band part " << std::endl;
249  theBand.print();
250  }
void print() const
Print matrix.
Definition: VMatrix.cc:189
VMatrix theMixed
Mixed part.
VSymMatrix theBorder
Border part.
void print() const
Print matrix.
Definition: VMatrix.cc:92
VMatrix theBand
Band part.
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 50 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().

51  {
52  numSize = nSize;
53  numBorder = nBorder;
54  numCol = nSize - nBorder;
55  numBand = 0;
58  theBand.resize((nBand + 1), numCol);
59  }
void resize(const unsigned int nRows, const unsigned int nCols)
Resize Matrix.
Definition: VMatrix.cc:55
unsigned int numCol
Band matrix size.
void resize(const unsigned int nRows)
Resize symmetric matrix.
Definition: VMatrix.cc:175
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 210 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().

211  {
212 
213  // decompose band
214  decomposeBand();
215  // invert band
216  VMatrix inverseBand = invertBand();
217  if (numBorder > 0) { // need to use block matrix decomposition to solve
218  // solve for mixed part
219  const VMatrix auxMat = solveBand(theMixed);// = Xt
220  const VMatrix auxMatT = auxMat.transpose();// = X
221  // solve for border part
222  const VVector auxVec = aRightHandSide.getVec(numBorder)
223  - auxMat * aRightHandSide.getVec(numCol, numBorder);// = b1 - Xt*b2
224  VSymMatrix inverseBorder = theBorder - theMixed * auxMatT;
225  inverseBorder.invert();// = E
226  const VVector borderSolution = inverseBorder * auxVec;// = x1
227  // solve for band part
228  const VVector bandSolution = solveBand(
229  aRightHandSide.getVec(numCol, numBorder));// = x
230  aSolution.putVec(borderSolution);
231  aSolution.putVec(bandSolution - auxMatT * borderSolution, numBorder);// = x2
232  // parts of inverse
233  theBorder = inverseBorder;// E
234  theMixed = inverseBorder * auxMat;// E*Xt (-mixed part of inverse) !!!
235  theBand = inverseBand + bandOfAVAT(auxMatT, inverseBorder);// band(D^-1 + X*E*Xt)
236  } else {
237  aSolution.putVec(solveBand(aRightHandSide));
238  theBand = inverseBand;
239  }
240  }
void putVec(const VVector &aVector, unsigned int start=0)
Put part of vector.
Definition: VMatrix.cc:284
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 297 of file BorderedBandMatrix.cc.

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

Referenced by solveAndInvertBorderedBand().

297  {
298 
299  int nRow = theBand.getNumRows();
300  int nCol = theBand.getNumCols();
301  VVector aSolution(aRightHandSide);
302  for (int i = 0; i < nCol; ++i) // forward substitution
303  {
304  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
305  aSolution(j + i) -= theBand(j, i) * aSolution(i);
306  }
307  }
308  for (int i = nCol - 1; i >= 0; i--) // backward substitution
309  {
310  double rxw = theBand(0, i) * aSolution(i);
311  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
312  rxw -= theBand(j, i) * aSolution(j + i);
313  }
314  aSolution(i) = rxw;
315  }
316  return aSolution;
317  }
unsigned int getNumRows() const
Get number of rows.
Definition: VMatrix.cc:79
T min(T a, T b)
Definition: MathUtil.h:58
unsigned int getNumCols() const
Get number of columns.
Definition: VMatrix.cc:87
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 326 of file BorderedBandMatrix.cc.

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

326  {
327 
328  int nRow = theBand.getNumRows();
329  int nCol = theBand.getNumCols();
330  VMatrix aSolution(aRightHandSide);
331  for (unsigned int iBorder = 0; iBorder < numBorder; iBorder++) {
332  for (int i = 0; i < nCol; ++i) // forward substitution
333  {
334  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
335  aSolution(iBorder, j + i) -= theBand(j, i)
336  * aSolution(iBorder, i);
337  }
338  }
339  for (int i = nCol - 1; i >= 0; i--) // backward substitution
340  {
341  double rxw = theBand(0, i) * aSolution(iBorder, i);
342  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
343  rxw -= theBand(j, i) * aSolution(iBorder, j + i);
344  }
345  aSolution(iBorder, i) = rxw;
346  }
347  }
348  return aSolution;
349  }
unsigned int getNumRows() const
Get number of rows.
Definition: VMatrix.cc:79
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:87
VMatrix theBand
Band part.

Member Data Documentation

unsigned int gbl::BorderedBandMatrix::numBand
private

Band width.

Definition at line 96 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 97 of file BorderedBandMatrix.h.

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

unsigned int gbl::BorderedBandMatrix::numSize
private

Matrix size.

Definition at line 94 of file BorderedBandMatrix.h.

Referenced by resize().

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

Border part.

Definition at line 98 of file BorderedBandMatrix.h.

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

VMatrix gbl::BorderedBandMatrix::theMixed
private

Mixed part.

Definition at line 99 of file BorderedBandMatrix.h.

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