libMesh::SparseMatrix< T > Class Template Reference

#include <sparse_matrix.h>

Inheritance diagram for libMesh::SparseMatrix< T >:

List of all members.

Public Member Functions

 SparseMatrix ()
virtual ~SparseMatrix ()
virtual bool initialized () const
void attach_dof_map (const DofMap &dof_map)
virtual bool need_full_sparsity_pattern () const
virtual void update_sparsity_pattern (const SparsityPattern::Graph &)
virtual void init (const numeric_index_type m, const numeric_index_type n, const numeric_index_type m_l, const numeric_index_type n_l, const numeric_index_type nnz=30, const numeric_index_type noz=10)=0
virtual void init ()=0
virtual void clear ()=0
virtual void zero ()=0
virtual void zero_rows (std::vector< numeric_index_type > &rows, T diag_value=0.0)
virtual void close () const =0
virtual numeric_index_type m () const =0
virtual numeric_index_type n () const =0
virtual numeric_index_type row_start () const =0
virtual numeric_index_type row_stop () const =0
virtual void set (const numeric_index_type i, const numeric_index_type j, const T value)=0
virtual void add (const numeric_index_type i, const numeric_index_type j, const T value)=0
virtual void add_matrix (const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
virtual void add_matrix (const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &dof_indices)=0
virtual void add (const T, SparseMatrix< T > &)=0
virtual T operator() (const numeric_index_type i, const numeric_index_type j) const =0
virtual Real l1_norm () const =0
virtual Real linfty_norm () const =0
virtual bool closed () const =0
void print (std::ostream &os=libMesh::out, const bool sparse=false) const
virtual void print_personal (std::ostream &os=libMesh::out) const =0
virtual void print_matlab (const std::string name="NULL") const
virtual void create_submatrix (SparseMatrix< T > &submatrix, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) const
virtual void reinit_submatrix (SparseMatrix< T > &submatrix, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) const
void vector_mult (NumericVector< T > &dest, const NumericVector< T > &arg) const
void vector_mult_add (NumericVector< T > &dest, const NumericVector< T > &arg) const
virtual void get_diagonal (NumericVector< T > &dest) const =0
virtual void get_transpose (SparseMatrix< T > &dest) const =0
template<>
void print (std::ostream &os, const bool sparse) const

Static Public Member Functions

static AutoPtr< SparseMatrix< T > > build (const SolverPackage solver_package=libMesh::default_solver_package())
static std::string get_info ()
static void print_info (std::ostream &out=libMesh::out)
static unsigned int n_objects ()
static void enable_print_counter_info ()
static void disable_print_counter_info ()

Protected Types

typedef std::map< std::string,
std::pair< unsigned int,
unsigned int > > 
Counts

Protected Member Functions

virtual void _get_submatrix (SparseMatrix< T > &, const std::vector< numeric_index_type > &, const std::vector< numeric_index_type > &, const bool) const
void increment_constructor_count (const std::string &name)
void increment_destructor_count (const std::string &name)

Protected Attributes

DofMap const * _dof_map
bool _is_initialized

Static Protected Attributes

static Counts _counts
static Threads::atomic
< unsigned int > 
_n_objects
static Threads::spin_mutex _mutex
static bool _enable_print_counter = true

Friends

std::ostream & operator<< (std::ostream &os, const SparseMatrix< T > &m)

Detailed Description

template<typename T>
class libMesh::SparseMatrix< T >

Generic sparse matrix. This class contains pure virtual members that must be overloaded in derived classes. Using a common base class allows for uniform access to sparse matrices from various different solver packages in different formats.

Author:
Benjamin S. Kirk, 2003

Definition at line 64 of file sparse_matrix.h.


Member Typedef Documentation

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts [protected, inherited]

Data structure to log the information. The log is identified by the class name.

Definition at line 113 of file reference_counter.h.


Constructor & Destructor Documentation

template<typename T >
libMesh::SparseMatrix< T >::SparseMatrix (  )  [inline]

Constructor; initializes the matrix to be empty, without any structure, i.e. the matrix is not usable at all. This constructor is therefore only useful for matrices which are members of a class. All other matrices should be created at a point in the data flow where all necessary information is available.

You have to initialize the matrix before usage with init(...).

Definition at line 41 of file sparse_matrix.C.

00041                                :
00042   _dof_map(NULL),
00043   _is_initialized(false)
00044 {}

template<typename T >
libMesh::SparseMatrix< T >::~SparseMatrix (  )  [inline, virtual]

Destructor. Free all memory, but do not release the memory of the sparsity structure.

Definition at line 50 of file sparse_matrix.C.

00051 {}


Member Function Documentation

template<typename T>
virtual void libMesh::SparseMatrix< T >::_get_submatrix ( SparseMatrix< T > &  ,
const std::vector< numeric_index_type > &  ,
const std::vector< numeric_index_type > &  ,
const   bool 
) const [inline, protected, virtual]

Protected implementation of the create_submatrix and reinit_submatrix routines. Note that this function must be redefined in derived classes for it to work properly!

Definition at line 398 of file sparse_matrix.h.

Referenced by libMesh::SparseMatrix< Number >::create_submatrix(), and libMesh::SparseMatrix< Number >::reinit_submatrix().

00402   {
00403     libMesh::err << "Error! This function is not yet implemented in the base class!"
00404                   << std::endl;
00405     libmesh_error();
00406   }

template<typename T>
virtual void libMesh::SparseMatrix< T >::add ( const   T,
SparseMatrix< T > &   
) [pure virtual]

Add a Sparse matrix _X, scaled with _a, to this, stores the result in this: $\texttt{this} = \_a*\_X + \texttt{this} $.

template<typename T>
virtual void libMesh::SparseMatrix< T >::add ( const numeric_index_type  i,
const numeric_index_type  j,
const T  value 
) [pure virtual]

Add value to the element (i,j). Throws an error if the entry does not exist. Still, it is allowed to store zero values in non-existent fields.

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

Referenced by libMesh::NewmarkSystem::compute_matrix().

template<typename T>
virtual void libMesh::SparseMatrix< T >::add_matrix ( const DenseMatrix< T > &  dm,
const std::vector< numeric_index_type > &  dof_indices 
) [pure virtual]

Same, but assumes the row and column maps are the same. Thus the matrix dm must be square.

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

template<typename T>
virtual void libMesh::SparseMatrix< T >::add_matrix ( const DenseMatrix< T > &  dm,
const std::vector< numeric_index_type > &  rows,
const std::vector< numeric_index_type > &  cols 
) [pure virtual]

Add the full matrix dm to the Sparse matrix. This is useful for adding an element matrix at assembly time

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

template<typename T>
void libMesh::SparseMatrix< T >::attach_dof_map ( const DofMap dof_map  )  [inline]

Get a pointer to the DofMap to use.

Definition at line 107 of file sparse_matrix.h.

Referenced by libMesh::DofMap::attach_matrix(), libMesh::Problem_Interface::computeJacobian(), and DMJacobian_libMesh().

00108   { _dof_map = &dof_map; }

template<typename T >
AutoPtr< SparseMatrix< T > > libMesh::SparseMatrix< T >::build ( const SolverPackage  solver_package = libMesh::default_solver_package()  )  [inline, static]

Builds a SparseMatrix<T> using the linear solver package specified by solver_package

Definition at line 91 of file sparse_matrix.C.

References libMesh::err, libMesh::LASPACK_SOLVERS, libMeshEnums::PETSC_SOLVERS, and libMesh::TRILINOS_SOLVERS.

00092 {
00093   // Build the appropriate vector
00094   switch (solver_package)
00095     {
00096 
00097 
00098 #ifdef LIBMESH_HAVE_LASPACK
00099     case LASPACK_SOLVERS:
00100       {
00101         AutoPtr<SparseMatrix<T> > ap(new LaspackMatrix<T>);
00102         return ap;
00103       }
00104 #endif
00105 
00106 
00107 #ifdef LIBMESH_HAVE_PETSC
00108     case PETSC_SOLVERS:
00109       {
00110         AutoPtr<SparseMatrix<T> > ap(new PetscMatrix<T>);
00111         return ap;
00112       }
00113 #endif
00114 
00115 
00116 #ifdef LIBMESH_HAVE_TRILINOS
00117     case TRILINOS_SOLVERS:
00118       {
00119         AutoPtr<SparseMatrix<T> > ap(new EpetraMatrix<T>);
00120         return ap;
00121       }
00122 #endif
00123 
00124 
00125     default:
00126       libMesh::err << "ERROR:  Unrecognized solver package: "
00127                     << solver_package
00128                     << std::endl;
00129       libmesh_error();
00130     }
00131 
00132   AutoPtr<SparseMatrix<T> > ap(NULL);
00133   return ap;
00134 }

template<typename T>
virtual void libMesh::SparseMatrix< T >::clear (  )  [pure virtual]

Release all memory and return to a state just like after having called the default constructor.

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

Referenced by libMesh::PetscMatrix< T >::_get_submatrix(), libMesh::PetscMatrix< T >::get_transpose(), and libMesh::EigenSystem::reinit().

template<typename T>
virtual bool libMesh::SparseMatrix< T >::closed (  )  const [pure virtual]

see if Sparse matrix has been closed and fully assembled yet

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

template<typename T>
virtual void libMesh::SparseMatrix< T >::create_submatrix ( SparseMatrix< T > &  submatrix,
const std::vector< numeric_index_type > &  rows,
const std::vector< numeric_index_type > &  cols 
) const [inline, virtual]

This function creates a matrix called "submatrix" which is defined by the row and column indices given in the "rows" and "cols" entries. Currently this operation is only defined for the PetscMatrix type.

Definition at line 341 of file sparse_matrix.h.

Referenced by libMesh::CondensedEigenSystem::solve().

00344   {
00345     this->_get_submatrix(submatrix,
00346                          rows,
00347                          cols,
00348                          false); // false means DO NOT REUSE submatrix
00349   }

void libMesh::ReferenceCounter::disable_print_counter_info (  )  [static, inherited]

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

00107 {
00108   _enable_print_counter = false;
00109   return;
00110 }

void libMesh::ReferenceCounter::enable_print_counter_info (  )  [static, inherited]

Methods to enable/disable the reference counter output from print_info()

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

00101 {
00102   _enable_print_counter = true;
00103   return;
00104 }

template<typename T>
virtual void libMesh::SparseMatrix< T >::get_diagonal ( NumericVector< T > &  dest  )  const [pure virtual]

Copies the diagonal part of the matrix into dest.

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

std::string libMesh::ReferenceCounter::get_info (  )  [static, inherited]

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().

Referenced by libMesh::ReferenceCounter::print_info().

00048 {
00049 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
00050 
00051   std::ostringstream oss;
00052 
00053   oss << '\n'
00054       << " ---------------------------------------------------------------------------- \n"
00055       << "| Reference count information                                                |\n"
00056       << " ---------------------------------------------------------------------------- \n";
00057 
00058   for (Counts::iterator it = _counts.begin();
00059        it != _counts.end(); ++it)
00060     {
00061       const std::string name(it->first);
00062       const unsigned int creations    = it->second.first;
00063       const unsigned int destructions = it->second.second;
00064 
00065       oss << "| " << name << " reference count information:\n"
00066           << "|  Creations:    " << creations    << '\n'
00067           << "|  Destructions: " << destructions << '\n';
00068     }
00069 
00070   oss << " ---------------------------------------------------------------------------- \n";
00071 
00072   return oss.str();
00073 
00074 #else
00075 
00076   return "";
00077 
00078 #endif
00079 }

template<typename T>
virtual void libMesh::SparseMatrix< T >::get_transpose ( SparseMatrix< T > &  dest  )  const [pure virtual]

Copies the transpose of the matrix into dest, which may be *this.

Referenced by libMesh::LinearSolver< T >::adjoint_solve(), and libMesh::ImplicitSystem::weighted_sensitivity_adjoint_solve().

void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name  )  [inline, protected, inherited]

Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.

Definition at line 163 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

00164 {
00165   Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
00166   std::pair<unsigned int, unsigned int>& p = _counts[name];
00167 
00168   p.first++;
00169 }

void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name  )  [inline, protected, inherited]

Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.

Definition at line 176 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

00177 {
00178   Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
00179   std::pair<unsigned int, unsigned int>& p = _counts[name];
00180 
00181   p.second++;
00182 }

template<typename T>
virtual void libMesh::SparseMatrix< T >::init (  )  [pure virtual]

Initialize using sparsity structure computed by dof_map.

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

template<typename T>
virtual void libMesh::SparseMatrix< T >::init ( const numeric_index_type  m,
const numeric_index_type  n,
const numeric_index_type  m_l,
const numeric_index_type  n_l,
const numeric_index_type  nnz = 30,
const numeric_index_type  noz = 10 
) [pure virtual]

Initialize a Sparse matrix that is of global dimension $ m \times n $ with local dimensions $ m_l \times n_l $. nnz is the number of on-processor nonzeros per row (defaults to 30). noz is the number of on-processor nonzeros per row (defaults to 10).

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

Referenced by libMesh::EigenSystem::init_matrices(), and libMesh::EigenSystem::reinit().

template<typename T>
virtual bool libMesh::SparseMatrix< T >::initialized (  )  const [inline, virtual]
Returns:
true if the matrix has been initialized, false otherwise.

Definition at line 102 of file sparse_matrix.h.

Referenced by libMesh::PetscMatrix< T >::_get_submatrix(), libMesh::EpetraMatrix< T >::add(), libMesh::PetscMatrix< T >::add(), libMesh::LaspackMatrix< T >::add(), libMesh::EpetraMatrix< T >::add_matrix(), libMesh::PetscMatrix< T >::add_matrix(), libMesh::LaspackMatrix< T >::add_matrix(), libMesh::ImplicitSystem::assemble(), libMesh::EpetraMatrix< T >::clear(), libMesh::PetscMatrix< T >::clear(), libMesh::LaspackMatrix< T >::clear(), libMesh::EpetraMatrix< T >::closed(), libMesh::PetscMatrix< T >::closed(), libMesh::EpetraMatrix< T >::init(), libMesh::PetscMatrix< T >::init(), libMesh::LaspackMatrix< T >::init(), libMesh::ImplicitSystem::init_matrices(), libMesh::EpetraMatrix< T >::l1_norm(), libMesh::PetscMatrix< T >::l1_norm(), libMesh::EpetraMatrix< T >::linfty_norm(), libMesh::PetscMatrix< T >::linfty_norm(), libMesh::EpetraMatrix< T >::m(), libMesh::PetscMatrix< T >::m(), libMesh::LaspackMatrix< T >::m(), libMesh::EpetraMatrix< T >::n(), libMesh::PetscMatrix< T >::n(), libMesh::LaspackMatrix< T >::n(), libMesh::EpetraMatrix< T >::operator()(), libMesh::PetscMatrix< T >::operator()(), libMesh::LaspackMatrix< T >::operator()(), libMesh::SparseMatrix< T >::print(), libMesh::EpetraMatrix< T >::print_matlab(), libMesh::PetscMatrix< T >::print_matlab(), libMesh::EpetraMatrix< T >::print_personal(), libMesh::PetscMatrix< T >::print_personal(), libMesh::EpetraMatrix< T >::row_start(), libMesh::PetscMatrix< T >::row_start(), libMesh::EpetraMatrix< T >::row_stop(), libMesh::PetscMatrix< T >::row_stop(), libMesh::EpetraMatrix< T >::set(), libMesh::PetscMatrix< T >::set(), libMesh::LaspackMatrix< T >::set(), libMesh::EpetraMatrix< T >::update_sparsity_pattern(), libMesh::LaspackMatrix< T >::update_sparsity_pattern(), libMesh::EpetraMatrix< T >::zero(), libMesh::PetscMatrix< T >::zero(), and libMesh::PetscMatrix< T >::zero_rows().

00102 { return _is_initialized; }

template<typename T>
virtual Real libMesh::SparseMatrix< T >::l1_norm (  )  const [pure virtual]

Return the l1-norm of the matrix, that is $|M|_1=max_{all columns j}\sum_{all rows i} |M_ij|$, (max. sum of columns). This is the natural matrix norm that is compatible to the l1-norm for vectors, i.e. $|Mv|_1\leq |M|_1 |v|_1$.

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

Referenced by libMesh::FEMSystem::assembly().

template<typename T>
virtual Real libMesh::SparseMatrix< T >::linfty_norm (  )  const [pure virtual]

Return the linfty-norm of the matrix, that is $|M|_\infty=max_{all rows i}\sum_{all columns j} |M_ij|$, (max. sum of rows). This is the natural matrix norm that is compatible to the linfty-norm of vectors, i.e. $|Mv|_\infty \leq |M|_\infty |v|_\infty$.

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

template<typename T>
virtual numeric_index_type libMesh::SparseMatrix< T >::m (  )  const [pure virtual]
template<typename T>
virtual numeric_index_type libMesh::SparseMatrix< T >::n (  )  const [pure virtual]
static unsigned int libMesh::ReferenceCounter::n_objects (  )  [inline, static, inherited]

Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 79 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

00080   { return _n_objects; }

template<typename T>
virtual bool libMesh::SparseMatrix< T >::need_full_sparsity_pattern (  )  const [inline, virtual]

returns true if this sparse matrix format needs to be fed the graph of the sparse matrix. This is true in the case of the LaspackMatrix, but not for the PetscMatrix. In the case where the full graph is not required we can efficiently approximate it to provide a good estimate of the required size of the sparse matrix.

Reimplemented in libMesh::LaspackMatrix< T >, and libMesh::EpetraMatrix< T >.

Definition at line 117 of file sparse_matrix.h.

Referenced by libMesh::DofMap::attach_matrix().

00118   { return false; }

template<typename T>
virtual T libMesh::SparseMatrix< T >::operator() ( const numeric_index_type  i,
const numeric_index_type  j 
) const [pure virtual]

Return the value of the entry (i,j). This may be an expensive operation, and you should always be careful where you call this function.

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

template<>
void libMesh::SparseMatrix< Complex >::print ( std::ostream &  os,
const bool  sparse 
) const [inline]

Definition at line 57 of file sparse_matrix.C.

References libMesh::SparseMatrix< T >::m(), and libMesh::SparseMatrix< T >::n().

00058 {
00059   // std::complex<>::operator<<() is defined, but use this form
00060 
00061   if(sparse)
00062     {
00063       libmesh_not_implemented();
00064     }
00065 
00066   os << "Real part:" << std::endl;
00067   for (numeric_index_type i=0; i<this->m(); i++)
00068     {
00069       for (numeric_index_type j=0; j<this->n(); j++)
00070         os << std::setw(8) << (*this)(i,j).real() << " ";
00071       os << std::endl;
00072     }
00073 
00074   os << std::endl << "Imaginary part:" << std::endl;
00075   for (numeric_index_type i=0; i<this->m(); i++)
00076     {
00077       for (numeric_index_type j=0; j<this->n(); j++)
00078         os << std::setw(8) << (*this)(i,j).imag() << " ";
00079       os << std::endl;
00080     }
00081 }

template<typename T >
void libMesh::SparseMatrix< T >::print ( std::ostream &  os = libMesh::out,
const bool  sparse = false 
) const [inline]

Print the contents of the matrix to the screen in a uniform style, regardless of matrix/solver package being used.

Definition at line 168 of file sparse_matrix.C.

References libMesh::SparseMatrix< T >::_dof_map, libMesh::CommWorld, libMesh::DofMap::end_dof(), libMesh::DofMap::first_dof(), libMesh::SparseMatrix< T >::initialized(), libMesh::SparseMatrix< T >::m(), libMesh::SparseMatrix< T >::n(), libMesh::n_processors(), libMesh::processor_id(), libMesh::Parallel::Communicator::receive(), and libMesh::Parallel::Communicator::send().

Referenced by libMesh::LaspackMatrix< T >::print_personal().

00169 {
00170   parallel_only();
00171 
00172   libmesh_assert (this->initialized());
00173 
00174   if(!this->_dof_map)
00175   {
00176     os << std::endl << "Error!  Trying to print a matrix with no dof_map set!" << std::endl << std::endl;
00177     libmesh_error();
00178   }
00179 
00180   // We'll print the matrix from processor 0 to make sure
00181   // it's serialized properly
00182   if (libMesh::processor_id() == 0)
00183     {
00184       libmesh_assert_equal_to (this->_dof_map->first_dof(), 0);
00185       for (numeric_index_type i=this->_dof_map->first_dof();
00186            i!=this->_dof_map->end_dof(); ++i)
00187         {
00188           if(sparse)
00189             {
00190               for (numeric_index_type j=0; j<this->n(); j++)
00191                 {
00192                   T c = (*this)(i,j);
00193                   if (c != static_cast<T>(0.0))
00194                     {
00195                       os << i << " " << j << " " << c << std::endl;
00196                     }
00197                 }
00198             }
00199           else
00200             {
00201               for (numeric_index_type j=0; j<this->n(); j++)
00202                 os << (*this)(i,j) << " ";
00203               os << std::endl;
00204             }
00205         }
00206 
00207       std::vector<numeric_index_type> ibuf, jbuf;
00208       std::vector<T> cbuf;
00209       numeric_index_type currenti = this->_dof_map->end_dof();
00210       for (processor_id_type p=1; p < libMesh::n_processors(); ++p)
00211         {
00212           CommWorld.receive(p, ibuf);
00213           CommWorld.receive(p, jbuf);
00214           CommWorld.receive(p, cbuf);
00215           libmesh_assert_equal_to (ibuf.size(), jbuf.size());
00216           libmesh_assert_equal_to (ibuf.size(), cbuf.size());
00217 
00218           if (ibuf.empty())
00219             continue;
00220           libmesh_assert_greater_equal (ibuf.front(), currenti);
00221           libmesh_assert_greater_equal (ibuf.back(), ibuf.front());
00222 
00223           std::size_t currentb = 0;
00224           for (;currenti <= ibuf.back(); ++currenti)
00225             {
00226               if(sparse)
00227                 {
00228                   for (numeric_index_type j=0; j<this->n(); j++)
00229                     {
00230                       if (currentb < ibuf.size() &&
00231                           ibuf[currentb] == currenti &&
00232                           jbuf[currentb] == j)
00233                         {
00234                           os << currenti << " " << j << " " << cbuf[currentb] << std::endl;
00235                           currentb++;
00236                         }
00237                     }
00238                 }
00239               else
00240                 {
00241                   for (numeric_index_type j=0; j<this->n(); j++)
00242                     {
00243                       if (currentb < ibuf.size() &&
00244                           ibuf[currentb] == currenti &&
00245                           jbuf[currentb] == j)
00246                         {
00247                           os << cbuf[currentb] << " ";
00248                           currentb++;
00249                         }
00250                       else
00251                         os << static_cast<T>(0.0) << " ";
00252                     }
00253                   os << std::endl;
00254                 }
00255             }
00256         }
00257       if(!sparse)
00258         {
00259           for (; currenti != this->m(); ++currenti)
00260             {
00261               for (numeric_index_type j=0; j<this->n(); j++)
00262                 os << static_cast<T>(0.0) << " ";
00263               os << std::endl;
00264             }
00265         }
00266     }
00267   else
00268     {
00269       std::vector<numeric_index_type> ibuf, jbuf;
00270       std::vector<T> cbuf;
00271 
00272       // We'll assume each processor has access to entire
00273       // matrix rows, so (*this)(i,j) is valid if i is a local index.
00274       for (numeric_index_type i=this->_dof_map->first_dof();
00275            i!=this->_dof_map->end_dof(); ++i)
00276         {
00277           for (numeric_index_type j=0; j<this->n(); j++)
00278             {
00279               T c = (*this)(i,j);
00280               if (c != static_cast<T>(0.0))
00281                 {
00282                   ibuf.push_back(i);
00283                   jbuf.push_back(j);
00284                   cbuf.push_back(c);
00285                 }
00286             }
00287         }
00288       CommWorld.send(0,ibuf);
00289       CommWorld.send(0,jbuf);
00290       CommWorld.send(0,cbuf);
00291     }
00292 }

void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out  )  [static, inherited]

Prints the reference information, by default to libMesh::out.

Definition at line 88 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().

00089 {
00090   if( _enable_print_counter ) out_stream << ReferenceCounter::get_info();
00091 }

template<typename T>
virtual void libMesh::SparseMatrix< T >::print_matlab ( const std::string  name = "NULL"  )  const [inline, virtual]

Print the contents of the matrix in Matlab's sparse matrix format. Optionally prints the matrix to the file named name. If name is not specified it is dumped to the screen.

Reimplemented in libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

Definition at line 329 of file sparse_matrix.h.

00330   {
00331     libMesh::err << "ERROR: Not Implemented in base class yet!" << std::endl;
00332     libMesh::err << "ERROR writing MATLAB file " << name << std::endl;
00333     libmesh_error();
00334   }

template<typename T>
virtual void libMesh::SparseMatrix< T >::print_personal ( std::ostream &  os = libMesh::out  )  const [pure virtual]

Print the contents of the matrix to the screen in a package-personalized style, if available.

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

template<typename T>
virtual void libMesh::SparseMatrix< T >::reinit_submatrix ( SparseMatrix< T > &  submatrix,
const std::vector< numeric_index_type > &  rows,
const std::vector< numeric_index_type > &  cols 
) const [inline, virtual]

This function is similar to the one above, but it allows you to reuse the existing sparsity pattern of "submatrix" instead of reallocating it again. This should hopefully be more efficient if you are frequently extracting submatrices of the same size.

Definition at line 357 of file sparse_matrix.h.

00360   {
00361     this->_get_submatrix(submatrix,
00362                          rows,
00363                          cols,
00364                          true); // true means REUSE submatrix
00365   }

template<typename T>
virtual numeric_index_type libMesh::SparseMatrix< T >::row_start (  )  const [pure virtual]

return row_start, the index of the first matrix row stored on this processor

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

template<typename T>
virtual numeric_index_type libMesh::SparseMatrix< T >::row_stop (  )  const [pure virtual]

return row_stop, the index of the last matrix row (+1) stored on this processor

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

template<typename T>
virtual void libMesh::SparseMatrix< T >::set ( const numeric_index_type  i,
const numeric_index_type  j,
const T  value 
) [pure virtual]

Set the element (i,j) to value. Throws an error if the entry does not exist. Still, it is allowed to store zero values in non-existent fields.

Implemented in libMesh::LaspackMatrix< T >, libMesh::PetscMatrix< T >, and libMesh::EpetraMatrix< T >.

template<typename T>
virtual void libMesh::SparseMatrix< T >::update_sparsity_pattern ( const SparsityPattern::Graph  )  [inline, virtual]

Updates the matrix sparsity pattern. When your SparseMatrix<T> implementation does not need this data simply do not overload this method.

Reimplemented in libMesh::LaspackMatrix< T >, and libMesh::EpetraMatrix< T >.

Definition at line 125 of file sparse_matrix.h.

Referenced by libMesh::DofMap::attach_matrix().

00125 {}

template<typename T>
void libMesh::SparseMatrix< T >::vector_mult ( NumericVector< T > &  dest,
const NumericVector< T > &  arg 
) const [inline]

Multiplies the matrix with arg and stores the result in dest.

Definition at line 138 of file sparse_matrix.C.

References libMesh::SparseMatrix< T >::vector_mult_add(), and libMesh::NumericVector< T >::zero().

00140 {
00141   dest.zero();
00142   this->vector_mult_add(dest,arg);
00143 }

template<typename T>
void libMesh::SparseMatrix< T >::vector_mult_add ( NumericVector< T > &  dest,
const NumericVector< T > &  arg 
) const [inline]

Multiplies the matrix with arg and adds the result to dest.

Definition at line 148 of file sparse_matrix.C.

References libMesh::NumericVector< T >::add_vector().

Referenced by libMesh::SparseMatrix< T >::vector_mult(), and libMesh::ImplicitSystem::weighted_sensitivity_adjoint_solve().

00150 {
00151   /* This functionality is actually implemented in the \p
00152      NumericVector class.  */
00153   dest.add_vector(arg,*this);
00154 }

template<typename T>
void libMesh::SparseMatrix< T >::zero_rows ( std::vector< numeric_index_type > &  rows,
diag_value = 0.0 
) [inline, virtual]

Set all row entries to 0 then puts diag_value in the diagonal entry

Reimplemented in libMesh::PetscMatrix< T >.

Definition at line 159 of file sparse_matrix.C.

00160 {
00161   /* This functionality isn't implemented or stubbed in every subclass yet */
00162   libmesh_not_implemented();
00163 }


Friends And Related Function Documentation

template<typename T>
std::ostream& operator<< ( std::ostream &  os,
const SparseMatrix< T > &  m 
) [friend]

Same as the print method above, but allows you to print to a stream in the standard syntax.

template <typename u>=""> friend std::ostream& operator << (std::ostream& os, const SparseMatrix<U>& m);

Obscure C++ note 1: the above syntax, which does not require any prior declaration of operator<<, declares *any* instantiation of SparseMatrix<X> is friend to *any* instantiation of operator<<(ostream&, SparseMatrix<Y>&). It would not happen in practice, but in principle it means that SparseMatrix<Complex> would be friend to operator<<(ostream&, SparseMatrix<Real>).

Obscure C++ note 2: The form below, which requires a previous declaration of the operator<<(stream&, SparseMatrix<T>&) function (see top of this file), means that any instantiation of SparseMatrix<T> is friend to the specialization operator<<(ostream&, SparseMatrix<T>&), but e.g. SparseMatrix<U> is *not* friend to the same function. So this is slightly different to the form above...

This method seems to be the "preferred" technique, see http://www.parashift.com/c++-faq-lite/template-friends.html


Member Data Documentation

bool libMesh::ReferenceCounter::_enable_print_counter = true [static, protected, inherited]

Flag to control whether reference count information is printed when print_info is called.

Definition at line 137 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 131 of file reference_counter.h.

Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects [static, protected, inherited]

The number of objects. Print the reference count information when the number returns to 0.

Definition at line 126 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().


The documentation for this class was generated from the following files:

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:39 UTC

Hosted By:
SourceForge.net Logo