21 #include "libmesh/libmesh_config.h"
23 #ifdef LIBMESH_HAVE_EIGEN
25 #include "libmesh/eigen_sparse_vector.h"
26 #include "libmesh/eigen_sparse_matrix.h"
27 #include "libmesh/dense_matrix.h"
28 #include "libmesh/dof_map.h"
29 #include "libmesh/sparsity_pattern.h"
51 libmesh_assert_equal_to (m_in, m_l);
52 libmesh_assert_equal_to (n_in, n_l);
53 libmesh_assert_equal_to (m_in, n_in);
54 libmesh_assert_greater (nnz, 0);
56 _mat.resize(m_in, n_in);
57 _mat.reserve(Eigen::Matrix<numeric_index_type, Eigen::Dynamic, 1>::Constant(m_in,nnz));
72 libmesh_assert(this->_dof_map);
89 libmesh_assert_equal_to (n_rows, n_cols);
90 libmesh_assert_equal_to (m_l, n_rows);
91 libmesh_assert_equal_to (n_l, n_cols);
93 const std::vector<numeric_index_type> & n_nz = this->_sp->get_n_nz();
98 const std::vector<numeric_index_type> & n_oz = this->_sp->get_n_oz();
102 libmesh_assert_equal_to (n_nz.size(), n_l);
103 libmesh_assert_equal_to (n_oz.size(), n_l);
111 _mat.resize(n_rows,n_cols);
116 libmesh_assert_equal_to (n_rows, this->m());
117 libmesh_assert_equal_to (n_cols, this->n());
122 template <
typename T>
124 const std::vector<numeric_index_type> & rows,
125 const std::vector<numeric_index_type> & cols)
129 unsigned int n_rows = cast_int<unsigned int>(rows.size());
130 unsigned int n_cols = cast_int<unsigned int>(cols.size());
131 libmesh_assert_equal_to (dm.
m(), n_rows);
132 libmesh_assert_equal_to (dm.
n(), n_cols);
135 for (
unsigned int i=0; i<n_rows; i++)
136 for (
unsigned int j=0; j<n_cols; j++)
137 this->add(rows[i],cols[j],dm(i,j));
142 template <
typename T>
147 dest.
_vec = _mat.diagonal();
152 template <
typename T>
157 dest.
_mat = _mat.transpose();
162 template <
typename T>
171 template <
typename T>
182 template <
typename T>
191 const std::vector<numeric_index_type> & n_nz = this->_sp->get_n_nz();
198 template <
typename T>
203 auto ret = std::make_unique<EigenSparseMatrix<T>>(*this);
208 return std::unique_ptr<SparseMatrix<T>>(ret.release());
213 template <
typename T>
216 return std::make_unique<EigenSparseMatrix<T>>(*this);
221 template <
typename T>
226 return cast_int<numeric_index_type>(_mat.rows());
231 template <
typename T>
236 return cast_int<numeric_index_type>(_mat.cols());
241 template <
typename T>
249 template <
typename T>
257 template <
typename T>
263 libmesh_assert_less (i, this->m());
264 libmesh_assert_less (j, this->n());
266 _mat.coeffRef(i,j) = value;
271 template <
typename T>
277 libmesh_assert_less (i, this->m());
278 libmesh_assert_less (j, this->n());
280 _mat.coeffRef(i,j) += value;
285 template <
typename T>
287 const std::vector<numeric_index_type> & dof_indices)
289 this->add_matrix (dm, dof_indices, dof_indices);
294 template <
typename T>
298 libmesh_assert_equal_to (this->m(), X_in.
m());
299 libmesh_assert_equal_to (this->n(), X_in.
n());
302 cast_ref<const EigenSparseMatrix<T> &> (X_in);
310 template <
typename T>
315 libmesh_assert_less (i, this->m());
316 libmesh_assert_less (j, this->n());
318 return _mat.coeff(i,j);
323 template <
typename T>
330 std::vector<Real> abs_col_sums(this->n());
334 for (
auto row : make_range(this->m()))
336 EigenSM::InnerIterator it(_mat, row);
338 abs_col_sums[it.col()] +=
std::abs(it.value());
341 return *(std::max_element(abs_col_sums.begin(), abs_col_sums.end()));
346 template <
typename T>
349 Real max_abs_row_sum = 0.;
353 for (
auto row : make_range(this->m()))
355 Real current_abs_row_sum = 0.;
356 EigenSM::InnerIterator it(_mat, row);
358 current_abs_row_sum +=
std::abs(it.value());
360 max_abs_row_sum = std::max(max_abs_row_sum, current_abs_row_sum);
363 return max_abs_row_sum;
374 #endif // #ifdef LIBMESH_HAVE_EIGEN
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value) override
设置元素 (i,j) 为 value .
EigenSparseMatrix 类包装了来自 Eigen 库的稀疏矩阵对象。
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
使用指定的大小初始化 SparseMatrix。
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) override
将 value 添加到元素 (i,j) .
virtual numeric_index_type m() const override
unsigned int n() const
返回矩阵的列维度。
virtual numeric_index_type row_stop() const override
DataType _vec
Actual Eigen::SparseVector<> we are wrapping.
virtual std::unique_ptr< SparseMatrix< T > > clone() const override
DataType _mat
实际的 Eigen::SparseMatrix<> 对象。
virtual void get_transpose(SparseMatrix< T > &dest) const override
将矩阵的转置复制到 dest 中,dest 可能是 *this。
ADRealEigenVector< T, D, asd > abs(const ADRealEigenVector< T, D, asd > &)
计算自动微分实数向量的绝对值。
这是一个通用的稀疏矩阵类。该类包含了必须在派生类中覆盖的纯虚拟成员。 使用一个公共的基类允许从不同的求解器包中以不同的格式统一访问稀疏矩阵。
virtual void zero() override
将所有条目设置为 0。
virtual numeric_index_type n() const override
virtual void get_diagonal(NumericVector< T > &dest) const override
复制矩阵的对角线部分到 dest。
dof_id_type numeric_index_type
bool _is_initialized
Flag that tells if init() has been called.
unsigned int m() const
返回矩阵的行维度。
EigenSparseMatrix(const Parallel::Communicator &comm)
构造函数;初始化矩阵为空,没有任何结构,即矩阵无法使用。
virtual numeric_index_type m() const =0
virtual T operator()(const numeric_index_type i, const numeric_index_type j) const override
virtual numeric_index_type row_start() const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real linfty_norm() const override
获取矩阵的 -范数,即最大行和:
virtual void clear() override
恢复 SparseMatrix<T> 到原始状态。
bool initialized()
Checks that library initialization has been done.
virtual std::unique_ptr< SparseMatrix< T > > zero_clone() const override
virtual Real l1_norm() const override
获取矩阵的 -范数,即最大列和:
定义用于有限元类型计算的密集矩阵。 用于在求和成全局矩阵之前存储单元刚度矩阵。所有被覆盖的虚函数都记录在dense_matrix_base.h中。
This class provides a nice interface to the Eigen C++-based data structures for serial vectors...
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
将完整矩阵 dm 添加到 SparseMatrix。这对于在装配时添加元素矩阵很有用。