GeneralBrokenLines  trunk_rev46
 All Classes Files Functions Variables Typedefs Pages
BorderedBandMatrix.cpp
Go to the documentation of this file.
1 /*
2  * BorderedBandMatrix.cpp
3  *
4  * Created on: Aug 14, 2011
5  * Author: kleinwrt
6  */
7 
8 #include "BorderedBandMatrix.h"
9 
12 }
13 
15 }
16 
18 
23 void BorderedBandMatrix::resize(unsigned int nSize, unsigned int nBorder,
24  unsigned int nBand) {
25  numSize = nSize;
26  numBorder = nBorder;
27  numCol = nSize - nBorder;
28  numBand = 0;
31  theBand.resize((nBand + 1), numCol);
32 }
33 
35 
44  const std::vector<unsigned int>* anIndex,
45  const std::vector<double>* aVector) {
46  int nBorder = numBorder;
47  for (unsigned int i = 0; i < anIndex->size(); ++i) {
48  int iIndex = (*anIndex)[i] - 1; // anIndex has to be sorted
49  for (unsigned int j = 0; j <= i; ++j) {
50  int jIndex = (*anIndex)[j] - 1;
51  if (iIndex < nBorder) {
52  theBorder(iIndex, jIndex) += (*aVector)[i] * aWeight
53  * (*aVector)[j];
54  } else if (jIndex < nBorder) {
55  theMixed(jIndex, iIndex - nBorder) += (*aVector)[i] * aWeight
56  * (*aVector)[j];
57  } else {
58  unsigned int nBand = iIndex - jIndex;
59  theBand(nBand, jIndex - nBorder) += (*aVector)[i] * aWeight
60  * (*aVector)[j];
61  numBand = std::max(numBand, nBand); // update band width
62  }
63  }
64  }
65 }
66 
68 
73  const std::vector<unsigned int> anIndex) const {
74 
75  TMatrixDSym aMatrix(anIndex.size());
76  int nBorder = numBorder;
77  for (unsigned int i = 0; i < anIndex.size(); ++i) {
78  int iIndex = anIndex[i] - 1; // anIndex has to be sorted
79  for (unsigned int j = 0; j <= i; ++j) {
80  int jIndex = anIndex[j] - 1;
81  if (iIndex < nBorder) {
82  aMatrix(i, j) = theBorder(iIndex, jIndex); // border part of inverse
83  } else if (jIndex < nBorder) {
84  aMatrix(i, j) = -theMixed(jIndex, iIndex - nBorder); // mixed part of inverse
85  } else {
86  unsigned int nBand = iIndex - jIndex;
87  aMatrix(i, j) = theBand(nBand, jIndex - nBorder); // band part of inverse
88  }
89  aMatrix(j, i) = aMatrix(i, j);
90  }
91  }
92  return aMatrix;
93 }
94 
96 
123  const VVector &aRightHandSide, VVector &aSolution) {
124 
125  // decompose band
126  decomposeBand();
127  // invert band
128  VMatrix inverseBand = invertBand();
129  if (numBorder > 0) { // need to use block matrix decomposition to solve
130  // solve for mixed part
131  const VMatrix auxMat = solveBand(theMixed); // = Xt
132  const VMatrix auxMatT = auxMat.transpose(); // = X
133  // solve for border part
134  const VVector auxVec = aRightHandSide.getVec(numBorder)
135  - auxMat * aRightHandSide.getVec(numCol, numBorder); // = b1 - Xt*b2
136  VSymMatrix inverseBorder = theBorder - theMixed * auxMatT;
137  inverseBorder.invert(); // = E
138  const VVector borderSolution = inverseBorder * auxVec; // = x1
139  // solve for band part
140  const VVector bandSolution = solveBand(
141  aRightHandSide.getVec(numCol, numBorder)); // = x
142  aSolution.putVec(borderSolution);
143  aSolution.putVec(bandSolution - auxMatT * borderSolution, numBorder); // = x2
144  // parts of inverse
145  theBorder = inverseBorder; // E
146  theMixed = inverseBorder * auxMat; // E*Xt (-mixed part of inverse) !!!
147  theBand = inverseBand + bandOfAVAT(auxMatT, inverseBorder); // band(D^-1 + X*E*Xt)
148  } else {
149  aSolution.putVec(solveBand(aRightHandSide));
150  theBand = inverseBand;
151  }
152 }
153 
156  std::cout << "Border part " << std::endl;
157  theBorder.print();
158  std::cout << "Mixed part " << std::endl;
159  theMixed.print();
160  std::cout << "Band part " << std::endl;
161  theBand.print();
162 }
163 
164 /*============================================================================
165  from Dbandmatrix.F (MillePede-II by V. Blobel, Univ. Hamburg)
166  ============================================================================*/
168 
175 
176  int nRow = numBand + 1;
177  int nCol = numCol;
178  VVector auxVec(nCol);
179  for (int i = 0; i < nCol; ++i) {
180  auxVec(i) = theBand(0, i) * 16.0; // save diagonal elements
181  }
182  for (int i = 0; i < nCol; ++i) {
183  if ((theBand(0, i) + auxVec(i)) != theBand(0, i)) {
184  theBand(0, i) = 1.0 / theBand(0, i);
185  if (theBand(0, i) < 0.) {
186  throw 3; // not positive definite
187  }
188  } else {
189  theBand(0, i) = 0.0;
190  throw 2; // singular
191  }
192  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
193  double rxw = theBand(j, i) * theBand(0, i);
194  for (int k = 0; k < std::min(nRow, nCol - i) - j; ++k) {
195  theBand(k, i + j) -= theBand(k + j, i) * rxw;
196  }
197  theBand(j, i) = rxw;
198  }
199  }
200 }
201 
203 
209 VVector BorderedBandMatrix::solveBand(const VVector &aRightHandSide) const {
210 
211  int nRow = theBand.getNumRows();
212  int nCol = theBand.getNumCols();
213  VVector aSolution(aRightHandSide);
214  for (int i = 0; i < nCol; ++i) // forward substitution
215  {
216  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
217  aSolution(j + i) -= theBand(j, i) * aSolution(i);
218  }
219  }
220  for (int i = nCol - 1; i >= 0; i--) // backward substitution
221  {
222  double rxw = theBand(0, i) * aSolution(i);
223  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
224  rxw -= theBand(j, i) * aSolution(j + i);
225  }
226  aSolution(i) = rxw;
227  }
228  return aSolution;
229 }
230 
232 
238 VMatrix BorderedBandMatrix::solveBand(const VMatrix &aRightHandSide) const {
239 
240  int nRow = theBand.getNumRows();
241  int nCol = theBand.getNumCols();
242  VMatrix aSolution(aRightHandSide);
243  for (unsigned int iBorder = 0; iBorder < numBorder; iBorder++) {
244  for (int i = 0; i < nCol; ++i) // forward substitution
245  {
246  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
247  aSolution(iBorder, j + i) -= theBand(j, i)
248  * aSolution(iBorder, i);
249  }
250  }
251  for (int i = nCol - 1; i >= 0; i--) // backward substitution
252  {
253  double rxw = theBand(0, i) * aSolution(iBorder, i);
254  for (int j = 1; j < std::min(nRow, nCol - i); ++j) {
255  rxw -= theBand(j, i) * aSolution(iBorder, j + i);
256  }
257  aSolution(iBorder, i) = rxw;
258  }
259  }
260  return aSolution;
261 }
262 
264 
268 
269  int nRow = numBand + 1;
270  int nCol = numCol;
271  VMatrix inverseBand(nRow, nCol);
272 
273  for (int i = nCol - 1; i >= 0; i--) {
274  double rxw = theBand(0, i);
275  for (int j = i; j >= std::max(0, i - nRow + 1); j--) {
276  for (int k = j + 1; k < std::min(nCol, j + nRow); ++k) {
277  rxw -= inverseBand(abs(i - k), std::min(i, k))
278  * theBand(k - j, j);
279  }
280  inverseBand(i - j, j) = rxw;
281  rxw = 0.;
282  }
283  }
284  return inverseBand;
285 }
286 
288 
292  const VSymMatrix &aSymArray) const {
293  int nBand = numBand;
294  int nCol = numCol;
295  int nBorder = numBorder;
296  double sum;
297  VMatrix aBand((nBand + 1), nCol);
298  for (int i = 0; i < nCol; ++i) {
299  for (int j = std::max(0, i - nBand); j <= i; ++j) {
300  sum = 0.;
301  for (int l = 0; l < nBorder; ++l) { // diagonal
302  sum += anArray(i, l) * aSymArray(l, l) * anArray(j, l);
303  for (int k = 0; k < l; ++k) { // off diagonal
304  sum += anArray(i, l) * aSymArray(l, k) * anArray(j, k)
305  + anArray(i, k) * aSymArray(l, k) * anArray(j, l);
306  }
307  }
308  aBand(i - j, j) = sum;
309  }
310  }
311  return aBand;
312 }