laspack_matrix.h
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 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 #ifndef LIBMESH_LASPACK_MATRIX_H 00021 #define LIBMESH_LASPACK_MATRIX_H 00022 00023 #include "libmesh/libmesh_config.h" 00024 00025 #ifdef LIBMESH_HAVE_LASPACK 00026 00027 // Local includes 00028 #include "libmesh/sparse_matrix.h" 00029 00030 // Laspack includes 00031 #include <qmatrix.h> 00032 00033 // C++ includes 00034 #include <algorithm> 00035 #include <cstddef> 00036 00037 namespace libMesh 00038 { 00039 00040 00041 00042 // Forward declarations 00043 template <typename T> class DenseMatrix; 00044 template <typename T> class LaspackVector; 00045 template <typename T> class LaspackLinearSolver; 00046 00047 00048 00058 template <typename T> 00059 class LaspackMatrix : public SparseMatrix<T> 00060 { 00061 00062 public: 00078 LaspackMatrix (); 00079 00085 ~LaspackMatrix (); 00086 00090 bool need_full_sparsity_pattern() const 00091 { return true; } 00092 00098 void update_sparsity_pattern (const SparsityPattern::Graph &); 00099 00108 void init (const numeric_index_type m, 00109 const numeric_index_type n, 00110 const numeric_index_type m_l, 00111 const numeric_index_type n_l, 00112 const numeric_index_type nnz=30, 00113 const numeric_index_type noz=10); 00114 00118 void init (); 00119 00126 void clear (); 00127 00131 void zero (); 00132 00138 void close () const { const_cast<LaspackMatrix<T>*>(this)->_closed = true; } 00139 00144 numeric_index_type m () const; 00145 00150 numeric_index_type n () const; 00151 00156 numeric_index_type row_start () const; 00157 00162 numeric_index_type row_stop () const; 00163 00170 void set (const numeric_index_type i, 00171 const numeric_index_type j, 00172 const T value); 00173 00182 void add (const numeric_index_type i, 00183 const numeric_index_type j, 00184 const T value); 00185 00193 void add_matrix (const DenseMatrix<T> &dm, 00194 const std::vector<numeric_index_type> &rows, 00195 const std::vector<numeric_index_type> &cols); 00196 00201 void add_matrix (const DenseMatrix<T> &dm, 00202 const std::vector<numeric_index_type> &dof_indices); 00203 00211 void add (const T a, SparseMatrix<T> &X); 00212 00220 T operator () (const numeric_index_type i, 00221 const numeric_index_type j) const; 00222 00233 Real l1_norm () const { libmesh_error(); return 0.; } 00234 00246 Real linfty_norm () const { libmesh_error(); return 0.; } 00247 00252 bool closed() const { return _closed; } 00253 00258 void print_personal(std::ostream& os=libMesh::out) const { this->print(os); } 00259 00263 virtual void get_diagonal (NumericVector<T>& dest) const; 00264 00269 virtual void get_transpose (SparseMatrix<T>& dest) const; 00270 00271 private: 00272 00277 numeric_index_type pos (const numeric_index_type i, 00278 const numeric_index_type j) const; 00279 00283 QMatrix _QMat; 00284 00288 std::vector<numeric_index_type> _csr; 00289 00294 std::vector<std::vector<numeric_index_type>::const_iterator> _row_start; 00295 00299 bool _closed; 00300 00304 friend class LaspackVector<T>; 00305 friend class LaspackLinearSolver<T>; 00306 }; 00307 00308 } // namespace libMesh 00309 00310 #endif // #ifdef LIBMESH_HAVE_LASPACK 00311 #endif // #ifdef LIBMESH_LASPACK_MATRIX_H
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:47 UTC
Hosted By: