eigen_sparse_matrix.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 // C++ includes
21 
22 // Local includes
23 #include "libmesh/libmesh_config.h"
24 
25 #ifdef LIBMESH_HAVE_EIGEN
26 
29 #include "libmesh/dense_matrix.h"
30 #include "libmesh/dof_map.h"
32 
33 namespace libMesh
34 {
35 
36 
37 //-----------------------------------------------------------------------
38 // EigenSparseMatrix members
39 template <typename T>
41  const numeric_index_type n_in,
42  const numeric_index_type libmesh_dbg_var(m_l),
43  const numeric_index_type libmesh_dbg_var(n_l),
44  const numeric_index_type nnz,
45  const numeric_index_type,
46  const numeric_index_type)
47 {
48  // noz ignored... only used for multiple processors!
49  libmesh_assert_equal_to (m_in, m_l);
50  libmesh_assert_equal_to (n_in, n_l);
51  libmesh_assert_equal_to (m_in, n_in);
52  libmesh_assert_greater (nnz, 0);
53 
54  _mat.resize(m_in, n_in);
55  _mat.reserve(Eigen::Matrix<numeric_index_type, Eigen::Dynamic, 1>::Constant(m_in,nnz));
56 
57  this->_is_initialized = true;
58 }
59 
60 
61 
62 template <typename T>
64 {
65  // Ignore calls on initialized objects
66  if (this->initialized())
67  return;
68 
69  // We need the DofMap for this!
70  libmesh_assert(this->_dof_map);
71 
72  // Clear initialized matrices
73  if (this->initialized())
74  this->clear();
75 
76  const numeric_index_type n_rows = this->_dof_map->n_dofs();
77  const numeric_index_type n_cols = n_rows;
78 
79 #ifndef NDEBUG
80  // The following variables are only used for assertions,
81  // so avoid declaring them when asserts are inactive.
82  const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(0);
83  const numeric_index_type m_l = n_l;
84 #endif
85 
86  // Laspack Matrices only work for uniprocessor cases
87  libmesh_assert_equal_to (n_rows, n_cols);
88  libmesh_assert_equal_to (m_l, n_rows);
89  libmesh_assert_equal_to (n_l, n_cols);
90 
91  const std::vector<numeric_index_type> & n_nz = this->_dof_map->get_n_nz();
92 
93 #ifndef NDEBUG
94  // The following variables are only used for assertions,
95  // so avoid declaring them when asserts are inactive.
96  const std::vector<numeric_index_type> & n_oz = this->_dof_map->get_n_oz();
97 #endif
98 
99  // Make sure the sparsity pattern isn't empty
100  libmesh_assert_equal_to (n_nz.size(), n_l);
101  libmesh_assert_equal_to (n_oz.size(), n_l);
102 
103  if (n_rows==0)
104  {
105  _mat.resize(0,0);
106  return;
107  }
108 
109  _mat.resize(n_rows,n_cols);
110  _mat.reserve(n_nz);
111 
112  this->_is_initialized = true;
113 
114  libmesh_assert_equal_to (n_rows, this->m());
115  libmesh_assert_equal_to (n_cols, this->n());
116 }
117 
118 
119 
120 template <typename T>
122  const std::vector<numeric_index_type> & rows,
123  const std::vector<numeric_index_type> & cols)
124 
125 {
126  libmesh_assert (this->initialized());
127  unsigned int n_rows = cast_int<unsigned int>(rows.size());
128  unsigned int n_cols = cast_int<unsigned int>(cols.size());
129  libmesh_assert_equal_to (dm.m(), n_rows);
130  libmesh_assert_equal_to (dm.n(), n_cols);
131 
132 
133  for (unsigned int i=0; i<n_rows; i++)
134  for (unsigned int j=0; j<n_cols; j++)
135  this->add(rows[i],cols[j],dm(i,j));
136 }
137 
138 
139 
140 template <typename T>
142 {
143  EigenSparseVector<T> & dest = cast_ref<EigenSparseVector<T> &>(dest_in);
144 
145  dest._vec = _mat.diagonal();
146 }
147 
148 
149 
150 template <typename T>
152 {
153  EigenSparseMatrix<T> & dest = cast_ref<EigenSparseMatrix<T> &>(dest_in);
154 
155  dest._mat = _mat.transpose();
156 }
157 
158 
159 
160 template <typename T>
162  SparseMatrix<T>(comm_in),
163  _closed (false)
164 {
165 }
166 
167 
168 
169 template <typename T>
171 {
172  _mat.resize(0,0);
173 
174  _closed = false;
175  this->_is_initialized = false;
176 }
177 
178 
179 
180 template <typename T>
182 {
183  _mat.setZero();
184 }
185 
186 
187 
188 template <typename T>
190 {
191  libmesh_assert (this->initialized());
192 
193  return cast_int<numeric_index_type>(_mat.rows());
194 }
195 
196 
197 
198 template <typename T>
200 {
201  libmesh_assert (this->initialized());
202 
203  return cast_int<numeric_index_type>(_mat.cols());
204 }
205 
206 
207 
208 template <typename T>
210 {
211  return 0;
212 }
213 
214 
215 
216 template <typename T>
218 {
219  return this->m();
220 }
221 
222 
223 
224 template <typename T>
226  const numeric_index_type j,
227  const T value)
228 {
229  libmesh_assert (this->initialized());
230  libmesh_assert_less (i, this->m());
231  libmesh_assert_less (j, this->n());
232 
233  _mat.coeffRef(i,j) = value;
234 }
235 
236 
237 
238 template <typename T>
240  const numeric_index_type j,
241  const T value)
242 {
243  libmesh_assert (this->initialized());
244  libmesh_assert_less (i, this->m());
245  libmesh_assert_less (j, this->n());
246 
247  _mat.coeffRef(i,j) += value;
248 }
249 
250 
251 
252 template <typename T>
254  const std::vector<numeric_index_type> & dof_indices)
255 {
256  this->add_matrix (dm, dof_indices, dof_indices);
257 }
258 
259 
260 
261 template <typename T>
262 void EigenSparseMatrix<T>::add (const T a_in, const SparseMatrix<T> & X_in)
263 {
264  libmesh_assert (this->initialized());
265  libmesh_assert_equal_to (this->m(), X_in.m());
266  libmesh_assert_equal_to (this->n(), X_in.n());
267 
268  const EigenSparseMatrix<T> & X =
269  cast_ref<const EigenSparseMatrix<T> &> (X_in);
270 
271  _mat += X._mat*a_in;
272 }
273 
274 
275 
276 
277 template <typename T>
279  const numeric_index_type j) const
280 {
281  libmesh_assert (this->initialized());
282  libmesh_assert_less (i, this->m());
283  libmesh_assert_less (j, this->n());
284 
285  return _mat.coeff(i,j);
286 }
287 
288 
289 
290 template <typename T>
292 {
293  // There does not seem to be a straightforward way to iterate over
294  // the columns of an EigenSparseMatrix. So we use some extra
295  // storage and keep track of the column sums while going over the
296  // row entries...
297  std::vector<Real> abs_col_sums(this->n());
298 
299  // For a row-major Eigen SparseMatrix like we're using, the
300  // InnerIterator iterates over the non-zero entries of rows.
301  for (unsigned row=0; row<this->m(); ++row)
302  {
303  EigenSM::InnerIterator it(_mat, row);
304  for (; it; ++it)
305  abs_col_sums[it.col()] += std::abs(it.value());
306  }
307 
308  return *(std::max_element(abs_col_sums.begin(), abs_col_sums.end()));
309 }
310 
311 
312 
313 template <typename T>
315 {
316  Real max_abs_row_sum = 0.;
317 
318  // For a row-major Eigen SparseMatrix like we're using, the
319  // InnerIterator iterates over the non-zero entries of rows.
320  for (unsigned row=0; row<this->m(); ++row)
321  {
322  Real current_abs_row_sum = 0.;
323  EigenSM::InnerIterator it(_mat, row);
324  for (; it; ++it)
325  current_abs_row_sum += std::abs(it.value());
326 
327  max_abs_row_sum = std::max(max_abs_row_sum, current_abs_row_sum);
328  }
329 
330  return max_abs_row_sum;
331 }
332 
333 
334 //------------------------------------------------------------------
335 // Explicit instantiations
336 template class EigenSparseMatrix<Number>;
337 
338 } // namespace libMesh
339 
340 
341 #endif // #ifdef LIBMESH_HAVE_EIGEN
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value) override
virtual void init() override
double abs(double a)
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) override
virtual numeric_index_type m() const override
virtual numeric_index_type row_stop() const override
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: diff_context.h:40
unsigned int m() const
long double max(long double a, double b)
virtual void get_transpose(SparseMatrix< T > &dest) const override
virtual void zero() override
virtual numeric_index_type n() const override
virtual void get_diagonal(NumericVector< T > &dest) const override
dof_id_type numeric_index_type
Definition: id_types.h:92
EigenSparseMatrix(const Parallel::Communicator &comm)
virtual numeric_index_type m() const =0
virtual T operator()(const numeric_index_type i, const numeric_index_type j) const override
virtual numeric_index_type row_start() const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real linfty_norm() const override
virtual void clear() override
static const bool value
Definition: xdr_io.C:109
bool initialized()
Definition: libmesh.C:258
virtual Real l1_norm() const override
unsigned int n() const
A matrix object used for finite element assembly and numerics.
Definition: dense_matrix.h:54
virtual numeric_index_type n() const =0
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) override