libMesh::DenseSubMatrix< T > Class Template Reference

#include <dense_submatrix.h>

Inheritance diagram for libMesh::DenseSubMatrix< T >:

Public Member Functions

 DenseSubMatrix (DenseMatrix< T > &new_parent, const unsigned int ioff=0, const unsigned int joff=0, const unsigned int m=0, const unsigned int n=0)
 
 DenseSubMatrix (DenseSubMatrix &&)=default
 
 DenseSubMatrix (const DenseSubMatrix &)=default
 
DenseSubMatrixoperator= (const DenseSubMatrix &)=default
 
DenseSubMatrixoperator= (DenseSubMatrix &&)=default
 
virtual ~DenseSubMatrix ()=default
 
DenseMatrix< T > & parent ()
 
virtual void zero () override
 
operator() (const unsigned int i, const unsigned int j) const
 
T & operator() (const unsigned int i, const unsigned int j)
 
virtual T el (const unsigned int i, const unsigned int j) const override
 
virtual T & el (const unsigned int i, const unsigned int j) override
 
virtual void left_multiply (const DenseMatrixBase< T > &M2) override
 
virtual void right_multiply (const DenseMatrixBase< T > &M3) override
 
void reposition (const unsigned int ioff, const unsigned int joff, const unsigned int new_m, const unsigned int new_n)
 
unsigned int i_off () const
 
unsigned int j_off () const
 
void condense (const unsigned int i, const unsigned int j, const T val, DenseSubVector< T > &rhs)
 
unsigned int m () const
 
unsigned int n () const
 
void print (std::ostream &os=libMesh::out) const
 
void print_scientific (std::ostream &os, unsigned precision=8) const
 
template<typename T2 , typename T3 >
boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type add (const T2 factor, const DenseMatrixBase< T3 > &mat)
 

Protected Member Functions

void condense (const unsigned int i, const unsigned int j, const T val, DenseVectorBase< T > &rhs)
 

Static Protected Member Functions

static void multiply (DenseMatrixBase< T > &M1, const DenseMatrixBase< T > &M2, const DenseMatrixBase< T > &M3)
 

Protected Attributes

unsigned int _m
 
unsigned int _n
 

Private Attributes

DenseMatrix< T > & _parent_matrix
 
unsigned int _i_off
 
unsigned int _j_off
 

Detailed Description

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

Defines a dense submatrix for use in Finite Element-type computations. Useful for storing element stiffness matrices before summation into a global matrix, particularly when you have systems of equations. All overridden virtual functions are documented in dense_matrix_base.h.

Author
Benjamin S. Kirk
Date
2003

Definition at line 44 of file dense_submatrix.h.

Constructor & Destructor Documentation

◆ DenseSubMatrix() [1/3]

template<typename T >
libMesh::DenseSubMatrix< T >::DenseSubMatrix ( DenseMatrix< T > &  new_parent,
const unsigned int  ioff = 0,
const unsigned int  joff = 0,
const unsigned int  m = 0,
const unsigned int  n = 0 
)
inline

Constructor. Creates a dense submatrix of the matrix parent. The submatrix has dimensions $(m \times n)$, and the $(0,0)$ entry of the submatrix is located at the $(ioff,joff)$ location in the parent matrix.

Definition at line 158 of file dense_submatrix.h.

References libMesh::DenseSubMatrix< T >::reposition().

162  :
163  DenseMatrixBase<T>(new_m,new_n),
164  _parent_matrix(new_parent)
165 {
166  this->reposition (ioff, joff, new_m, new_n);
167 }
void reposition(const unsigned int ioff, const unsigned int joff, const unsigned int new_m, const unsigned int new_n)
DenseMatrix< T > & _parent_matrix

◆ DenseSubMatrix() [2/3]

template<typename T>
libMesh::DenseSubMatrix< T >::DenseSubMatrix ( DenseSubMatrix< T > &&  )
default

The 5 special functions can be defaulted for this class, as it does not manage any memory itself.

◆ DenseSubMatrix() [3/3]

template<typename T>
libMesh::DenseSubMatrix< T >::DenseSubMatrix ( const DenseSubMatrix< T > &  )
default

◆ ~DenseSubMatrix()

template<typename T>
virtual libMesh::DenseSubMatrix< T >::~DenseSubMatrix ( )
virtualdefault

Member Function Documentation

◆ add()

template<typename T >
template<typename T2 , typename T3 >
boostcopy::enable_if_c< ScalarTraits< T2 >::value, void >::type libMesh::DenseMatrixBase< T >::add ( const T2  factor,
const DenseMatrixBase< T3 > &  mat 
)
inlineinherited

Adds factor to every element in the matrix. This should only work if T += T2 * T3 is valid C++ and if T2 is scalar. Return type is void

Definition at line 188 of file dense_matrix_base.h.

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

190 {
191  libmesh_assert_equal_to (this->m(), mat.m());
192  libmesh_assert_equal_to (this->n(), mat.n());
193 
194  for (unsigned int j=0; j<this->n(); j++)
195  for (unsigned int i=0; i<this->m(); i++)
196  this->el(i,j) += factor*mat.el(i,j);
197 }
unsigned int m() const
virtual T el(const unsigned int i, const unsigned int j) const =0
unsigned int n() const

◆ condense() [1/2]

template<typename T>
void libMesh::DenseSubMatrix< T >::condense ( const unsigned int  i,
const unsigned int  j,
const T  val,
DenseSubVector< T > &  rhs 
)
inline

Condense-out the (i,j) entry of the matrix, forcing it to take on the value val. This is useful in numerical simulations for applying boundary conditions. Preserves the symmetry of the matrix.

Definition at line 125 of file dense_submatrix.h.

References libMesh::DenseSubMatrix< T >::i_off(), libMesh::DenseSubMatrix< T >::j_off(), libMesh::DenseSubVector< T >::parent(), and libMesh::DenseSubMatrix< T >::parent().

129  {
130  this->parent().condense(this->i_off()+i,
131  this->j_off()+j,
132  val, rhs.parent());
133  }
unsigned int j_off() const
DenseMatrix< T > & parent()
unsigned int i_off() const

◆ condense() [2/2]

template<typename T >
void libMesh::DenseMatrixBase< T >::condense ( const unsigned int  i,
const unsigned int  j,
const T  val,
DenseVectorBase< T > &  rhs 
)
protectedinherited

Condense-out the (i,j) entry of the matrix, forcing it to take on the value val. This is useful in numerical simulations for applying boundary conditions. Preserves the symmetry of the matrix.

Definition at line 58 of file dense_matrix_base.C.

References libMesh::DenseVectorBase< T >::el(), and libMesh::DenseVectorBase< T >::size().

Referenced by libMesh::DenseMatrix< Number >::condense().

62 {
63  libmesh_assert_equal_to (this->_m, rhs.size());
64  libmesh_assert_equal_to (iv, jv);
65 
66 
67  // move the known value into the RHS
68  // and zero the column
69  for (unsigned int i=0; i<this->m(); i++)
70  {
71  rhs.el(i) -= this->el(i,jv)*val;
72  this->el(i,jv) = 0.;
73  }
74 
75  // zero the row
76  for (unsigned int j=0; j<this->n(); j++)
77  this->el(iv,j) = 0.;
78 
79  this->el(iv,jv) = 1.;
80  rhs.el(iv) = val;
81 
82 }
unsigned int m() const
virtual T el(const unsigned int i, const unsigned int j) const =0
unsigned int n() const

◆ el() [1/2]

template<typename T>
virtual T libMesh::DenseSubMatrix< T >::el ( const unsigned int  i,
const unsigned int  j 
) const
inlineoverridevirtual
Returns
The (i,j) element of the matrix. Since internal data representations may differ, you must redefine this function.

Implements libMesh::DenseMatrixBase< T >.

Definition at line 89 of file dense_submatrix.h.

91  { return (*this)(i,j); }

◆ el() [2/2]

template<typename T>
virtual T& libMesh::DenseSubMatrix< T >::el ( const unsigned int  i,
const unsigned int  j 
)
inlineoverridevirtual
Returns
The (i,j) element of the matrix as a writable reference. Since internal data representations may differ, you must redefine this function.

Implements libMesh::DenseMatrixBase< T >.

Definition at line 93 of file dense_submatrix.h.

95  { return (*this)(i,j); }

◆ i_off()

template<typename T>
unsigned int libMesh::DenseSubMatrix< T >::i_off ( ) const
inline
Returns
The row offset into the parent matrix.

Definition at line 112 of file dense_submatrix.h.

References libMesh::DenseSubMatrix< T >::_i_off.

Referenced by libMesh::DenseSubMatrix< T >::condense().

112 { return _i_off; }

◆ j_off()

template<typename T>
unsigned int libMesh::DenseSubMatrix< T >::j_off ( ) const
inline
Returns
The column offset into the parent matrix.

Definition at line 117 of file dense_submatrix.h.

References libMesh::DenseSubMatrix< T >::_j_off.

Referenced by libMesh::DenseSubMatrix< T >::condense().

117 { return _j_off; }

◆ left_multiply()

template<typename T >
void libMesh::DenseSubMatrix< T >::left_multiply ( const DenseMatrixBase< T > &  M2)
overridevirtual

Performs the operation: (*this) <- M2 * (*this)

Implements libMesh::DenseMatrixBase< T >.

Definition at line 31 of file dense_submatrix.C.

32 {
33  // (*this) <- M2 * M3
34  // Where:
35  // (*this) = (m x n),
36  // M2 = (m x p),
37  // M3 = (p x n)
38 
39  // M3 is a simply a copy of *this
40  DenseSubMatrix<T> M3(*this);
41 
42  // Call the multiply function in the base class
43  this->multiply(*this, M2, M3);
44 }
static void multiply(DenseMatrixBase< T > &M1, const DenseMatrixBase< T > &M2, const DenseMatrixBase< T > &M3)

◆ m()

template<typename T>
unsigned int libMesh::DenseMatrixBase< T >::m ( ) const
inlineinherited
Returns
The row-dimension of the matrix.

Definition at line 102 of file dense_matrix_base.h.

References libMesh::DenseMatrixBase< T >::_m.

Referenced by libMesh::DenseMatrix< Number >::_multiply_blas(), libMesh::DenseMatrix< Number >::_svd_solve_lapack(), libMesh::DenseMatrixBase< T >::add(), libMesh::DenseMatrix< Number >::add(), libMesh::PetscMatrix< T >::add_block_matrix(), libMesh::SparseMatrix< ValOut >::add_block_matrix(), libMesh::EigenSparseMatrix< T >::add_matrix(), libMesh::LaspackMatrix< T >::add_matrix(), libMesh::EpetraMatrix< T >::add_matrix(), libMesh::PetscMatrix< T >::add_matrix(), libMesh::DofMap::build_constraint_matrix(), libMesh::DofMap::build_constraint_matrix_and_vector(), libMesh::DofMap::constrain_element_dyad_matrix(), libMesh::DofMap::constrain_element_matrix(), libMesh::DofMap::constrain_element_matrix_and_vector(), libMesh::DofMap::constrain_element_vector(), libMesh::DofMap::extract_local_vector(), libMesh::DenseMatrix< Number >::get_transpose(), libMesh::DofMap::heterogenously_constrain_element_matrix_and_vector(), libMesh::DofMap::heterogenously_constrain_element_vector(), libMesh::DenseMatrix< Number >::left_multiply(), libMesh::DenseMatrix< Number >::left_multiply_transpose(), libMesh::DofMap::max_constraint_error(), libMesh::DenseMatrixBase< T >::multiply(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::DenseMatrix< Number >::operator=(), libMesh::DenseMatrix< Number >::right_multiply(), and libMesh::DenseMatrix< Number >::right_multiply_transpose().

102 { return _m; }

◆ multiply()

template<typename T >
void libMesh::DenseMatrixBase< T >::multiply ( DenseMatrixBase< T > &  M1,
const DenseMatrixBase< T > &  M2,
const DenseMatrixBase< T > &  M3 
)
staticprotectedinherited

Helper function - Performs the computation M1 = M2 * M3 where: M1 = (m x n) M2 = (m x p) M3 = (p x n)

Definition at line 31 of file dense_matrix_base.C.

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

34 {
35  // Assertions to make sure we have been
36  // passed matrices of the correct dimension.
37  libmesh_assert_equal_to (M1.m(), M2.m());
38  libmesh_assert_equal_to (M1.n(), M3.n());
39  libmesh_assert_equal_to (M2.n(), M3.m());
40 
41  const unsigned int m_s = M2.m();
42  const unsigned int p_s = M2.n();
43  const unsigned int n_s = M1.n();
44 
45  // Do it this way because there is a
46  // decent chance (at least for constraint matrices)
47  // that M3(k,j) = 0. when right-multiplying.
48  for (unsigned int k=0; k<p_s; k++)
49  for (unsigned int j=0; j<n_s; j++)
50  if (M3.el(k,j) != 0.)
51  for (unsigned int i=0; i<m_s; i++)
52  M1.el(i,j) += M2.el(i,k) * M3.el(k,j);
53 }

◆ n()

template<typename T>
unsigned int libMesh::DenseMatrixBase< T >::n ( ) const
inlineinherited
Returns
The column-dimension of the matrix.

Definition at line 107 of file dense_matrix_base.h.

References libMesh::DenseMatrixBase< T >::_n.

Referenced by libMesh::DenseMatrix< Number >::_multiply_blas(), libMesh::DenseMatrix< Number >::_svd_solve_lapack(), libMesh::DenseMatrixBase< T >::add(), libMesh::DenseMatrix< Number >::add(), libMesh::PetscMatrix< T >::add_block_matrix(), libMesh::SparseMatrix< ValOut >::add_block_matrix(), libMesh::EigenSparseMatrix< T >::add_matrix(), libMesh::LaspackMatrix< T >::add_matrix(), libMesh::EpetraMatrix< T >::add_matrix(), libMesh::PetscMatrix< T >::add_matrix(), libMesh::DofMap::build_constraint_matrix(), libMesh::DofMap::build_constraint_matrix_and_vector(), libMesh::DofMap::constrain_element_dyad_matrix(), libMesh::DofMap::constrain_element_matrix(), libMesh::DofMap::constrain_element_matrix_and_vector(), libMesh::DofMap::constrain_element_vector(), libMesh::DofMap::extract_local_vector(), libMesh::DenseMatrix< Number >::get_transpose(), libMesh::DofMap::heterogenously_constrain_element_matrix_and_vector(), libMesh::DofMap::heterogenously_constrain_element_vector(), libMesh::DenseMatrix< Number >::left_multiply(), libMesh::DenseMatrix< Number >::left_multiply_transpose(), libMesh::DofMap::max_constraint_error(), libMesh::DenseMatrixBase< T >::multiply(), libMesh::WeightedPatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::PatchRecoveryErrorEstimator::EstimateError::operator()(), libMesh::DenseMatrix< Number >::operator=(), libMesh::DenseMatrix< Number >::right_multiply(), and libMesh::DenseMatrix< Number >::right_multiply_transpose().

107 { return _n; }

◆ operator()() [1/2]

template<typename T >
T libMesh::DenseSubMatrix< T >::operator() ( const unsigned int  i,
const unsigned int  j 
) const
inline
Returns
The (i,j) element of the submatrix.

Definition at line 203 of file dense_submatrix.h.

205 {
206  libmesh_assert_less (i, this->m());
207  libmesh_assert_less (j, this->n());
208  libmesh_assert_less (i + this->i_off(), _parent_matrix.m());
209  libmesh_assert_less (j + this->j_off(), _parent_matrix.n());
210 
211  return _parent_matrix (i + this->i_off(),
212  j + this->j_off());
213 }
unsigned int j_off() const
unsigned int m() const
unsigned int n() const
DenseMatrix< T > & _parent_matrix
unsigned int i_off() const

◆ operator()() [2/2]

template<typename T >
T & libMesh::DenseSubMatrix< T >::operator() ( const unsigned int  i,
const unsigned int  j 
)
inline
Returns
The (i,j) element of the submatrix as a writable reference.

Definition at line 218 of file dense_submatrix.h.

220 {
221  libmesh_assert_less (i, this->m());
222  libmesh_assert_less (j, this->n());
223  libmesh_assert_less (i + this->i_off(), _parent_matrix.m());
224  libmesh_assert_less (j + this->j_off(), _parent_matrix.n());
225 
226  return _parent_matrix (i + this->i_off(),
227  j + this->j_off());
228 }
unsigned int j_off() const
unsigned int m() const
unsigned int n() const
DenseMatrix< T > & _parent_matrix
unsigned int i_off() const

◆ operator=() [1/2]

template<typename T>
DenseSubMatrix& libMesh::DenseSubMatrix< T >::operator= ( const DenseSubMatrix< T > &  )
default

◆ operator=() [2/2]

template<typename T>
DenseSubMatrix& libMesh::DenseSubMatrix< T >::operator= ( DenseSubMatrix< T > &&  )
default

◆ parent()

template<typename T>
DenseMatrix<T>& libMesh::DenseSubMatrix< T >::parent ( )
inline
Returns
A reference to the parent matrix.

Definition at line 73 of file dense_submatrix.h.

References libMesh::DenseSubMatrix< T >::_parent_matrix.

Referenced by libMesh::DenseSubMatrix< T >::condense().

73 { return _parent_matrix; }
DenseMatrix< T > & _parent_matrix

◆ print()

template<typename T >
void libMesh::DenseMatrixBase< T >::print ( std::ostream &  os = libMesh::out) const
inherited

Pretty-print the matrix, by default to libMesh::out.

Definition at line 110 of file dense_matrix_base.C.

111 {
112  for (unsigned int i=0; i<this->m(); i++)
113  {
114  for (unsigned int j=0; j<this->n(); j++)
115  os << std::setw(8)
116  << this->el(i,j) << " ";
117 
118  os << std::endl;
119  }
120 
121  return;
122 }
unsigned int m() const
virtual T el(const unsigned int i, const unsigned int j) const =0
unsigned int n() const

◆ print_scientific()

template<typename T >
void libMesh::DenseMatrixBase< T >::print_scientific ( std::ostream &  os,
unsigned  precision = 8 
) const
inherited

Prints the matrix entries with more decimal places in scientific notation.

Definition at line 86 of file dense_matrix_base.C.

87 {
88  // save the initial format flags
89  std::ios_base::fmtflags os_flags = os.flags();
90 
91  // Print the matrix entries.
92  for (unsigned int i=0; i<this->m(); i++)
93  {
94  for (unsigned int j=0; j<this->n(); j++)
95  os << std::setw(15)
96  << std::scientific
97  << std::setprecision(precision)
98  << this->el(i,j) << " ";
99 
100  os << std::endl;
101  }
102 
103  // reset the original format flags
104  os.flags(os_flags);
105 }
unsigned int m() const
virtual T el(const unsigned int i, const unsigned int j) const =0
unsigned int n() const

◆ reposition()

template<typename T >
void libMesh::DenseSubMatrix< T >::reposition ( const unsigned int  ioff,
const unsigned int  joff,
const unsigned int  new_m,
const unsigned int  new_n 
)
inline

Changes the location of the submatrix in the parent matrix.

Definition at line 172 of file dense_submatrix.h.

Referenced by libMesh::DenseSubMatrix< T >::DenseSubMatrix().

176 {
177  _i_off = ioff;
178  _j_off = joff;
179  this->_m = new_m;
180  this->_n = new_n;
181 
182  // Make sure we still fit in the parent matrix.
183  libmesh_assert_less_equal ((this->i_off() + this->m()), _parent_matrix.m());
184  libmesh_assert_less_equal ((this->j_off() + this->n()), _parent_matrix.n());
185 }
unsigned int j_off() const
unsigned int m() const
unsigned int n() const
DenseMatrix< T > & _parent_matrix
unsigned int i_off() const

◆ right_multiply()

template<typename T >
void libMesh::DenseSubMatrix< T >::right_multiply ( const DenseMatrixBase< T > &  M3)
overridevirtual

Performs the operation: (*this) <- (*this) * M3

Implements libMesh::DenseMatrixBase< T >.

Definition at line 49 of file dense_submatrix.C.

50 {
51  // (*this) <- M2 * M3
52  // Where:
53  // (*this) = (m x n),
54  // M2 = (m x p),
55  // M3 = (p x n)
56 
57  // M2 is simply a copy of *this
58  DenseSubMatrix<T> M2(*this);
59 
60  // Call the multiply function in the base class
61  this->multiply(*this, M2, M3);
62 }
static void multiply(DenseMatrixBase< T > &M1, const DenseMatrixBase< T > &M2, const DenseMatrixBase< T > &M3)

◆ zero()

template<typename T >
void libMesh::DenseSubMatrix< T >::zero ( )
inlineoverridevirtual

Set every element in the matrix to 0. You must redefine what you mean by zeroing the matrix since it depends on how your values are stored.

Implements libMesh::DenseMatrixBase< T >.

Definition at line 191 of file dense_submatrix.h.

192 {
193  for (unsigned int i=0; i<this->m(); i++)
194  for (unsigned int j=0; j<this->n(); j++)
195  _parent_matrix(i + this->i_off(),
196  j + this->j_off()) = 0.;
197 }
unsigned int j_off() const
unsigned int m() const
unsigned int n() const
DenseMatrix< T > & _parent_matrix
unsigned int i_off() const

Member Data Documentation

◆ _i_off

template<typename T>
unsigned int libMesh::DenseSubMatrix< T >::_i_off
private

The row offset into the parent matrix.

Definition at line 145 of file dense_submatrix.h.

Referenced by libMesh::DenseSubMatrix< T >::i_off().

◆ _j_off

template<typename T>
unsigned int libMesh::DenseSubMatrix< T >::_j_off
private

The column offset into the parent matrix.

Definition at line 150 of file dense_submatrix.h.

Referenced by libMesh::DenseSubMatrix< T >::j_off().

◆ _m

template<typename T>
unsigned int libMesh::DenseMatrixBase< T >::_m
protectedinherited

The row dimension.

Definition at line 169 of file dense_matrix_base.h.

Referenced by libMesh::DenseMatrixBase< T >::m(), and libMesh::DenseMatrix< Number >::swap().

◆ _n

template<typename T>
unsigned int libMesh::DenseMatrixBase< T >::_n
protectedinherited

The column dimension.

Definition at line 174 of file dense_matrix_base.h.

Referenced by libMesh::DenseMatrixBase< T >::n(), and libMesh::DenseMatrix< Number >::swap().

◆ _parent_matrix

template<typename T>
DenseMatrix<T>& libMesh::DenseSubMatrix< T >::_parent_matrix
private

The parent matrix that contains this submatrix.

Definition at line 140 of file dense_submatrix.h.

Referenced by libMesh::DenseSubMatrix< T >::parent().


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