Go to the documentation of this file.
20 #ifndef LIBMESH_NUMERIC_VECTOR_H
21 #define LIBMESH_NUMERIC_VECTOR_H
24 #include "libmesh/libmesh_common.h"
25 #include "libmesh/enum_parallel_type.h"
26 #include "libmesh/id_types.h"
27 #include "libmesh/int_range.h"
28 #include "libmesh/reference_counted_object.h"
29 #include "libmesh/libmesh.h"
30 #include "libmesh/parallel_object.h"
31 #include "libmesh/dense_subvector.h"
32 #include "libmesh/dense_vector.h"
34 #ifdef LIBMESH_FORWARD_DECLARE_ENUMS
40 #include "libmesh/enum_solver_package.h"
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;
147 static std::unique_ptr<NumericVector<T>>
148 build(
const Parallel::Communicator &
comm,
177 virtual void close () = 0;
182 virtual void clear ();
188 virtual void zero () = 0;
196 virtual std::unique_ptr<NumericVector<T>>
zero_clone ()
const = 0;
203 virtual std::unique_ptr<NumericVector<T>>
clone ()
const = 0;
217 const bool fast =
false,
224 const bool fast =
false,
233 const std::vector<numeric_index_type> & ghost,
234 const bool fast =
false,
242 const bool fast =
false) = 0;
273 virtual T
sum()
const = 0;
360 virtual void get(
const std::vector<numeric_index_type> & index,
369 void get(
const std::vector<numeric_index_type> & index,
370 std::vector<T> & values)
const;
450 virtual void add (
const T s) = 0;
473 const std::vector<numeric_index_type> & dof_indices);
481 const std::vector<numeric_index_type> & dof_indices);
489 const std::vector<numeric_index_type> & dof_indices);
497 const std::vector<numeric_index_type> & dof_indices);
526 virtual void insert (
const T * v,
527 const std::vector<numeric_index_type> & dof_indices);
532 void insert (
const std::vector<T> & v,
533 const std::vector<numeric_index_type> & dof_indices);
539 const std::vector<numeric_index_type> & dof_indices);
545 const std::vector<numeric_index_type> & dof_indices);
551 const std::vector<numeric_index_type> & dof_indices);
556 virtual void scale (
const T factor) = 0;
561 virtual void abs() = 0;
575 virtual void localize (std::vector<T> & v_local)
const = 0;
588 const std::vector<numeric_index_type> & send_list)
const = 0;
619 virtual void localize (std::vector<T> & v_local,
620 const std::vector<numeric_index_type> & indices)
const = 0;
628 const std::vector<numeric_index_type> & send_list) = 0;
698 libmesh_not_implemented();
708 const std::vector<numeric_index_type> &)
const
710 libmesh_not_implemented();
745 template <
typename T>
758 template <
typename T>
768 libmesh_not_implemented();
774 template <
typename T>
785 libmesh_not_implemented();
791 template <
typename T>
796 const std::vector<numeric_index_type> & ,
803 libmesh_not_implemented();
809 template <
typename T>
819 template <
typename T>
824 const std::size_t num = index.size();
825 for (std::size_t i=0; i<num; i++)
827 values[i] = (*this)(index[i]);
833 template <
typename T>
836 std::vector<T> & values)
const
838 const std::size_t num = index.size();
843 this->
get(index, values.data());
848 template <
typename T>
851 const std::vector<numeric_index_type> & dof_indices)
855 this->add_vector(v.data(), dof_indices);
860 template <
typename T>
863 const std::vector<numeric_index_type> & dof_indices)
867 this->add_vector(&v(0), dof_indices);
872 template <
typename T>
875 const std::vector<numeric_index_type> & dof_indices)
879 this->insert(v.data(), dof_indices);
884 template <
typename T>
887 const std::vector<numeric_index_type> & dof_indices)
891 this->insert(&v(0), dof_indices);
896 template <
typename T>
899 const std::vector<numeric_index_type> & dof_indices)
903 this->insert(&v(0), dof_indices);
916 os <<
"Size\tglobal = " << this->size()
917 <<
"\t\tlocal = " << this->local_size() << std::endl;
920 os <<
"#\tReal part\t\tImaginary part" << std::endl;
923 << (*this)(i).
real() <<
"\t\t"
924 << (*this)(i).
imag() << std::endl;
929 template <
typename T>
934 os <<
"Size\tglobal = " << this->size()
935 <<
"\t\tlocal = " << this->local_size() << std::endl;
937 os <<
"#\tValue" << std::endl;
939 os << i <<
"\t" << (*this)(i) << std::endl;
950 std::vector<Complex> v(this->size());
954 if (this->processor_id())
957 os <<
"Size\tglobal = " << this->size() << std::endl;
958 os <<
"#\tReal part\t\tImaginary part" << std::endl;
961 << v[i].real() <<
"\t\t"
962 << v[i].imag() << std::endl;
966 template <
typename T>
972 std::vector<T> v(this->size());
976 if (this->processor_id())
979 os <<
"Size\tglobal = " << this->size() << std::endl;
980 os <<
"#\tValue" << std::endl;
982 os << i <<
"\t" << v[i] << std::endl;
987 template <
typename T>
1001 #ifdef LIBMESH_DEFAULT_QUADRUPLE_PRECISION
1002 namespace boost {
namespace multiprecision {
namespace detail {
1003 template <
typename T,
typename To>
1004 struct is_lossy_conversion<
libMesh::NumericVector<T>, To> {
1012 #endif // LIBMESH_NUMERIC_VECTOR_H
Generic shell matrix, i.e.
virtual void zero()=0
Set all entries to zero.
SolverPackage
Defines an enum for various linear solver packages.
virtual void add(const numeric_index_type i, const T value)=0
Adds value to each entry of the vector.
virtual T operator()(const numeric_index_type i) const =0
virtual void reciprocal()=0
Computes the component-wise reciprocal, .
virtual void print_matlab(const std::string &="") const
Print the contents of the vector in Matlab's sparse matrix format.
virtual void swap(NumericVector< T > &v)
Swaps the contents of this with v.
virtual bool empty() const override
virtual Real subset_linfty_norm(const std::set< numeric_index_type > &indices) const
virtual numeric_index_type last_local_index() const =0
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...
The libMesh namespace provides an interface to certain functionality in the library.
virtual NumericVector< T > & operator+=(const NumericVector< T > &v)=0
Adds v to *this, .
virtual void close()=0
Calls the NumericVector's internal assembly routines, ensuring that the values are consistent across ...
virtual void init(const numeric_index_type n, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC)=0
Change the dimension of the vector to n.
static const Real TOLERANCE
const Parallel::Communicator & comm() const
SolverPackage default_solver_package()
virtual std::unique_ptr< NumericVector< T > > clone() const =0
bool _is_initialized
true once init() has been called.
virtual Real max() const =0
bool _is_closed
Flag which tracks whether the vector's values are consistent on all processors after insertion or add...
NumericVector< T > & operator*=(const T a)
Scales the vector by a, .
virtual ~NumericVector()=default
While this class doesn't manage any memory, the derived class might and users may be deleting through...
virtual void pointwise_mult(const NumericVector< T > &vec1, const NumericVector< T > &vec2)=0
Computes (summation not implied) i.e.
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.
virtual bool closed() const
NumericVector(const Parallel::Communicator &comm_in, const ParallelType ptype=AUTOMATIC)
Dummy-Constructor.
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...
virtual numeric_index_type size() const =0
virtual Real min() const =0
Provides a uniform interface to vector storage schemes for different linear algebra libraries.
ParallelType _type
Type of vector.
virtual NumericVector< T > & operator-=(const NumericVector< T > &v)=0
Subtracts v from *this, .
The IntRange templated class is intended to make it easy to loop over integers which are indices of a...
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].
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.
virtual Real subset_l2_norm(const std::set< numeric_index_type > &indices) const
virtual unsigned int size() const override
virtual bool initialized() const
virtual T el(const numeric_index_type i) const
uint8_t processor_id_type
virtual int local_relative_compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
virtual bool empty() const override
virtual void add_vector_transpose(const NumericVector< T > &v, const SparseMatrix< T > &A)=0
Computes , i.e.
virtual numeric_index_type local_size() const =0
virtual NumericVector< T > & operator=(const NumericVector< T > &v)=0
This looks like a copy assignment operator, but note that, unlike normal copy assignment operators,...
static PetscErrorCode Mat * A
virtual void print_global(std::ostream &os=libMesh::out) const
Prints the global contents of the vector, by default to libMesh::out.
virtual unsigned int size() const override
virtual Real l1_norm() const =0
virtual void clear()
Restores the NumericVector<T> to a pristine state.
dof_id_type numeric_index_type
Defines a dense subvector for use in finite element computations.
virtual Real linfty_norm() const =0
bool initialized()
Checks that library initialization has been done.
void swap(Iterator &lhs, Iterator &rhs)
swap, used to implement op=
const Elem & get(const ElemType type_in)
friend std::ostream & operator<<(std::ostream &os, const NumericVector< T > &v)
Same as above but allows you to use stream syntax.
bool _is_initialized
Flag that tells if init() has been called.
virtual void set(const numeric_index_type i, const T value)=0
Sets v(i) = value.
virtual void localize_to_one(std::vector< T > &v_local, const processor_id_type proc_id=0) const =0
Creates a local copy of the global vector in v_local only on processor proc_id.
virtual Real l2_norm() const =0
virtual void scale(const T factor)=0
Scale each element of the vector by the given factor.
virtual void create_subvector(NumericVector< T > &, const std::vector< numeric_index_type > &) const
Fills in subvector from this vector using the indices in rows.
virtual void print(std::ostream &os=libMesh::out) const
Prints the local contents of the vector, by default to libMesh::out.
ParallelType
Defines an enum for parallel data structure types.
virtual void get(const std::vector< numeric_index_type > &index, T *values) const
Access multiple components at once.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
An object whose state is distributed along a set of processors.
virtual void conjugate()=0
Negates the imaginary component of each entry in the vector.
virtual Real subset_l1_norm(const std::set< numeric_index_type > &indices) const
NumericVector< T > & operator/=(const T a)
Scales the vector by 1/a, .
virtual int global_relative_compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
virtual std::unique_ptr< NumericVector< T > > zero_clone() const =0
ParallelType type() const
boost::multiprecision::float128 imag(const boost::multiprecision::float128)
void ErrorVector unsigned int
virtual T dot(const NumericVector< T > &v) const =0
Defines a dense vector for use in Finite Element-type computations.
virtual void abs()=0
Sets for each entry in the vector.
virtual numeric_index_type first_local_index() const =0
virtual int compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
boost::multiprecision::float128 real(const boost::multiprecision::float128 in)