libmesh解析
本工作只是尝试解析原libmesh的代码,供学习使用
 全部  命名空间 文件 函数 变量 类型定义 枚举 枚举值 友元 
sparse_matrix.h
浏览该文件的文档.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2023 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_SPARSE_MATRIX_H
21 #define LIBMESH_SPARSE_MATRIX_H
22 
23 
24 // Local includes
25 #include "libmesh/libmesh.h"
26 #include "libmesh/libmesh_common.h"
27 #include "libmesh/reference_counted_object.h"
28 #include "libmesh/id_types.h"
29 #include "libmesh/parallel_object.h"
30 #include "libmesh/enum_parallel_type.h" // PARALLEL
31 #include "libmesh/enum_matrix_build_type.h" // AUTOMATIC
32 
33 // C++ includes
34 #include <cstddef>
35 #include <iomanip>
36 #include <vector>
37 #include <memory>
38 
39 namespace libMesh
40 {
41 
42 // forward declarations
43 template <typename T> class SparseMatrix;
44 template <typename T> class DenseMatrix;
45 class DofMap;
46 namespace SparsityPattern {
47  class Build;
48  class Graph;
49 }
50 template <typename T> class NumericVector;
51 
52 // 在定义这个模板辅助函数之前,必须先声明它。
53 template <typename T>
54 std::ostream & operator << (std::ostream & os, const SparseMatrix<T> & m);
55 
56 
64 template <typename T>
65 class SparseMatrix : public ReferenceCountedObject<SparseMatrix<T>>,
66  public ParallelObject
67 {
68 public:
79  explicit
80  SparseMatrix (const Parallel::Communicator & comm);
81 
85  SparseMatrix (SparseMatrix &&) = default;
86  SparseMatrix (const SparseMatrix &) = default;
87  SparseMatrix & operator= (const SparseMatrix &) = default;
88  SparseMatrix & operator= (SparseMatrix &&) = default;
89  virtual ~SparseMatrix () = default;
90 
94  static std::unique_ptr<SparseMatrix<T>>
95  build(const Parallel::Communicator & comm,
96  const SolverPackage solver_package = libMesh::default_solver_package(),
97  const MatrixBuildType matrix_build_type = MatrixBuildType::AUTOMATIC);
98 
102  virtual bool initialized() const { return _is_initialized; }
103 
110  void attach_dof_map (const DofMap & dof_map);
111 
119 
126  virtual bool need_full_sparsity_pattern() const
127  { return false; }
128 
134 
146  virtual void init (const numeric_index_type m,
147  const numeric_index_type n,
148  const numeric_index_type m_l,
149  const numeric_index_type n_l,
150  const numeric_index_type nnz=30,
151  const numeric_index_type noz=10,
152  const numeric_index_type blocksize=1) = 0;
153 
161  virtual void init (ParallelType type = PARALLEL) = 0;
162 
166  virtual void clear () = 0;
167 
171  virtual void zero () = 0;
172 
178  virtual std::unique_ptr<SparseMatrix<T>> zero_clone () const = 0;
179 
185  virtual std::unique_ptr<SparseMatrix<T>> clone () const = 0;
186 
190  virtual void zero_rows (std::vector<numeric_index_type> & rows, T diag_value = 0.0);
191 
195  virtual void close () = 0;
196 
202  virtual void flush () { close(); }
203 
207  virtual numeric_index_type m () const = 0;
208 
212  virtual numeric_index_type local_m () const { return row_stop() - row_start(); }
213 
217  virtual numeric_index_type n () const = 0;
218 
222  virtual numeric_index_type row_start () const = 0;
223 
227  virtual numeric_index_type row_stop () const = 0;
228 
237  virtual void set (const numeric_index_type i,
238  const numeric_index_type j,
239  const T value) = 0;
240 
249  virtual void add (const numeric_index_type i,
250  const numeric_index_type j,
251  const T value) = 0;
252 
256  virtual void add_matrix (const DenseMatrix<T> & dm,
257  const std::vector<numeric_index_type> & rows,
258  const std::vector<numeric_index_type> & cols) = 0;
259 
263  virtual void add_matrix (const DenseMatrix<T> & dm,
264  const std::vector<numeric_index_type> & dof_indices) = 0;
265 
270  virtual void add_block_matrix (const DenseMatrix<T> & dm,
271  const std::vector<numeric_index_type> & brows,
272  const std::vector<numeric_index_type> & bcols);
273 
278  virtual void add_block_matrix (const DenseMatrix<T> & dm,
279  const std::vector<numeric_index_type> & dof_indices)
280  { this->add_block_matrix (dm, dof_indices, dof_indices); }
281 
285  virtual void add (const T a, const SparseMatrix<T> & X) = 0;
286 
290  virtual void matrix_matrix_mult (SparseMatrix<T> & /*X*/, SparseMatrix<T> & /*Y*/, bool /*reuse*/)
291  { libmesh_not_implemented(); }
292 
297  virtual void add_sparse_matrix (const SparseMatrix<T> & /*spm*/,
298  const std::map<numeric_index_type, numeric_index_type> & /*rows*/,
299  const std::map<numeric_index_type, numeric_index_type> & /*cols*/,
300  const T /*scalar*/)
301  { libmesh_not_implemented(); }
302 
312  virtual T operator () (const numeric_index_type i,
313  const numeric_index_type j) const = 0;
314 
315 
324  virtual Real l1_norm () const = 0;
325 
335  virtual Real linfty_norm () const = 0;
336 
340  virtual bool closed() const = 0;
341 
345  void print(std::ostream & os=libMesh::out, const bool sparse=false) const;
346 
365  friend std::ostream & operator << <>(std::ostream & os, const SparseMatrix<T> & m);
366 
370  virtual void print_personal(std::ostream & os=libMesh::out) const = 0;
371 
376  virtual void print_matlab(const std::string & /*name*/ = "") const
377  {
378  libmesh_not_implemented();
379  }
380 
388  virtual void create_submatrix(SparseMatrix<T> & submatrix,
389  const std::vector<numeric_index_type> & rows,
390  const std::vector<numeric_index_type> & cols) const
391  {
392  this->_get_submatrix(submatrix,
393  rows,
394  cols,
395  false); // false means DO NOT REUSE submatrix
396  }
397 
404  virtual void create_submatrix_nosort(SparseMatrix<T> & /*submatrix*/,
405  const std::vector<numeric_index_type> & /*rows*/,
406  const std::vector<numeric_index_type> & /*cols*/) const
407  {
408  libmesh_not_implemented();
409  }
410 
415  virtual void reinit_submatrix(SparseMatrix<T> & submatrix,
416  const std::vector<numeric_index_type> & rows,
417  const std::vector<numeric_index_type> & cols) const
418  {
419  this->_get_submatrix(submatrix,
420  rows,
421  cols,
422  true); // true 表示重用 submatrix
423  }
424 
433  void vector_mult (NumericVector<T> & dest,
434  const NumericVector<T> & arg) const;
435 
439  void vector_mult_add (NumericVector<T> & dest,
440  const NumericVector<T> & arg) const;
441 
449  virtual void get_diagonal (NumericVector<T> & dest) const = 0;
450 
458  virtual void get_transpose (SparseMatrix<T> & dest) const = 0;
459 
469  virtual void get_row(numeric_index_type /*i*/,
470  std::vector<numeric_index_type> & /*indices*/,
471  std::vector<T> & /*values*/) const
472  {
473  libmesh_not_implemented();
474  }
475 
476 protected:
477 
483  virtual void _get_submatrix(SparseMatrix<T> & /*submatrix*/,
484  const std::vector<numeric_index_type> & /*rows*/,
485  const std::vector<numeric_index_type> & /*cols*/,
486  const bool /*reuse_submatrix*/) const
487  {
488  libmesh_not_implemented();
489  }
490 
494  DofMap const * _dof_map;
495 
501 
506 };
507 
508 
509 
510 //-----------------------------------------------------------------------
511 // SparseMatrix 内联成员
512 
513 // 对于 SGI MIPSpro,此实现必须出现在 print() 成员的完全特化之后。
514 //
515 // 通常,更容易在头文件中定义这些 friend 函数。
516 template <typename T>
517 std::ostream & operator << (std::ostream & os, const SparseMatrix<T> & m)
518 {
519  m.print(os);
520  return os;
521 }
522 
523 
524 } // namespace libMesh
525 
526 
527 #endif // LIBMESH_SPARSE_MATRIX_H
virtual void add_block_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &brows, const std::vector< numeric_index_type > &bcols)
将完整矩阵 dm 添加到 SparseMatrix。这对于在装配时添加元素矩阵很有用。 矩阵被假定为块矩阵,brow、bcol 对应于*块*行和列索引。
Definition: sparse_matrix.C:77
virtual void print_personal(std::ostream &os=libMesh::out) const =0
以个性化的风格(如果可用)将矩阵的内容打印到屏幕上。
virtual void create_submatrix_nosort(SparseMatrix< T > &, const std::vector< numeric_index_type > &, const std::vector< numeric_index_type > &) const
与上述函数类似,此函数创建由 rows 和 cols 向量中的索引定义的 submatrix。 注意:rows 和 cols 可以是无序的; 如果索引已排序,则对于效率更高的操作,请使用上面的函数; r...
SparseMatrix(const Parallel::Communicator &comm)
构造函数;初始化矩阵为空,没有任何结构,即矩阵无法使用。
Definition: sparse_matrix.C:48
This helper class can be called on multiple threads to compute the sparsity pattern (or graph) of the...
virtual ~SparseMatrix()=default
virtual void print_matlab(const std::string &="") const
以 Matlab 的稀疏矩阵格式打印矩阵的内容。可选择将矩阵打印到名为 name 的文件中。 如果未指定 name ,则将其转储到屏幕上。
void vector_mult(NumericVector< T > &dest, const NumericVector< T > &arg) const
将矩阵乘以 NumericVector arg 并将结果存储在 NumericVector dest 中。
virtual void get_row(numeric_index_type, std::vector< numeric_index_type > &, std::vector< T > &) const
从矩阵获取一行。
SparsityPattern::Build const * _sp
与此对象关联的 sparsity pattern。在需要时, 应查询其入口计数(或使用 need_full_sparsity_pattern,模式)。
void attach_sparsity_pattern(const SparsityPattern::Build &sp)
设置要使用的稀疏性模式的指针。在矩阵需要比系统中的大(或更小以提高效率)的模式, 或者在 DofMap 未计算系统稀疏性模式的情况下使用。
Definition: sparse_matrix.C:68
提供了不同线性代数库的向量存储方案的统一接口。
Definition: dof_map.h:67
void vector_mult_add(NumericVector< T > &dest, const NumericVector< T > &arg) const
将矩阵乘以 NumericVector arg ,并将结果添加到 NumericVector dest 中。
virtual void matrix_matrix_mult(SparseMatrix< T > &, SparseMatrix< T > &, bool)
计算 Y = A*X,其中 X 为矩阵。
virtual numeric_index_type local_m() const
获取此进程拥有的行数。
static std::unique_ptr< SparseMatrix< T > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package(), const MatrixBuildType matrix_build_type=MatrixBuildType::AUTOMATIC)
使用由 solver_package 指定的线性求解器包构建一个 SparseMatrix&lt;T&gt;。
virtual numeric_index_type row_stop() const =0
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value)=0
设置元素 (i,j) 为 value .
virtual void clear()=0
恢复 SparseMatrix&lt;T&gt; 到原始状态。
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value)=0
将 value 添加到元素 (i,j) .
virtual void add_sparse_matrix(const SparseMatrix< T > &, const std::map< numeric_index_type, numeric_index_type > &, const std::map< numeric_index_type, numeric_index_type > &, const T)
将 scalar* spm 添加到此矩阵的行和列中: A(rows[i], cols[j]) += scalar * spm(i,j)
virtual bool initialized() const
这是一个通用的稀疏矩阵类。该类包含了必须在派生类中覆盖的纯虚拟成员。 使用一个公共的基类允许从不同的求解器包中以不同的格式统一访问稀疏矩阵。
Definition: dof_map.h:66
SolverPackage default_solver_package()
Definition: libmesh.C:967
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:169
virtual void flush()
对于 PETSc 矩阵,此函数类似于 close,但不收缩内存。 当我们希望在 ADD_VALUES 和 INSERT_VALUES 之间切换时,这很有用。 在使用矩阵之前应该调用 close。 ...
virtual void create_submatrix(SparseMatrix< T > &submatrix, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) const
此函数创建一个名为 &quot;submatrix&quot; 的矩阵,其定义由 &quot;rows&quot; 和 &quot;cols&quot; 条目中的行和列索引给定。 目前,此操作仅对 PetscMatrix 类型定义。 注意:rows 和 cols...
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols)=0
将完整矩阵 dm 添加到 SparseMatrix。这对于在装配时添加元素矩阵很有用。
virtual void zero()=0
将所有条目设置为 0。
dof_id_type numeric_index_type
Definition: id_types.h:99
bool _is_initialized
标志,指示矩阵是否已初始化。
virtual std::unique_ptr< SparseMatrix< T > > zero_clone() const =0
virtual numeric_index_type m() const =0
virtual bool need_full_sparsity_pattern() const
virtual void zero_rows(std::vector< numeric_index_type > &rows, T diag_value=0.0)
将所有行条目设置为 0,然后将 diag_value 放在对角线条目。
virtual void reinit_submatrix(SparseMatrix< T > &submatrix, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) const
此函数与上述函数类似,但允许您重用 &quot;submatrix&quot; 的现有稀疏性模式,而不是再次分配它。 如果经常提取相同大小的子矩阵,这应该更有效。
virtual void get_diagonal(NumericVector< T > &dest) const =0
复制矩阵的对角线部分到 dest。
OStreamProxy out
virtual bool closed() const =0
void attach_dof_map(const DofMap &dof_map)
设置要使用的 DofMap 的指针。如果不使用单独的稀疏性模式, 则使用来自 DofMap 的模式。
Definition: sparse_matrix.C:58
DofMap const * _dof_map
与此对象关联的 DofMap 对象。可以查询其在处理器上的自由度计数。
virtual void close()=0
调用 SparseMatrix 的内部装配例程,确保值在处理器之间保持一致。
virtual void update_sparsity_pattern(const SparsityPattern::Graph &)
更新矩阵的稀疏性模式。当您的 SparseMatrix&lt;T&gt; 实现不需要此数据时, 只需不覆盖此方法。
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual numeric_index_type row_start() const =0
virtual void get_transpose(SparseMatrix< T > &dest) const =0
将矩阵的转置复制到 dest 中,dest 可能是 *this。
void print(std::ostream &os=libMesh::out, const bool sparse=false) const
将矩阵的内容以统一的样式打印到屏幕上,而不考虑正在使用的矩阵/求解器包。
virtual Real linfty_norm() const =0
获取矩阵的 -范数,即最大行和:
SparseMatrix & operator=(const SparseMatrix &)=default
virtual void _get_submatrix(SparseMatrix< T > &, const std::vector< numeric_index_type > &, const std::vector< numeric_index_type > &, const bool) const
创建子矩阵和 reinit_submatrix 例程的受保护实现。
virtual Real l1_norm() const =0
获取矩阵的 -范数,即最大列和:
virtual std::unique_ptr< SparseMatrix< T > > clone() const =0
定义用于有限元类型计算的密集矩阵。 用于在求和成全局矩阵之前存储单元刚度矩阵。所有被覆盖的虚函数都记录在dense_matrix_base.h中。
Definition: dof_map.h:65
virtual T operator()(const numeric_index_type i, const numeric_index_type j) const =0
virtual numeric_index_type n() const =0
virtual void add_block_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &dof_indices)
与 add_block_matrix() 相同,但假定行和列映射相同。 因此矩阵 dm 必须是方阵。
virtual void init(const numeric_index_type m, const numeric_index_type n, const numeric_index_type m_l, const numeric_index_type n_l, const numeric_index_type nnz=30, const numeric_index_type noz=10, const numeric_index_type blocksize=1)=0
使用指定的大小初始化 SparseMatrix。