![]() |
GeneralBrokenLines
V01-11-00
|
00001 /* 00002 * BorderedBandMatrix.cpp 00003 * 00004 * Created on: Aug 14, 2011 00005 * Author: kleinwrt 00006 */ 00007 00008 #include "BorderedBandMatrix.h" 00009 00011 BorderedBandMatrix::BorderedBandMatrix() { 00012 00013 } 00014 00015 BorderedBandMatrix::~BorderedBandMatrix() { 00016 // TODO Auto-generated destructor stub 00017 } 00018 00020 00025 void BorderedBandMatrix::resize(unsigned int nSize, unsigned int nBorder, 00026 unsigned int nBand) { 00027 numSize = nSize; 00028 numBorder = nBorder; 00029 numCol = nSize - nBorder; 00030 numBand = 0; 00031 theBorder.resize(numBorder); 00032 theMixed.resize(numBorder, numCol); 00033 theBand.resize((nBand + 1), numCol); 00034 } 00035 00037 00045 void BorderedBandMatrix::addBlockMatrix(double aWeight, 00046 const std::vector<unsigned int>* anIndex, 00047 const std::vector<double>* aVector) { 00048 int nBorder = numBorder; 00049 for (unsigned int i = 0; i < anIndex->size(); i++) { 00050 int iIndex = (*anIndex)[i] - 1; // anIndex has to be sorted 00051 for (unsigned int j = 0; j <= i; j++) { 00052 int jIndex = (*anIndex)[j] - 1; 00053 if (iIndex < nBorder) { 00054 theBorder(iIndex, jIndex) += (*aVector)[i] * aWeight 00055 * (*aVector)[j]; 00056 } else if (jIndex < nBorder) { 00057 theMixed(jIndex, iIndex - nBorder) += (*aVector)[i] * aWeight 00058 * (*aVector)[j]; 00059 } else { 00060 unsigned int nBand = iIndex - jIndex; 00061 theBand(nBand, jIndex - nBorder) += (*aVector)[i] * aWeight 00062 * (*aVector)[j]; 00063 numBand = std::max(numBand, nBand); // update band width 00064 } 00065 } 00066 } 00067 } 00068 00070 00074 TMatrixDSym BorderedBandMatrix::getBlockMatrix( 00075 std::vector<unsigned int> anIndex) const { 00076 00077 TMatrixDSym aMatrix(anIndex.size()); 00078 int nBorder = numBorder; 00079 for (unsigned int i = 0; i < anIndex.size(); i++) { 00080 int iIndex = anIndex[i] - 1; // anIndex has to be sorted 00081 for (unsigned int j = 0; j <= i; j++) { 00082 int jIndex = anIndex[j] - 1; 00083 if (iIndex < nBorder) { 00084 aMatrix(i, j) = theBorder(iIndex, jIndex); // border part of inverse 00085 } else if (jIndex < nBorder) { 00086 aMatrix(i, j) = -theMixed(jIndex, iIndex - nBorder); // mixed part of inverse 00087 } else { 00088 unsigned int nBand = iIndex - jIndex; 00089 aMatrix(i, j) = theBand(nBand, jIndex - nBorder); // band part of inverse 00090 } 00091 aMatrix(j, i) = aMatrix(i, j); 00092 } 00093 } 00094 return aMatrix; 00095 } 00096 00098 00124 void BorderedBandMatrix::solveAndInvertBorderedBand( 00125 const VVector &aRightHandSide, VVector &aSolution) { 00126 00127 // decompose band 00128 decomposeBand(); // TODO: check for positive definiteness 00129 // invert band 00130 VMatrix inverseBand = invertBand(); 00131 if (numBorder > 0) { // need to use block matrix decomposition to solve 00132 // solve for mixed part 00133 VMatrix auxMat = solveBand(theMixed); // = Xt 00134 VMatrix auxMatT = auxMat.transpose(); // = X 00135 // solve for border part 00136 VVector auxVec = aRightHandSide.getVec(numBorder) 00137 - auxMat * aRightHandSide.getVec(numCol, numBorder); // = b1 - Xt*b2 00138 VSymMatrix inverseBorder = theBorder - theMixed * auxMatT; 00139 inverseBorder.invert(); // = E 00140 VVector borderSolution = inverseBorder * auxVec; // = x1 00141 // solve for band part 00142 VVector bandSolution = solveBand( 00143 aRightHandSide.getVec(numCol, numBorder)); // = x 00144 aSolution.putVec(borderSolution); 00145 aSolution.putVec(bandSolution - auxMatT * borderSolution, numBorder); // = x2 00146 // parts of inverse 00147 theBorder = inverseBorder; // E 00148 theMixed = inverseBorder * auxMat; // E*Xt (-mixed part of inverse) !!! 00149 theBand = inverseBand + bandOfAVAT(auxMatT, inverseBorder); // band(D^-1 + X*E*Xt) 00150 } else { 00151 aSolution.putVec(solveBand(aRightHandSide)); 00152 theBand = inverseBand; 00153 } 00154 } 00155 00157 void BorderedBandMatrix::printMatrix() const { 00158 std::cout << "Border part " << std::endl; 00159 theBorder.print(); 00160 std::cout << "Mixed part " << std::endl; 00161 theMixed.print(); 00162 std::cout << "Band part " << std::endl; 00163 theBand.print(); 00164 } 00165 00166 /*============================================================================ 00167 from Dbandmatrix.F (MillePede-II by V. Blobel, Univ. Hamburg) 00168 ============================================================================*/ 00170 00176 void BorderedBandMatrix::decomposeBand() { 00177 00178 int nRow = numBand + 1; 00179 int nCol = numCol; 00180 VVector auxVec(nCol); 00181 for (int i = 0; i < nCol; i++) { 00182 auxVec(i) = theBand(0, i) * 16.0; // save diagonal elements 00183 } 00184 for (int i = 0; i < nCol; i++) { 00185 if ((theBand(0, i) + auxVec(i)) != theBand(0, i)) { 00186 theBand(0, i) = 1.0 / theBand(0, i); 00187 if (theBand(0, i) < 0.) { 00188 throw 3; // not positive definite 00189 } 00190 } else { 00191 theBand(0, i) = 0.0; 00192 throw 2; // singular 00193 } 00194 for (int j = 1; j < std::min(nRow, nCol - i); j++) { 00195 double rxw = theBand(j, i) * theBand(0, i); 00196 for (int k = 0; k < std::min(nRow, nCol - i) - j; k++) { 00197 theBand(k, i + j) -= theBand(k + j, i) * rxw; 00198 } 00199 theBand(j, i) = rxw; 00200 } 00201 } 00202 } 00203 00205 00211 VVector BorderedBandMatrix::solveBand(const VVector &aRightHandSide) const { 00212 00213 int nRow = theBand.getNumRows(); 00214 int nCol = theBand.getNumCols(); 00215 VVector aSolution(aRightHandSide); 00216 for (int i = 0; i < nCol; i++) // forward substitution 00217 { 00218 for (int j = 1; j < std::min(nRow, nCol - i); j++) { 00219 aSolution(j + i) -= theBand(j, i) * aSolution(i); 00220 } 00221 } 00222 for (int i = nCol - 1; i >= 0; i--) // backward substitution 00223 { 00224 double rxw = theBand(0, i) * aSolution(i); 00225 for (int j = 1; j < std::min(nRow, nCol - i); j++) { 00226 rxw -= theBand(j, i) * aSolution(j + i); 00227 } 00228 aSolution(i) = rxw; 00229 } 00230 return aSolution; 00231 } 00232 00234 00240 VMatrix BorderedBandMatrix::solveBand(const VMatrix &aRightHandSide) const { 00241 00242 int nRow = theBand.getNumRows(); 00243 int nCol = theBand.getNumCols(); 00244 VMatrix aSolution(aRightHandSide); 00245 for (unsigned int iBorder = 0; iBorder < numBorder; iBorder++) { 00246 for (int i = 0; i < nCol; i++) // forward substitution 00247 { 00248 for (int j = 1; j < std::min(nRow, nCol - i); j++) { 00249 aSolution(iBorder, j + i) -= theBand(j, i) 00250 * aSolution(iBorder, i); 00251 } 00252 } 00253 for (int i = nCol - 1; i >= 0; i--) // backward substitution 00254 { 00255 double rxw = theBand(0, i) * aSolution(iBorder, i); 00256 for (int j = 1; j < std::min(nRow, nCol - i); j++) { 00257 rxw -= theBand(j, i) * aSolution(iBorder, j + i); 00258 } 00259 aSolution(iBorder, i) = rxw; 00260 } 00261 } 00262 return aSolution; 00263 } 00264 00266 00269 VMatrix BorderedBandMatrix::invertBand() { 00270 00271 int nRow = numBand + 1; 00272 int nCol = numCol; 00273 VMatrix inverseBand(nRow, nCol); 00274 00275 for (int i = nCol - 1; i >= 0; i--) { 00276 double rxw = theBand(0, i); 00277 for (int j = i; j >= std::max(0, i - nRow + 1); j--) { 00278 for (int k = j + 1; k < std::min(nCol, j + nRow); k++) { 00279 rxw -= inverseBand(abs(i - k), std::min(i, k)) 00280 * theBand(k - j, j); 00281 } 00282 inverseBand(i - j, j) = rxw; 00283 rxw = 0.; 00284 } 00285 } 00286 return inverseBand; 00287 } 00288 00290 00293 VMatrix BorderedBandMatrix::bandOfAVAT(const VMatrix &anArray, 00294 const VSymMatrix &aSymArray) const { 00295 int nBand = numBand; 00296 int nCol = numCol; 00297 int nBorder = numBorder; 00298 double sum; 00299 VMatrix aBand((nBand + 1), nCol); 00300 for (int i = 0; i < nCol; i++) { 00301 for (int j = std::max(0, i - nBand); j <= i; j++) { 00302 sum = 0.; 00303 for (int l = 0; l < nBorder; l++) { // diagonal 00304 sum += anArray(i, l) * aSymArray(l, l) * anArray(j, l); 00305 for (int k = 0; k < l; k++) { // off diagonal 00306 sum += anArray(i, l) * aSymArray(l, k) * anArray(j, k) 00307 + anArray(i, k) * aSymArray(l, k) * anArray(j, l); 00308 } 00309 } 00310 aBand(i - j, j) = sum; 00311 } 00312 } 00313 return aBand; 00314 }
1.7.6.1