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"
35 #include "libmesh/enum_solver_package.h"
36 #include "libmesh/int_range.h"
48 std::unique_ptr<NumericVector<T>>
52 switch (solver_package)
55 #ifdef LIBMESH_HAVE_LASPACK
57 return libmesh_make_unique<LaspackVector<T>>(comm,
AUTOMATIC);
60 #ifdef LIBMESH_HAVE_PETSC
62 return libmesh_make_unique<PetscVector<T>>(comm,
AUTOMATIC);
65 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA
67 return libmesh_make_unique<EpetraVector<T>>(comm,
AUTOMATIC);
70 #ifdef LIBMESH_HAVE_EIGEN
72 return libmesh_make_unique<EigenSparseVector<T>>(comm,
AUTOMATIC);
76 return libmesh_make_unique<DistributedVector<T>>(comm,
AUTOMATIC);
84 const std::vector<numeric_index_type> & dof_indices)
87 this->set (dof_indices[i], v[i]);
94 const std::vector<numeric_index_type> & dof_indices)
96 libmesh_assert_equal_to (V.
size(), dof_indices.size());
99 this->set (dof_indices[i], V(i));
104 template <
typename T>
106 const Real threshold)
const
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());
113 int first_different_i = std::numeric_limits<int>::max();
118 if (
std::abs((*
this)(i) - other_vector(i)) > threshold)
119 first_different_i = i;
123 while (first_different_i==std::numeric_limits<int>::max()
124 && i<last_local_index());
127 this->comm().min(first_different_i);
129 if (first_different_i == std::numeric_limits<int>::max())
132 return first_different_i;
136 template <
typename T>
138 const Real threshold)
const
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());
145 int first_different_i = std::numeric_limits<int>::max();
150 if (
std::abs((*
this)(i) - other_vector(i)) > threshold *
152 first_different_i = i;
156 while (first_different_i==std::numeric_limits<int>::max()
157 && i<last_local_index());
160 this->comm().min(first_different_i);
162 if (first_different_i == std::numeric_limits<int>::max())
165 return first_different_i;
169 template <
typename T>
171 const Real threshold)
const
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());
178 int first_different_i = std::numeric_limits<int>::max();
181 const Real my_norm = this->linfty_norm();
183 const Real abs_threshold = std::max(my_norm, other_norm) * threshold;
187 if (
std::abs((*
this)(i) - other_vector(i) ) > abs_threshold)
188 first_different_i = i;
192 while (first_different_i==std::numeric_limits<int>::max()
193 && i<last_local_index());
196 this->comm().min(first_different_i);
198 if (first_different_i == std::numeric_limits<int>::max())
201 return first_different_i;
318 for (
const auto & index : indices)
321 this->comm().sum(
norm);
333 for (
const auto & index : indices)
336 this->comm().sum(
norm);
348 for (
const auto & index : indices)
355 this->comm().max(
norm);
362 template <
typename T>
364 const std::vector<numeric_index_type> & dof_indices)
367 this->add (dof_indices[i], v[i]);
372 template <
typename T>
374 const std::vector<numeric_index_type> & dof_indices)
376 const std::size_t n = dof_indices.size();
377 libmesh_assert_equal_to(v.
size(), n);
379 this->add (dof_indices[i], v(i));
384 template <
typename T>