libmesh解析
本工作只是尝试解析原libmesh的代码,供学习使用
 全部  命名空间 文件 函数 变量 类型定义 枚举 枚举值 友元 
sparse_matrix.C
浏览该文件的文档.
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 // Local Includes
21 #include "libmesh/dof_map.h"
22 #include "libmesh/dense_matrix.h"
23 #include "libmesh/diagonal_matrix.h"
24 #include "libmesh/laspack_matrix.h"
25 #include "libmesh/eigen_sparse_matrix.h"
26 #include "libmesh/parallel.h"
27 #include "libmesh/petsc_matrix.h"
28 #include "libmesh/sparse_matrix.h"
29 #include "libmesh/trilinos_epetra_matrix.h"
30 #include "libmesh/numeric_vector.h"
31 #include "libmesh/enum_solver_package.h"
32 
33 
34 // C++ includes
35 #include <memory>
36 
37 
38 namespace libMesh
39 {
40 
41 
42 //------------------------------------------------------------------
43 // SparseMatrix Methods
44 
45 
46 // Constructor
47 template <typename T>
48 SparseMatrix<T>::SparseMatrix (const Parallel::Communicator & comm_in) :
49  ParallelObject(comm_in),
50  _dof_map(nullptr),
51  _sp(nullptr),
52  _is_initialized(false)
53 {}
54 
55 
56 
57 template <typename T>
59 {
60  _dof_map = &dof_map;
61  if (!_sp)
62  _sp = dof_map.get_sparsity_pattern();
63 }
64 
65 
66 
67 template <typename T>
69 {
70  _sp = &sp;
71 }
72 
73 
74 
75 // default implementation is to fall back to non-blocked method
76 template <typename T>
78  const std::vector<numeric_index_type> & brows,
79  const std::vector<numeric_index_type> & bcols)
80 {
81  libmesh_assert_equal_to (dm.m() / brows.size(), dm.n() / bcols.size());
82 
83  const numeric_index_type blocksize = cast_int<numeric_index_type>
84  (dm.m() / brows.size());
85 
86  libmesh_assert_equal_to (dm.m()%blocksize, 0);
87  libmesh_assert_equal_to (dm.n()%blocksize, 0);
88 
89  std::vector<numeric_index_type> rows, cols;
90 
91  rows.reserve(blocksize*brows.size());
92  cols.reserve(blocksize*bcols.size());
93 
94  for (auto & row : brows)
95  {
96  numeric_index_type i = row * blocksize;
97 
98  for (unsigned int v=0; v<blocksize; v++)
99  rows.push_back(i++);
100  }
101 
102  for (auto & col : bcols)
103  {
104  numeric_index_type j = col * blocksize;
105 
106  for (unsigned int v=0; v<blocksize; v++)
107  cols.push_back(j++);
108  }
109 
110  this->add_matrix (dm, rows, cols);
111 }
112 
113 
114 
115 // Full specialization of print method for Complex datatypes
116 template <>
117 void SparseMatrix<Complex>::print(std::ostream & os, const bool sparse) const
118 {
119  // std::complex<>::operator<<() is defined, but use this form
120 
121  if (sparse)
122  {
123  libmesh_not_implemented();
124  }
125 
126  os << "Real part:" << std::endl;
127  for (auto i : make_range(this->m()))
128  {
129  for (auto j : make_range(this->n()))
130  os << std::setw(8) << (*this)(i,j).real() << " ";
131  os << std::endl;
132  }
133 
134  os << std::endl << "Imaginary part:" << std::endl;
135  for (auto i : make_range(this->m()))
136  {
137  for (auto j : make_range(this->n()))
138  os << std::setw(8) << (*this)(i,j).imag() << " ";
139  os << std::endl;
140  }
141 }
142 
143 
144 
145 
146 
147 
148 // Full specialization for Real datatypes
149 template <typename T>
150 std::unique_ptr<SparseMatrix<T>>
151 SparseMatrix<T>::build(const Parallel::Communicator & comm,
152  const SolverPackage solver_package,
153  const MatrixBuildType matrix_build_type /* = AUTOMATIC */)
154 {
155  // Avoid unused parameter warnings when no solver packages are enabled.
156  libmesh_ignore(comm);
157 
158  if (matrix_build_type == MatrixBuildType::DIAGONAL)
159  return std::make_unique<DiagonalMatrix<T>>(comm);
160 
161  // Build the appropriate vector
162  switch (solver_package)
163  {
164 
165 #ifdef LIBMESH_HAVE_LASPACK
166  case LASPACK_SOLVERS:
167  return std::make_unique<LaspackMatrix<T>>(comm);
168 #endif
169 
170 
171 #ifdef LIBMESH_HAVE_PETSC
172  case PETSC_SOLVERS:
173  return std::make_unique<PetscMatrix<T>>(comm);
174 #endif
175 
176 
177 #ifdef LIBMESH_TRILINOS_HAVE_EPETRA
178  case TRILINOS_SOLVERS:
179  return std::make_unique<EpetraMatrix<T>>(comm);
180 #endif
181 
182 
183 #ifdef LIBMESH_HAVE_EIGEN
184  case EIGEN_SOLVERS:
185  return std::make_unique<EigenSparseMatrix<T>>(comm);
186 #endif
187 
188  default:
189  libmesh_error_msg("ERROR: Unrecognized solver package: " << solver_package);
190  }
191 }
192 
193 
194 template <typename T>
196  const NumericVector<T> & arg) const
197 {
198  dest.zero();
199  this->vector_mult_add(dest,arg);
200 }
201 
202 
203 
204 template <typename T>
206  const NumericVector<T> & arg) const
207 {
208  /* This functionality is actually implemented in the \p
209  NumericVector class. */
210  dest.add_vector(arg,*this);
211 }
212 
213 
214 
215 template <typename T>
216 void SparseMatrix<T>::zero_rows (std::vector<numeric_index_type> &, T)
217 {
218  /* This functionality isn't implemented or stubbed in every subclass yet */
219  libmesh_not_implemented();
220 }
221 
222 
223 
224 template <typename T>
225 void SparseMatrix<T>::print(std::ostream & os, const bool sparse) const
226 {
227  parallel_object_only();
228 
229  libmesh_assert (this->initialized());
230 
231  libmesh_error_msg_if(!this->_dof_map, "Error! Trying to print a matrix with no dof_map set!");
232 
233  // We'll print the matrix from processor 0 to make sure
234  // it's serialized properly
235  if (this->processor_id() == 0)
236  {
237  libmesh_assert_equal_to (this->_dof_map->first_dof(), 0);
238  for (numeric_index_type i=this->_dof_map->first_dof();
239  i!=this->_dof_map->end_dof(); ++i)
240  {
241  if (sparse)
242  {
243  for (auto j : make_range(this->n()))
244  {
245  T c = (*this)(i,j);
246  if (c != static_cast<T>(0.0))
247  {
248  os << i << " " << j << " " << c << std::endl;
249  }
250  }
251  }
252  else
253  {
254  for (auto j : make_range(this->n()))
255  os << (*this)(i,j) << " ";
256  os << std::endl;
257  }
258  }
259 
260  std::vector<numeric_index_type> ibuf, jbuf;
261  std::vector<T> cbuf;
262  numeric_index_type currenti = this->_dof_map->end_dof();
263  for (auto p : IntRange<processor_id_type>(1, this->n_processors()))
264  {
265  this->comm().receive(p, ibuf);
266  this->comm().receive(p, jbuf);
267  this->comm().receive(p, cbuf);
268  libmesh_assert_equal_to (ibuf.size(), jbuf.size());
269  libmesh_assert_equal_to (ibuf.size(), cbuf.size());
270 
271  if (ibuf.empty())
272  continue;
273  libmesh_assert_greater_equal (ibuf.front(), currenti);
274  libmesh_assert_greater_equal (ibuf.back(), ibuf.front());
275 
276  std::size_t currentb = 0;
277  for (;currenti <= ibuf.back(); ++currenti)
278  {
279  if (sparse)
280  {
281  for (numeric_index_type j=0; j<this->n(); j++)
282  {
283  if (currentb < ibuf.size() &&
284  ibuf[currentb] == currenti &&
285  jbuf[currentb] == j)
286  {
287  os << currenti << " " << j << " " << cbuf[currentb] << std::endl;
288  currentb++;
289  }
290  }
291  }
292  else
293  {
294  for (auto j : make_range(this->n()))
295  {
296  if (currentb < ibuf.size() &&
297  ibuf[currentb] == currenti &&
298  jbuf[currentb] == j)
299  {
300  os << cbuf[currentb] << " ";
301  currentb++;
302  }
303  else
304  os << static_cast<T>(0.0) << " ";
305  }
306  os << std::endl;
307  }
308  }
309  }
310  if (!sparse)
311  {
312  for (; currenti != this->m(); ++currenti)
313  {
314  for (numeric_index_type j=0; j<this->n(); j++)
315  os << static_cast<T>(0.0) << " ";
316  os << std::endl;
317  }
318  }
319  }
320  else
321  {
322  std::vector<numeric_index_type> ibuf, jbuf;
323  std::vector<T> cbuf;
324 
325  // We'll assume each processor has access to entire
326  // matrix rows, so (*this)(i,j) is valid if i is a local index.
327  for (numeric_index_type i=this->_dof_map->first_dof();
328  i!=this->_dof_map->end_dof(); ++i)
329  {
330  for (auto j : make_range(this->n()))
331  {
332  T c = (*this)(i,j);
333  if (c != static_cast<T>(0.0))
334  {
335  ibuf.push_back(i);
336  jbuf.push_back(j);
337  cbuf.push_back(c);
338  }
339  }
340  }
341  this->comm().send(0,ibuf);
342  this->comm().send(0,jbuf);
343  this->comm().send(0,cbuf);
344  }
345 }
346 
347 
348 
349 //------------------------------------------------------------------
350 // Explicit instantiations
351 template class LIBMESH_EXPORT SparseMatrix<Number>;
352 
353 } // namespace libMesh
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
SparseMatrix(const Parallel::Communicator &comm)
构造函数;初始化矩阵为空,没有任何结构,即矩阵无法使用。
Definition: sparse_matrix.C:48
unsigned int n() const
返回矩阵的列维度。
This helper class can be called on multiple threads to compute the sparsity pattern (or graph) of the...
boost::multiprecision::float128 real(const boost::multiprecision::float128 in)
EIGEN_SOLVERS
Definition: libmesh.C:249
TRILINOS_SOLVERS
Definition: libmesh.C:247
virtual void add_vector(const T *v, const std::vector< numeric_index_type > &dof_indices)
计算 ,其中 v 是一个指针, 每个 dof_indices[i] 指定了要添加的值 v[i] 的位置。 这应该在子类中进行重写以提高效率。请注意,此方法的库实现是线程安全的。 ...
void vector_mult(NumericVector< T > &dest, const NumericVector< T > &arg) const
将矩阵乘以 NumericVector arg 并将结果存储在 NumericVector dest 中。
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 中。
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 void zero()=0
将所有条目设置为零。等同于 v = 0,但更明显且更快。
这是一个通用的稀疏矩阵类。该类包含了必须在派生类中覆盖的纯虚拟成员。 使用一个公共的基类允许从不同的求解器包中以不同的格式统一访问稀疏矩阵。
Definition: dof_map.h:66
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:169
void libmesh_ignore(const Args &...)
dof_id_type numeric_index_type
Definition: id_types.h:99
bool _is_initialized
Flag that tells if init() has been called.
Definition: libmesh.C:242
unsigned int m() const
返回矩阵的行维度。
LASPACK_SOLVERS
Definition: libmesh.C:251
const SparsityPattern::Build * get_sparsity_pattern() const
Definition: dof_map.h:548
virtual void zero_rows(std::vector< numeric_index_type > &rows, T diag_value=0.0)
将所有行条目设置为 0,然后将 diag_value 放在对角线条目。
void attach_dof_map(const DofMap &dof_map)
设置要使用的 DofMap 的指针。如果不使用单独的稀疏性模式, 则使用来自 DofMap 的模式。
Definition: sparse_matrix.C:58
void print(std::ostream &os=libMesh::out, const bool sparse=false) const
将矩阵的内容以统一的样式打印到屏幕上,而不考虑正在使用的矩阵/求解器包。
bool initialized()
Checks that library initialization has been done.
Definition: libmesh.C:261
boost::multiprecision::float128 imag(const boost::multiprecision::float128)
定义用于有限元类型计算的密集矩阵。 用于在求和成全局矩阵之前存储单元刚度矩阵。所有被覆盖的虚函数都记录在dense_matrix_base.h中。
Definition: dof_map.h:65