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"
46 template <
typename T>
class NumericVector;
47 template <
typename T>
class DenseVector;
48 template <
typename T>
class DenseSubVector;
49 template <
typename T>
class SparseMatrix;
50 template <
typename T>
class ShellMatrix;
51 enum SolverPackage : int;
66 class NumericVector :
public ReferenceCountedObject<NumericVector<T>>,
79 const ParallelType ptype = AUTOMATIC);
91 const ParallelType ptype = AUTOMATIC);
104 const ParallelType ptype = AUTOMATIC);
118 const std::vector<numeric_index_type> & ghost,
119 const ParallelType ptype = AUTOMATIC);
128 virtual NumericVector<T> &
operator= (
const NumericVector<T> & v) = 0;
150 static std::unique_ptr<NumericVector<T>>
151 build(
const Parallel::Communicator & comm,
185 virtual void close () = 0;
190 virtual void clear ();
195 virtual void zero () = 0;
204 virtual std::unique_ptr<NumericVector<T>>
zero_clone ()
const = 0;
213 virtual std::unique_ptr<NumericVector<T>>
clone ()
const = 0;
227 const bool fast =
false,
228 const ParallelType ptype = AUTOMATIC) = 0;
238 const bool fast =
false,
239 const ParallelType ptype = AUTOMATIC) = 0;
252 const std::vector<numeric_index_type> & ghost,
253 const bool fast =
false,
254 const ParallelType ptype = AUTOMATIC) = 0;
263 const bool fast =
false) = 0;
302 virtual T
sum()
const = 0;
425 virtual void get(
const std::vector<numeric_index_type> & index,
435 void get(
const std::vector<numeric_index_type> & index,
436 std::vector<T> & values)
const;
528 virtual void add (
const T s) = 0;
558 const std::vector<numeric_index_type> & dof_indices);
570 const std::vector<numeric_index_type> & dof_indices);
580 const std::vector<numeric_index_type> & dof_indices);
590 const std::vector<numeric_index_type> & dof_indices);
628 virtual void insert (
const T * v,
629 const std::vector<numeric_index_type> & dof_indices);
638 const std::vector<numeric_index_type> & dof_indices);
648 const std::vector<numeric_index_type> & dof_indices);
658 const std::vector<numeric_index_type> & dof_indices);
665 virtual void scale (
const T factor) = 0;
670 virtual void abs() = 0;
688 virtual void localize (std::vector<T> & v_local)
const = 0;
704 const std::vector<numeric_index_type> & send_list)
const = 0;
738 virtual void localize (std::vector<T> & v_local,
739 const std::vector<numeric_index_type> & indices)
const = 0;
750 const std::vector<numeric_index_type> & send_list) = 0;
847 friend std::ostream & operator << (std::ostream & os, const NumericVector<T> & v)
861 libmesh_not_implemented();
872 const std::vector<numeric_index_type> & rows)
const
874 libmesh_not_implemented();
936 template <
typename T>
939 const ParallelType ptype) :
940 ParallelObject(comm_in),
949 template <
typename T>
953 const ParallelType ptype) :
954 ParallelObject(comm_in),
959 libmesh_not_implemented();
965 template <
typename T>
970 const ParallelType ptype) :
971 ParallelObject(comm_in),
976 libmesh_not_implemented();
982 template <
typename T>
987 const std::vector<numeric_index_type> & ,
988 const ParallelType ptype) :
989 ParallelObject(comm_in),
994 libmesh_not_implemented();
1000 template <
typename T>
1010 template <
typename T>
1015 const std::size_t num = index.size();
1016 for (std::size_t i=0; i<num; i++)
1018 values[i] = (*this)(index[i]);
1024 template <
typename T>
1027 std::vector<T> & values)
const
1029 const std::size_t num = index.size();
1034 this->
get(index, values.data());
1039 template <
typename T>
1042 const std::vector<numeric_index_type> & dof_indices)
1044 libmesh_assert(v.size() == dof_indices.size());
1046 this->add_vector(v.data(), dof_indices);
1051 template <
typename T>
1054 const std::vector<numeric_index_type> & dof_indices)
1056 libmesh_assert(v.
size() == dof_indices.size());
1058 this->add_vector(&v(0), dof_indices);
1063 template <
typename T>
1066 const std::vector<numeric_index_type> & dof_indices)
1068 libmesh_assert(v.size() == dof_indices.size());
1070 this->insert(v.data(), dof_indices);
1075 template <
typename T>
1078 const std::vector<numeric_index_type> & dof_indices)
1080 libmesh_assert(v.
size() == dof_indices.size());
1082 this->insert(&v(0), dof_indices);
1087 template <
typename T>
1090 const std::vector<numeric_index_type> & dof_indices)
1092 libmesh_assert(v.
size() == dof_indices.size());
1094 this->insert(&v(0), dof_indices);
1106 os <<
"Size\tglobal = " << this->size()
1107 <<
"\t\tlocal = " << this->local_size() << std::endl;
1110 os <<
"#\tReal part\t\tImaginary part" << std::endl;
1111 for (
auto i : index_range(*
this))
1113 << (*this)(i).
real() <<
"\t\t"
1114 << (*this)(i).
imag() << std::endl;
1119 template <
typename T>
1124 os <<
"Size\tglobal = " << this->size()
1125 <<
"\t\tlocal = " << this->local_size() << std::endl;
1127 os <<
"#\tValue" << std::endl;
1128 for (
auto i : index_range(*
this))
1129 os << i <<
"\t" << (*this)(i) << std::endl;
1140 std::vector<Complex> v(this->size());
1144 if (this->processor_id())
1147 os <<
"Size\tglobal = " << this->size() << std::endl;
1148 os <<
"#\tReal part\t\tImaginary part" << std::endl;
1149 for (
auto i : make_range(v.size()))
1151 << v[i].real() <<
"\t\t"
1152 << v[i].imag() << std::endl;
1156 template <
typename T>
1162 std::vector<T> v(this->size());
1166 if (this->processor_id())
1169 os <<
"Size\tglobal = " << this->size() << std::endl;
1170 os <<
"#\tValue" << std::endl;
1171 for (
auto i : make_range(v.size()))
1172 os << i <<
"\t" << v[i] << std::endl;
1177 template <
typename T>
1181 std::swap(_is_closed, v._is_closed);
1183 std::swap(_type, v._type);
1191 #ifdef LIBMESH_DEFAULT_QUADRUPLE_PRECISION
1192 namespace boost {
namespace multiprecision {
namespace detail {
1193 template <
typename T,
typename To>
1194 struct is_lossy_conversion<libMesh::NumericVector<T>, To> {
1196 static const bool value = type::value;
1202 #endif // LIBMESH_NUMERIC_VECTOR_H
virtual bool closed() const
检查向量是否已经关闭并准备好进行计算。
virtual void insert(const T *v, const std::vector< numeric_index_type > &dof_indices)
将 v 的条目插入到 *this 中,位置由 dof_indices 指定。请注意,此方法的库实现是线程安全的。
std::mutex _numeric_vector_mutex
用于执行线程安全操作的互斥锁。
virtual bool empty() const overridefinal
virtual int compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
比较 this 与 other_vector 的等效性,(在给定的 threshold 内) 如果等效则返回 -1 ,或者返回第一个索引,其中 abs(a[i]-b[i]) 超过阈值。 ...
boost::multiprecision::float128 real(const boost::multiprecision::float128 in)
virtual Real subset_linfty_norm(const std::set< numeric_index_type > &indices) const
获取指定条目的向量的最大绝对值,即指定条目的 -范数。
virtual void create_subvector(NumericVector< T > &subvector, const std::vector< numeric_index_type > &rows) const
使用 rows 中的索引从该向量中填充 subvector。类似于 SparseMatrix 类的 create_submatrix() 函数, 当前仅对 PetscVectors 实现了此功能。 ...
static constexpr Real TOLERANCE
virtual T sum() const =0
获取向量中所有值的总和。
virtual void add_vector_transpose(const NumericVector< T > &v, const SparseMatrix< T > &A)=0
计算 , 即将矩阵 A 的转置与 NumericVector v 的乘积添加到 this。
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)
计算 ,其中 v 是一个指针, 每个 dof_indices[i] 指定了要添加的值 v[i] 的位置。 这应该在子类中进行重写以提高效率。请注意,此方法的库实现是线程安全的。 ...
virtual T el(const numeric_index_type i) const
获取向量的第 i 个条目。
virtual NumericVector< T > & operator=(const NumericVector< T > &v)=0
一个复制赋值运算符
virtual bool empty() const overridefinal
ParallelType & type()
获取向量的类型。
bool readable() const
检查该向量是否能够用于全局操作。
virtual std::unique_ptr< NumericVector< T > > clone() const =0
返回一个包装了此向量副本的智能指针。
bool _is_initialized
在调用 init() 后设置为 true。
virtual bool initialized() const
检查向量是否已经初始化。
virtual T dot(const NumericVector< T > &v) const =0
计算 ,即 (*this) 与向量 v 的点积。
virtual void init(const numeric_index_type n, const numeric_index_type n_local, const bool fast=false, const ParallelType ptype=AUTOMATIC)=0
更改向量的维度为 n 。如果可能的话 ,该向量的保留内存保持不变。 如果 n==0 ,所有内存都将被释放。因此,如果要调整向量的大小并释放不需要的内存, 必须首先调用 init(0) ,然后调用 ini...
uint8_t processor_id_type
virtual void zero()=0
将所有条目设置为零。等同于 v = 0,但更明显且更快。
定义了一个用于有限元计算的稠密子向量。 在将元素载荷向量累加到全局向量之前存储这些载荷向量时特别有用,尤其是在存在方程组的情况下。 所有重写的虚拟函数在 dense_vector_base.h 中有文档说明。
这是一个通用的稀疏矩阵类。该类包含了必须在派生类中覆盖的纯虚拟成员。 使用一个公共的基类允许从不同的求解器包中以不同的格式统一访问稀疏矩阵。
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
在处理器 proc_id 上创建全局向量的本地副本。 默认情况下,数据发送到处理器 0。此方法对于从一个处理器输出数据非常有用。
virtual Real l2_norm() const =0
获取向量的 -范数,即条目平方和的平方根。
virtual NumericVector< T > & operator+=(const NumericVector< T > &v)=0
将向量加上 v , 。等价于 u.add(1, v)。
virtual void print_matlab(const std::string &="") const
以Matlab的稀疏矩阵格式打印向量内容。可选择将向量打印到名为 name 的文件中。 如果未指定 name,则内容将被打印到屏幕上。
virtual ~NumericVector()=default
虽然这个类不管理任何内存,但派生类可能会管理内存,用户可能会通过指向这个基类的指针删除内存。
dof_id_type numeric_index_type
bool _is_initialized
Flag that tells if init() has been called.
NumericVector(const Parallel::Communicator &comm_in, const ParallelType ptype=AUTOMATIC)
虚拟构造函数。维度=0。
NumericVector< T > & operator/=(const T a)
将向量缩放为 1/a , 。等价于 u.scale(1.
NumericVector< T > & operator*=(const T a)
将向量缩放为 a , 。等价于 u.scale(a)。
virtual Real subset_l1_norm(const std::set< numeric_index_type > &indices) const
获取指定条目的向量的 -范数,即指定条目的绝对值之和。
virtual int local_relative_compare(const NumericVector< T > &other_vector, const Real threshold=TOLERANCE) const
比较该向量与另一个向量的局部相对差异。
virtual void print_global(std::ostream &os=libMesh::out) const
打印全局向量的内容,默认输出到 libMesh::out 流。
virtual Real max() const =0
获取向量中的最大值,或者在复数情况下获取最大的实部。
virtual Real min() const =0
获取向量中的最小值,或者在复数情况下获取最小的实部。
virtual Real l1_norm() const =0
获取向量的 -范数,即条目的绝对值之和。
virtual void close()=0
调用 NumericVector 的内部组装函数,确保值在处理器之间一致。
static std::unique_ptr< NumericVector< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
构建一个 NumericVector 对象。
virtual unsigned int size() const overridefinal
virtual void pointwise_mult(const NumericVector< T > &vec1, const NumericVector< T > &vec2)=0
比较该向量与另一个向量的全局相对差异。
virtual void abs()=0
设置 ,对向量中的每个条目进行绝对值操作。
virtual numeric_index_type first_local_index() const =0
获取实际存储在该处理器上的第一个向量元素的索引。
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual NumericVector< T > & operator-=(const NumericVector< T > &v)=0
将 v 从 *this 减去, 。等价于 u.add(-1, v)。
virtual void swap(NumericVector< T > &v)
交换该向量的内容与向量 v 的内容。子类应提供足够的间接性以使此操作成为 O(1) 的头部交换操作。
virtual void get(const std::vector< numeric_index_type > &index, T *values) const
一次访问多个组件。 values 将 *不会* 重新分配空间;它应该已经具有足够的空间。 默认实现对每个索引调用 operator() ,但某些实现可能在此处提供更快的方法。 ...
virtual numeric_index_type local_size() const =0
获取向量的本地大小,即 index_stop - index_start。
virtual std::size_t max_allowed_id() const =0
返回 NumericVector 可以包含的最大条目数(在所有处理器上)。
ParallelType type() const
获取向量的类型。
virtual void print(std::ostream &os=libMesh::out) const
打印本地向量的内容,默认输出到 libMesh::out 流。
bool initialized()
Checks that library initialization has been done.
virtual void clear()
将 NumericVector<T> 恢复到原始状态。
virtual void set(const numeric_index_type i, const T value)=0
设置 v(i) = value 。 请注意,此方法的库实现是线程安全的, 例如,将在写入向量之前锁定 _numeric_vector_mutex 。
定义用于有限元计算的稠密向量类。该类基本上是为了补充 DenseMatrix 类而设计的。 它相对于 std::vector 具有额外的功能,使其在有限元中特别有用,特别是对于方程组。 所有重写的虚拟函...
boost::multiprecision::float128 imag(const boost::multiprecision::float128)
virtual void reciprocal()=0
计算每个向量条目的分量倒数, 。
通用的Shell矩阵,即一个仅定义其对向量的作用的矩阵。此类包含必须在派生类中重写的纯虚拟成员。
virtual void add(const numeric_index_type i, const T value)=0
将 value 添加到由 i 指定的向量条目。 请注意,此方法的库实现是线程安全的, 例如,将在写入向量之前锁定 _numeric_vector_mutex 。
virtual T operator()(const numeric_index_type i) const =0
获取向量的第 i 个条目的副本。
virtual unsigned int size() const overridefinal
virtual numeric_index_type last_local_index() const =0
获取实际存储在该处理器上的最后一个向量元素的索引+1。
bool compatible(const NumericVector< T > &v) const
检查该向量和向量 v 是否能够一起用于全局操作。
virtual Real linfty_norm() const =0
获取向量的 -范数,即向量的最大绝对值。
virtual Real subset_l2_norm(const std::set< numeric_index_type > &indices) const
获取指定条目的向量的 -范数,即指定条目平方和的平方根。
virtual void pointwise_divide(const NumericVector< T > &vec1, const NumericVector< T > &vec2)=0
计算该向量与另一个向量的逐点除法。
bool _is_closed
用于跟踪向量的值在在一些或全部处理器上进行插入或添加值操作后是否在所有处理器上保持一致的标志。
virtual void conjugate()=0
反转向量中每个条目的虚部。
Real l2_norm_diff(const NumericVector< T > &other_vec) const
获取 -范数的向量差值 , 其中 是 this。
virtual void localize(std::vector< T > &v_local) const =0
创建全局向量的副本并存储在本地向量 v_local 中。