libMesh::LaspackMatrix< T > Class Template Reference

#include <laspack_matrix.h>

Inheritance diagram for libMesh::LaspackMatrix< T >:

List of all members.

Public Member Functions

 LaspackMatrix ()
 ~LaspackMatrix ()
bool need_full_sparsity_pattern () const
void update_sparsity_pattern (const SparsityPattern::Graph &)
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)
void init ()
void clear ()
void zero ()
void close () const
numeric_index_type m () const
numeric_index_type n () const
numeric_index_type row_start () const
numeric_index_type row_stop () const
void set (const numeric_index_type i, const numeric_index_type j, const T value)
void add (const numeric_index_type i, const numeric_index_type j, const T value)
void add_matrix (const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)
void add_matrix (const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &dof_indices)
void add (const T a, SparseMatrix< T > &X)
operator() (const numeric_index_type i, const numeric_index_type j) const
Real l1_norm () const
Real linfty_norm () const
bool closed () const
void print_personal (std::ostream &os=libMesh::out) const
virtual void get_diagonal (NumericVector< T > &dest) const
virtual void get_transpose (SparseMatrix< T > &dest) const
virtual bool initialized () const
void attach_dof_map (const DofMap &dof_map)
virtual void zero_rows (std::vector< numeric_index_type > &rows, T diag_value=0.0)
virtual void add (const T, SparseMatrix< T > &)=0
void print (std::ostream &os=libMesh::out, const bool sparse=false) const
template<>
void print (std::ostream &os, const bool sparse) const
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_transpose (SparseMatrix< T > &dest) const =0

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

Private Member Functions

numeric_index_type pos (const numeric_index_type i, const numeric_index_type j) const

Private Attributes

QMatrix _QMat
std::vector< numeric_index_type_csr
std::vector< std::vector
< numeric_index_type >
::const_iterator > 
_row_start
bool _closed

Friends

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

Detailed Description

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

The LaspackMatrix class wraps a QMatrix object from the Laspack library. Currently Laspack only supports real datatypes, so this class is a full specialization of SparseMatrix<T> with T = Real.

Author:
Benjamin S. Kirk, 2003

Definition at line 59 of file laspack_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::LaspackMatrix< T >::LaspackMatrix (  )  [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 238 of file laspack_matrix.C.

00238                                  :
00239   _closed (false)
00240 {
00241 }

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

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

Definition at line 246 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::clear().

00247 {
00248   this->clear ();
00249 }


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, inherited]

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, inherited]

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

template<typename T >
void libMesh::LaspackMatrix< T >::add ( const T  a,
SparseMatrix< T > &  X 
) [inline]

Add a Sparse matrix X, scaled with a, to this, stores the result in this: $\texttt{this} += a*X $. LASPACK does not provide a true axpy for matrices, so a hand-coded version with hopefully acceptable performance is provided.

Definition at line 383 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_QMat, libMesh::LaspackMatrix< T >::_row_start, libMesh::SparseMatrix< T >::initialized(), libMesh::SparseMatrix< T >::m(), libMesh::LaspackMatrix< T >::m(), libMesh::SparseMatrix< T >::n(), and libMesh::LaspackMatrix< T >::n().

00384 {
00385   libmesh_assert (this->initialized());
00386   libmesh_assert_equal_to (this->m(), X_in.m());
00387   libmesh_assert_equal_to (this->n(), X_in.n());
00388 
00389   LaspackMatrix<T>* X = libmesh_cast_ptr<LaspackMatrix<T>*> (&X_in);
00390   _LPNumber         a = static_cast<_LPNumber>          (a_in);
00391 
00392   libmesh_assert(X);
00393 
00394   // loops taken from LaspackMatrix<T>::zero ()
00395 
00396   const numeric_index_type n_rows = this->m();
00397 
00398   for (numeric_index_type row=0; row<n_rows; row++)
00399     {
00400       const std::vector<numeric_index_type>::const_iterator
00401         r_start = _row_start[row];
00402 
00403       const numeric_index_type len = (_row_start[row+1] - _row_start[row]);
00404 
00405       // Make sure we agree on the row length
00406       libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1));
00407       // compare matrix sparsity structures
00408       libmesh_assert_equal_to (len, Q_GetLen(&(X->_QMat), row+1));
00409 
00410 
00411       for (numeric_index_type l=0; l<len; l++)
00412         {
00413           const numeric_index_type j = *(r_start + l);
00414 
00415           // Make sure the data structures are working
00416           libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l));
00417 
00418           const _LPNumber value = a * Q_GetEl(const_cast<QMatrix*>(&(X->_QMat)), row+1, j+1);
00419           Q_AddVal   (&_QMat, row+1, l, value);
00420         }
00421     }
00422 }

template<typename T >
void libMesh::LaspackMatrix< T >::add ( const numeric_index_type  i,
const numeric_index_type  j,
const T  value 
) [inline, 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.

Implements libMesh::SparseMatrix< T >.

Definition at line 355 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_QMat, libMesh::LaspackMatrix< T >::_row_start, libMesh::SparseMatrix< T >::initialized(), libMesh::LaspackMatrix< T >::m(), libMesh::LaspackMatrix< T >::n(), and libMesh::LaspackMatrix< T >::pos().

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

00358 {
00359   libmesh_assert (this->initialized());
00360   libmesh_assert_less (i, this->m());
00361   libmesh_assert_less (j, this->n());
00362 
00363   const numeric_index_type position = this->pos(i,j);
00364 
00365   // Sanity check
00366   libmesh_assert_equal_to (*(_row_start[i]+position), j);
00367 
00368   Q_AddVal (&_QMat, i+1, position, value);
00369 }

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

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

Implements libMesh::SparseMatrix< T >.

Definition at line 374 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::add_matrix().

00376 {
00377   this->add_matrix (dm, dof_indices, dof_indices);
00378 }

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

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

Implements libMesh::SparseMatrix< T >.

Definition at line 202 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::add(), libMesh::SparseMatrix< T >::initialized(), libMesh::DenseMatrixBase< T >::m(), and libMesh::DenseMatrixBase< T >::n().

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

00206 {
00207   libmesh_assert (this->initialized());
00208   unsigned int n_rows = libmesh_cast_int<unsigned int>(rows.size());
00209   unsigned int n_cols = libmesh_cast_int<unsigned int>(cols.size());
00210   libmesh_assert_equal_to (dm.m(), n_rows);
00211   libmesh_assert_equal_to (dm.n(), n_cols);
00212 
00213 
00214   for (unsigned int i=0; i<n_rows; i++)
00215     for (unsigned int j=0; j<n_cols; j++)
00216       this->add(rows[i],cols[j],dm(i,j));
00217 }

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

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, inherited]

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 >
void libMesh::LaspackMatrix< T >::clear (  )  [inline, virtual]

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

Implements libMesh::SparseMatrix< T >.

Definition at line 254 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_closed, libMesh::LaspackMatrix< T >::_csr, libMesh::SparseMatrix< T >::_is_initialized, libMesh::LaspackMatrix< T >::_QMat, libMesh::LaspackMatrix< T >::_row_start, and libMesh::SparseMatrix< T >::initialized().

Referenced by libMesh::LaspackMatrix< T >::init(), libMesh::LaspackMatrix< T >::update_sparsity_pattern(), and libMesh::LaspackMatrix< T >::~LaspackMatrix().

00255 {
00256   if (this->initialized())
00257     {
00258       Q_Destr(&_QMat);
00259     }
00260 
00261   _csr.clear();
00262   _row_start.clear();
00263   _closed = false;
00264   this->_is_initialized = false;
00265 }

template<typename T>
void libMesh::LaspackMatrix< T >::close (  )  const [inline, virtual]

Close the matrix. Dummy routine. After calling this method closed() is true and the matrix can be used in computations.

Implements libMesh::SparseMatrix< T >.

Definition at line 138 of file laspack_matrix.h.

References libMesh::LaspackMatrix< T >::_closed.

Referenced by libMesh::LaspackLinearSolver< T >::adjoint_solve(), and libMesh::LaspackLinearSolver< T >::solve().

00138 { const_cast<LaspackMatrix<T>*>(this)->_closed = true; }

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

see if Laspack matrix has been closed and fully assembled yet

Implements libMesh::SparseMatrix< T >.

Definition at line 252 of file laspack_matrix.h.

References libMesh::LaspackMatrix< T >::_closed.

00252 { return _closed; }

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, inherited]

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 >
void libMesh::LaspackMatrix< T >::get_diagonal ( NumericVector< T > &  dest  )  const [inline, virtual]

Copies the diagonal part of the matrix into dest.

Implements libMesh::SparseMatrix< T >.

Definition at line 222 of file laspack_matrix.C.

00223 {
00224   libmesh_not_implemented();
00225 }

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, inherited]

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().

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

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

Definition at line 230 of file laspack_matrix.C.

00231 {
00232   libmesh_not_implemented();
00233 }

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 >
void libMesh::LaspackMatrix< T >::init (  )  [inline, virtual]

Initialize using sparsity structure computed by dof_map.

Implements libMesh::SparseMatrix< T >.

Definition at line 151 of file laspack_matrix.C.

References libMesh::SparseMatrix< T >::_dof_map, libMesh::SparseMatrix< T >::_is_initialized, libMesh::LaspackMatrix< T >::_QMat, libMesh::LaspackMatrix< T >::clear(), libMesh::DofMap::get_n_nz(), libMesh::DofMap::get_n_oz(), libMesh::SparseMatrix< T >::initialized(), libMesh::LaspackMatrix< T >::m(), libMesh::DofMap::n_dofs(), and libMesh::DofMap::n_dofs_on_processor().

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

00152 {
00153   // Ignore calls on initialized objects
00154   if (this->initialized())
00155     return;
00156 
00157   // We need the DofMap for this!
00158   libmesh_assert(this->_dof_map);
00159 
00160   // Clear intialized matrices
00161   if (this->initialized())
00162     this->clear();
00163 
00164   const numeric_index_type n_rows   = this->_dof_map->n_dofs();
00165 #ifndef NDEBUG
00166   // The following variables are only used for assertions,
00167   // so avoid declaring them when asserts are inactive.
00168   const numeric_index_type n_cols   = n_rows;
00169   const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(0);
00170   const numeric_index_type m_l = n_l;
00171 #endif
00172 
00173   // Laspack Matrices only work for uniprocessor cases
00174   libmesh_assert_equal_to (n_rows, n_cols);
00175   libmesh_assert_equal_to (m_l, n_rows);
00176   libmesh_assert_equal_to (n_l, n_cols);
00177 
00178 #ifndef NDEBUG
00179   // The following variables are only used for assertions,
00180   // so avoid declaring them when asserts are inactive.
00181   const std::vector<numeric_index_type>& n_nz = this->_dof_map->get_n_nz();
00182   const std::vector<numeric_index_type>& n_oz = this->_dof_map->get_n_oz();
00183 #endif
00184 
00185   // Make sure the sparsity pattern isn't empty
00186   libmesh_assert_equal_to (n_nz.size(), n_l);
00187   libmesh_assert_equal_to (n_oz.size(), n_l);
00188 
00189   if (n_rows==0)
00190     return;
00191 
00192   Q_Constr(&_QMat, const_cast<char*>("Mat"), n_rows, _LPFalse, Rowws, Normal, _LPTrue);
00193 
00194   this->_is_initialized = true;
00195 
00196   libmesh_assert_equal_to (n_rows, this->m());
00197 }

template<typename T>
void libMesh::LaspackMatrix< 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 
) [virtual]

Initialize a Laspack 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 30).

Implements libMesh::SparseMatrix< T >.

template<typename T>
virtual bool libMesh::SparseMatrix< T >::initialized (  )  const [inline, virtual, inherited]
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>
Real libMesh::LaspackMatrix< T >::l1_norm (  )  const [inline, 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$.

Implements libMesh::SparseMatrix< T >.

Definition at line 233 of file laspack_matrix.h.

00233 { libmesh_error(); return 0.; }

template<typename T>
Real libMesh::LaspackMatrix< T >::linfty_norm (  )  const [inline, 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$.

Implements libMesh::SparseMatrix< T >.

Definition at line 246 of file laspack_matrix.h.

00246 { libmesh_error(); return 0.; }

template<typename T >
numeric_index_type libMesh::LaspackMatrix< T >::m (  )  const [inline, virtual]
template<typename T >
numeric_index_type libMesh::LaspackMatrix< T >::n (  )  const [inline, virtual]
Returns:
n, the column-dimension of the matrix where the marix is $ M \times N $.

Implements libMesh::SparseMatrix< T >.

Definition at line 309 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_QMat, and libMesh::SparseMatrix< T >::initialized().

Referenced by libMesh::LaspackMatrix< T >::add(), libMesh::LaspackMatrix< T >::operator()(), libMesh::LaspackMatrix< T >::pos(), and libMesh::LaspackMatrix< T >::set().

00310 {
00311   libmesh_assert (this->initialized());
00312 
00313   return static_cast<numeric_index_type>(Q_GetDim(const_cast<QMatrix*>(&_QMat)));
00314 }

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>
bool libMesh::LaspackMatrix< T >::need_full_sparsity_pattern (  )  const [inline, virtual]

The LaspackMatrix needs the full sparsity pattern.

Reimplemented from libMesh::SparseMatrix< T >.

Definition at line 90 of file laspack_matrix.h.

00091   { return true; }

template<typename T >
T libMesh::LaspackMatrix< T >::operator() ( const numeric_index_type  i,
const numeric_index_type  j 
) const [inline, 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.

Implements libMesh::SparseMatrix< T >.

Definition at line 428 of file laspack_matrix.C.

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

00430 {
00431   libmesh_assert (this->initialized());
00432   libmesh_assert_less (i, this->m());
00433   libmesh_assert_less (j, this->n());
00434 
00435   return Q_GetEl (const_cast<QMatrix*>(&_QMat), i+1, j+1);
00436 }

template<typename T >
numeric_index_type libMesh::LaspackMatrix< T >::pos ( const numeric_index_type  i,
const numeric_index_type  j 
) const [inline, private]
Returns:
the position in the compressed row storage scheme of the $ (i,j) $ element.

Definition at line 441 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_csr, libMesh::LaspackMatrix< T >::_row_start, libMesh::LaspackMatrix< T >::m(), and libMesh::LaspackMatrix< T >::n().

Referenced by libMesh::LaspackMatrix< T >::add(), libMesh::LaspackMatrix< T >::set(), and libMesh::LaspackMatrix< T >::update_sparsity_pattern().

00443 {
00444   libmesh_assert_less (i, this->m());
00445   libmesh_assert_less (j, this->n());
00446   libmesh_assert_less (i+1, _row_start.size());
00447   libmesh_assert (_row_start.back() == _csr.end());
00448 
00449   // note this requires the _csr to be
00450   std::pair<std::vector<numeric_index_type>::const_iterator,
00451             std::vector<numeric_index_type>::const_iterator> p =
00452     std::equal_range (_row_start[i],
00453                       _row_start[i+1],
00454                       j);
00455 
00456   // Make sure the row contains the element j
00457   libmesh_assert (p.first != p.second);
00458 
00459   // Make sure the values match
00460   libmesh_assert (*p.first == j);
00461 
00462   // Return the position in the compressed row
00463   return std::distance (_row_start[i], p.first);
00464 }

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

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, inherited]

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, inherited]

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>
void libMesh::LaspackMatrix< T >::print_personal ( std::ostream &  os = libMesh::out  )  const [inline, virtual]

Print the contents of the matrix, by default to libMesh::out. Currently identical to print().

Implements libMesh::SparseMatrix< T >.

Definition at line 258 of file laspack_matrix.h.

References libMesh::SparseMatrix< T >::print().

00258 { this->print(os); }

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, inherited]

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 >
numeric_index_type libMesh::LaspackMatrix< T >::row_start (  )  const [inline, virtual]

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

Implements libMesh::SparseMatrix< T >.

Definition at line 319 of file laspack_matrix.C.

00320 {
00321   return 0;
00322 }

template<typename T >
numeric_index_type libMesh::LaspackMatrix< T >::row_stop (  )  const [inline, virtual]

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

Implements libMesh::SparseMatrix< T >.

Definition at line 327 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::m().

00328 {
00329   return this->m();
00330 }

template<typename T >
void libMesh::LaspackMatrix< T >::set ( const numeric_index_type  i,
const numeric_index_type  j,
const T  value 
) [inline, 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.

Implements libMesh::SparseMatrix< T >.

Definition at line 335 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_QMat, libMesh::LaspackMatrix< T >::_row_start, libMesh::SparseMatrix< T >::initialized(), libMesh::LaspackMatrix< T >::m(), libMesh::LaspackMatrix< T >::n(), and libMesh::LaspackMatrix< T >::pos().

00338 {
00339   libmesh_assert (this->initialized());
00340   libmesh_assert_less (i, this->m());
00341   libmesh_assert_less (j, this->n());
00342 
00343   const numeric_index_type position = this->pos(i,j);
00344 
00345   // Sanity check
00346   libmesh_assert_equal_to (*(_row_start[i]+position), j);
00347   libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, i+1, position));
00348 
00349   Q_SetEntry (&_QMat, i+1, position, j+1, value);
00350 }

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

Updates the matrix sparsity pattern. This will tell the underlying matrix storage scheme how to map the $ (i,j) $ elements.

Reimplemented from libMesh::SparseMatrix< T >.

Definition at line 39 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_csr, libMesh::SparseMatrix< T >::_dof_map, libMesh::LaspackMatrix< T >::_QMat, libMesh::LaspackMatrix< T >::_row_start, libMesh::LaspackMatrix< T >::clear(), libMesh::LaspackMatrix< T >::init(), libMesh::SparseMatrix< T >::initialized(), libMesh::LaspackMatrix< T >::m(), and libMesh::LaspackMatrix< T >::pos().

00040 {
00041   // clear data, start over
00042   this->clear ();
00043 
00044   // big trouble if this fails!
00045   libmesh_assert(this->_dof_map);
00046 
00047   const numeric_index_type n_rows = sparsity_pattern.size();
00048 
00049   // Initialize the _row_start data structure,
00050   // allocate storage for the _csr array
00051   {
00052     std::size_t size = 0;
00053 
00054     for (numeric_index_type row=0; row<n_rows; row++)
00055       size += sparsity_pattern[row].size();
00056 
00057     _csr.resize       (size);
00058     _row_start.reserve(n_rows + 1);
00059   }
00060 
00061 
00062   // Initize the _csr data structure.
00063   {
00064     std::vector<numeric_index_type>::iterator position = _csr.begin();
00065 
00066     _row_start.push_back (position);
00067 
00068     for (numeric_index_type row=0; row<n_rows; row++)
00069       {
00070         // insert the row indices
00071         for (SparsityPattern::Row::const_iterator col = sparsity_pattern[row].begin();
00072              col != sparsity_pattern[row].end(); ++col)
00073           {
00074             libmesh_assert (position != _csr.end());
00075             *position = *col;
00076             ++position;
00077           }
00078 
00079         _row_start.push_back (position);
00080       }
00081   }
00082 
00083 
00084   // Initialize the matrix
00085   libmesh_assert (!this->initialized());
00086   this->init ();
00087   libmesh_assert (this->initialized());
00088   //libMesh::out << "n_rows=" << n_rows << std::endl;
00089   //libMesh::out << "m()=" << m() << std::endl;
00090   libmesh_assert_equal_to (n_rows, this->m());
00091 
00092   // Tell the matrix about its structure.  Initialize it
00093   // to zero.
00094   for (numeric_index_type i=0; i<n_rows; i++)
00095     {
00096       const std::vector<numeric_index_type>::const_iterator
00097         rs = _row_start[i];
00098 
00099       const numeric_index_type length = _row_start[i+1] - rs;
00100 
00101       Q_SetLen (&_QMat, i+1, length);
00102 
00103       for (numeric_index_type l=0; l<length; l++)
00104         {
00105           const numeric_index_type j = *(rs+l);
00106 
00107           // sanity check
00108           //libMesh::out << "m()=" << m() << std::endl;
00109           //libMesh::out << "(i,j,l) = (" << i
00110           //              << "," << j
00111           //              << "," << l
00112           //              << ")" << std::endl;
00113           //libMesh::out << "pos(i,j)=" << pos(i,j)
00114           //              << std::endl;
00115           libmesh_assert_equal_to (this->pos(i,j), l);
00116           Q_SetEntry (&_QMat, i+1, l, j+1, 0.);
00117         }
00118     }
00119 
00120   // That's it!
00121   //libmesh_here();
00122 }

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

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, inherited]

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::LaspackMatrix< T >::zero (  )  [inline, virtual]

Set all entries to 0.

Implements libMesh::SparseMatrix< T >.

Definition at line 270 of file laspack_matrix.C.

References libMesh::LaspackMatrix< T >::_QMat, libMesh::LaspackMatrix< T >::_row_start, and libMesh::LaspackMatrix< T >::m().

00271 {
00272   const numeric_index_type n_rows = this->m();
00273 
00274   for (numeric_index_type row=0; row<n_rows; row++)
00275     {
00276       const std::vector<numeric_index_type>::const_iterator
00277         r_start = _row_start[row];
00278 
00279       const numeric_index_type len = (_row_start[row+1] - _row_start[row]);
00280 
00281       // Make sure we agree on the row length
00282       libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1));
00283 
00284       for (numeric_index_type l=0; l<len; l++)
00285         {
00286           const numeric_index_type j = *(r_start + l);
00287 
00288           // Make sure the data structures are working
00289           libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l));
00290 
00291           Q_SetEntry (&_QMat, row+1, l, j+1, 0.);
00292         }
00293     }
00294 }

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

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>
friend class LaspackLinearSolver< T > [friend]

Definition at line 305 of file laspack_matrix.h.

template<typename T>
friend class LaspackVector< T > [friend]

Make other Laspack datatypes friends

Definition at line 304 of file laspack_matrix.h.

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

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

template<typename T>
bool libMesh::LaspackMatrix< T >::_closed [private]

Flag indicating if the matrix has been closed yet.

Definition at line 299 of file laspack_matrix.h.

Referenced by libMesh::LaspackMatrix< T >::clear(), libMesh::LaspackMatrix< T >::close(), and libMesh::LaspackMatrix< T >::closed().

template<typename T>
std::vector<numeric_index_type> libMesh::LaspackMatrix< T >::_csr [private]
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().

template<typename T>
std::vector<std::vector<numeric_index_type>::const_iterator> libMesh::LaspackMatrix< T >::_row_start [private]

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

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

Hosted By:
SourceForge.net Logo