linearAlgebra.cpp

Go to the documentation of this file.
00001 
00013 #include <dlrCommon/exception.h>
00014 #include <dlrLinearAlgebra/linearAlgebra.h>
00015 #include <dlrLinearAlgebra/clapack.h>
00016 #include <dlrNumeric/utilities.h>
00017 
00018 
00019 namespace dlr {
00020 
00025   namespace linearAlgebra {
00026 
00027     void
00028     choleskyFactorization(const Array2D<Float64>& inputArray,
00029                           Array2D<Float64>& kArray,
00030                           bool isUpperTriangular)
00031     {
00032       // Argument checking.
00033       if(inputArray.size() == 0) {
00034         DLR_THROW3(ValueException,
00035                    "choleskyFactorization()",
00036                    "Argument inputArray cannot have zero size.");
00037       }
00038       if(inputArray.rows() != inputArray.columns()) {
00039         DLR_THROW3(ValueException,
00040                    "choleskyFactorization()",
00041                    "Argument inputArray must be square.");
00042       }
00043     
00044       // Transpose A to match LAPACK's convention.  Since inputArray
00045       // is symmetric, we don't really need to transpose!  Also, we
00046       // only need to copy the upper/lower triangular part, depending
00047       // on the value of isUpperTriangular.  For convenience, we copy
00048       // the whole matrix, though.
00049       size_t dimension = inputArray.rows();
00050       Array2D<Float64> aColumnMajor(dimension, dimension);
00051       aColumnMajor = 0.0;
00052       {
00053         Array2D<Float64>::const_iterator inPtr = inputArray.begin();
00054         Array2D<Float64>::iterator outPtr = aColumnMajor.begin();
00055         if(isUpperTriangular) {
00056           for(size_t index0 = 0; index0 < dimension; ++index0) {
00057             inPtr += index0;
00058             outPtr += index0;
00059             for(size_t index1 = index0; index1 < dimension; ++index1) {
00060               *(outPtr++) = *(inPtr++);
00061             }
00062           }
00063         } else {
00064           for(size_t index0 = 0; index0 < dimension; ++index0) {
00065             for(size_t index1 = 0; index1 <= index0; ++index1) {
00066               *(outPtr++) = *(inPtr++);
00067             }
00068             inPtr += (dimension - index0 - 1);
00069             outPtr += (dimension - index0 - 1);
00070           }
00071         }        
00072       }
00073 
00074 
00075       // Set up arguments for the LAPACK call.
00076 
00077       // Set UPLO opposite of what you might expect because LAPACK is
00078       // column major.
00079       char UPLO = isUpperTriangular ? 'L' : 'U';
00080       Int32 N = static_cast<Int32>(dimension);
00081       Int32 LDA = static_cast<Int32>(dimension);
00082       Int32 INFO;
00083 
00084       // Dispatch to the LAPACK routine that actually does the
00085       // factorization.
00086       dpotrf_(&UPLO, &N, aColumnMajor.data(), &LDA, &INFO);
00087 
00088       // Check for errors.
00089       if(INFO < 0L) {
00090         std::ostringstream message;
00091         message << "Call to dpotrf_ returns " << INFO
00092                 << ".  Something is wrong.";
00093         DLR_THROW3(ValueException, "choleskyFactorization()",
00094                    message.str().c_str());
00095       } else if(INFO > 0L) {
00096         std::ostringstream message;
00097         DLR_THROW3(ValueException, "choleskyFactorization()",
00098                    "Input matrix is not positive definite.");
00099       }
00100 
00101       // Recover the result.
00102       kArray = aColumnMajor;
00103     }
00104 
00105     
00106     Float64
00107     determinant(const Array2D<Float64>& A)
00108     {
00109       // First argument checking.
00110       if(A.columns() != A.rows()) {
00111         DLR_THROW3(ValueException,
00112                    "determinant(const Array2D<Float64>&)",
00113                    "Input array is not square.");
00114       }
00115 
00116       // In this routine, we take advantage of the fact that the
00117       // determinant of a matrix is related to the product of the
00118       // diagonal elements of its LU factorization.
00119     
00120       // Start by computing the LU factorization of A.
00121       Array2D<Float64> AColumnMajor = A.transpose();
00122       Int32 M = static_cast<Int32>(A.rows());
00123       Int32 N = static_cast<Int32>(A.columns());
00124       Int32 LDA = static_cast<Int32>(A.rows());
00125       Array1D<Int32> IPIV(M); // Really should be min(M, N), but A is square.
00126       Int32 info;
00127       dgetrf_(&M, &N, AColumnMajor.data(), &LDA, IPIV.data(), &info);
00128 
00129       // Check for errors.
00130       if(info != 0L) {
00131         std::ostringstream message;
00132         message << "Call to dgetrf_ returns " << info
00133                 << ".  Something is wrong.";
00134         DLR_THROW3(ValueException,
00135                    "determinant(const Array2D<Float64>&)",
00136                    message.str().c_str());
00137       }
00138 
00139       // Compute the product of the diagonal elements, and find out if
00140       // the determinant is equal to + or - the product of the diagonal
00141       // element. Do this second thing by counting how many row-swaps
00142       // were conducteed.
00143       Float64 determinant = 1.0;
00144       size_t swapCount = 0;
00145       for(size_t index = 0; index < AColumnMajor.rows(); ++index) {
00146         determinant *= AColumnMajor(index, index);
00147         if(IPIV[index] != static_cast<Int32>(index + 1)) {
00148           // Note(xxx): check that this is right.
00149           ++swapCount;
00150         }
00151       }
00152 
00153       // Finally, correct the sign
00154       if((swapCount % 2) == 1) {
00155         determinant *= -1.0;
00156       }
00157       return determinant;
00158     }
00159 
00160   
00161     Array1D<Float64>
00162     eigenvaluesSymmetric(const Array2D<Float64>& inputArray)
00163     {
00164       // Argument checking.
00165       if(inputArray.size() == 0) {
00166         DLR_THROW3(ValueException,
00167                    "eigenvaluesSymmetric()",
00168                    "Argument inputArray cannot have zero size.");
00169       }
00170       if(inputArray.rows() != inputArray.columns()) {
00171         DLR_THROW3(ValueException,
00172                    "eigenvaluesSymmetric()",
00173                    "Argument inputArray must be square.");
00174       }
00175     
00176       // Transpose A to match LAPACK's convention.  Since inputArray is
00177       // symmetric, we don't really need to transpose!  Also, we only
00178       // need to copy the upper triangular part.
00179       size_t dimension = inputArray.rows();
00180       Array2D<Float64> aColumnMajor(dimension, dimension);
00181       {
00182         Array2D<Float64>::const_iterator inPtr = inputArray.begin();
00183         Array2D<Float64>::iterator outPtr = aColumnMajor.begin();
00184         for(size_t index0 = 0; index0 < dimension; ++index0) {
00185           inPtr += index0;
00186           outPtr += index0;
00187           for(size_t index1 = index0; index1 < dimension; ++index1) {
00188             *(outPtr++) = *(inPtr++);
00189           }
00190         }
00191       }
00192 
00193       // Set up storage for return values.
00194       Array1D<Float64> eigenvalues(dimension);
00195 
00196       // Set up arguments for the LAPACK call.
00197       char JOBZ = 'N';  // Compute eigenvalues only.
00198       char UPLO = 'L';  // Get input from upper triangle of A.
00199       Int32 N = static_cast<Int32>(inputArray.rows());
00200       Int32 LDA = N;
00201       Float64 WORK;
00202       Int32 LWORK = -1;
00203       Int32 INFO;
00204 
00205       // Call once to request optimal workspace size.
00206       dsyev_(&JOBZ,  &UPLO, &N, aColumnMajor.data(), &LDA,
00207              eigenvalues.data(), &WORK, &LWORK, &INFO);
00208 
00209       // Check for errors.
00210       if(INFO != 0L) {
00211         std::ostringstream message;
00212         message << "First call to dsyev_ returns " << INFO
00213                 << ".  Something is wrong.";
00214         DLR_THROW3(ValueException,
00215                    "eigenvectorsSymmetric()",
00216                    message.str().c_str());
00217       }
00218     
00219       // Resize workspace.
00220       LWORK = static_cast<Int32>(WORK);
00221       Array1D<Float64> doubleWorkSpace(static_cast<size_t>(LWORK));
00222 
00223       // Call again to really compute the eigenvectors.
00224       dsyev_(&JOBZ,  &UPLO, &N, aColumnMajor.data(), &LDA,
00225              eigenvalues.data(), doubleWorkSpace.data(), &LWORK, &INFO);
00226 
00227       // Check for errors.
00228       if(INFO != 0L) {
00229         std::ostringstream message;
00230         message << "First call to dsyev_ returns " << INFO
00231                 << ".  Something is wrong.";
00232         DLR_THROW3(ValueException,
00233                    "eigenvectorsSymmetric()",
00234                    message.str().c_str());
00235       }
00236 
00237       // Our convention differs from LAPACK's about the order of eigenvalues.
00238       std::reverse(eigenvalues.begin(), eigenvalues.end());
00239       return eigenvalues;
00240     }
00241 
00242 
00243     void
00244     eigenvectorsSymmetric(const Array2D<Float64>& inputArray,
00245                           Array1D<Float64>& eigenvalues,
00246                           Array2D<Float64>& eigenvectors)
00247     {
00248       // Argument checking.
00249       if(inputArray.size() == 0) {
00250         DLR_THROW3(ValueException,
00251                    "eigenvectorsSymmetric()",
00252                    "Argument inputArray cannot have zero size.");
00253       }
00254       if(inputArray.rows() != inputArray.columns()) {
00255         DLR_THROW3(ValueException,
00256                    "eigenvectorsSymmetric()",
00257                    "Argument inputArray must be square.");
00258       }
00259     
00260       // Transpose A to match LAPACK's convention.  Since inputArray is
00261       // symmetric, we don't really need to transpose!  Also, we only
00262       // need to copy the upper triangular part.
00263       size_t dimension = inputArray.rows();
00264       Array2D<Float64> aColumnMajor(dimension, dimension);
00265       {
00266         Array2D<Float64>::const_iterator inPtr = inputArray.begin();
00267         Array2D<Float64>::iterator outPtr = aColumnMajor.begin();
00268         for(size_t index0 = 0; index0 < dimension; ++index0) {
00269           inPtr += index0;
00270           outPtr += index0;
00271           for(size_t index1 = index0; index1 < dimension; ++index1) {
00272             *(outPtr++) = *(inPtr++);
00273           }
00274         }
00275       }
00276 
00277       // Set up storage for return values.
00278       Array1D<Float64> eigenvaluesTmp(dimension);
00279 
00280       // Set up arguments for the LAPACK call.
00281       char JOBZ = 'V';  // Compute eigenvalues and eigenvectors.
00282       char UPLO = 'L';  // Get input from lower triangle of A.
00283       Int32 N = static_cast<Int32>(dimension);
00284       Int32 LDA = static_cast<Int32>(dimension);
00285       Float64 WORK;
00286       Int32 LWORK = -1;
00287       Int32 INFO;
00288 
00289       // Call once to request optimal workspace size.
00290       dsyev_(&JOBZ,  &UPLO, &N, aColumnMajor.data(), &LDA,
00291              eigenvaluesTmp.data(), &WORK, &LWORK, &INFO);
00292 
00293       // Check for errors.
00294       if(INFO != 0L) {
00295         std::ostringstream message;
00296         message << "First call to dsyev_ returns " << INFO
00297                 << ".  Something is wrong.";
00298         DLR_THROW3(ValueException,
00299                    "eigenvectorsSymmetric()",
00300                    message.str().c_str());
00301       }
00302     
00303       // Resize workspace.
00304       LWORK = static_cast<Int32>(WORK);
00305       Array1D<Float64> doubleWorkSpace(static_cast<size_t>(LWORK));
00306 
00307       // Call again to really compute the eigenvectors.
00308       dsyev_(&JOBZ,  &UPLO, &N, aColumnMajor.data(), &LDA,
00309              eigenvaluesTmp.data(), doubleWorkSpace.data(), &LWORK, &INFO);
00310 
00311       // Check for errors.
00312       if(INFO != 0L) {
00313         std::ostringstream message;
00314         message << "First call to dsyev_ returns " << INFO
00315                 << ".  Something is wrong.";
00316         DLR_THROW3(ValueException,
00317                    "eigenvectorsSymmetric()",
00318                    message.str().c_str());
00319       }
00320 
00321       // Recover the result.  Eigenvectors are left in aColumnMajor, but
00322       // unfortunately they're transposed and in the reverse order from
00323       // our convention.  How awkward.
00324 
00325       // Check the sizes of the return references.
00326       if(eigenvalues.size() != dimension) {
00327         eigenvalues.reinit(dimension);
00328       }
00329       if(eigenvectors.rows() != dimension
00330          || eigenvectors.columns() != dimension) {
00331         eigenvectors.reinit(dimension, dimension);
00332       }
00333 
00334       // Copy eigenvalues in reverse order.
00335       {
00336         Array1D<Float64>::iterator outPtr = eigenvalues.end() - 1;
00337         Array1D<Float64>::iterator inPtr0 = eigenvaluesTmp.begin();
00338         Array1D<Float64>::iterator inPtr1 = eigenvaluesTmp.end();
00339         while(inPtr0 != inPtr1) {
00340           *outPtr = *inPtr0;
00341           ++inPtr0;
00342           --outPtr;
00343         }
00344       }
00345 
00346       // Transpose and reverse and copy, oh my!
00347       for(size_t index0 = 0; index0 < dimension; ++index0) {
00348         Float64* outPtr = eigenvectors.data(index0);
00349         Float64* inPtr = aColumnMajor.data((dimension - index0 - 1) * dimension);
00350         for(size_t index1 = 0; index1 < dimension; ++index1) {
00351           *outPtr = *inPtr;
00352           ++inPtr;
00353           outPtr += dimension;
00354         }
00355       }
00356     }
00357 
00358 
00359     Array2D<Float64>
00360     inverse(const Array2D<Float64>& A)
00361     {
00362       // First argument checking.
00363       if(A.columns() != A.rows()) {
00364         DLR_THROW3(ValueException,
00365                    "inverse(const Array2D<Float64>&)",
00366                    "Input array is not square.");
00367       }
00368       // Now set up some linear equations to solve.
00369       Array2D<Float64> AInverse =
00370         identity(A.rows(), A.rows(), type_tag<Float64>());
00371       Array2D<Float64> AA = A.transpose();
00372       // And solve for the inverse matrix.
00373       linearSolveInPlace(AA, AInverse); //Modifies AInverse.
00374       return AInverse;
00375     }
00376 
00377 
00378     // This function computes the best linear fit between the two input
00379     // arrays.
00380     std::pair<Float64, Float64>
00381     linearFit(const Array1D<Float64>& array0,
00382               const Array1D<Float64>& array1)
00383     {
00384       // We're looking for constants a and b which most nearly (in the
00385       // least squares sense) satisfy the equation
00386       //
00387       //   a * array0 + b = array1
00388       //
00389       // Which can be rewritten
00390       //
00391       //   [array0[0], 1]   [a] = [array1[0]]
00392       //   [array0[1], 1] * [b]   [array1[1]]
00393       //   [array0[2], 1]         [array1[2]]
00394       //   ...                    ...
00395       // 
00396       // Solving this using the Moore-Penrose pseudoinverse gives
00397       //
00398       //                                             -1
00399       //   [a]  =  [dot(array0, array0), sum(array0)]  * [dot(array0, array1)]
00400       //   [b]     [sum(array0),         N          ]    [sum(array1)        ]
00401 
00402       // First some argument checking.
00403       if(array0.size() != array1.size()) {
00404         std::ostringstream message;
00405         message << "Arguments array0 and array1 must have the same size, "
00406                 << "but are of size " << array0.size()
00407                 << " and " << array1.size() << " respectively." << std::endl;
00408         DLR_THROW3(ValueException,
00409                    "linearFit(const Array1D<Float64>&, const Array1D<Float64>&)",
00410                    message.str().c_str());                 
00411       }
00412       if(array0.size() == 0) {
00413         DLR_THROW3(ValueException,
00414                    "linearFit(const Array1D<Float64>&, const Array1D<Float64>&)",
00415                    "Arguments cannot have zero size.");
00416       }
00417 
00418       // Do linear regression.
00419       Array2D<Float64> AMatrix(2, 2);
00420       AMatrix(0, 0) = dot(array0, array0);
00421       AMatrix(0, 1) = sum(array0);
00422       AMatrix(1, 0) = sum(array0);
00423       AMatrix(1, 1) = array0.size();
00424 
00425       Array1D<Float64> bVector(2);
00426       bVector(0) = dot(array0, array1);
00427       bVector(1) = sum(array1);
00428     
00429       Array2D<Float64> AInverse = inverse(AMatrix);
00430 
00431       Array1D<Float64> result = matrixMultiply(AInverse, bVector);
00432       return std::make_pair(result(0), result(1));
00433     }
00434 
00435   
00436     // This function solves the system of equations A*x = b, where A and
00437     // b are known Array2D<double> instances.
00438     Array1D<Float64>
00439     linearLeastSquares(const Array2D<Float64>& A, const Array1D<Float64>& b)
00440     {
00441       // First some argument checking.
00442       if(A.size() == 0) {
00443         DLR_THROW3(
00444           ValueException,
00445           "linearLeastSquares(const Array2D<Float64>&, const Array1D<Float64>&)",
00446           "Input array A must have nonzero size.");
00447       }
00448       if(A.rows() != b.size()) {
00449         DLR_THROW3(
00450           ValueException,
00451           "linearLeastSquares(const Array2D<Float64>&, const Array1D<Float64>&)",
00452           "The number of rows in input array A must be the same as the number "
00453           "of elements in b.");
00454       }
00455 
00456       // This two-line implementation is slow.
00457       // Array2D<Float64> APInv = pseudoinverse(A);
00458       // return matrixMultiply(APInv, b);
00459 
00460       // Set up scalar arguments for the LAPACK routine.
00461       char trans = 'N';
00462       Int32 rows = static_cast<Int32>(A.rows());
00463       Int32 columns = static_cast<Int32>(A.columns());
00464       Int32 nrhs = static_cast<Int32>(1);
00465       Int32 ldb = static_cast<Int32>(std::max(rows, columns));
00466       Float64 temporaryWorkspace;
00467       Int32 lwork = -1;  // Request info on optimal workspace size.
00468       Int32 info;
00469 
00470       // Set up array arguments for the LAPACK routine.
00471       Array2D<Float64> AColumnMajor = A.transpose();
00472       Array1D<Float64> bCopy(ldb);
00473       std::copy(b.begin(), b.end(), bCopy.begin());
00474     
00475       // Now invoke the LAPACK routine to find the optimal workspace
00476       // size.
00477       dgels_(&trans, &rows, &columns, &nrhs, AColumnMajor.data(), &rows,
00478              bCopy.data(), &ldb, &temporaryWorkspace, &lwork, &info);
00479 
00480       // Check for errors.
00481       if(info != 0L) {
00482         std::ostringstream message;
00483         message << "First call to dgels_ returns " << info
00484                 << ".  Something is wrong.";
00485         DLR_THROW3(ValueException,
00486                    "linearLeastSquares()",
00487                    message.str().c_str());
00488       }
00489     
00490       // Resize workspace.
00491       lwork = static_cast<Int32>(temporaryWorkspace);
00492       Array1D<Float64> doubleWorkspace(static_cast<size_t>(lwork));
00493 
00494       // Call again to solve the system of equations.
00495       dgels_(&trans, &rows, &columns, &nrhs, AColumnMajor.data(), &rows,
00496              bCopy.data(), &ldb, doubleWorkspace.data(), &lwork, &info);
00497     
00498       // Check for errors.
00499       if(info != 0L) {
00500         std::ostringstream message;
00501         message << "Second call to dgels_ returns " << info
00502                 << ".  Something is wrong.";
00503         DLR_THROW3(ValueException,
00504                    "linearLeastSquares()",
00505                    message.str().c_str());
00506       }
00507 
00508       if(bCopy.size() == A.columns()) {
00509         return bCopy;
00510       } else {
00511         Array1D<Float64> returnValue(A.columns());
00512         std::copy(bCopy.begin(), bCopy.begin() + A.columns(),
00513                   returnValue.begin());
00514         return returnValue;
00515       }
00516     }
00517 
00518 
00519     // WARNING:  linearSolveInPlace() destructively modifies 
00520     // both arguments!
00521     //
00522     // This function solves the system of equations A*x = b, where A is
00523     // a known matrix, and b is a known vector.
00524     void
00525     linearSolveInPlace(Array2D<Float64>& A, Array1D<Float64>& b)
00526     {
00527       Array2D<Float64> bMatrix(b.size(), 1, b.data());
00528       linearSolveInPlace(A, bMatrix);
00529     }
00530 
00531 
00532     // This function is identical to linearSolveInPlace(Array2D<Float64>,
00533     // Array1D<Float64>&), except that b (and therefore x) is not
00534     // constrained to be a vector.
00535     void
00536     linearSolveInPlace(Array2D<Float64> &A, Array2D<Float64> &b)
00537     {
00538       // First some argument checking.
00539       if(A.rows() != b.rows()) {
00540         DLR_THROW3(ValueException,
00541                    "linearSolveInPlace(Array2D<Float64>&, Array2D<Float64>&)",
00542                    "Input arrays A and b must have the same number of rows.");
00543       }
00544       if(A.rows() != A.columns()) {
00545         DLR_THROW3(ValueException,
00546                    "linearSolveInPlace(Array2D<Float64>&, Array2D<Float64>&)",
00547                    "Input array A must be square.");
00548       }
00549 
00550       // Grr.  Have to transpose to match lapack.
00551       Array2D<Float64> AColumnMajor = A.transpose();
00552     
00553       // Now invoke the LAPACK routine.
00554       Int32 rows = static_cast<Int32>(A.rows());
00555       Int32 xColumns = static_cast<Int32>(b.columns());
00556       Int32 *iPiv = new Int32[rows];
00557       Int32 info;
00558       dgesv_(&rows, &xColumns, AColumnMajor.data(), &rows, iPiv, b.data(), &rows,
00559              &info );
00560       // Clean up
00561       delete[] iPiv;
00562       // And do some error checking.
00563       if(info != 0L) {
00564         std::ostringstream message;
00565         message << "Call to dgesv_ returns " << info << ".  Something is wrong."
00566                 << "  Perhaps the the input equations are poorly conditioned "
00567                 << "or perhaps there is no solution.";
00568         DLR_THROW3(ValueException,
00569                    "linearSolveInPlace(Array2D<Float64>&, Array2D<Float64>&)",
00570                    message.str().c_str());
00571       }
00572     }
00573 
00574 
00575     // This function solves the system of equations A*x = b, where A is
00576     // a known tridiagonal matrix and b is a known vector.
00577     Array1D<Float64>
00578     linearSolveTridiagonal(const Array1D<Float64>& subDiagonal,
00579                            const Array1D<Float64>& centerDiagonal,
00580                            const Array1D<Float64>& superDiagonal,
00581                            const Array1D<Float64>& bVector)
00582     {
00583       // First some argument checking.
00584       if(subDiagonal.size() != superDiagonal.size()) {
00585         DLR_THROW3(ValueException,
00586                    "linearSolveTridiagonal()",
00587                    "Input arguments subDiagonal and superDiagonal "
00588                    "must have the same size.");
00589       }
00590       if(centerDiagonal.size() != bVector.size()) {
00591         DLR_THROW3(ValueException,
00592                    "linearSolveTridiagonal()",
00593                    "Input arguments centerDiagonal and bVector "
00594                    "must have the same size.");
00595       }
00596       if(centerDiagonal.size() != subDiagonal.size() + 1) {
00597         DLR_THROW3(ValueException,
00598                    "linearSolveTridiagonal()",
00599                    "Input argument centerDiagonal must have one more element "
00600                    "than input argument subDiagonal.");
00601       }
00602 
00603     
00604       // Now invoke the LAPACK routine.
00605       Int32 N = static_cast<Int32>(centerDiagonal.size());
00606       Int32 NRHS = 1;
00607       Array1D<Float64> subDiagonalCopy = subDiagonal.copy();
00608       Array1D<Float64> centerDiagonalCopy = centerDiagonal.copy();
00609       Array1D<Float64> superDiagonalCopy = superDiagonal.copy();
00610       Array1D<Float64> xVector = bVector.copy();
00611       Int32 INFO;
00612 
00613       dgtsv_(&N, &NRHS, subDiagonalCopy.data(), centerDiagonalCopy.data(),
00614              superDiagonalCopy.data(), xVector.data(), &N, &INFO);
00615 
00616       // And do some error checking.
00617       if(INFO != 0L) {
00618         std::ostringstream message;
00619         message << "Call to dgtsv_ returns " << INFO << ".  Something is wrong."
00620                 << "  Perhaps the the input equations are poorly conditioned "
00621                 << "or perhaps there is no solution.";
00622         DLR_THROW3(ValueException,
00623                    "linearSolveTridiagonal()",
00624                    message.str().c_str());
00625       }
00626 
00627       return xVector;
00628     }
00629     
00630 
00631     // This function accepts an Array2D<Float64> instance having at least
00632     // as many rows as columns, and returns the Moore-Penrose
00633     // pseudoinverse.
00634     Array2D<Float64>
00635     pseudoinverse(const Array2D<Float64>& A)
00636     {
00637       Array2D<Float64> ATranspose = A.transpose();
00638       Array2D<Float64> ATA = matrixMultiply(ATranspose, A);
00639       return matrixMultiply(inverse(ATA), ATranspose);
00640     }
00641 
00642 
00643     // This function computes the singular value decomposition of a
00644     // matrix.
00645     void
00646     singularValueDecomposition(const Array2D<Float64>& inputArray,
00647                                Array2D<Float64>& uArray,
00648                                Array1D<Float64>& sigmaArray,
00649                                Array2D<Float64>& vTransposeArray)
00650     {
00651       // Argument checking.
00652       if(inputArray.size() == 0) {
00653         DLR_THROW3(ValueException,
00654                    "singularValueDecomposition()",
00655                    "Argument inputArray cannot have zero size.");
00656       }
00657     
00658       // We do not transpose A. Instead, we let LAPACK operate on the
00659       // untransposed matrix as if it were column major and then swap U
00660       // and VT at the end.  We make a copy here since the contents of
00661       // whatever array we pass to LAPACK will be destroyed.
00662       Array2D<Float64> aColumnMajor = inputArray.copy();
00663 
00664       // Set up storage for return values.
00665       size_t numberOfSingularValues =
00666         std::min(inputArray.rows(), inputArray.columns());
00667 
00668       if((vTransposeArray.rows() != numberOfSingularValues)
00669          || (vTransposeArray.columns() != inputArray.columns())) {
00670         vTransposeArray.reinit(numberOfSingularValues, inputArray.columns());
00671       }
00672       // Shallow copy.
00673       Array2D<Float64> uColumnMajor = vTransposeArray;
00674 
00675       if((uArray.rows() != inputArray.rows())
00676          || (uArray.columns() != numberOfSingularValues)) {
00677         uArray.reinit(inputArray.rows(), numberOfSingularValues);
00678       }
00679       // Shallow copy.
00680       Array2D<Float64> vTColumnMajor = uArray;
00681 
00682       if(sigmaArray.size() != numberOfSingularValues) {
00683         sigmaArray.reinit(numberOfSingularValues);
00684       }
00685       Array1D<Int32> integerWorkSpace(8 * numberOfSingularValues);
00686     
00687       // Set up arguments for the LAPACK call.
00688       char JOBZ = 'S';  // Compute only min(M, N) columns of U and rows of VT.
00689       Int32 M = static_cast<Int32>(inputArray.columns());
00690       Int32 N = static_cast<Int32>(inputArray.rows());
00691       Int32 LDA = M;
00692       Int32 LDU = M;
00693       Int32 LDVT = static_cast<Int32>(numberOfSingularValues);
00694       Float64 WORK;
00695       Int32 LWORK = -1;  // Request info on optimal workspace size.
00696       Int32 INFO;
00697 
00698       // Call once to request optimal workspace size.
00699       dgesdd_(&JOBZ, &M, &N, aColumnMajor.data(), &LDA,
00700               sigmaArray.data(), uColumnMajor.data(), &LDU,
00701               vTColumnMajor.data(), &LDVT,
00702               &WORK, &LWORK,
00703               integerWorkSpace.data(), &INFO);
00704 
00705       // Check for errors.
00706       if(INFO != 0L) {
00707         std::ostringstream message;
00708         message << "First call to dgesdd_ returns " << INFO
00709                 << ".  Something is wrong.";
00710         DLR_THROW3(ValueException,
00711                    "singularValueDecomposition()",
00712                    message.str().c_str());
00713       }
00714     
00715       // Resize workspace.
00716       LWORK = static_cast<Int32>(WORK);
00717       Array1D<Float64> doubleWorkSpace(static_cast<size_t>(LWORK));
00718 
00719       // Call again to do the SVD.
00720       dgesdd_(&JOBZ, &M, &N, aColumnMajor.data(), &LDA,
00721               sigmaArray.data(), uColumnMajor.data(), &LDU,
00722               vTColumnMajor.data(), &LDVT,
00723               doubleWorkSpace.data(), &LWORK,
00724               integerWorkSpace.data(), &INFO);
00725 
00726       // Check for errors.
00727       if(INFO != 0L) {
00728         std::ostringstream message;
00729         message << "Second call to dgesdd_ returns " << INFO
00730                 << ".  Something is wrong.";
00731         DLR_THROW3(ValueException,
00732                    "singularValueDecomposition()",
00733                    message.str().c_str());
00734       }
00735     } 
00736 
00737 
00738 
00739     void
00740     singularValueDecompositionSimple(const Array2D<Float64>& inputArray,
00741                                      Array2D<Float64>& uArray,
00742                                      Array1D<Float64>& sigmaArray,
00743                                      Array2D<Float64>& vTransposeArray)
00744     {
00745       // Argument checking.
00746       if(inputArray.size() == 0) {
00747         DLR_THROW3(ValueException,
00748                    "singularValueDecomposition()",
00749                    "Argument inputArray cannot have zero size.");
00750       }
00751     
00752       // We do not transpose A. Instead, we let LAPACK operate on the
00753       // untransposed matrix as if it were column major and then swap U
00754       // and VT at the end.  We make a copy here since the contents of
00755       // whatever array we pass to LAPACK will be destroyed.
00756       Array2D<Float64> aColumnMajor = inputArray.copy();
00757 
00758       // Set up storage for return values.
00759       size_t numberOfSingularValues =
00760         std::min(inputArray.rows(), inputArray.columns());
00761 
00762       // Make sure V is appropriately sized.
00763       if((vTransposeArray.rows() != numberOfSingularValues)
00764          || (vTransposeArray.columns() != inputArray.columns())) {
00765         vTransposeArray.reinit(numberOfSingularValues, inputArray.columns());
00766       }
00767 
00768       // Shallow copy.
00769       Array2D<Float64> uColumnMajor = vTransposeArray;
00770 
00771       // Make sure U is appropriately sized.
00772       if((uArray.rows() != inputArray.rows())
00773          || (uArray.columns() != numberOfSingularValues)) {
00774         uArray.reinit(inputArray.rows(), numberOfSingularValues);
00775       }
00776     
00777       // Shallow copy.
00778       Array2D<Float64> vTColumnMajor = uArray;
00779 
00780       if(sigmaArray.size() != numberOfSingularValues) {
00781         sigmaArray.reinit(numberOfSingularValues);
00782       }
00783 
00784       // Integer workspace.
00785       Array1D<Int32> integerWorkspace(8 * numberOfSingularValues);
00786 
00787       // Set up arguments for the LAPACK call.
00788       char JOBZ = 'S';  // Compute only min(M, N) columns of U.
00789       Int32 M = static_cast<Int32>(inputArray.columns());
00790       Int32 N = static_cast<Int32>(inputArray.rows());
00791       Int32 LDA = M;
00792       Int32 LDU = M;
00793       Int32 LDVT = static_cast<Int32>(numberOfSingularValues);
00794       Float64 WORK;
00795       Int32 LWORK = -1;  // Request info on optimal workspace size.
00796       Int32 INFO;
00797 
00798       // Call once to request optimal workspace size.
00799       dgesdd_(&JOBZ, &M, &N, aColumnMajor.data(), &LDA,
00800               sigmaArray.data(), uColumnMajor.data(), &LDU,
00801               vTColumnMajor.data(), &LDVT,
00802               &WORK, &LWORK, integerWorkspace.data(), &INFO);
00803 
00804       // Check for errors.
00805       if(INFO != 0L) {
00806         std::ostringstream message;
00807         message << "First call to dgesdd_ returns " << INFO
00808                 << ".  Something is wrong.";
00809         DLR_THROW3(ValueException,
00810                    "singularValueDecomposition()",
00811                    message.str().c_str());
00812       }
00813     
00814       // Resize workspace.
00815       LWORK = static_cast<Int32>(WORK);
00816       Array1D<Float64> doubleWorkSpace(static_cast<size_t>(LWORK));
00817 
00818       // Call again to do the SVD.
00819       dgesdd_(&JOBZ, &M, &N, aColumnMajor.data(), &LDA,
00820               sigmaArray.data(), uColumnMajor.data(), &LDU,
00821               vTColumnMajor.data(), &LDVT,
00822               doubleWorkSpace.data(), &LWORK, integerWorkspace.data(), &INFO);
00823 
00824       // Check for errors.
00825       if(INFO != 0L) {
00826         std::ostringstream message;
00827         message << "Second call to dgesdd_ returns " << INFO
00828                 << ".  Something is wrong.";
00829         DLR_THROW3(ValueException,
00830                    "singularValueDecomposition()",
00831                    message.str().c_str());
00832       }
00833     } 
00834 
00835   
00836 //   // This function computes the singular value decomposition of a
00837 //   // matrix, but is less efficient than the routine above: It
00838 //   // explicitly transposes the matrix before passing it to LAPACK.
00839 //   void
00840 //   singularValueDecomposition(
00841 //     const Array2D<double>& inputArray,
00842 //     Array2D<double>& uArray,
00843 //     Array1D<double>& sigmaArray,
00844 //     Array2D<double>& vTransposeArray)
00845 //   {
00846 //     // Argument checking.
00847 //     if(inputArray.size() == 0) {
00848 //       DLR_THROW3(ValueException,
00849 //                  "singularValueDecomposition()",
00850 //                  "Argument inputArray cannot have zero size.");
00851 //     }
00852     
00853 //     // Transpose A to match LAPACK's convention.
00854 //     Array2D<double> aColumnMajor = inputArray.transpose();
00855 
00856 //     // Set up storage for return values.
00857 //     size_t numberOfSingularValues =
00858 //       std::min(inputArray.rows(), inputArray.columns());
00859 //     Array2D<double> uColumnMajor(numberOfSingularValues, inputArray.rows());
00860 //     Array2D<double> vTColumnMajor(
00861 //       inputArray.columns(), numberOfSingularValues);
00862 //     if(sigmaArray.size() != numberOfSingularValues) {
00863 //       sigmaArray.reinit(numberOfSingularValues);
00864 //     }
00865 //     Array1D<Int32> integerWorkSpace(8 * numberOfSingularValues);
00866     
00867 //     // Set up arguments for the LAPACK call.
00868 //     char JOBZ = 'S';  // Compute only min(M, N) columns of U and rows of VT.
00869 //     Int32 M = static_cast<Int32>(inputArray.rows());
00870 //     Int32 N = static_cast<Int32>(inputArray.columns());
00871 //     Int32 LDA = M;
00872 //     Int32 LDU = M;
00873 //     Int32 LDVT = static_cast<Int32>(numberOfSingularValues);
00874 //     double WORK;
00875 //     Int32 LWORK = -1;  // Request info on optimal workspace size.
00876 //     Int32 INFO;
00877 
00878 //     // Call once to request optimal workspace size.
00879 //     dgesdd_(&JOBZ, &M, &N, aColumnMajor.data(), &LDA,
00880 //             sigmaArray.data(), uColumnMajor.data(), &LDU,
00881 //             vTColumnMajor.data(), &LDVT,
00882 //             &WORK, &LWORK,
00883 //             integerWorkSpace.data(), &INFO);
00884 
00885 //     // Check for errors.
00886 //     if(INFO != 0L) {
00887 //       std::ostringstream message;
00888 //       message << "First call to dgesdd_ returns " << INFO
00889 //               << ".  Something is wrong.";
00890 //       DLR_THROW3(ValueException,
00891 //                  "singularValueDecomposition()",
00892 //                  message.str().c_str());
00893 //     }
00894     
00895 //     // Resize workspace.
00896 //     LWORK = static_cast<Int32>(WORK);
00897 //     Array1D<double> doubleWorkSpace(static_cast<size_t>(LWORK));
00898 
00899 //     // Call again to do the SVD.
00900 //     dgesdd_(&JOBZ, &M, &N, aColumnMajor.data(), &LDA,
00901 //             sigmaArray.data(), uColumnMajor.data(), &LDU,
00902 //             vTColumnMajor.data(), &LDVT,
00903 //             doubleWorkSpace.data(), &LWORK,
00904 //             integerWorkSpace.data(), &INFO);
00905 
00906 //     // Check for errors.
00907 //     if(INFO != 0L) {
00908 //       std::ostringstream message;
00909 //       message << "Second call to dgesdd_ returns " << INFO
00910 //               << ".  Something is wrong.";
00911 //       DLR_THROW3(ValueException,
00912 //                  "singularValueDecomposition()",
00913 //                  message.str().c_str());
00914 //     }
00915 
00916 //     // Recover the result.
00917 //     uArray = uColumnMajor.transpose();
00918 //     vTransposeArray = vTColumnMajor.transpose();
00919 //   } 
00920   
00921 
00922     // This function computes the singular values a matrix without
00923     // computing the associated U and V matrices.
00924     Array1D<Float64>
00925     singularValues(const Array2D<Float64>& inputArray)
00926     {
00927       // Argument checking.
00928       if(inputArray.size() == 0) {
00929         DLR_THROW3(ValueException,
00930                    "singularValues()",
00931                    "Argument inputArray cannot have zero size.");
00932       }
00933     
00934       // Transpose A to match LAPACK's convention.
00935       //
00936       // Actually, it's slightly quicker not to transpose, and singular
00937       // values remain the same, so we just copy.
00938       // 
00939       // Array2D<Float64> aColumnMajor = inputArray.transpose();
00940       Array2D<Float64> aColumnMajor = inputArray.copy();
00941 
00942       // Set up storage for return values.
00943       size_t numberOfSingularValues =
00944         std::min(inputArray.rows(), inputArray.columns());
00945       Array1D<Float64> sigmaArray(numberOfSingularValues);
00946       Array1D<Int32> integerWorkSpace(8 * numberOfSingularValues);
00947     
00948       // Set up arguments for the LAPACK call.
00949       char JOBZ = 'N';  // Compute no columns of U and no rows of VT.
00950       // Since we didn't transpose, we have to reverse rows & columns.
00951       // Int32 M = static_cast<Int32>(inputArray.rows());
00952       // Int32 N = static_cast<Int32>(inputArray.columns());
00953       Int32 M = static_cast<Int32>(inputArray.columns());
00954       Int32 N = static_cast<Int32>(inputArray.rows());
00955       Int32 LDA = M;
00956       Float64 U;         // U is not referenced by the LAPACK call.
00957       Int32 LDU = 1;
00958       Float64 VT;        // VT is not referenced by the LAPACK call.
00959       Int32 LDVT = 1;
00960       Float64 WORK;
00961       Int32 LWORK = -1;  // Request info on optimal workspace size.
00962       Int32 INFO;
00963 
00964       // Call once to request optimal workspace size.
00965       dgesdd_(&JOBZ, &M, &N, aColumnMajor.data(), &LDA,
00966               sigmaArray.data(), &U, &LDU, &VT, &LDVT,
00967               &WORK, &LWORK,
00968               integerWorkSpace.data(), &INFO);
00969 
00970       // Check for errors.
00971       if(INFO != 0L) {
00972         std::ostringstream message;
00973         message << "First call to dgesdd_ returns " << INFO
00974                 << ".  Something is wrong.";
00975         DLR_THROW3(ValueException,
00976                    "singularValues()",
00977                    message.str().c_str());
00978       }
00979     
00980       // Resize workspace.
00981       LWORK = static_cast<Int32>(WORK);
00982 
00983       // Sanity check to correct a bug in LAPACK3.
00984       Int32 minimumLWORK =
00985         3 * std::min(M,N) + std::max(std::max(M, N), 6 * std::min(M, N));
00986       if(LWORK < minimumLWORK) {
00987         // std::cerr << "Warning: changing LWORK from " << LWORK << " to "
00988         //           << minimumLWORK << "." << std::endl;
00989         LWORK = minimumLWORK;
00990       }
00991       Array1D<Float64> doubleWorkSpace(static_cast<size_t>(LWORK));
00992       
00993       // Call again to do the SVD.
00994       dgesdd_(&JOBZ, &M, &N, aColumnMajor.data(), &LDA,
00995               sigmaArray.data(), &U, &LDU, &VT, &LDVT,
00996               doubleWorkSpace.data(), &LWORK,
00997               integerWorkSpace.data(), &INFO);
00998 
00999       // Check for errors.
01000       if(INFO != 0L) {
01001         std::ostringstream message;
01002         message << "First call to dgesdd_ returns " << INFO
01003                 << ".  Something is wrong.";
01004         DLR_THROW3(ValueException,
01005                    "singularValues()",
01006                    message.str().c_str());
01007       }
01008 
01009       // Return the result.
01010       return sigmaArray;
01011     }
01012   
01013   } // namespace linearAlgebra
01014   
01015 } // namespace dlr

Generated on Tue Jan 6 22:59:48 2009 for dlrLinearAlgebra Utility Library by  doxygen 1.5.6