$extrastylesheet
dense_matrix.C
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 
00019 // C++ Includes
00020 #include <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
00021 #include <cmath> // for sqrt
00022 
00023 // Local Includes
00024 #include "libmesh/dense_matrix.h"
00025 #include "libmesh/dense_vector.h"
00026 #include "libmesh/libmesh.h"
00027 
00028 namespace libMesh
00029 {
00030 
00031 
00032 
00033 // ------------------------------------------------------------
00034 // Dense Matrix member functions
00035 
00036 template<typename T>
00037 void DenseMatrix<T>::left_multiply (const DenseMatrixBase<T>& M2)
00038 {
00039   if (this->use_blas_lapack)
00040     this->_multiply_blas(M2, LEFT_MULTIPLY);
00041   else
00042     {
00043       // (*this) <- M2 * (*this)
00044       // Where:
00045       // (*this) = (m x n),
00046       // M2      = (m x p),
00047       // M3      = (p x n)
00048 
00049       // M3 is a copy of *this before it gets resize()d
00050       DenseMatrix<T> M3(*this);
00051 
00052       // Resize *this so that the result can fit
00053       this->resize (M2.m(), M3.n());
00054 
00055       // Call the multiply function in the base class
00056       this->multiply(*this, M2, M3);
00057     }
00058 }
00059 
00060 
00061 
00062 template<typename T>
00063 template<typename T2>
00064 void DenseMatrix<T>::left_multiply (const DenseMatrixBase<T2>& M2)
00065 {
00066   // (*this) <- M2 * (*this)
00067   // Where:
00068   // (*this) = (m x n),
00069   // M2      = (m x p),
00070   // M3      = (p x n)
00071 
00072   // M3 is a copy of *this before it gets resize()d
00073   DenseMatrix<T> M3(*this);
00074 
00075   // Resize *this so that the result can fit
00076   this->resize (M2.m(), M3.n());
00077 
00078   // Call the multiply function in the base class
00079   this->multiply(*this, M2, M3);
00080 }
00081 
00082 
00083 
00084 template<typename T>
00085 void DenseMatrix<T>::left_multiply_transpose(const DenseMatrix<T>& A)
00086 {
00087   if (this->use_blas_lapack)
00088     this->_multiply_blas(A, LEFT_MULTIPLY_TRANSPOSE);
00089   else
00090     {
00091       //Check to see if we are doing (A^T)*A
00092       if (this == &A)
00093         {
00094           //libmesh_here();
00095           DenseMatrix<T> B(*this);
00096 
00097           // Simple but inefficient way
00098           // return this->left_multiply_transpose(B);
00099 
00100           // More efficient, but more code way
00101           // If A is mxn, the result will be a square matrix of Size n x n.
00102           const unsigned int n_rows = A.m();
00103           const unsigned int n_cols = A.n();
00104 
00105           // resize() *this and also zero out all entries.
00106           this->resize(n_cols,n_cols);
00107 
00108           // Compute the lower-triangular part
00109           for (unsigned int i=0; i<n_cols; ++i)
00110             for (unsigned int j=0; j<=i; ++j)
00111               for (unsigned int k=0; k<n_rows; ++k) // inner products are over n_rows
00112                 (*this)(i,j) += B(k,i)*B(k,j);
00113 
00114           // Copy lower-triangular part into upper-triangular part
00115           for (unsigned int i=0; i<n_cols; ++i)
00116             for (unsigned int j=i+1; j<n_cols; ++j)
00117               (*this)(i,j) = (*this)(j,i);
00118         }
00119 
00120       else
00121         {
00122           DenseMatrix<T> B(*this);
00123 
00124           this->resize (A.n(), B.n());
00125 
00126           libmesh_assert_equal_to (A.m(), B.m());
00127           libmesh_assert_equal_to (this->m(), A.n());
00128           libmesh_assert_equal_to (this->n(), B.n());
00129 
00130           const unsigned int m_s = A.n();
00131           const unsigned int p_s = A.m();
00132           const unsigned int n_s = this->n();
00133 
00134           // Do it this way because there is a
00135           // decent chance (at least for constraint matrices)
00136           // that A.transpose(i,k) = 0.
00137           for (unsigned int i=0; i<m_s; i++)
00138             for (unsigned int k=0; k<p_s; k++)
00139               if (A.transpose(i,k) != 0.)
00140                 for (unsigned int j=0; j<n_s; j++)
00141                   (*this)(i,j) += A.transpose(i,k)*B(k,j);
00142         }
00143     }
00144 
00145 }
00146 
00147 
00148 
00149 template<typename T>
00150 template<typename T2>
00151 void DenseMatrix<T>::left_multiply_transpose(const DenseMatrix<T2>& A)
00152 {
00153   //Check to see if we are doing (A^T)*A
00154   if (this == &A)
00155     {
00156       //libmesh_here();
00157       DenseMatrix<T> B(*this);
00158 
00159       // Simple but inefficient way
00160       // return this->left_multiply_transpose(B);
00161 
00162       // More efficient, but more code way
00163       // If A is mxn, the result will be a square matrix of Size n x n.
00164       const unsigned int n_rows = A.m();
00165       const unsigned int n_cols = A.n();
00166 
00167       // resize() *this and also zero out all entries.
00168       this->resize(n_cols,n_cols);
00169 
00170       // Compute the lower-triangular part
00171       for (unsigned int i=0; i<n_cols; ++i)
00172         for (unsigned int j=0; j<=i; ++j)
00173           for (unsigned int k=0; k<n_rows; ++k) // inner products are over n_rows
00174             (*this)(i,j) += B(k,i)*B(k,j);
00175 
00176       // Copy lower-triangular part into upper-triangular part
00177       for (unsigned int i=0; i<n_cols; ++i)
00178         for (unsigned int j=i+1; j<n_cols; ++j)
00179           (*this)(i,j) = (*this)(j,i);
00180     }
00181 
00182   else
00183     {
00184       DenseMatrix<T> B(*this);
00185 
00186       this->resize (A.n(), B.n());
00187 
00188       libmesh_assert_equal_to (A.m(), B.m());
00189       libmesh_assert_equal_to (this->m(), A.n());
00190       libmesh_assert_equal_to (this->n(), B.n());
00191 
00192       const unsigned int m_s = A.n();
00193       const unsigned int p_s = A.m();
00194       const unsigned int n_s = this->n();
00195 
00196       // Do it this way because there is a
00197       // decent chance (at least for constraint matrices)
00198       // that A.transpose(i,k) = 0.
00199       for (unsigned int i=0; i<m_s; i++)
00200         for (unsigned int k=0; k<p_s; k++)
00201           if (A.transpose(i,k) != 0.)
00202             for (unsigned int j=0; j<n_s; j++)
00203               (*this)(i,j) += A.transpose(i,k)*B(k,j);
00204     }
00205 }
00206 
00207 
00208 
00209 template<typename T>
00210 void DenseMatrix<T>::right_multiply (const DenseMatrixBase<T>& M3)
00211 {
00212   if (this->use_blas_lapack)
00213     this->_multiply_blas(M3, RIGHT_MULTIPLY);
00214   else
00215     {
00216       // (*this) <- M3 * (*this)
00217       // Where:
00218       // (*this) = (m x n),
00219       // M2      = (m x p),
00220       // M3      = (p x n)
00221 
00222       // M2 is a copy of *this before it gets resize()d
00223       DenseMatrix<T> M2(*this);
00224 
00225       // Resize *this so that the result can fit
00226       this->resize (M2.m(), M3.n());
00227 
00228       this->multiply(*this, M2, M3);
00229     }
00230 }
00231 
00232 
00233 
00234 template<typename T>
00235 template<typename T2>
00236 void DenseMatrix<T>::right_multiply (const DenseMatrixBase<T2>& M3)
00237 {
00238   // (*this) <- M3 * (*this)
00239   // Where:
00240   // (*this) = (m x n),
00241   // M2      = (m x p),
00242   // M3      = (p x n)
00243 
00244   // M2 is a copy of *this before it gets resize()d
00245   DenseMatrix<T> M2(*this);
00246 
00247   // Resize *this so that the result can fit
00248   this->resize (M2.m(), M3.n());
00249 
00250   this->multiply(*this, M2, M3);
00251 }
00252 
00253 
00254 
00255 
00256 template<typename T>
00257 void DenseMatrix<T>::right_multiply_transpose (const DenseMatrix<T>& B)
00258 {
00259   if (this->use_blas_lapack)
00260     this->_multiply_blas(B, RIGHT_MULTIPLY_TRANSPOSE);
00261   else
00262     {
00263       //Check to see if we are doing B*(B^T)
00264       if (this == &B)
00265         {
00266           //libmesh_here();
00267           DenseMatrix<T> A(*this);
00268 
00269           // Simple but inefficient way
00270           // return this->right_multiply_transpose(A);
00271 
00272           // More efficient, more code
00273           // If B is mxn, the result will be a square matrix of Size m x m.
00274           const unsigned int n_rows = B.m();
00275           const unsigned int n_cols = B.n();
00276 
00277           // resize() *this and also zero out all entries.
00278           this->resize(n_rows,n_rows);
00279 
00280           // Compute the lower-triangular part
00281           for (unsigned int i=0; i<n_rows; ++i)
00282             for (unsigned int j=0; j<=i; ++j)
00283               for (unsigned int k=0; k<n_cols; ++k) // inner products are over n_cols
00284                 (*this)(i,j) += A(i,k)*A(j,k);
00285 
00286           // Copy lower-triangular part into upper-triangular part
00287           for (unsigned int i=0; i<n_rows; ++i)
00288             for (unsigned int j=i+1; j<n_rows; ++j)
00289               (*this)(i,j) = (*this)(j,i);
00290         }
00291 
00292       else
00293         {
00294           DenseMatrix<T> A(*this);
00295 
00296           this->resize (A.m(), B.m());
00297 
00298           libmesh_assert_equal_to (A.n(), B.n());
00299           libmesh_assert_equal_to (this->m(), A.m());
00300           libmesh_assert_equal_to (this->n(), B.m());
00301 
00302           const unsigned int m_s = A.m();
00303           const unsigned int p_s = A.n();
00304           const unsigned int n_s = this->n();
00305 
00306           // Do it this way because there is a
00307           // decent chance (at least for constraint matrices)
00308           // that B.transpose(k,j) = 0.
00309           for (unsigned int j=0; j<n_s; j++)
00310             for (unsigned int k=0; k<p_s; k++)
00311               if (B.transpose(k,j) != 0.)
00312                 for (unsigned int i=0; i<m_s; i++)
00313                   (*this)(i,j) += A(i,k)*B.transpose(k,j);
00314         }
00315     }
00316 }
00317 
00318 
00319 
00320 template<typename T>
00321 template<typename T2>
00322 void DenseMatrix<T>::right_multiply_transpose (const DenseMatrix<T2>& B)
00323 {
00324   //Check to see if we are doing B*(B^T)
00325   if (this == &B)
00326     {
00327       //libmesh_here();
00328       DenseMatrix<T> A(*this);
00329 
00330       // Simple but inefficient way
00331       // return this->right_multiply_transpose(A);
00332 
00333       // More efficient, more code
00334       // If B is mxn, the result will be a square matrix of Size m x m.
00335       const unsigned int n_rows = B.m();
00336       const unsigned int n_cols = B.n();
00337 
00338       // resize() *this and also zero out all entries.
00339       this->resize(n_rows,n_rows);
00340 
00341       // Compute the lower-triangular part
00342       for (unsigned int i=0; i<n_rows; ++i)
00343         for (unsigned int j=0; j<=i; ++j)
00344           for (unsigned int k=0; k<n_cols; ++k) // inner products are over n_cols
00345             (*this)(i,j) += A(i,k)*A(j,k);
00346 
00347       // Copy lower-triangular part into upper-triangular part
00348       for (unsigned int i=0; i<n_rows; ++i)
00349         for (unsigned int j=i+1; j<n_rows; ++j)
00350           (*this)(i,j) = (*this)(j,i);
00351     }
00352 
00353   else
00354     {
00355       DenseMatrix<T> A(*this);
00356 
00357       this->resize (A.m(), B.m());
00358 
00359       libmesh_assert_equal_to (A.n(), B.n());
00360       libmesh_assert_equal_to (this->m(), A.m());
00361       libmesh_assert_equal_to (this->n(), B.m());
00362 
00363       const unsigned int m_s = A.m();
00364       const unsigned int p_s = A.n();
00365       const unsigned int n_s = this->n();
00366 
00367       // Do it this way because there is a
00368       // decent chance (at least for constraint matrices)
00369       // that B.transpose(k,j) = 0.
00370       for (unsigned int j=0; j<n_s; j++)
00371         for (unsigned int k=0; k<p_s; k++)
00372           if (B.transpose(k,j) != 0.)
00373             for (unsigned int i=0; i<m_s; i++)
00374               (*this)(i,j) += A(i,k)*B.transpose(k,j);
00375     }
00376 }
00377 
00378 
00379 
00380 
00381 template<typename T>
00382 void DenseMatrix<T>::vector_mult (DenseVector<T>& dest,
00383                                   const DenseVector<T>& arg) const
00384 {
00385   // Make sure the input sizes are compatible
00386   libmesh_assert_equal_to (this->n(), arg.size());
00387 
00388   // Resize and clear dest.
00389   // Note: DenseVector::resize() also zeros the vector.
00390   dest.resize(this->m());
00391 
00392   // Short-circuit if the matrix is empty
00393   if(this->m() == 0 || this->n() == 0)
00394     return;
00395 
00396   if (this->use_blas_lapack)
00397     this->_matvec_blas(1., 0., dest, arg);
00398   else
00399     {
00400       const unsigned int n_rows = this->m();
00401       const unsigned int n_cols = this->n();
00402 
00403       for(unsigned int i=0; i<n_rows; i++)
00404         for(unsigned int j=0; j<n_cols; j++)
00405           dest(i) += (*this)(i,j)*arg(j);
00406     }
00407 }
00408 
00409 
00410 
00411 template<typename T>
00412 template<typename T2>
00413 void DenseMatrix<T>::vector_mult (DenseVector<typename CompareTypes<T,T2>::supertype>& dest,
00414                                   const DenseVector<T2>& arg) const
00415 {
00416   // Make sure the input sizes are compatible
00417   libmesh_assert_equal_to (this->n(), arg.size());
00418 
00419   // Resize and clear dest.
00420   // Note: DenseVector::resize() also zeros the vector.
00421   dest.resize(this->m());
00422 
00423   // Short-circuit if the matrix is empty
00424   if(this->m() == 0 || this->n() == 0)
00425     return;
00426 
00427   const unsigned int n_rows = this->m();
00428   const unsigned int n_cols = this->n();
00429 
00430   for(unsigned int i=0; i<n_rows; i++)
00431     for(unsigned int j=0; j<n_cols; j++)
00432       dest(i) += (*this)(i,j)*arg(j);
00433 }
00434 
00435 
00436 
00437 template<typename T>
00438 void DenseMatrix<T>::vector_mult_transpose (DenseVector<T>& dest,
00439                                             const DenseVector<T>& arg) const
00440 {
00441   // Make sure the input sizes are compatible
00442   libmesh_assert_equal_to (this->m(), arg.size());
00443 
00444   // Resize and clear dest.
00445   // Note: DenseVector::resize() also zeros the vector.
00446   dest.resize(this->n());
00447 
00448   // Short-circuit if the matrix is empty
00449   if(this->m() == 0)
00450     return;
00451 
00452   if (this->use_blas_lapack)
00453     {
00454       this->_matvec_blas(1., 0., dest, arg, /*trans=*/true);
00455     }
00456   else
00457     {
00458       const unsigned int n_rows = this->m();
00459       const unsigned int n_cols = this->n();
00460 
00461       // WORKS
00462       // for(unsigned int j=0; j<n_cols; j++)
00463       //   for(unsigned int i=0; i<n_rows; i++)
00464       //     dest(j) += (*this)(i,j)*arg(i);
00465 
00466       // ALSO WORKS, (i,j) just swapped
00467       for(unsigned int i=0; i<n_cols; i++)
00468         for(unsigned int j=0; j<n_rows; j++)
00469           dest(i) += (*this)(j,i)*arg(j);
00470     }
00471 }
00472 
00473 
00474 
00475 template<typename T>
00476 template<typename T2>
00477 void DenseMatrix<T>::vector_mult_transpose (DenseVector<typename CompareTypes<T,T2>::supertype>& dest,
00478                                             const DenseVector<T2>& arg) const
00479 {
00480   // Make sure the input sizes are compatible
00481   libmesh_assert_equal_to (this->m(), arg.size());
00482 
00483   // Resize and clear dest.
00484   // Note: DenseVector::resize() also zeros the vector.
00485   dest.resize(this->n());
00486 
00487   // Short-circuit if the matrix is empty
00488   if(this->m() == 0)
00489     return;
00490 
00491   const unsigned int n_rows = this->m();
00492   const unsigned int n_cols = this->n();
00493 
00494   // WORKS
00495   // for(unsigned int j=0; j<n_cols; j++)
00496   //   for(unsigned int i=0; i<n_rows; i++)
00497   //     dest(j) += (*this)(i,j)*arg(i);
00498 
00499   // ALSO WORKS, (i,j) just swapped
00500   for(unsigned int i=0; i<n_cols; i++)
00501     for(unsigned int j=0; j<n_rows; j++)
00502       dest(i) += (*this)(j,i)*arg(j);
00503 }
00504 
00505 
00506 
00507 template<typename T>
00508 void DenseMatrix<T>::vector_mult_add (DenseVector<T>& dest,
00509                                       const T factor,
00510                                       const DenseVector<T>& arg) const
00511 {
00512   // Short-circuit if the matrix is empty
00513   if(this->m() == 0)
00514     {
00515       dest.resize(0);
00516       return;
00517     }
00518 
00519   if (this->use_blas_lapack)
00520     this->_matvec_blas(factor, 1., dest, arg);
00521   else
00522     {
00523       DenseVector<T> temp(arg.size());
00524       this->vector_mult(temp, arg);
00525       dest.add(factor, temp);
00526     }
00527 }
00528 
00529 
00530 
00531 template<typename T>
00532 template<typename T2, typename T3>
00533 void DenseMatrix<T>::vector_mult_add (DenseVector<typename CompareTypes<T, typename CompareTypes<T2,T3>::supertype>::supertype>& dest,
00534                                       const T2 factor,
00535                                       const DenseVector<T3>& arg) const
00536 {
00537   // Short-circuit if the matrix is empty
00538   if (this->m() == 0)
00539     {
00540       dest.resize(0);
00541       return;
00542     }
00543 
00544   DenseVector<typename CompareTypes<T,T3>::supertype>
00545     temp(arg.size());
00546   this->vector_mult(temp, arg);
00547   dest.add(factor, temp);
00548 }
00549 
00550 
00551 
00552 template<typename T>
00553 void DenseMatrix<T>::get_principal_submatrix (unsigned int sub_m,
00554                                               unsigned int sub_n,
00555                                               DenseMatrix<T>& dest) const
00556 {
00557   libmesh_assert( (sub_m <= this->m()) && (sub_n <= this->n()) );
00558 
00559   dest.resize(sub_m, sub_n);
00560   for(unsigned int i=0; i<sub_m; i++)
00561     for(unsigned int j=0; j<sub_n; j++)
00562       dest(i,j) = (*this)(i,j);
00563 }
00564 
00565 
00566 
00567 template<typename T>
00568 void DenseMatrix<T>::get_principal_submatrix (unsigned int sub_m, DenseMatrix<T>& dest) const
00569 {
00570   get_principal_submatrix(sub_m, sub_m, dest);
00571 }
00572 
00573 
00574 
00575 template<typename T>
00576 void DenseMatrix<T>::get_transpose (DenseMatrix<T>& dest) const
00577 {
00578   dest.resize(this->n(), this->m());
00579 
00580   for (unsigned int i=0; i<dest.m(); i++)
00581     for (unsigned int j=0; j<dest.n(); j++)
00582       dest(i,j) = (*this)(j,i);
00583 }
00584 
00585 
00586 
00587 
00588 template<typename T>
00589 void DenseMatrix<T>::lu_solve (const DenseVector<T>& b,
00590                                DenseVector<T>& x)
00591 {
00592   // Check to be sure that the matrix is square before attempting
00593   // an LU-solve.  In general, one can compute the LU factorization of
00594   // a non-square matrix, but:
00595   //
00596   // Overdetermined systems (m>n) have a solution only if enough of
00597   // the equations are linearly-dependent.
00598   //
00599   // Underdetermined systems (m<n) typically have infinitely many
00600   // solutions.
00601   //
00602   // We don't want to deal with either of these ambiguous cases here...
00603   libmesh_assert_equal_to (this->m(), this->n());
00604 
00605   switch(this->_decomposition_type)
00606     {
00607     case NONE:
00608       {
00609         if (this->use_blas_lapack)
00610           this->_lu_decompose_lapack();
00611         else
00612           this->_lu_decompose ();
00613         break;
00614       }
00615 
00616     case LU_BLAS_LAPACK:
00617       {
00618         // Already factored, just need to call back_substitute.
00619         if (this->use_blas_lapack)
00620           break;
00621       }
00622 
00623     case LU:
00624       {
00625         // Already factored, just need to call back_substitute.
00626         if ( !(this->use_blas_lapack) )
00627           break;
00628       }
00629 
00630     default:
00631       libmesh_error_msg("Error! This matrix already has a different decomposition...");
00632     }
00633 
00634   if (this->use_blas_lapack)
00635     this->_lu_back_substitute_lapack (b, x);
00636   else
00637     this->_lu_back_substitute (b, x);
00638 }
00639 
00640 
00641 
00642 
00643 
00644 
00645 template<typename T>
00646 void DenseMatrix<T>::_lu_back_substitute (const DenseVector<T>& b,
00647                                           DenseVector<T>& x ) const
00648 {
00649   const unsigned int
00650     n_cols = this->n();
00651 
00652   libmesh_assert_equal_to (this->m(), n_cols);
00653   libmesh_assert_equal_to (this->m(), b.size());
00654 
00655   x.resize (n_cols);
00656 
00657   // A convenient reference to *this
00658   const DenseMatrix<T>& A = *this;
00659 
00660   // Temporary vector storage.  We use this instead of
00661   // modifying the RHS.
00662   DenseVector<T> z = b;
00663 
00664   // Lower-triangular "top to bottom" solve step, taking into account pivots
00665   for (unsigned int i=0; i<n_cols; ++i)
00666     {
00667       // Swap
00668       if (_pivots[i] != static_cast<int>(i))
00669         std::swap( z(i), z(_pivots[i]) );
00670 
00671       x(i) = z(i);
00672 
00673       for (unsigned int j=0; j<i; ++j)
00674         x(i) -= A(i,j)*x(j);
00675 
00676       x(i) /= A(i,i);
00677     }
00678 
00679   // Upper-triangular "bottom to top" solve step
00680   const unsigned int last_row = n_cols-1;
00681 
00682   for (int i=last_row; i>=0; --i)
00683     {
00684       for (int j=i+1; j<static_cast<int>(n_cols); ++j)
00685         x(i) -= A(i,j)*x(j);
00686     }
00687 }
00688 
00689 
00690 
00691 
00692 
00693 
00694 
00695 
00696 template<typename T>
00697 void DenseMatrix<T>::_lu_decompose ()
00698 {
00699   // If this function was called, there better not be any
00700   // previous decomposition of the matrix.
00701   libmesh_assert_equal_to (this->_decomposition_type, NONE);
00702 
00703   // Get the matrix size and make sure it is square
00704   const unsigned int
00705     n_rows = this->m();
00706 
00707   // A convenient reference to *this
00708   DenseMatrix<T>& A = *this;
00709 
00710   _pivots.resize(n_rows);
00711 
00712   for (unsigned int i=0; i<n_rows; ++i)
00713     {
00714       // Find the pivot row by searching down the i'th column
00715       _pivots[i] = i;
00716 
00717       // std::abs(complex) must return a Real!
00718       Real the_max = std::abs( A(i,i) );
00719       for (unsigned int j=i+1; j<n_rows; ++j)
00720         {
00721           Real candidate_max = std::abs( A(j,i) );
00722           if (the_max < candidate_max)
00723             {
00724               the_max = candidate_max;
00725               _pivots[i] = j;
00726             }
00727         }
00728 
00729       // libMesh::out << "the_max=" << the_max << " found at row " << _pivots[i] << std::endl;
00730 
00731       // If the max was found in a different row, interchange rows.
00732       // Here we interchange the *entire* row, in Gaussian elimination
00733       // you would only interchange the subrows A(i,j) and A(p(i),j), for j>i
00734       if (_pivots[i] != static_cast<int>(i))
00735         {
00736           for (unsigned int j=0; j<n_rows; ++j)
00737             std::swap( A(i,j), A(_pivots[i], j) );
00738         }
00739 
00740 
00741       // If the max abs entry found is zero, the matrix is singular
00742       if (A(i,i) == libMesh::zero)
00743         libmesh_error_msg("Matrix A is singular!");
00744 
00745       // Scale upper triangle entries of row i by the diagonal entry
00746       // Note: don't scale the diagonal entry itself!
00747       const T diag_inv = 1. / A(i,i);
00748       for (unsigned int j=i+1; j<n_rows; ++j)
00749         A(i,j) *= diag_inv;
00750 
00751       // Update the remaining sub-matrix A[i+1:m][i+1:m]
00752       // by subtracting off (the diagonal-scaled)
00753       // upper-triangular part of row i, scaled by the
00754       // i'th column entry of each row.  In terms of
00755       // row operations, this is:
00756       // for each r > i
00757       //   SubRow(r) = SubRow(r) - A(r,i)*SubRow(i)
00758       //
00759       // If we were scaling the i'th column as well, like
00760       // in Gaussian elimination, this would 'zero' the
00761       // entry in the i'th column.
00762       for (unsigned int row=i+1; row<n_rows; ++row)
00763         for (unsigned int col=i+1; col<n_rows; ++col)
00764           A(row,col) -= A(row,i) * A(i,col);
00765 
00766     } // end i loop
00767 
00768   // Set the flag for LU decomposition
00769   this->_decomposition_type = LU;
00770 }
00771 
00772 
00773 
00774 template<typename T>
00775 void DenseMatrix<T>::svd (DenseVector<T>& sigma)
00776 {
00777   // We use the LAPACK svd implementation
00778   _svd_lapack(sigma);
00779 }
00780 
00781 
00782 template<typename T>
00783 void DenseMatrix<T>::svd (DenseVector<T>& sigma, DenseMatrix<T>& U, DenseMatrix<T>& VT)
00784 {
00785   // We use the LAPACK svd implementation
00786   _svd_lapack(sigma, U, VT);
00787 }
00788 
00789 
00790 
00791 template<typename T>
00792 void DenseMatrix<T>::evd (DenseVector<T>& lambda_real, DenseVector<T>& lambda_imag)
00793 {
00794   // We use the LAPACK eigenvalue problem implementation
00795   _evd_lapack(lambda_real, lambda_imag);
00796 }
00797 
00798 
00799 
00800 template<typename T>
00801 T DenseMatrix<T>::det ()
00802 {
00803   switch(this->_decomposition_type)
00804     {
00805     case NONE:
00806       {
00807         // First LU decompose the matrix.
00808         // Note that the lu_decompose routine will check to see if the
00809         // matrix is square so we don't worry about it.
00810         if (this->use_blas_lapack)
00811           this->_lu_decompose_lapack();
00812         else
00813           this->_lu_decompose ();
00814       }
00815     case LU:
00816     case LU_BLAS_LAPACK:
00817       {
00818         // Already decomposed, don't do anything
00819         break;
00820       }
00821     default:
00822       libmesh_error_msg("Error! Can't compute the determinant under the current decomposition.");
00823     }
00824 
00825   // A variable to keep track of the running product of diagonal terms.
00826   T determinant = 1.;
00827 
00828   // Loop over diagonal terms, computing the product.  In practice,
00829   // be careful because this value could easily become too large to
00830   // fit in a double or float.  To be safe, one should keep track of
00831   // the power (of 10) of the determinant in a separate variable
00832   // and maintain an order 1 value for the determinant itself.
00833   unsigned int n_interchanges = 0;
00834   for (unsigned int i=0; i<this->m(); i++)
00835     {
00836       if (this->_decomposition_type==LU)
00837         if (_pivots[i] != static_cast<int>(i))
00838           n_interchanges++;
00839 
00840       // Lapack pivots are 1-based!
00841       if (this->_decomposition_type==LU_BLAS_LAPACK)
00842         if (_pivots[i] != static_cast<int>(i+1))
00843           n_interchanges++;
00844 
00845       determinant *= (*this)(i,i);
00846     }
00847 
00848   // Compute sign of determinant, depends on number of row interchanges!
00849   // The sign should be (-1)^{n}, where n is the number of interchanges.
00850   Real sign = n_interchanges % 2 == 0 ? 1. : -1.;
00851 
00852   return sign*determinant;
00853 }
00854 
00855 
00856 
00857 // The cholesky solve function first decomposes the matrix
00858 // with cholesky_decompose and then uses the cholesky_back_substitute
00859 // routine to find the solution x.
00860 template <typename T>
00861 template <typename T2>
00862 void DenseMatrix<T>::cholesky_solve (const DenseVector<T2>& b,
00863                                      DenseVector<T2>& x)
00864 {
00865   // Check for a previous decomposition
00866   switch(this->_decomposition_type)
00867     {
00868     case NONE:
00869       {
00870         this->_cholesky_decompose ();
00871         break;
00872       }
00873 
00874     case CHOLESKY:
00875       {
00876         // Already factored, just need to call back_substitute.
00877         break;
00878       }
00879 
00880     default:
00881       libmesh_error_msg("Error! This matrix already has a different decomposition...");
00882     }
00883 
00884   // Perform back substitution
00885   this->_cholesky_back_substitute (b, x);
00886 }
00887 
00888 
00889 
00890 
00891 // This algorithm is based on the Cholesky decomposition in
00892 // the Numerical Recipes in C book.
00893 template<typename T>
00894 void DenseMatrix<T>::_cholesky_decompose ()
00895 {
00896   // If we called this function, there better not be any
00897   // previous decomposition of the matrix.
00898   libmesh_assert_equal_to (this->_decomposition_type, NONE);
00899 
00900   // Shorthand notation for number of rows and columns.
00901   const unsigned int
00902     n_rows = this->m(),
00903     n_cols = this->n();
00904 
00905   // Just to be really sure...
00906   libmesh_assert_equal_to (n_rows, n_cols);
00907 
00908   // A convenient reference to *this
00909   DenseMatrix<T>& A = *this;
00910 
00911   for (unsigned int i=0; i<n_rows; ++i)
00912     {
00913       for (unsigned int j=i; j<n_cols; ++j)
00914         {
00915           for (unsigned int k=0; k<i; ++k)
00916             A(i,j) -= A(i,k) * A(j,k);
00917 
00918           if (i == j)
00919             {
00920 #ifndef LIBMESH_USE_COMPLEX_NUMBERS
00921               if (A(i,j) <= 0.0)
00922                 libmesh_error_msg("Error! Can only use Cholesky decomposition with symmetric positive definite matrices.");
00923 #endif
00924 
00925               A(i,i) = std::sqrt(A(i,j));
00926             }
00927           else
00928             A(j,i) = A(i,j) / A(i,i);
00929         }
00930     }
00931 
00932   // Set the flag for CHOLESKY decomposition
00933   this->_decomposition_type = CHOLESKY;
00934 }
00935 
00936 
00937 
00938 template <typename T>
00939 template <typename T2>
00940 void DenseMatrix<T>::_cholesky_back_substitute (const DenseVector<T2>& b,
00941                                                 DenseVector<T2>& x) const
00942 {
00943   // Shorthand notation for number of rows and columns.
00944   const unsigned int
00945     n_rows = this->m(),
00946     n_cols = this->n();
00947 
00948   // Just to be really sure...
00949   libmesh_assert_equal_to (n_rows, n_cols);
00950 
00951   // A convenient reference to *this
00952   const DenseMatrix<T>& A = *this;
00953 
00954   // Now compute the solution to Ax =b using the factorization.
00955   x.resize(n_rows);
00956 
00957   // Solve for Ly=b
00958   for (unsigned int i=0; i<n_cols; ++i)
00959     {
00960       T2 temp = b(i);
00961 
00962       for (unsigned int k=0; k<i; ++k)
00963         temp -= A(i,k)*x(k);
00964 
00965       x(i) = temp / A(i,i);
00966     }
00967 
00968   // Solve for L^T x = y
00969   for (unsigned int i=0; i<n_cols; ++i)
00970     {
00971       const unsigned int ib = (n_cols-1)-i;
00972 
00973       for (unsigned int k=(ib+1); k<n_cols; ++k)
00974         x(ib) -= A(k,ib) * x(k);
00975 
00976       x(ib) /= A(ib,ib);
00977     }
00978 }
00979 
00980 
00981 
00982 
00983 
00984 
00985 
00986 
00987 // This routine is commented out since it is not really a memory
00988 // efficient implementation.  Also, you don't *need* the inverse
00989 // for anything, instead just use lu_solve to solve Ax=b.
00990 // template<typename T>
00991 // void DenseMatrix<T>::inverse ()
00992 // {
00993 //   // First LU decompose the matrix
00994 //   // Note that the lu_decompose routine will check to see if the
00995 //   // matrix is square so we don't worry about it.
00996 //   if (!this->_lu_decomposed)
00997 //     this->_lu_decompose();
00998 
00999 //   // A unit vector which will be used as a rhs
01000 //   // to pick off a single value each time.
01001 //   DenseVector<T> e;
01002 //   e.resize(this->m());
01003 
01004 //   // An empty vector which will be used to hold the solution
01005 //   // to the back substitutions.
01006 //   DenseVector<T> x;
01007 //   x.resize(this->m());
01008 
01009 //   // An empty dense matrix to store the resulting inverse
01010 //   // temporarily until we can overwrite A.
01011 //   DenseMatrix<T> inv;
01012 //   inv.resize(this->m(), this->n());
01013 
01014 //   // Resize the passed in matrix to hold the inverse
01015 //   inv.resize(this->m(), this->n());
01016 
01017 //   for (unsigned int j=0; j<this->n(); ++j)
01018 //     {
01019 //       e.zero();
01020 //       e(j) = 1.;
01021 //       this->_lu_back_substitute(e, x, false);
01022 //       for (unsigned int i=0; i<this->n(); ++i)
01023 // inv(i,j) = x(i);
01024 //     }
01025 
01026 //   // Now overwrite all the entries
01027 //   *this = inv;
01028 // }
01029 
01030 
01031 //--------------------------------------------------------------
01032 // Explicit instantiations
01033 #define LIBMESH_VMA_INSTANTIATE(T1,T2,T3)               \
01034   template void DenseMatrix<T1>::vector_mult_add        \
01035   (DenseVector<                                         \
01036    CompareTypes<T1,                                     \
01037    CompareTypes<T2,T3>::supertype>::supertype>& dest,   \
01038    const T2 factor,                                     \
01039    const DenseVector<T3>& arg) const
01040 
01041 template class DenseMatrix<Real>;
01042 template void DenseMatrix<Real>::cholesky_solve(const DenseVector<Real>&, DenseVector<Real>&);
01043 template void DenseMatrix<Real>::_cholesky_back_substitute(const DenseVector<Real>&, DenseVector<Real>&) const;
01044 template void DenseMatrix<Real>::cholesky_solve(const DenseVector<Complex>&, DenseVector<Complex>&);
01045 template void DenseMatrix<Real>::_cholesky_back_substitute(const DenseVector<Complex>&, DenseVector<Complex>&) const;
01046 LIBMESH_VMA_INSTANTIATE(Real,int,Real);
01047 #ifndef LIBMESH_DEFAULT_SINGLE_PRECISION
01048 LIBMESH_VMA_INSTANTIATE(Real,float,Real);
01049 #endif
01050 #ifndef LIBMESH_DEFAULT_DOUBLE_PRECISION
01051 LIBMESH_VMA_INSTANTIATE(Real,double,Real);
01052 #endif
01053 
01054 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
01055 template class DenseMatrix<Complex>;
01056 template void DenseMatrix<Complex>::cholesky_solve(const DenseVector<Complex>&,DenseVector<Complex>&);
01057 template void DenseMatrix<Complex>::_cholesky_back_substitute(const DenseVector<Complex>&, DenseVector<Complex>&) const;
01058 template void DenseMatrix<Real>::vector_mult (DenseVector<CompareTypes<Real,Complex>::supertype>& dest,
01059                                               const DenseVector<Complex>& arg) const;
01060 template void DenseMatrix<Real>::vector_mult_transpose (DenseVector<CompareTypes<Real,Complex>::supertype>& dest,
01061                                                         const DenseVector<Complex>& arg) const;
01062 LIBMESH_VMA_INSTANTIATE(Real,int,Complex);
01063 LIBMESH_VMA_INSTANTIATE(Complex,int,Complex);
01064 LIBMESH_VMA_INSTANTIATE(Complex,int,Real);
01065 
01066 // complex<int> and complex<float_foo> don't interact well
01067 //LIBMESH_VMA_INSTANTIATE(Real,std::complex<int>,Complex);
01068 //LIBMESH_VMA_INSTANTIATE(Complex,std::complex<int>,Complex);
01069 //LIBMESH_VMA_INSTANTIATE(Complex,std::complex<int>,Real);
01070 
01071 LIBMESH_VMA_INSTANTIATE(Real,float,Complex);
01072 LIBMESH_VMA_INSTANTIATE(Complex,float,Complex);
01073 LIBMESH_VMA_INSTANTIATE(Complex,float,Real);
01074 LIBMESH_VMA_INSTANTIATE(Real,std::complex<float>,Complex);
01075 #ifndef LIBMESH_DEFAULT_SINGLE_PRECISION
01076 LIBMESH_VMA_INSTANTIATE(Complex,std::complex<float>,Complex);
01077 #endif
01078 LIBMESH_VMA_INSTANTIATE(Complex,std::complex<float>,Real);
01079 
01080 LIBMESH_VMA_INSTANTIATE(Real,double,Complex);
01081 LIBMESH_VMA_INSTANTIATE(Complex,double,Complex);
01082 LIBMESH_VMA_INSTANTIATE(Complex,double,Real);
01083 LIBMESH_VMA_INSTANTIATE(Real,std::complex<double>,Complex);
01084 #ifndef LIBMESH_DEFAULT_DOUBLE_PRECISION
01085 LIBMESH_VMA_INSTANTIATE(Complex,std::complex<double>,Complex);
01086 #endif
01087 LIBMESH_VMA_INSTANTIATE(Complex,std::complex<double>,Real);
01088 #endif
01089 
01090 } // namespace libMesh