libMesh
numeric_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 <cstdlib> // *must* precede <cmath> for proper std:abs() on PGI, Sun Studio CC
22 #include <cmath> // for std::abs
23 #include <limits>
24 
25 // Local Includes
26 #include "libmesh/numeric_vector.h"
27 #include "libmesh/distributed_vector.h"
28 #include "libmesh/laspack_vector.h"
29 #include "libmesh/eigen_sparse_vector.h"
30 #include "libmesh/petsc_vector.h"
31 #include "libmesh/trilinos_epetra_vector.h"
32 #include "libmesh/shell_matrix.h"
33 #include "libmesh/tensor_tools.h"
34 #include "libmesh/auto_ptr.h" // libmesh_make_unique
35 #include "libmesh/enum_solver_package.h"
36 #include "libmesh/int_range.h"
37 
38 namespace libMesh
39 {
40 
41 
42 
43 //------------------------------------------------------------------
44 // NumericVector methods
45 
46 // Full specialization for Real datatypes
47 template <typename T>
48 std::unique_ptr<NumericVector<T>>
49 NumericVector<T>::build(const Parallel::Communicator & comm, const SolverPackage solver_package)
50 {
51  // Build the appropriate vector
52  switch (solver_package)
53  {
54 
55 #ifdef LIBMESH_HAVE_LASPACK
56  case LASPACK_SOLVERS:
57  return libmesh_make_unique<LaspackVector<T>>(comm, AUTOMATIC);
58 #endif
59 
60 #ifdef LIBMESH_HAVE_PETSC
61  case PETSC_SOLVERS:
62  return libmesh_make_unique<PetscVector<T>>(comm, AUTOMATIC);
63 #endif
64 
65 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA
66  case TRILINOS_SOLVERS:
67  return libmesh_make_unique<EpetraVector<T>>(comm, AUTOMATIC);
68 #endif
69 
70 #ifdef LIBMESH_HAVE_EIGEN
71  case EIGEN_SOLVERS:
72  return libmesh_make_unique<EigenSparseVector<T>>(comm, AUTOMATIC);
73 #endif
74 
75  default:
76  return libmesh_make_unique<DistributedVector<T>>(comm, AUTOMATIC);
77  }
78 }
79 
80 
81 
82 template <typename T>
83 void NumericVector<T>::insert (const T * v,
84  const std::vector<numeric_index_type> & dof_indices)
85 {
86  for (auto i : index_range(dof_indices))
87  this->set (dof_indices[i], v[i]);
88 }
89 
90 
91 
92 template <typename T>
94  const std::vector<numeric_index_type> & dof_indices)
95 {
96  libmesh_assert_equal_to (V.size(), dof_indices.size());
97 
98  for (auto i : index_range(dof_indices))
99  this->set (dof_indices[i], V(i));
100 }
101 
102 
103 
104 template <typename T>
105 int NumericVector<T>::compare (const NumericVector<T> & other_vector,
106  const Real threshold) const
107 {
108  libmesh_assert (this->initialized());
109  libmesh_assert (other_vector.initialized());
110  libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
111  libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
112 
113  int first_different_i = std::numeric_limits<int>::max();
114  numeric_index_type i = first_local_index();
115 
116  do
117  {
118  if (std::abs((*this)(i) - other_vector(i)) > threshold)
119  first_different_i = i;
120  else
121  i++;
122  }
123  while (first_different_i==std::numeric_limits<int>::max()
124  && i<last_local_index());
125 
126  // Find the correct first differing index in parallel
127  this->comm().min(first_different_i);
128 
129  if (first_different_i == std::numeric_limits<int>::max())
130  return -1;
131 
132  return first_different_i;
133 }
134 
135 
136 template <typename T>
138  const Real threshold) const
139 {
140  libmesh_assert (this->initialized());
141  libmesh_assert (other_vector.initialized());
142  libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
143  libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
144 
145  int first_different_i = std::numeric_limits<int>::max();
146  numeric_index_type i = first_local_index();
147 
148  do
149  {
150  if (std::abs((*this)(i) - other_vector(i)) > threshold *
151  std::max(std::abs((*this)(i)), std::abs(other_vector(i))))
152  first_different_i = i;
153  else
154  i++;
155  }
156  while (first_different_i==std::numeric_limits<int>::max()
157  && i<last_local_index());
158 
159  // Find the correct first differing index in parallel
160  this->comm().min(first_different_i);
161 
162  if (first_different_i == std::numeric_limits<int>::max())
163  return -1;
164 
165  return first_different_i;
166 }
167 
168 
169 template <typename T>
171  const Real threshold) const
172 {
173  libmesh_assert (this->initialized());
174  libmesh_assert (other_vector.initialized());
175  libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
176  libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
177 
178  int first_different_i = std::numeric_limits<int>::max();
179  numeric_index_type i = first_local_index();
180 
181  const Real my_norm = this->linfty_norm();
182  const Real other_norm = other_vector.linfty_norm();
183  const Real abs_threshold = std::max(my_norm, other_norm) * threshold;
184 
185  do
186  {
187  if (std::abs((*this)(i) - other_vector(i) ) > abs_threshold)
188  first_different_i = i;
189  else
190  i++;
191  }
192  while (first_different_i==std::numeric_limits<int>::max()
193  && i<last_local_index());
194 
195  // Find the correct first differing index in parallel
196  this->comm().min(first_different_i);
197 
198  if (first_different_i == std::numeric_limits<int>::max())
199  return -1;
200 
201  return first_different_i;
202 }
203 
204 /*
205 // Full specialization for float datatypes (DistributedVector wants this)
206 
207 template <>
208 int NumericVector<float>::compare (const NumericVector<float> & other_vector,
209 const Real threshold) const
210 {
211 libmesh_assert (this->initialized());
212 libmesh_assert (other_vector.initialized());
213 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
214 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
215 
216 int rvalue = -1;
217 numeric_index_type i = first_local_index();
218 
219 do
220 {
221 if (std::abs((*this)(i) - other_vector(i) ) > threshold)
222 rvalue = i;
223 else
224 i++;
225 }
226 while (rvalue==-1 && i<last_local_index());
227 
228 return rvalue;
229 }
230 
231 // Full specialization for double datatypes
232 template <>
233 int NumericVector<double>::compare (const NumericVector<double> & other_vector,
234 const Real threshold) const
235 {
236 libmesh_assert (this->initialized());
237 libmesh_assert (other_vector.initialized());
238 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
239 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
240 
241 int rvalue = -1;
242 numeric_index_type i = first_local_index();
243 
244 do
245 {
246 if (std::abs((*this)(i) - other_vector(i) ) > threshold)
247 rvalue = i;
248 else
249 i++;
250 }
251 while (rvalue==-1 && i<last_local_index());
252 
253 return rvalue;
254 }
255 
256 #ifdef LIBMESH_DEFAULT_TRIPLE_PRECISION
257 // Full specialization for long double datatypes
258 template <>
259 int NumericVector<long double>::compare (const NumericVector<long double> & other_vector,
260 const Real threshold) const
261 {
262 libmesh_assert (this->initialized());
263 libmesh_assert (other_vector.initialized());
264 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
265 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
266 
267 int rvalue = -1;
268 numeric_index_type i = first_local_index();
269 
270 do
271 {
272 if (std::abs((*this)(i) - other_vector(i) ) > threshold)
273 rvalue = i;
274 else
275 i++;
276 }
277 while (rvalue==-1 && i<last_local_index());
278 
279 return rvalue;
280 }
281 #endif
282 
283 
284 // Full specialization for Complex datatypes
285 template <>
286 int NumericVector<Complex>::compare (const NumericVector<Complex> & other_vector,
287 const Real threshold) const
288 {
289 libmesh_assert (this->initialized());
290 libmesh_assert (other_vector.initialized());
291 libmesh_assert_equal_to (this->first_local_index(), other_vector.first_local_index());
292 libmesh_assert_equal_to (this->last_local_index(), other_vector.last_local_index());
293 
294 int rvalue = -1;
295 numeric_index_type i = first_local_index();
296 
297 do
298 {
299 if ((std::abs((*this)(i).real() - other_vector(i).real()) > threshold) || (std::abs((*this)(i).imag() - other_vector(i).imag()) > threshold))
300 rvalue = i;
301 else
302 i++;
303 }
304 while (rvalue==-1 && i<this->last_local_index());
305 
306 return rvalue;
307 }
308 */
309 
310 
311 template <class T>
312 Real NumericVector<T>::subset_l1_norm (const std::set<numeric_index_type> & indices) const
313 {
314  const NumericVector<T> & v = *this;
315 
316  Real norm = 0;
317 
318  for (const auto & index : indices)
319  norm += std::abs(v(index));
320 
321  this->comm().sum(norm);
322 
323  return norm;
324 }
325 
326 template <class T>
327 Real NumericVector<T>::subset_l2_norm (const std::set<numeric_index_type> & indices) const
328 {
329  const NumericVector<T> & v = *this;
330 
331  Real norm = 0;
332 
333  for (const auto & index : indices)
334  norm += TensorTools::norm_sq(v(index));
335 
336  this->comm().sum(norm);
337 
338  return std::sqrt(norm);
339 }
340 
341 template <class T>
342 Real NumericVector<T>::subset_linfty_norm (const std::set<numeric_index_type> & indices) const
343 {
344  const NumericVector<T> & v = *this;
345 
346  Real norm = 0;
347 
348  for (const auto & index : indices)
349  {
350  Real value = std::abs(v(index));
351  if (value > norm)
352  norm = value;
353  }
354 
355  this->comm().max(norm);
356 
357  return norm;
358 }
359 
360 
361 
362 template <typename T>
364  const std::vector<numeric_index_type> & dof_indices)
365 {
366  for (auto i : index_range(dof_indices))
367  this->add (dof_indices[i], v[i]);
368 }
369 
370 
371 
372 template <typename T>
374  const std::vector<numeric_index_type> & dof_indices)
375 {
376  const std::size_t n = dof_indices.size();
377  libmesh_assert_equal_to(v.size(), n);
378  for (numeric_index_type i=0; i != n; i++)
379  this->add (dof_indices[i], v(i));
380 }
381 
382 
383 
384 template <typename T>
386  const ShellMatrix<T> & a)
387 {
388  a.vector_mult_add(*this,v);
389 }
390 
391 
392 
393 //------------------------------------------------------------------
394 // Explicit instantiations
395 template class NumericVector<Number>;
396 
397 } // namespace libMesh
libMesh::ShellMatrix
Generic shell matrix, i.e.
Definition: eigen_preconditioner.h:36
libMesh::SolverPackage
SolverPackage
Defines an enum for various linear solver packages.
Definition: enum_solver_package.h:34
libMesh::PETSC_SOLVERS
Definition: enum_solver_package.h:36
libMesh::NumericVector::subset_linfty_norm
virtual Real subset_linfty_norm(const std::set< numeric_index_type > &indices) const
Definition: numeric_vector.C:342
libMesh::NumericVector::last_local_index
virtual numeric_index_type last_local_index() const =0
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
std::sqrt
MetaPhysicL::DualNumber< T, D > sqrt(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::NumericVector::insert
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices)
Inserts the entries of v in *this at the locations specified by v.
Definition: numeric_vector.C:83
libMesh::NumericVector::build
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
Builds a NumericVector on the processors in communicator comm using the linear solver package specifi...
Definition: numeric_vector.C:49
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::NumericVector::add_vector
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Computes , where v is a pointer and each dof_indices[i] specifies where to add value v[i].
Definition: numeric_vector.C:363
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::NumericVector::subset_l2_norm
virtual Real subset_l2_norm(const std::set< numeric_index_type > &indices) const
Definition: numeric_vector.C:327
libMesh::NumericVector::initialized
virtual bool initialized() const
Definition: numeric_vector.h:155
libMesh::NumericVector::local_relative_compare
virtual int local_relative_compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
Definition: numeric_vector.C:137
libMesh::LASPACK_SOLVERS
Definition: enum_solver_package.h:38
libMesh::numeric_index_type
dof_id_type numeric_index_type
Definition: id_types.h:99
libMesh::NumericVector::linfty_norm
virtual Real linfty_norm() const =0
libMesh::TRILINOS_SOLVERS
Definition: enum_solver_package.h:37
libMesh::initialized
bool initialized()
Checks that library initialization has been done.
Definition: libmesh.C:265
libMesh::TensorTools::norm_sq
T norm_sq(std::complex< T > a)
Definition: tensor_tools.h:85
value
static const bool value
Definition: xdr_io.C:56
std::norm
MetaPhysicL::DualNumber< T, D > norm(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::AUTOMATIC
Definition: enum_parallel_type.h:34
libMesh::ShellMatrix::vector_mult_add
virtual void vector_mult_add(NumericVector< T > &dest, const NumericVector< T > &arg) const =0
Multiplies the matrix with arg and adds the result to dest.
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::NumericVector::subset_l1_norm
virtual Real subset_l1_norm(const std::set< numeric_index_type > &indices) const
Definition: numeric_vector.C:312
libMesh::NumericVector::global_relative_compare
virtual int global_relative_compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
Definition: numeric_vector.C:170
libMesh::EIGEN_SOLVERS
Definition: enum_solver_package.h:40
libMesh::NumericVector::first_local_index
virtual numeric_index_type first_local_index() const =0
libMesh::NumericVector::compare
virtual int compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
Definition: numeric_vector.C:105