libMesh
eigen_sparse_vector.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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 #include <algorithm> // for std::min
22 #include <limits>
23 
24 // Local Includes
25 #include "libmesh/dense_subvector.h"
26 #include "libmesh/dense_vector.h"
27 #include "libmesh/eigen_sparse_vector.h"
28 #include "libmesh/eigen_sparse_matrix.h"
29 #include "libmesh/int_range.h"
30 
31 #ifdef LIBMESH_HAVE_EIGEN
32 
33 namespace libMesh
34 {
35 
36 template <typename T>
38 {
39  libmesh_assert (this->closed());
40  libmesh_assert (this->initialized());
41 
42  return _vec.sum();
43 }
44 
45 
46 
47 template <typename T>
49 {
50  libmesh_assert (this->closed());
51  libmesh_assert (this->initialized());
52 
53  return _vec.lpNorm<1>();
54 }
55 
56 
57 
58 template <typename T>
60 {
61  libmesh_assert (this->closed());
62  libmesh_assert (this->initialized());
63 
64  return _vec.lpNorm<2>();
65 }
66 
67 
68 
69 template <typename T>
71 {
72  libmesh_assert (this->closed());
73  libmesh_assert (this->initialized());
74 
75  return _vec.lpNorm<Eigen::Infinity>();
76 }
77 
78 
79 
80 template <typename T>
82 {
83  libmesh_assert (this->closed());
84 
85  const EigenSparseVector<T> & v = cast_ref<const EigenSparseVector<T> &>(v_in);
86 
87  _vec += v._vec;
88 
89  return *this;
90 }
91 
92 
93 
94 
95 template <typename T>
97 {
98  libmesh_assert (this->closed());
99 
100  const EigenSparseVector<T> & v = cast_ref<const EigenSparseVector<T> &>(v_in);
101 
102  _vec -= v._vec;
103 
104  return *this;
105 }
106 
107 
108 
109 template <typename T>
111 {
112  libmesh_assert (this->closed());
113  libmesh_assert_equal_to(size(), v_in.size());
114 
115  const EigenSparseVector<T> & v = cast_ref<const EigenSparseVector<T> &>(v_in);
116 
117  _vec = _vec.cwiseProduct(v._vec);
118 
119  return *this;
120 }
121 
122 
123 
124 template <typename T>
126 {
127  libmesh_assert (this->closed());
128  libmesh_assert_equal_to(size(), v_in.size());
129 
130  const EigenSparseVector<T> & v = cast_ref<const EigenSparseVector<T> &>(v_in);
131 
132  _vec = _vec.cwiseQuotient(v._vec);
133 
134  return *this;
135 }
136 
137 
138 
139 
140 template <typename T>
142 {
143 #ifndef NDEBUG
144  const numeric_index_type n = this->size();
145 
146  for (numeric_index_type i=0; i<n; i++)
147  // Don't divide by zero!
148  libmesh_assert_not_equal_to ((*this)(i), T(0));
149 #endif
150 
151  _vec = _vec.cwiseInverse();
152 }
153 
154 
155 
156 template <typename T>
158 {
159  _vec = _vec.conjugate();
160 }
161 
162 
163 
164 template <typename T>
166 {
167  _vec += EigenSV::Constant(this->size(), v);
168 
169  this->_is_closed = false;
170 }
171 
172 
173 
174 
175 template <typename T>
177 {
178  libmesh_assert (this->initialized());
179 
180  const EigenSparseVector<T> & v = cast_ref<const EigenSparseVector<T> &>(v_in);
181 
182  _vec += v._vec;
183 }
184 
185 
186 
187 template <typename T>
188 void EigenSparseVector<T>::add (const T a, const NumericVector<T> & v_in)
189 {
190  libmesh_assert (this->initialized());
191 
192  const EigenSparseVector<T> & v = cast_ref<const EigenSparseVector<T> &>(v_in);
193 
194  _vec += v._vec*a;
195 }
196 
197 
198 
199 template <typename T>
201  const SparseMatrix<T> & mat_in)
202 {
203  // Make sure the data passed in are really in Eigen types
204  const EigenSparseVector<T> * e_vec = cast_ptr<const EigenSparseVector<T> *>(&vec_in);
205  const EigenSparseMatrix<T> * mat = cast_ptr<const EigenSparseMatrix<T> *>(&mat_in);
206 
207  libmesh_assert(e_vec);
208  libmesh_assert(mat);
209 
210  _vec += mat->_mat*e_vec->_vec;
211 }
212 
213 
214 
215 template <typename T>
217  const SparseMatrix<T> & mat_in)
218 {
219  // Make sure the data passed in are really in Eigen types
220  const EigenSparseVector<T> * e_vec = cast_ptr<const EigenSparseVector<T> *>(&vec_in);
221  const EigenSparseMatrix<T> * mat = cast_ptr<const EigenSparseMatrix<T> *>(&mat_in);
222 
223  libmesh_assert(e_vec);
224  libmesh_assert(mat);
225 
226  _vec += mat->_mat.transpose()*e_vec->_vec;
227 }
228 
229 
230 
231 template <typename T>
232 void EigenSparseVector<T>::scale (const T factor)
233 {
234  libmesh_assert (this->initialized());
235 
236  _vec *= factor;
237 }
238 
239 
240 
241 template <typename T>
243 {
244  libmesh_assert (this->initialized());
245 
246  const numeric_index_type n = this->size();
247 
248  for (numeric_index_type i=0; i!=n; ++i)
249  this->set(i,std::abs((*this)(i)));
250 }
251 
252 
253 
254 template <typename T>
256 {
257  libmesh_assert (this->initialized());
258 
259  // Make sure the NumericVector passed in is really a EigenSparseVector
260  const EigenSparseVector<T> * v = cast_ptr<const EigenSparseVector<T> *>(&v_in);
261  libmesh_assert(v);
262 
263  return _vec.dot(v->_vec);
264 }
265 
266 
267 
268 template <typename T>
271 {
272  libmesh_assert (this->initialized());
273  libmesh_assert (this->closed());
274 
275  _vec.fill(s);
276 
277  return *this;
278 }
279 
280 
281 
282 template <typename T>
285 {
286  // Make sure the NumericVector passed in is really a EigenSparseVector
287  const EigenSparseVector<T> * v =
288  cast_ptr<const EigenSparseVector<T> *>(&v_in);
289 
290  libmesh_assert(v);
291 
292  *this = *v;
293 
294  return *this;
295 }
296 
297 
298 
299 template <typename T>
302 {
303  libmesh_assert (this->initialized());
304  libmesh_assert (v.closed());
305  libmesh_assert_equal_to (this->size(), v.size());
306 
307  _vec = v._vec;
308 
309  this->_is_closed = true;
310 
311  return *this;
312 }
313 
314 
315 
316 template <typename T>
318 EigenSparseVector<T>::operator = (const std::vector<T> & v)
319 {
324  if (this->size() == v.size())
325  for (auto i : index_range(v))
326  this->set (i, v[i]);
327 
328  else
329  libmesh_error_msg("this->size() = " << this->size() << " must be equal to v.size() = " << v.size());
330 
331  return *this;
332 }
333 
334 
335 template <typename T>
337 {
338  // Make sure the NumericVector passed in is really a EigenSparseVector
339  EigenSparseVector<T> * v_local =
340  cast_ptr<EigenSparseVector<T> *>(&v_local_in);
341 
342  libmesh_assert(v_local);
343 
344  *v_local = *this;
345 }
346 
347 
348 
349 template <typename T>
351  const std::vector<numeric_index_type> & libmesh_dbg_var(send_list)) const
352 {
353  // Make sure the NumericVector passed in is really a EigenSparseVector
354  EigenSparseVector<T> * v_local =
355  cast_ptr<EigenSparseVector<T> *>(&v_local_in);
356 
357  libmesh_assert(v_local);
358  libmesh_assert_less_equal (send_list.size(), v_local->size());
359 
360  *v_local = *this;
361 }
362 
363 
364 
365 template <typename T>
366 void EigenSparseVector<T>::localize (std::vector<T> & v_local,
367  const std::vector<numeric_index_type> & indices) const
368 {
369  // EigenSparseVectors are serial, so we can just copy values
370  v_local.resize(indices.size());
371 
372  for (auto i : index_range(v_local))
373  v_local[i] = (*this)(indices[i]);
374 }
375 
376 
377 
378 template <typename T>
379 void EigenSparseVector<T>::localize (const numeric_index_type libmesh_dbg_var(first_local_idx),
380  const numeric_index_type libmesh_dbg_var(last_local_idx),
381  const std::vector<numeric_index_type> & libmesh_dbg_var(send_list))
382 {
383  libmesh_assert_equal_to (first_local_idx, 0);
384  libmesh_assert_equal_to (last_local_idx+1, this->size());
385 
386  libmesh_assert_less_equal (send_list.size(), this->size());
387 
388  this->_is_closed = true;
389 }
390 
391 
392 
393 template <typename T>
394 void EigenSparseVector<T>::localize (std::vector<T> & v_local) const
395 
396 {
397  v_local.resize(this->size());
398 
399  for (auto i : index_range(v_local))
400  v_local[i] = (*this)(i);
401 }
402 
403 
404 
405 template <typename T>
406 void EigenSparseVector<T>::localize_to_one (std::vector<T> & v_local,
407  const processor_id_type libmesh_dbg_var(pid)) const
408 {
409  libmesh_assert_equal_to (pid, 0);
410 
411  this->localize (v_local);
412 }
413 
414 
415 
416 template <typename T>
418  const NumericVector<T> & /*vec2*/)
419 {
420  libmesh_not_implemented();
421 }
422 
423 
424 
425 template <typename T>
427 {
428  libmesh_assert (this->initialized());
429  if (!this->size())
430  return -std::numeric_limits<Real>::max();
431 
432 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
433  Real the_max = libmesh_real((*this)(0));
434 
435  const numeric_index_type n = this->size();
436 
437  for (numeric_index_type i=1; i<n; i++)
438  the_max = std::max (the_max, libmesh_real((*this)(i)));
439 
440  return the_max;
441 #else
442  return libmesh_real(_vec.maxCoeff());
443 #endif
444 }
445 
446 
447 
448 template <typename T>
450 {
451  libmesh_assert (this->initialized());
452  if (!this->size())
453  return std::numeric_limits<Real>::max();
454 
455 #ifdef LIBMESH_USE_COMPLEX_NUMBERS
456  Real the_min = libmesh_real((*this)(0));
457 
458  const numeric_index_type n = this->size();
459 
460  for (numeric_index_type i=1; i<n; i++)
461  the_min = std::min (the_min, libmesh_real((*this)(i)));
462 
463  return the_min;
464 #else
465  return libmesh_real(_vec.minCoeff());
466 #endif
467 }
468 
469 
470 //------------------------------------------------------------------
471 // Explicit instantiations
472 template class EigenSparseVector<Number>;
473 
474 } // namespace libMesh
475 
476 
477 #endif // #ifdef LIBMESH_HAVE_EIGEN
libMesh::EigenSparseVector::max
virtual Real max() const override
Definition: eigen_sparse_vector.C:426
libMesh::EigenSparseVector::operator=
EigenSparseVector< T > & operator=(const EigenSparseVector< T > &v)
Copy assignment operator.
Definition: eigen_sparse_vector.C:301
libMesh::EigenSparseVector::operator*=
virtual NumericVector< T > & operator*=(const NumericVector< T > &v_in) override
Computes the component-wise multiplication of this vector's entries by another's, .
Definition: eigen_sparse_vector.C:110
libMesh::EigenSparseVector::size
virtual numeric_index_type size() const override
Definition: eigen_sparse_vector.h:431
libMesh::EigenSparseMatrix
The EigenSparseMatrix class wraps a sparse matrix object from the Eigen library.
Definition: eigen_sparse_matrix.h:54
libMesh::EigenSparseVector::abs
virtual void abs() override
Sets for each entry in the vector.
Definition: eigen_sparse_vector.C:242
libMesh::libmesh_real
T libmesh_real(T a)
Definition: libmesh_common.h:166
libMesh::EigenSparseVector::add_vector_transpose
virtual void add_vector_transpose(const NumericVector< T > &v, const SparseMatrix< T > &A) override
Computes , i.e.
Definition: eigen_sparse_vector.C:216
libMesh::EigenSparseVector::scale
virtual void scale(const T factor) override
Scale each element of the vector by the given factor.
Definition: eigen_sparse_vector.C:232
libMesh::index_range
IntRange< std::size_t > index_range(const std::vector< T > &vec)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
Definition: int_range.h:106
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::EigenSparseVector::localize_to_one
virtual void localize_to_one(std::vector< T > &v_local, const processor_id_type proc_id=0) const override
Creates a local copy of the global vector in v_local only on processor proc_id.
Definition: eigen_sparse_vector.C:406
libMesh::EigenSparseVector::l1_norm
virtual Real l1_norm() const override
Definition: eigen_sparse_vector.C:48
libMesh::EigenSparseVector::add
virtual void add(const numeric_index_type i, const T value) override
Adds value to each entry of the vector.
Definition: eigen_sparse_vector.h:489
libMesh::SparseMatrix
Generic sparse matrix.
Definition: vector_fe_ex5.C:45
libMesh::NumericVector::closed
virtual bool closed() const
Definition: numeric_vector.h:171
libMesh::NumericVector::size
virtual numeric_index_type size() const =0
libMesh::NumericVector
Provides a uniform interface to vector storage schemes for different linear algebra libraries.
Definition: vector_fe_ex5.C:43
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::EigenSparseVector::dot
virtual T dot(const NumericVector< T > &v) const override
Definition: eigen_sparse_vector.C:255
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::EigenSparseVector::operator+=
virtual NumericVector< T > & operator+=(const NumericVector< T > &v) override
Adds v to *this, .
Definition: eigen_sparse_vector.C:81
libMesh::EigenSparseVector
This class provides a nice interface to the Eigen C++-based data structures for serial vectors.
Definition: eigen_sparse_matrix.h:42
libMesh::processor_id_type
uint8_t processor_id_type
Definition: id_types.h:104
libMesh::EigenSparseVector::localize
virtual void localize(std::vector< T > &v_local) const override
Creates a copy of the global vector in the local vector v_local.
Definition: eigen_sparse_vector.C:394
libMesh::EigenSparseVector::linfty_norm
virtual Real linfty_norm() const override
Definition: eigen_sparse_vector.C:70
libMesh::EigenSparseVector::conjugate
virtual void conjugate() override
Negates the imaginary component of each entry in the vector.
Definition: eigen_sparse_vector.C:157
libMesh::EigenSparseVector::pointwise_mult
virtual void pointwise_mult(const NumericVector< T > &vec1, const NumericVector< T > &vec2) override
Computes (summation not implied) i.e.
Definition: eigen_sparse_vector.C:417
libMesh::EigenSparseVector::_vec
DataType _vec
Actual Eigen::SparseVector<> we are wrapping.
Definition: eigen_sparse_vector.h:242
libMesh::EigenSparseVector::sum
virtual T sum() const override
Definition: eigen_sparse_vector.C:37
libMesh::numeric_index_type
dof_id_type numeric_index_type
Definition: id_types.h:99
libMesh::initialized
bool initialized()
Checks that library initialization has been done.
Definition: libmesh.C:265
libMesh::EigenSparseVector::add_vector
virtual void add_vector(const NumericVector< T > &v, const SparseMatrix< T > &A) override
Computes , i.e.
Definition: eigen_sparse_vector.C:200
libMesh::EigenSparseVector::reciprocal
virtual void reciprocal() override
Computes the component-wise reciprocal, .
Definition: eigen_sparse_vector.C:141
libMesh::EigenSparseVector::min
virtual Real min() const override
Definition: eigen_sparse_vector.C:449
libMesh::EigenSparseVector::operator-=
virtual NumericVector< T > & operator-=(const NumericVector< T > &v) override
Subtracts v from *this, .
Definition: eigen_sparse_vector.C:96
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::EigenSparseVector::l2_norm
virtual Real l2_norm() const override
Definition: eigen_sparse_vector.C:59
libMesh::closed
bool closed()
Checks that the library has been closed.
Definition: libmesh.C:272
libMesh::EigenSparseVector::operator/=
virtual NumericVector< T > & operator/=(const NumericVector< T > &v_in) override
Computes the component-wise division of this vector's entries by another's, .
Definition: eigen_sparse_vector.C:125