20 #ifndef LIBMESH_NUMERIC_VECTOR_H 21 #define LIBMESH_NUMERIC_VECTOR_H 34 #ifdef LIBMESH_FORWARD_DECLARE_ENUMS 54 template <
typename T>
class NumericVector;
55 template <
typename T>
class DenseVector;
56 template <
typename T>
class DenseSubVector;
57 template <
typename T>
class SparseMatrix;
58 template <
typename T>
class ShellMatrix;
75 class NumericVector :
public ReferenceCountedObject<NumericVector<T>>,
112 const std::vector<numeric_index_type> & ghost,
126 virtual NumericVector<T> &
operator= (
const NumericVector<T> & v) = 0;
141 static std::unique_ptr<NumericVector<T>>
142 build(
const Parallel::Communicator &
comm,
171 virtual void close () = 0;
176 virtual void clear ();
182 virtual void zero () = 0;
190 virtual std::unique_ptr<NumericVector<T>>
zero_clone ()
const = 0;
197 virtual std::unique_ptr<NumericVector<T>>
clone ()
const = 0;
211 const bool fast =
false,
218 const bool fast =
false,
227 const std::vector<numeric_index_type> & ghost,
228 const bool fast =
false,
236 const bool fast =
false) = 0;
267 virtual T
sum()
const = 0;
354 virtual void get(
const std::vector<numeric_index_type> & index,
363 void get(
const std::vector<numeric_index_type> & index,
364 std::vector<T> & values)
const;
435 virtual void add (
const T s) = 0;
458 const std::vector<numeric_index_type> & dof_indices);
466 const std::vector<numeric_index_type> & dof_indices);
474 const std::vector<numeric_index_type> & dof_indices);
482 const std::vector<numeric_index_type> & dof_indices);
511 virtual void insert (
const T * v,
512 const std::vector<numeric_index_type> & dof_indices);
517 void insert (
const std::vector<T> & v,
518 const std::vector<numeric_index_type> & dof_indices);
524 const std::vector<numeric_index_type> & dof_indices);
530 const std::vector<numeric_index_type> & dof_indices);
536 const std::vector<numeric_index_type> & dof_indices);
541 virtual void scale (
const T factor) = 0;
546 virtual void abs() = 0;
560 virtual void localize (std::vector<T> & v_local)
const = 0;
573 const std::vector<numeric_index_type> & send_list)
const = 0;
604 virtual void localize (std::vector<T> & v_local,
605 const std::vector<numeric_index_type> & indices)
const = 0;
613 const std::vector<numeric_index_type> & send_list) = 0;
670 friend std::ostream & operator << (std::ostream & os, const NumericVector<T> & v)
683 libmesh_not_implemented();
693 const std::vector<numeric_index_type> &)
const 695 libmesh_not_implemented();
730 template <
typename T>
743 template <
typename T>
753 libmesh_not_implemented();
759 template <
typename T>
770 libmesh_not_implemented();
776 template <
typename T>
781 const std::vector<numeric_index_type> & ,
788 libmesh_not_implemented();
794 template <
typename T>
804 template <
typename T>
809 const std::size_t num = index.size();
810 for (std::size_t i=0; i<num; i++)
812 values[i] = (*this)(index[i]);
818 template <
typename T>
821 std::vector<T> & values)
const 823 const std::size_t num = index.size();
828 this->
get(index, values.data());
833 template <
typename T>
836 const std::vector<numeric_index_type> & dof_indices)
838 libmesh_assert(v.size() == dof_indices.size());
840 this->add_vector(v.data(), dof_indices);
845 template <
typename T>
848 const std::vector<numeric_index_type> & dof_indices)
850 libmesh_assert(v.
size() == dof_indices.size());
852 this->add_vector(&v(0), dof_indices);
857 template <
typename T>
860 const std::vector<numeric_index_type> & dof_indices)
862 libmesh_assert(v.size() == dof_indices.size());
864 this->insert(v.data(), dof_indices);
869 template <
typename T>
872 const std::vector<numeric_index_type> & dof_indices)
874 libmesh_assert(v.
size() == dof_indices.size());
876 this->insert(&v(0), dof_indices);
881 template <
typename T>
884 const std::vector<numeric_index_type> & dof_indices)
886 libmesh_assert(v.
size() == dof_indices.size());
888 this->insert(&v(0), dof_indices);
901 os <<
"Size\tglobal = " << this->size()
902 <<
"\t\tlocal = " << this->local_size() << std::endl;
905 os <<
"#\tReal part\t\tImaginary part" << std::endl;
908 << (*
this)(i).real() <<
"\t\t" 909 << (*this)(i).imag() << std::endl;
914 template <
typename T>
919 os <<
"Size\tglobal = " << this->size()
920 <<
"\t\tlocal = " << this->local_size() << std::endl;
922 os <<
"#\tValue" << std::endl;
924 os << i <<
"\t" << (*
this)(i) << std::endl;
935 std::vector<Complex> v(this->size());
939 if (this->processor_id())
942 os <<
"Size\tglobal = " << this->size() << std::endl;
943 os <<
"#\tReal part\t\tImaginary part" << std::endl;
946 << v[i].real() <<
"\t\t" 947 << v[i].imag() << std::endl;
951 template <
typename T>
957 std::vector<T> v(this->size());
961 if (this->processor_id())
964 os <<
"Size\tglobal = " << this->size() << std::endl;
965 os <<
"#\tValue" << std::endl;
967 os << i <<
"\t" << v[i] << std::endl;
972 template <
typename T>
985 #endif // LIBMESH_NUMERIC_VECTOR_H virtual unsigned int size() const override
virtual void create_subvector(NumericVector< T > &, const std::vector< numeric_index_type > &) const
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 void print_global(std::ostream &os=libMesh::out) const
virtual T el(const numeric_index_type i) const
virtual void add_vector_transpose(const NumericVector< T > &v, const SparseMatrix< T > &A)=0
virtual bool initialized() const
virtual numeric_index_type size() const =0
virtual std::unique_ptr< NumericVector< T > > zero_clone() const =0
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
virtual NumericVector< T > & operator=(const NumericVector< T > &v)=0
uint8_t processor_id_type
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
const Parallel::Communicator & comm() const
virtual std::unique_ptr< NumericVector< T > > clone() const =0
static const Real TOLERANCE
virtual Real subset_linfty_norm(const std::set< numeric_index_type > &indices) const
virtual T dot(const NumericVector< T > &v) const =0
virtual void init(const numeric_index_type n, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC)=0
virtual Real subset_l1_norm(const std::set< numeric_index_type > &indices) const
virtual void print_matlab(const std::string &="") const
virtual bool empty() const override
SolverPackage default_solver_package()
virtual void scale(const T factor)=0
virtual void localize_to_one(std::vector< T > &v_local, const processor_id_type proc_id=0) const =0
virtual Real l2_norm() const =0
virtual NumericVector< T > & operator+=(const NumericVector< T > &v)=0
dof_id_type numeric_index_type
NumericVector(const Parallel::Communicator &comm_in, const ParallelType ptype=AUTOMATIC)
NumericVector< T > & operator/=(const T a)
NumericVector< T > & operator*=(const T a)
virtual int global_relative_compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
virtual void print(std::ostream &os=libMesh::out) const
virtual unsigned int size() const override
virtual Real max() const =0
virtual Real min() const =0
An object whose state is distributed along a set of processors.
virtual Real l1_norm() const =0
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
virtual void pointwise_mult(const NumericVector< T > &vec1, const NumericVector< T > &vec2)=0
virtual bool closed() const
virtual numeric_index_type first_local_index() const =0
ParallelType type() const
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 NumericVector< T > & operator-=(const NumericVector< T > &v)=0
virtual void swap(NumericVector< T > &v)
virtual int compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
static PetscErrorCode Mat * A
virtual bool empty() const override
void swap(Iterator &lhs, Iterator &rhs)
virtual numeric_index_type local_size() const =0
virtual void reciprocal()=0
virtual void add(const numeric_index_type i, const T value)=0
virtual T operator()(const numeric_index_type i) const =0
OStreamProxy out(std::cout)
virtual numeric_index_type last_local_index() const =0
virtual Real linfty_norm() const =0
const Elem & get(const ElemType type_in)
virtual void conjugate()=0
virtual void localize(std::vector< T > &v_local) const =0