$extrastylesheet
trilinos_epetra_matrix.h
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 #ifndef LIBMESH_TRILINOS_EPETRA_MATRIX_H
00021 #define LIBMESH_TRILINOS_EPETRA_MATRIX_H
00022 
00023 #include "libmesh/libmesh_common.h"
00024 
00025 #ifdef LIBMESH_HAVE_TRILINOS
00026 
00027 // Trilinos includes
00028 #include <Epetra_FECrsMatrix.h>
00029 #include <Epetra_Map.h>
00030 #include <EpetraExt_MatrixMatrix.h>
00031 #include <Epetra_MpiComm.h>
00032 
00033 // Local includes
00034 #include "libmesh/sparse_matrix.h"
00035 
00036 // C++ includes
00037 #include <algorithm>
00038 #include <cstddef>
00039 
00040 namespace libMesh
00041 {
00042 
00043 // Forward Declarations
00044 template <typename T> class DenseMatrix;
00045 
00046 
00047 
00056 template <typename T>
00057 class EpetraMatrix : public SparseMatrix<T>
00058 {
00059 public:
00075   EpetraMatrix (const Parallel::Communicator &comm
00076                 LIBMESH_CAN_DEFAULT_TO_COMMWORLD);
00077 
00085   EpetraMatrix (Epetra_FECrsMatrix * m,
00086                 const Parallel::Communicator &comm
00087                 LIBMESH_CAN_DEFAULT_TO_COMMWORLD);
00088 
00094   virtual ~EpetraMatrix ();
00095 
00099   bool need_full_sparsity_pattern () const
00100   { return true; }
00101 
00107   void update_sparsity_pattern (const SparsityPattern::Graph &);
00108 
00119   void init (const numeric_index_type m,
00120              const numeric_index_type n,
00121              const numeric_index_type m_l,
00122              const numeric_index_type n_l,
00123              const numeric_index_type nnz=30,
00124              const numeric_index_type noz=10,
00125              const numeric_index_type blocksize=1);
00126 
00130   void init ();
00131 
00138   void clear ();
00139 
00144   void zero ();
00145 
00151   void close () const;
00152 
00157   numeric_index_type m () const;
00158 
00163   numeric_index_type n () const;
00164 
00169   numeric_index_type row_start () const;
00170 
00175   numeric_index_type row_stop () const;
00176 
00183   void set (const numeric_index_type i,
00184             const numeric_index_type j,
00185             const T value);
00186 
00195   void add (const numeric_index_type i,
00196             const numeric_index_type j,
00197             const T value);
00198 
00206   void add_matrix (const DenseMatrix<T> &dm,
00207                    const std::vector<numeric_index_type> &rows,
00208                    const std::vector<numeric_index_type> &cols);
00209 
00214   void add_matrix (const DenseMatrix<T> &dm,
00215                    const std::vector<numeric_index_type> &dof_indices);
00216 
00226   void add (const T a, SparseMatrix<T> &X);
00227 
00235   T operator () (const numeric_index_type i,
00236                  const numeric_index_type j) const;
00237 
00249   Real l1_norm () const;
00250 
00263   Real linfty_norm () const;
00264 
00269   bool closed() const;
00270 
00274   void print_personal(std::ostream& os=libMesh::out) const;
00275 
00279   void get_diagonal (NumericVector<T>& dest) const;
00280 
00285   virtual void get_transpose (SparseMatrix<T>& dest) const;
00286 
00290   void swap (EpetraMatrix<T> &);
00291 
00297   Epetra_FECrsMatrix * mat () { libmesh_assert(_mat); return _mat; }
00298 
00299   const Epetra_FECrsMatrix * mat () const { libmesh_assert(_mat); return _mat; }
00300 
00301 
00302 private:
00303 
00308   Epetra_FECrsMatrix * _mat;
00309 
00313   Epetra_Map * _map;
00314 
00318   Epetra_CrsGraph * _graph;
00319 
00324   bool _destroy_mat_on_exit;
00325 
00330   bool _use_transpose;
00331 };
00332 
00333 } // namespace libMesh
00334 
00335 #endif // #ifdef LIBMESH_HAVE_TRILINOS
00336 #endif // LIBMESH_TRILINOS_EPETRA_MATRIX_H