21 #ifndef LIBMESH_EIGEN_SPARSE_VECTOR_H 22 #define LIBMESH_EIGEN_SPARSE_VECTOR_H 28 #ifdef LIBMESH_HAVE_EIGEN 38 template <
typename T>
class EigenSparseMatrix;
39 template <
typename T>
class EigenSparseLinearSolver;
40 template <
typename T>
class SparseMatrix;
51 class EigenSparseVector final :
public NumericVector<T>
87 const std::vector<numeric_index_type> & ghost,
95 EigenSparseVector<T> &
operator= (
const EigenSparseVector<T> & v);
111 virtual void close ()
override;
113 virtual void clear ()
override;
115 virtual void zero ()
override;
117 virtual std::unique_ptr<NumericVector<T>>
zero_clone ()
const override;
119 virtual std::unique_ptr<NumericVector<T>>
clone ()
const override;
123 const bool fast=
false,
127 const bool fast=
false,
132 const std::vector<numeric_index_type> & ghost,
133 const bool fast =
false,
137 const bool fast =
false)
override;
145 virtual Real min ()
const override;
147 virtual Real max ()
const override;
149 virtual T
sum ()
const override;
181 virtual void add (
const T s)
override;
199 virtual void scale (
const T factor)
override;
201 virtual void abs()
override;
205 virtual void localize (std::vector<T> & v_local)
const override;
210 const std::vector<numeric_index_type> & send_list)
const override;
212 virtual void localize (std::vector<T> & v_local,
213 const std::vector<numeric_index_type> & indices)
const override;
217 const std::vector<numeric_index_type> & send_list)
override;
253 template <
typename T>
264 template <
typename T>
271 this->
init(n, n,
false, ptype);
276 template <
typename T>
284 this->
init(n, n_local,
false, ptype);
289 template <
typename T>
294 const std::vector<numeric_index_type> & ghost,
298 this->
init(N, n_local, ghost,
false, ptype);
303 template <
typename T>
313 libmesh_error_msg(
"Error: EigenSparseVectors can only be used in serial!");
325 this->_is_closed =
true;
337 template <
typename T>
343 this->
init(n,n,fast,ptype);
347 template <
typename T>
351 const std::vector<numeric_index_type> & libmesh_dbg_var(ghost),
355 libmesh_assert(ghost.empty());
356 this->
init(n,n_local,fast,ptype);
372 template <
typename T>
379 this->_is_closed =
true;
385 template <
typename T>
393 this->_is_closed =
false;
399 template <
typename T>
inline 403 libmesh_assert (this->
closed());
410 template <
typename T>
415 cloned_vector->
init(*
this);
416 return std::unique_ptr<NumericVector<T>>(cloned_vector);
421 template <
typename T>
426 cloned_vector->
init(*
this,
true);
427 *cloned_vector = *
this;
428 return std::unique_ptr<NumericVector<T>>(cloned_vector);
433 template <
typename T>
444 template <
typename T>
455 template <
typename T>
466 template <
typename T>
477 template <
typename T>
482 libmesh_assert_less (i, this->size());
487 this->_is_closed =
false;
493 template <
typename T>
498 libmesh_assert_less (i, this->size());
503 this->_is_closed =
false;
509 template <
typename T>
514 libmesh_assert ( ((i >= this->first_local_index()) &&
515 (i < this->last_local_index())) );
522 template <
typename T>
530 std::swap (this->_is_closed, v._is_closed);
539 #endif // #ifdef LIBMESH_HAVE_EIGEN 540 #endif // LIBMESH_EIGEN_SPARSE_VECTOR_H
virtual Real linfty_norm() const override
virtual void localize(std::vector< T > &v_local) const override
virtual void conjugate() override
virtual void pointwise_mult(const NumericVector< T > &vec1, const NumericVector< T > &vec2) override
virtual T sum() const override
virtual void add_vector(const NumericVector< T > &v, const SparseMatrix< T > &A) override
virtual void reciprocal() override
virtual numeric_index_type size() const =0
virtual NumericVector< T > & operator-=(const NumericVector< T > &v) override
virtual Real min() const override
uint8_t processor_id_type
virtual std::unique_ptr< NumericVector< T > > clone() const override
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
virtual numeric_index_type last_local_index() const override
virtual void close() override
virtual void init(const numeric_index_type n, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC)=0
virtual void swap(NumericVector< T > &v) override
virtual std::unique_ptr< NumericVector< T > > zero_clone() const override
virtual Real l2_norm() const override
dof_id_type numeric_index_type
virtual void clear() override
virtual Real max() const override
virtual NumericVector< T > & operator/=(const NumericVector< T > &v_in) override
Eigen::Matrix< Number, Eigen::Dynamic, 1 > EigenSV
virtual numeric_index_type size() const override
const DataType & vec() const
void init(triangulateio &t)
virtual numeric_index_type local_size() const override
EigenSparseVector< T > & operator=(const EigenSparseVector< T > &v)
virtual void init(const numeric_index_type N, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC) override
virtual void add_vector_transpose(const NumericVector< T > &v, const SparseMatrix< T > &A) override
virtual void abs() override
virtual void localize_to_one(std::vector< T > &v_local, const processor_id_type proc_id=0) const override
virtual void set(const numeric_index_type i, const T value) override
virtual void scale(const T factor) override
virtual numeric_index_type first_local_index() const override
ParallelType type() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static PetscErrorCode Mat * A
virtual Real l1_norm() const override
virtual void add(const numeric_index_type i, const T value) override
void swap(Iterator &lhs, Iterator &rhs)
virtual numeric_index_type local_size() const =0
virtual void zero() override
EigenSparseVector(const Parallel::Communicator &comm_in, const ParallelType=AUTOMATIC)
virtual T operator()(const numeric_index_type i) const override
virtual ~EigenSparseVector()=default
virtual T dot(const NumericVector< T > &v) const override
virtual NumericVector< T > & operator+=(const NumericVector< T > &v) override