libmesh解析
本工作只是尝试解析原libmesh的代码,供学习使用
 全部  命名空间 文件 函数 变量 类型定义 枚举 枚举值 友元 
laspack_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/libmesh_config.h"
22 
23 #ifdef LIBMESH_HAVE_LASPACK
24 
25 #include "libmesh/laspack_matrix.h"
26 #include "libmesh/dense_matrix.h"
27 #include "libmesh/dof_map.h"
28 #include "libmesh/sparsity_pattern.h"
29 
30 
31 // C++ Includes
32 #include <memory>
33 
34 
35 namespace libMesh
36 {
37 
38 
39 //-----------------------------------------------------------------------
40 // LaspackMatrix members
41 template <typename T>
43 {
44  // clear data, start over
45  this->clear ();
46 
47  // big trouble if this fails!
48  libmesh_assert(this->_dof_map);
49 
50  const numeric_index_type n_rows = sparsity_pattern.size();
51 
52  // Initialize the _row_start data structure,
53  // allocate storage for the _csr array
54  {
55  std::size_t size = 0;
56 
57  for (numeric_index_type row=0; row<n_rows; row++)
58  size += sparsity_pattern[row].size();
59 
60  _csr.resize (size);
61  _row_start.reserve(n_rows + 1);
62  }
63 
64 
65  // Initialize the _csr data structure.
66  {
67  std::vector<numeric_index_type>::iterator position = _csr.begin();
68 
69  _row_start.push_back (position);
70 
71  for (numeric_index_type row=0; row<n_rows; row++)
72  {
73  // insert the row indices
74  for (const auto & col : sparsity_pattern[row])
75  {
76  libmesh_assert (position != _csr.end());
77  *position = col;
78  ++position;
79  }
80 
81  _row_start.push_back (position);
82  }
83  }
84 
85 
86  // Initialize the matrix
87  libmesh_assert (!this->initialized());
88  this->init ();
89  libmesh_assert (this->initialized());
90  //libMesh::out << "n_rows=" << n_rows << std::endl;
91  //libMesh::out << "m()=" << m() << std::endl;
92  libmesh_assert_equal_to (n_rows, this->m());
93 
94  // Tell the matrix about its structure. Initialize it
95  // to zero.
96  for (numeric_index_type i=0; i<n_rows; i++)
97  {
98  auto rs = _row_start[i];
99 
100  const numeric_index_type length = _row_start[i+1] - rs;
101 
102  Q_SetLen (&_QMat, i+1, length);
103 
104  for (numeric_index_type l=0; l<length; l++)
105  {
106  const numeric_index_type j = *(rs+l);
107 
108  // sanity check
109  //libMesh::out << "m()=" << m() << std::endl;
110  //libMesh::out << "(i,j,l) = (" << i
111  // << "," << j
112  // << "," << l
113  // << ")" << std::endl;
114  //libMesh::out << "pos(i,j)=" << pos(i,j)
115  // << std::endl;
116  libmesh_assert_equal_to (this->pos(i,j), l);
117  Q_SetEntry (&_QMat, i+1, l, j+1, 0.);
118  }
119  }
120 
121  // That's it!
122  //libmesh_here();
123 }
124 
125 
126 
127 template <typename T>
128 void LaspackMatrix<T>::init (const numeric_index_type libmesh_dbg_var(m_in),
129  const numeric_index_type libmesh_dbg_var(n_in),
130  const numeric_index_type libmesh_dbg_var(m_l),
131  const numeric_index_type libmesh_dbg_var(n_l),
132  const numeric_index_type libmesh_dbg_var(nnz),
133  const numeric_index_type,
134  const numeric_index_type)
135 {
136  // noz ignored... only used for multiple processors!
137  libmesh_assert_equal_to (m_in, m_l);
138  libmesh_assert_equal_to (n_in, n_l);
139  libmesh_assert_equal_to (m_in, n_in);
140  libmesh_assert_greater (nnz, 0);
141 
142  libmesh_error_msg("ERROR: Only the init() member that uses the DofMap is implemented for Laspack matrices!");
143 
144  this->_is_initialized = true;
145 }
146 
147 
148 
149 template <typename T>
150 void LaspackMatrix<T>::init (const ParallelType)
151 {
152  // Ignore calls on initialized objects
153  if (this->initialized())
154  return;
155 
156  // We need the DofMap for this!
157  libmesh_assert(this->_dof_map);
158 
159  // Clear initialized matrices
160  if (this->initialized())
161  this->clear();
162 
163  const numeric_index_type n_rows = this->_dof_map->n_dofs();
164 #ifndef NDEBUG
165  // The following variables are only used for assertions,
166  // so avoid declaring them when asserts are inactive.
167  const numeric_index_type n_cols = n_rows;
168  const numeric_index_type n_l = this->_dof_map->n_dofs_on_processor(0);
169  const numeric_index_type m_l = n_l;
170 #endif
171 
172  // Laspack Matrices only work for uniprocessor cases
173  libmesh_assert_equal_to (n_rows, n_cols);
174  libmesh_assert_equal_to (m_l, n_rows);
175  libmesh_assert_equal_to (n_l, n_cols);
176 
177 #ifndef NDEBUG
178  // The following variables are only used for assertions,
179  // so avoid declaring them when asserts are inactive.
180  const std::vector<numeric_index_type> & n_nz = this->_sp->get_n_nz();
181  const std::vector<numeric_index_type> & n_oz = this->_sp->get_n_oz();
182 #endif
183 
184  // Make sure the sparsity pattern isn't empty
185  libmesh_assert_equal_to (n_nz.size(), n_l);
186  libmesh_assert_equal_to (n_oz.size(), n_l);
187 
188  if (n_rows==0)
189  return;
190 
191  Q_Constr(&_QMat, const_cast<char *>("Mat"), n_rows, _LPFalse, Rowws, Normal, _LPTrue);
192 
193  this->_is_initialized = true;
194 
195  libmesh_assert_equal_to (n_rows, this->m());
196 }
197 
198 
199 
200 template <typename T>
202  const std::vector<numeric_index_type> & rows,
203  const std::vector<numeric_index_type> & cols)
204 
205 {
206  libmesh_assert (this->initialized());
207  unsigned int n_rows = cast_int<unsigned int>(rows.size());
208  unsigned int n_cols = cast_int<unsigned int>(cols.size());
209  libmesh_assert_equal_to (dm.m(), n_rows);
210  libmesh_assert_equal_to (dm.n(), n_cols);
211 
212 
213  for (unsigned int i=0; i<n_rows; i++)
214  for (unsigned int j=0; j<n_cols; j++)
215  this->add(rows[i],cols[j],dm(i,j));
216 }
217 
218 
219 
220 template <typename T>
222 {
223  libmesh_not_implemented();
224 }
225 
226 
227 
228 template <typename T>
230 {
231  libmesh_not_implemented();
232 }
233 
234 
235 
236 template <typename T>
237 LaspackMatrix<T>::LaspackMatrix (const Parallel::Communicator & comm) :
238  SparseMatrix<T>(comm),
239  _closed (false)
240 {
241 }
242 
243 
244 
245 template <typename T>
247 {
248  this->clear ();
249 }
250 
251 
252 
253 template <typename T>
255 {
256  if (this->initialized())
257  {
258  Q_Destr(&_QMat);
259  }
260 
261  _csr.clear();
262  _row_start.clear();
263  _closed = false;
264  this->_is_initialized = false;
265 }
266 
267 
268 
269 template <typename T>
271 {
272  const numeric_index_type n_rows = this->m();
273 
274  for (numeric_index_type row=0; row<n_rows; row++)
275  {
276  auto r_start = _row_start[row];
277 
278  const numeric_index_type len = (_row_start[row+1] - _row_start[row]);
279 
280  // Make sure we agree on the row length
281  libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1));
282 
283  for (numeric_index_type l=0; l<len; l++)
284  {
285  const numeric_index_type j = *(r_start + l);
286 
287  // Make sure the data structures are working
288  libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l));
289 
290  Q_SetEntry (&_QMat, row+1, l, j+1, 0.);
291  }
292  }
293 
294  this->close();
295 }
296 
297 
298 
299 template <typename T>
300 std::unique_ptr<SparseMatrix<T>> LaspackMatrix<T>::zero_clone () const
301 {
302  // This function is marked as "not implemented" since it hasn't been
303  // tested, the code below might serve as a possible implementation.
304  libmesh_not_implemented();
305 
306  // Make empty copy with matching comm, initialize, zero, and return.
307  auto mat_copy = std::make_unique<LaspackMatrix<T>>(this->comm());
308  mat_copy->init();
309  mat_copy->zero();
310 
311  // Work around an issue on older compilers. We are able to simply
312  // "return mat_copy;" on newer compilers
313  return std::unique_ptr<SparseMatrix<T>>(mat_copy.release());
314 }
315 
316 
317 
318 template <typename T>
319 std::unique_ptr<SparseMatrix<T>> LaspackMatrix<T>::clone () const
320 {
321  // We don't currently have a faster implementation than making a
322  // zero clone and then filling in the values.
323  auto mat_copy = this->zero_clone();
324  mat_copy->add(1., *this);
325 
326  // Work around an issue on older compilers. We are able to simply
327  // "return mat_copy;" on newer compilers
328  return std::unique_ptr<SparseMatrix<T>>(mat_copy.release());
329 }
330 
331 template <typename T>
333 {
334  libmesh_assert (this->initialized());
335 
336  return static_cast<numeric_index_type>(Q_GetDim(const_cast<QMatrix*>(&_QMat)));
337 }
338 
339 
340 
341 template <typename T>
343 {
344  libmesh_assert (this->initialized());
345 
346  return static_cast<numeric_index_type>(Q_GetDim(const_cast<QMatrix*>(&_QMat)));
347 }
348 
349 
350 
351 template <typename T>
353 {
354  return 0;
355 }
356 
357 
358 
359 template <typename T>
361 {
362  return this->m();
363 }
364 
365 
366 
367 template <typename T>
369  const numeric_index_type j,
370  const T value)
371 {
372  libmesh_assert (this->initialized());
373  libmesh_assert_less (i, this->m());
374  libmesh_assert_less (j, this->n());
375 
376  const numeric_index_type position = this->pos(i,j);
377 
378  // Sanity check
379  libmesh_assert_equal_to (*(_row_start[i]+position), j);
380  libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, i+1, position));
381 
382  Q_SetEntry (&_QMat, i+1, position, j+1, value);
383 }
384 
385 
386 
387 template <typename T>
389  const numeric_index_type j,
390  const T value)
391 {
392  libmesh_assert (this->initialized());
393  libmesh_assert_less (i, this->m());
394  libmesh_assert_less (j, this->n());
395 
396  const numeric_index_type position = this->pos(i,j);
397 
398  // Sanity check
399  libmesh_assert_equal_to (*(_row_start[i]+position), j);
400 
401  Q_AddVal (&_QMat, i+1, position, value);
402 }
403 
404 
405 
406 template <typename T>
408  const std::vector<numeric_index_type> & dof_indices)
409 {
410  this->add_matrix (dm, dof_indices, dof_indices);
411 }
412 
413 
414 
415 template <typename T>
416 void LaspackMatrix<T>::add (const T a_in, const SparseMatrix<T> & X_in)
417 {
418  libmesh_assert (this->initialized());
419  libmesh_assert_equal_to (this->m(), X_in.m());
420  libmesh_assert_equal_to (this->n(), X_in.n());
421 
422  const LaspackMatrix<T> * X =
423  cast_ptr<const LaspackMatrix<T> *> (&X_in);
424 
425  _LPNumber a = static_cast<_LPNumber> (a_in);
426 
427  libmesh_assert(X);
428 
429  // loops taken from LaspackMatrix<T>::zero ()
430 
431  const numeric_index_type n_rows = this->m();
432 
433  for (numeric_index_type row=0; row<n_rows; row++)
434  {
435  auto r_start = _row_start[row];
436 
437  const numeric_index_type len = (_row_start[row+1] - _row_start[row]);
438 
439  // Make sure we agree on the row length
440  libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1));
441  // compare matrix sparsity structures
442  libmesh_assert_equal_to (len, Q_GetLen(&(X->_QMat), row+1));
443 
444 
445  for (numeric_index_type l=0; l<len; l++)
446  {
447  const numeric_index_type j = *(r_start + l);
448 
449  // Make sure the data structures are working
450  libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l));
451 
452  const _LPNumber value = a * Q_GetEl(const_cast<QMatrix*>(&(X->_QMat)), row+1, j+1);
453  Q_AddVal (&_QMat, row+1, l, value);
454  }
455  }
456 }
457 
458 
459 
460 
461 template <typename T>
463  const numeric_index_type j) const
464 {
465  libmesh_assert (this->initialized());
466  libmesh_assert_less (i, this->m());
467  libmesh_assert_less (j, this->n());
468 
469  return Q_GetEl (const_cast<QMatrix*>(&_QMat), i+1, j+1);
470 }
471 
472 
473 
474 template <typename T>
476  const numeric_index_type j) const
477 {
478  libmesh_assert_less (i, this->m());
479  libmesh_assert_less (j, this->n());
480  libmesh_assert_less (i+1, _row_start.size());
481  libmesh_assert (_row_start.back() == _csr.end());
482 
483  // note this requires the _csr to be sorted
484  auto p = std::equal_range (_row_start[i], _row_start[i+1], j);
485 
486  // Make sure the row contains the element j
487  libmesh_assert (p.first != p.second);
488 
489  // Make sure the values match
490  libmesh_assert (*p.first == j);
491 
492  // Return the position in the compressed row
493  return std::distance (_row_start[i], p.first);
494 }
495 
496 
497 
498 template <typename T>
500 {
501  libmesh_assert(this->initialized());
502 
503  this->_closed = true;
504 
505  // We've probably changed some entries so we need to tell LASPACK
506  // that cached data is now invalid.
507  *_QMat.DiagElAlloc = _LPFalse;
508  *_QMat.ElSorted = _LPFalse;
509  if (*_QMat.ILUExists)
510  {
511  *_QMat.ILUExists = _LPFalse;
512  Q_Destr(_QMat.ILU);
513  }
514 }
515 
516 
517 
518 //------------------------------------------------------------------
519 // Explicit instantiations
520 template class LIBMESH_EXPORT LaspackMatrix<Number>;
521 
522 } // namespace libMesh
523 
524 
525 #endif // #ifdef LIBMESH_HAVE_LASPACK
virtual std::unique_ptr< SparseMatrix< T > > zero_clone() const override
创建一个与LaspackMatrix对象具有相同属性但所有元素为零的新矩阵。 这个函数用于创建一个与原矩阵具有相同大小和属性的全零矩阵。
virtual void get_diagonal(NumericVector< T > &dest) const override
获取矩阵的对角线元素,并存储到指定的NumericVector对象中。 这个函数用于获取矩阵的对角线元素。
unsigned int n() const
返回矩阵的列维度。
提供了不同线性代数库的向量存储方案的统一接口。
Definition: dof_map.h:67
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) override
将 处的值增加 value。
这是一个通用的稀疏矩阵类。该类包含了必须在派生类中覆盖的纯虚拟成员。 使用一个公共的基类允许从不同的求解器包中以不同的格式统一访问稀疏矩阵。
Definition: dof_map.h:66
virtual std::unique_ptr< SparseMatrix< T > > clone() const override
克隆LaspackMatrix对象,创建一个具有相同属性和数据的新矩阵。 这个函数用于创建一个与原矩阵具有相同大小、属性和数据的新矩阵。
virtual void get_transpose(SparseMatrix< T > &dest) const override
获取矩阵的转置,并存储到指定的SparseMatrix对象中。 这个函数用于获取矩阵的转置矩阵。
virtual void close() override
关闭LaspackMatrix对象,标记矩阵已经完成构建,不再修改其结构。 这个函数用于标记矩阵的状态为已完成构建,之后不再修改其结构。
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
返回矩阵的行维度。
virtual numeric_index_type m() const =0
virtual T operator()(const numeric_index_type i, const numeric_index_type j) const override
获取LaspackMatrix对象的指定位置 处的值。 这个函数用于获取矩阵的特定元素的值。
virtual numeric_index_type row_start() const override
返回LaspackMatrix对象的行起始位置。 这个函数用于获取矩阵的行起始位置,通常用于稀疏矩阵的索引操作。
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) override
初始化LaspackMatrix对象的属性,包括矩阵的行数、列数、非零元素的预估数量等。这些参数将影响矩阵的内存分配和性能。
virtual void zero() override
将LaspackMatrix对象的所有元素置零。 这个函数用于将矩阵中的所有元素设置为零。
virtual void clear() override
清空LaspackMatrix对象,释放分配的内存和资源。 这个函数用于清除矩阵的数据和状态。
virtual numeric_index_type m() const override
返回LaspackMatrix对象的行数。 这个函数用于获取矩阵的行数。
virtual numeric_index_type n() const override
返回LaspackMatrix对象的列数。 这个函数用于获取矩阵的列数。
bool initialized()
Checks that library initialization has been done.
Definition: libmesh.C:261
virtual void update_sparsity_pattern(const SparsityPattern::Graph &) override
更新矩阵的稀疏性模式。这将告诉底层矩阵存储方案如何映射 元素。
virtual numeric_index_type row_stop() const override
返回LaspackMatrix对象的行结束位置。 这个函数用于获取矩阵的行结束位置,通常用于稀疏矩阵的索引操作。
LaspackMatrix类封装了Laspack库中的QMatrix对象。 目前,Laspack仅支持实数数据类型,因此这个类是对 SparseMatrix&lt;T&gt; 的全特化,其中 T = Real。 所...
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value) override
设置LaspackMatrix对象的指定位置 处的值为 value。
定义用于有限元类型计算的密集矩阵。 用于在求和成全局矩阵之前存储单元刚度矩阵。所有被覆盖的虚函数都记录在dense_matrix_base.h中。
Definition: dof_map.h:65
numeric_index_type pos(const numeric_index_type i, const numeric_index_type j) const
virtual numeric_index_type n() const =0
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) override
将一个DenseMatrix对象的元素加到LaspackMatrix对象的指定行和列。
LaspackMatrix(const Parallel::Communicator &comm)
构造函数;将矩阵初始化为空,没有任何结构,即矩阵无法使用。因此,此构造函数仅适用于类的成员矩阵。 所有其他矩阵应在所有必要信息都可用的数据流的某一点创建。