$extrastylesheet
sparse_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 
00020 // C++ includes
00021 
00022 // Local Includes
00023 #include "libmesh/dof_map.h"
00024 #include "libmesh/dense_matrix.h"
00025 #include "libmesh/laspack_matrix.h"
00026 #include "libmesh/eigen_sparse_matrix.h"
00027 #include "libmesh/parallel.h"
00028 #include "libmesh/petsc_matrix.h"
00029 #include "libmesh/sparse_matrix.h"
00030 #include "libmesh/trilinos_epetra_matrix.h"
00031 #include "libmesh/numeric_vector.h"
00032 
00033 namespace libMesh
00034 {
00035 
00036 
00037 //------------------------------------------------------------------
00038 // SparseMatrix Methods
00039 
00040 
00041 // Constructor
00042 template <typename T>
00043 SparseMatrix<T>::SparseMatrix (const Parallel::Communicator &comm_in) :
00044   ParallelObject(comm_in),
00045   _dof_map(NULL),
00046   _is_initialized(false)
00047 {}
00048 
00049 
00050 
00051 // Destructor
00052 template <typename T>
00053 SparseMatrix<T>::~SparseMatrix ()
00054 {}
00055 
00056 
00057 
00058 
00059 // default implementation is to fall back to non-blocked method
00060 template <typename T>
00061 void SparseMatrix<T>::add_block_matrix (const DenseMatrix<T> &dm,
00062                                         const std::vector<numeric_index_type> &brows,
00063                                         const std::vector<numeric_index_type> &bcols)
00064 {
00065   libmesh_assert_equal_to (dm.m() / brows.size(), dm.n() / bcols.size());
00066 
00067   const numeric_index_type blocksize = cast_int<numeric_index_type>
00068     (dm.m() / brows.size());
00069 
00070   libmesh_assert_equal_to (dm.m()%blocksize, 0);
00071   libmesh_assert_equal_to (dm.n()%blocksize, 0);
00072 
00073   std::vector<numeric_index_type> rows, cols;
00074 
00075   rows.reserve(blocksize*brows.size());
00076   cols.reserve(blocksize*bcols.size());
00077 
00078   for (unsigned int ib=0; ib<brows.size(); ib++)
00079     {
00080       numeric_index_type i=brows[ib]*blocksize;
00081 
00082       for (unsigned int v=0; v<blocksize; v++)
00083         rows.push_back(i++);
00084     }
00085 
00086   for (unsigned int jb=0; jb<bcols.size(); jb++)
00087     {
00088       numeric_index_type j=bcols[jb]*blocksize;
00089 
00090       for (unsigned int v=0; v<blocksize; v++)
00091         cols.push_back(j++);
00092     }
00093 
00094   this->add_matrix (dm, rows, cols);
00095 }
00096 
00097 
00098 
00099 // Full specialization of print method for Complex datatypes
00100 template <>
00101 void SparseMatrix<Complex>::print(std::ostream& os, const bool sparse) const
00102 {
00103   // std::complex<>::operator<<() is defined, but use this form
00104 
00105   if(sparse)
00106     {
00107       libmesh_not_implemented();
00108     }
00109 
00110   os << "Real part:" << std::endl;
00111   for (numeric_index_type i=0; i<this->m(); i++)
00112     {
00113       for (numeric_index_type j=0; j<this->n(); j++)
00114         os << std::setw(8) << (*this)(i,j).real() << " ";
00115       os << std::endl;
00116     }
00117 
00118   os << std::endl << "Imaginary part:" << std::endl;
00119   for (numeric_index_type i=0; i<this->m(); i++)
00120     {
00121       for (numeric_index_type j=0; j<this->n(); j++)
00122         os << std::setw(8) << (*this)(i,j).imag() << " ";
00123       os << std::endl;
00124     }
00125 }
00126 
00127 
00128 
00129 
00130 
00131 
00132 // Full specialization for Real datatypes
00133 template <typename T>
00134 UniquePtr<SparseMatrix<T> >
00135 SparseMatrix<T>::build(const Parallel::Communicator &comm,
00136                        const SolverPackage solver_package)
00137 {
00138   // Build the appropriate vector
00139   switch (solver_package)
00140     {
00141 
00142 #ifdef LIBMESH_HAVE_LASPACK
00143     case LASPACK_SOLVERS:
00144       return UniquePtr<SparseMatrix<T> >(new LaspackMatrix<T>(comm));
00145 #endif
00146 
00147 
00148 #ifdef LIBMESH_HAVE_PETSC
00149     case PETSC_SOLVERS:
00150       return UniquePtr<SparseMatrix<T> >(new PetscMatrix<T>(comm));
00151 #endif
00152 
00153 
00154 #ifdef LIBMESH_HAVE_TRILINOS
00155     case TRILINOS_SOLVERS:
00156       return UniquePtr<SparseMatrix<T> >(new EpetraMatrix<T>(comm));
00157 #endif
00158 
00159 
00160 #ifdef LIBMESH_HAVE_EIGEN
00161     case EIGEN_SOLVERS:
00162       return UniquePtr<SparseMatrix<T> >(new EigenSparseMatrix<T>(comm));
00163 #endif
00164 
00165     default:
00166       libmesh_error_msg("ERROR:  Unrecognized solver package: " << solver_package);
00167     }
00168 
00169   libmesh_error_msg("We'll never get here!");
00170   return UniquePtr<SparseMatrix<T> >();
00171 }
00172 
00173 
00174 template <typename T>
00175 void SparseMatrix<T>::vector_mult (NumericVector<T>& dest,
00176                                    const NumericVector<T>& arg) const
00177 {
00178   dest.zero();
00179   this->vector_mult_add(dest,arg);
00180 }
00181 
00182 
00183 
00184 template <typename T>
00185 void SparseMatrix<T>::vector_mult_add (NumericVector<T>& dest,
00186                                        const NumericVector<T>& arg) const
00187 {
00188   /* This functionality is actually implemented in the \p
00189      NumericVector class.  */
00190   dest.add_vector(arg,*this);
00191 }
00192 
00193 
00194 
00195 template <typename T>
00196 void SparseMatrix<T>::zero_rows (std::vector<numeric_index_type> &, T)
00197 {
00198   /* This functionality isn't implemented or stubbed in every subclass yet */
00199   libmesh_not_implemented();
00200 }
00201 
00202 
00203 
00204 template <typename T>
00205 void SparseMatrix<T>::print(std::ostream& os, const bool sparse) const
00206 {
00207   parallel_object_only();
00208 
00209   libmesh_assert (this->initialized());
00210 
00211   if(!this->_dof_map)
00212     libmesh_error_msg("Error!  Trying to print a matrix with no dof_map set!");
00213 
00214   // We'll print the matrix from processor 0 to make sure
00215   // it's serialized properly
00216   if (this->processor_id() == 0)
00217     {
00218       libmesh_assert_equal_to (this->_dof_map->first_dof(), 0);
00219       for (numeric_index_type i=this->_dof_map->first_dof();
00220            i!=this->_dof_map->end_dof(); ++i)
00221         {
00222           if(sparse)
00223             {
00224               for (numeric_index_type j=0; j<this->n(); j++)
00225                 {
00226                   T c = (*this)(i,j);
00227                   if (c != static_cast<T>(0.0))
00228                     {
00229                       os << i << " " << j << " " << c << std::endl;
00230                     }
00231                 }
00232             }
00233           else
00234             {
00235               for (numeric_index_type j=0; j<this->n(); j++)
00236                 os << (*this)(i,j) << " ";
00237               os << std::endl;
00238             }
00239         }
00240 
00241       std::vector<numeric_index_type> ibuf, jbuf;
00242       std::vector<T> cbuf;
00243       numeric_index_type currenti = this->_dof_map->end_dof();
00244       for (processor_id_type p=1; p < this->n_processors(); ++p)
00245         {
00246           this->comm().receive(p, ibuf);
00247           this->comm().receive(p, jbuf);
00248           this->comm().receive(p, cbuf);
00249           libmesh_assert_equal_to (ibuf.size(), jbuf.size());
00250           libmesh_assert_equal_to (ibuf.size(), cbuf.size());
00251 
00252           if (ibuf.empty())
00253             continue;
00254           libmesh_assert_greater_equal (ibuf.front(), currenti);
00255           libmesh_assert_greater_equal (ibuf.back(), ibuf.front());
00256 
00257           std::size_t currentb = 0;
00258           for (;currenti <= ibuf.back(); ++currenti)
00259             {
00260               if(sparse)
00261                 {
00262                   for (numeric_index_type j=0; j<this->n(); j++)
00263                     {
00264                       if (currentb < ibuf.size() &&
00265                           ibuf[currentb] == currenti &&
00266                           jbuf[currentb] == j)
00267                         {
00268                           os << currenti << " " << j << " " << cbuf[currentb] << std::endl;
00269                           currentb++;
00270                         }
00271                     }
00272                 }
00273               else
00274                 {
00275                   for (numeric_index_type j=0; j<this->n(); j++)
00276                     {
00277                       if (currentb < ibuf.size() &&
00278                           ibuf[currentb] == currenti &&
00279                           jbuf[currentb] == j)
00280                         {
00281                           os << cbuf[currentb] << " ";
00282                           currentb++;
00283                         }
00284                       else
00285                         os << static_cast<T>(0.0) << " ";
00286                     }
00287                   os << std::endl;
00288                 }
00289             }
00290         }
00291       if(!sparse)
00292         {
00293           for (; currenti != this->m(); ++currenti)
00294             {
00295               for (numeric_index_type j=0; j<this->n(); j++)
00296                 os << static_cast<T>(0.0) << " ";
00297               os << std::endl;
00298             }
00299         }
00300     }
00301   else
00302     {
00303       std::vector<numeric_index_type> ibuf, jbuf;
00304       std::vector<T> cbuf;
00305 
00306       // We'll assume each processor has access to entire
00307       // matrix rows, so (*this)(i,j) is valid if i is a local index.
00308       for (numeric_index_type i=this->_dof_map->first_dof();
00309            i!=this->_dof_map->end_dof(); ++i)
00310         {
00311           for (numeric_index_type j=0; j<this->n(); j++)
00312             {
00313               T c = (*this)(i,j);
00314               if (c != static_cast<T>(0.0))
00315                 {
00316                   ibuf.push_back(i);
00317                   jbuf.push_back(j);
00318                   cbuf.push_back(c);
00319                 }
00320             }
00321         }
00322       this->comm().send(0,ibuf);
00323       this->comm().send(0,jbuf);
00324       this->comm().send(0,cbuf);
00325     }
00326 }
00327 
00328 
00329 
00330 //------------------------------------------------------------------
00331 // Explicit instantiations
00332 template class SparseMatrix<Number>;
00333 
00334 } // namespace libMesh