47 std::unique_ptr<NumericVector<T>>
51 switch (solver_package)
54 #ifdef LIBMESH_HAVE_LASPACK 56 return libmesh_make_unique<LaspackVector<T>>(comm,
AUTOMATIC);
59 #ifdef LIBMESH_HAVE_PETSC 61 return libmesh_make_unique<PetscVector<T>>(comm,
AUTOMATIC);
64 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA 66 return libmesh_make_unique<EpetraVector<T>>(comm,
AUTOMATIC);
69 #ifdef LIBMESH_HAVE_EIGEN 71 return libmesh_make_unique<EigenSparseVector<T>>(comm,
AUTOMATIC);
75 return libmesh_make_unique<DistributedVector<T>>(comm,
AUTOMATIC);
83 const std::vector<numeric_index_type> & dof_indices)
86 this->
set (dof_indices[i], v[i]);
93 const std::vector<numeric_index_type> & dof_indices)
95 libmesh_assert_equal_to (V.
size(), dof_indices.size());
98 this->
set (dof_indices[i], V(i));
103 template <
typename T>
105 const Real threshold)
const 109 libmesh_assert_equal_to (this->first_local_index(), other_vector.
first_local_index());
110 libmesh_assert_equal_to (this->last_local_index(), other_vector.
last_local_index());
117 if (
std::abs((*
this)(i) - other_vector(i)) > threshold)
118 first_different_i = i;
123 && i<last_local_index());
126 this->comm().min(first_different_i);
131 return first_different_i;
135 template <
typename T>
137 const Real threshold)
const 141 libmesh_assert_equal_to (this->first_local_index(), other_vector.
first_local_index());
142 libmesh_assert_equal_to (this->last_local_index(), other_vector.
last_local_index());
149 if (
std::abs((*
this)(i) - other_vector(i)) > threshold *
151 first_different_i = i;
156 && i<last_local_index());
159 this->comm().min(first_different_i);
164 return first_different_i;
168 template <
typename T>
170 const Real threshold)
const 174 libmesh_assert_equal_to (this->first_local_index(), other_vector.
first_local_index());
175 libmesh_assert_equal_to (this->last_local_index(), other_vector.
last_local_index());
180 const Real my_norm = this->linfty_norm();
182 const Real abs_threshold =
std::max(my_norm, other_norm) * threshold;
186 if (
std::abs((*
this)(i) - other_vector(i) ) > abs_threshold)
187 first_different_i = i;
192 && i<last_local_index());
195 this->comm().min(first_different_i);
200 return first_different_i;
317 for (
const auto & index : indices)
320 this->comm().sum(norm);
332 for (
const auto & index : indices)
335 this->comm().sum(norm);
337 return std::sqrt(norm);
347 for (
const auto & index : indices)
354 this->comm().max(norm);
361 template <
typename T>
363 const std::vector<numeric_index_type> & dof_indices)
365 for (std::size_t i=0, n = dof_indices.size(); i != n; i++)
366 this->add (dof_indices[i], v[i]);
371 template <
typename T>
373 const std::vector<numeric_index_type> & dof_indices)
375 const std::size_t n = dof_indices.size();
376 libmesh_assert_equal_to(v.
size(), n);
378 this->add (dof_indices[i], v(i));
383 template <
typename T>
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices)
virtual Real subset_l2_norm(const std::set< numeric_index_type > &indices) const
virtual bool initialized() const
virtual numeric_index_type size() const =0
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
virtual Real subset_linfty_norm(const std::set< numeric_index_type > &indices) const
long double max(long double a, double b)
virtual Real subset_l1_norm(const std::set< numeric_index_type > &indices) const
virtual void vector_mult_add(NumericVector< T > &dest, const NumericVector< T > &arg) const =0
dof_id_type numeric_index_type
virtual int global_relative_compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
virtual numeric_index_type first_local_index() const =0
virtual int local_relative_compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual int compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
virtual numeric_index_type last_local_index() const =0
virtual Real linfty_norm() const =0