21 #include "libmesh/libmesh_config.h" 23 #ifdef LIBMESH_HAVE_LASPACK 25 #include "libmesh/laspack_matrix.h" 26 #include "libmesh/dense_matrix.h" 27 #include "libmesh/dof_map.h" 28 #include "libmesh/sparsity_pattern.h" 58 size += sparsity_pattern[row].size();
61 _row_start.reserve(n_rows + 1);
67 std::vector<numeric_index_type>::iterator position = _csr.begin();
69 _row_start.push_back (position);
74 for (
const auto & col : sparsity_pattern[row])
81 _row_start.push_back (position);
92 libmesh_assert_equal_to (n_rows, this->m());
98 auto rs = _row_start[i];
102 Q_SetLen (&_QMat, i+1, length);
116 libmesh_assert_equal_to (this->pos(i,j), l);
117 Q_SetEntry (&_QMat, i+1, l, j+1, 0.);
127 template <
typename T>
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);
144 libmesh_not_implemented_msg(
"Laspack does not support distributed matrices");
147 libmesh_not_implemented_msg(
"Laspack does not support rectangular matrices");
150 libmesh_warning(
"Using inefficient LaspackMatrix allocation via this init() method");
157 _csr.resize (m_in * m_in);
158 _row_start.reserve(m_in + 1);
163 std::vector<numeric_index_type>::iterator position = _csr.begin();
165 _row_start.push_back (position);
177 _row_start.push_back (position);
181 Q_Constr(&_QMat, const_cast<char *>(
"Mat"), m_in, _LPFalse, Rowws, Normal, _LPTrue);
189 auto rs = _row_start[i];
193 Q_SetLen (&_QMat, i+1, length);
199 libmesh_assert_equal_to (this->pos(i,j), l);
200 Q_SetEntry (&_QMat, i+1, l, j+1, 0.);
207 template <
typename T>
231 libmesh_assert_equal_to (n_rows, n_cols);
232 libmesh_assert_equal_to (m_l, n_rows);
233 libmesh_assert_equal_to (n_l, n_cols);
238 const std::vector<numeric_index_type> & n_nz = this->_sp->get_n_nz();
239 const std::vector<numeric_index_type> & n_oz = this->_sp->get_n_oz();
243 libmesh_assert_equal_to (n_nz.size(), n_l);
244 libmesh_assert_equal_to (n_oz.size(), n_l);
249 Q_Constr(&_QMat, const_cast<char *>(
"Mat"), n_rows, _LPFalse, Rowws, Normal, _LPTrue);
253 libmesh_assert_equal_to (n_rows, this->m());
258 template <
typename T>
260 const std::vector<numeric_index_type> & rows,
261 const std::vector<numeric_index_type> & cols)
265 unsigned int n_rows = cast_int<unsigned int>(rows.size());
266 unsigned int n_cols = cast_int<unsigned int>(cols.size());
267 libmesh_assert_equal_to (dm.
m(), n_rows);
268 libmesh_assert_equal_to (dm.
n(), n_cols);
271 for (
unsigned int i=0; i<n_rows; i++)
272 for (
unsigned int j=0; j<n_cols; j++)
273 this->add(rows[i],cols[j],dm(i,j));
278 template <
typename T>
284 std::vector<Real> abs_col_sums(this->n());
290 auto r_start = _row_start[row];
295 libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1));
302 libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l));
304 abs_col_sums[j] += std::abs(Q_GetVal (&_QMat, row+1, l));
308 return *(std::max_element(abs_col_sums.begin(), abs_col_sums.end()));
313 template <
typename T>
316 libmesh_not_implemented();
321 template <
typename T>
324 libmesh_not_implemented();
329 template <
typename T>
331 std::vector<numeric_index_type> & indices,
332 std::vector<T> & values)
const 340 libmesh_assert_equal_to (len, Q_GetLen(&_QMat, i+1));
342 auto r_start = _row_start[i];
349 libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, i+1, l));
351 indices.push_back(j);
353 values.push_back(Q_GetVal (&_QMat, i+1, l));
359 template <
typename T>
368 template <
typename T>
376 template <
typename T>
392 template <
typename T>
399 auto r_start = _row_start[row];
404 libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1));
411 libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l));
413 Q_SetEntry (&_QMat, row+1, l, j+1, 0.);
422 template <
typename T>
426 auto mat_copy = std::make_unique<LaspackMatrix<T>>(this->comm());
428 mat_copy->attach_dof_map(*this->_dof_map);
430 mat_copy->init(this->m(), this->n(), this->local_m(),
431 this->local_n(), this->local_n());
440 template <
typename T>
445 auto mat_copy = this->zero_clone();
446 mat_copy->add(1., *
this);
451 template <
typename T>
461 template <
typename T>
471 template <
typename T>
479 template <
typename T>
487 template <
typename T>
495 template <
typename T>
503 template <
typename T>
509 libmesh_assert_less (i, this->m());
510 libmesh_assert_less (j, this->n());
515 libmesh_assert_equal_to (*(_row_start[i]+position), j);
516 libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, i+1, position));
518 Q_SetEntry (&_QMat, i+1, position, j+1,
value);
523 template <
typename T>
529 libmesh_assert_less (i, this->m());
530 libmesh_assert_less (j, this->n());
535 libmesh_assert_equal_to (*(_row_start[i]+position), j);
537 Q_AddVal (&_QMat, i+1, position,
value);
542 template <
typename T>
544 const std::vector<numeric_index_type> & dof_indices)
546 this->add_matrix (dm, dof_indices, dof_indices);
551 template <
typename T>
555 libmesh_assert_equal_to (this->m(), X_in.
m());
556 libmesh_assert_equal_to (this->n(), X_in.
n());
559 cast_ptr<const LaspackMatrix<T> *> (&X_in);
561 _LPNumber a =
static_cast<_LPNumber
> (a_in);
571 auto r_start = _row_start[row];
576 libmesh_assert_equal_to (len, Q_GetLen(&_QMat, row+1));
578 libmesh_assert_equal_to (len, Q_GetLen(&(X->_QMat), row+1));
586 libmesh_assert_equal_to ((j+1), Q_GetPos (&_QMat, row+1, l));
588 const _LPNumber
value = a * Q_GetEl(const_cast<QMatrix*>(&(X->_QMat)), row+1, j+1);
589 Q_AddVal (&_QMat, row+1, l,
value);
597 template <
typename T>
602 libmesh_assert_less (i, this->m());
603 libmesh_assert_less (j, this->n());
605 return Q_GetEl (const_cast<QMatrix*>(&_QMat), i+1, j+1);
610 template <
typename T>
614 libmesh_assert_less (i, this->m());
615 libmesh_assert_less (j, this->n());
616 libmesh_assert_less (i+1, _row_start.size());
620 auto p = std::equal_range (_row_start[i], _row_start[i+1], j);
634 template <
typename T>
639 this->_closed =
true;
643 *_QMat.DiagElAlloc = _LPFalse;
644 *_QMat.ElSorted = _LPFalse;
645 if (*_QMat.ILUExists)
647 *_QMat.ILUExists = _LPFalse;
661 #endif // #ifdef LIBMESH_HAVE_LASPACK 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
Initialize SparseMatrix with the specified sizes.
virtual std::unique_ptr< SparseMatrix< T > > zero_clone() const override
virtual void get_diagonal(NumericVector< T > &dest) const override
Copies the diagonal part of the matrix into dest.
virtual numeric_index_type col_start() const override
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
The libMesh namespace provides an interface to certain functionality in the library.
Real distance(const Point &p)
virtual void add(const numeric_index_type i, const numeric_index_type j, const T value) override
Add value to the element (i,j).
virtual std::unique_ptr< SparseMatrix< T > > clone() const override
numeric_index_type pos(const numeric_index_type i, const numeric_index_type j) const
virtual void get_transpose(SparseMatrix< T > &dest) const override
Copies the transpose of the matrix into dest, which may be *this.
virtual void close() override
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
dof_id_type numeric_index_type
bool _is_initialized
Flag that tells if init() has been called.
virtual Real l1_norm() const override
virtual numeric_index_type m() const =0
virtual T operator()(const numeric_index_type i, const numeric_index_type j) const override
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary.
virtual numeric_index_type row_start() const override
virtual void zero() override
Set all entries to 0.
virtual void clear() override
Restores the SparseMatrix<T> to a pristine state.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual numeric_index_type m() const override
virtual numeric_index_type n() const override
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
bool initialized()
Checks that library initialization has been done.
virtual void update_sparsity_pattern(const SparsityPattern::Graph &) override
Updates the matrix sparsity pattern.
virtual numeric_index_type row_stop() const override
The LaspackMatrix class wraps a QMatrix object from the Laspack library.
virtual void set(const numeric_index_type i, const numeric_index_type j, const T value) override
Set the element (i,j) to value.
Defines a dense matrix for use in Finite Element-type computations.
virtual void get_row(numeric_index_type i, std::vector< numeric_index_type > &indices, std::vector< T > &values) const override
Get a row from the matrix.
virtual numeric_index_type n() const =0
virtual numeric_index_type col_stop() const override
virtual void add_matrix(const DenseMatrix< T > &dm, const std::vector< numeric_index_type > &rows, const std::vector< numeric_index_type > &cols) override
Add the full matrix dm to the SparseMatrix.
LaspackMatrix(const Parallel::Communicator &comm)
Constructor; initializes the matrix to be empty, without any structure, i.e.
ParallelType
Defines an enum for parallel data structure types.