16#ifndef dealii_sparse_matrix_ez_h
17#define dealii_sparse_matrix_ez_h
33template <
typename number>
35template <
typename number>
103template <
typename number>
171 static_cast<unsigned short>(-1);
193 const unsigned short index);
245 const unsigned short index);
323 const unsigned int default_increment = 1);
358 unsigned int default_increment = 1,
417 template <
class StreamType>
434 std::vector<size_type> &used_by_line,
435 const bool compute_by_line)
const;
462 const bool elide_zero_values =
true);
491 template <
typename number2>
493 add(
const std::vector<size_type> &indices,
495 const bool elide_zero_values =
true);
502 template <
typename number2>
504 add(
const std::vector<size_type> &row_indices,
505 const std::vector<size_type> &col_indices,
507 const bool elide_zero_values =
true);
518 template <
typename number2>
521 const std::vector<size_type> &col_indices,
522 const std::vector<number2> & values,
523 const bool elide_zero_values =
true);
534 template <
typename number2>
539 const number2 * values,
540 const bool elide_zero_values =
true,
541 const bool col_indices_are_sorted =
false);
564 template <
typename MatrixType>
566 copy_from(
const MatrixType &source,
const bool elide_zero_values =
true);
575 template <
typename MatrixType>
577 add(
const number factor,
const MatrixType &matrix);
610 template <
typename somenumber>
619 template <
typename somenumber>
627 template <
typename somenumber>
636 template <
typename somenumber>
659 template <
typename somenumber>
663 const number omega = 1.)
const;
668 template <
typename somenumber>
672 const number om = 1.,
673 const std::vector<std::size_t> &pos_right_of_diagonal =
674 std::vector<std::size_t>())
const;
680 template <
typename somenumber>
684 const number om = 1.)
const;
690 template <
typename somenumber>
694 const number om = 1.)
const;
704 template <
typename MatrixTypeA,
typename MatrixTypeB>
707 const MatrixTypeB &B,
773 const unsigned int precision = 3,
774 const bool scientific =
true,
775 const unsigned int width = 0,
776 const char * zero_string =
" ",
777 const double denominator = 1.)
const;
817 <<
"The entry with index (" << arg1 <<
',' << arg2
818 <<
") does not exist.");
823 <<
"An entry with index (" << arg1 <<
',' << arg2
824 <<
") cannot be allocated.");
852 template <
typename somenumber>
864 template <
typename somenumber>
869 somenumber * partial_sum)
const;
876 template <
typename somenumber>
882 somenumber * partial_sum)
const;
915template <
typename number>
917 const number & value)
924template <
typename number>
931template <
typename number>
935 , diagonal(invalid_diagonal)
940template <
typename number>
944 const unsigned short i)
951template <
typename number>
959template <
typename number>
963 return matrix->data[matrix->row_info[a_row].start + a_index].column;
967template <
typename number>
976template <
typename number>
980 return matrix->data[matrix->row_info[a_row].start + a_index].value;
984template <
typename number>
988 const unsigned short i)
1018template <
typename number>
1025 ++(accessor.a_index);
1029 if (accessor.a_index >= accessor.matrix->row_info[accessor.a_row].length)
1031 accessor.a_index = 0;
1038 while (accessor.a_row < accessor.matrix->m() &&
1039 accessor.matrix->row_info[accessor.a_row].length == 0);
1045template <
typename number>
1053template <
typename number>
1061template <
typename number>
1071template <
typename number>
1076 return !(*
this == other);
1080template <
typename number>
1092template <
typename number>
1100template <
typename number>
1108template <
typename number>
1120 if (entry->
column == col)
1130template <
typename number>
1135 return t->
locate(row, col);
1139template <
typename number>
1156 while (i <
end &&
data[i].column < col)
1160 if (i !=
end &&
data[i].column == col)
1199 Entry temp = *entry;
1226 std::swap(
data[j], temp);
1237template <
typename number>
1242 const bool elide_zero_values)
1249 if (elide_zero_values && value == 0.)
1252 if (entry !=
nullptr)
1258 entry->
value = value;
1264template <
typename number>
1280 entry->
value += value;
1284template <
typename number>
1285template <
typename number2>
1289 const bool elide_zero_values)
1292 for (
size_type i = 0; i < indices.size(); ++i)
1293 for (
size_type j = 0; j < indices.size(); ++j)
1294 if ((full_matrix(i, j) != 0) || (elide_zero_values ==
false))
1295 add(indices[i], indices[j], full_matrix(i, j));
1300template <
typename number>
1301template <
typename number2>
1304 const std::vector<size_type> &col_indices,
1306 const bool elide_zero_values)
1309 for (
size_type i = 0; i < row_indices.size(); ++i)
1310 for (
size_type j = 0; j < col_indices.size(); ++j)
1311 if ((full_matrix(i, j) != 0) || (elide_zero_values ==
false))
1312 add(row_indices[i], col_indices[j], full_matrix(i, j));
1317template <
typename number>
1318template <
typename number2>
1321 const std::vector<size_type> &col_indices,
1322 const std::vector<number2> & values,
1323 const bool elide_zero_values)
1326 for (
size_type j = 0; j < col_indices.size(); ++j)
1327 if ((values[j] != 0) || (elide_zero_values ==
false))
1328 add(row, col_indices[j], values[j]);
1333template <
typename number>
1334template <
typename number2>
1339 const number2 * values,
1340 const bool elide_zero_values,
1345 if ((
std::abs(values[j]) != 0) || (elide_zero_values ==
false))
1346 add(row, col_indices[j], values[j]);
1351template <
typename number>
1357 return entry->
value;
1363template <
typename number>
1369 return entry->
value;
1375template <
typename number>
1383template <
typename number>
1390template <
typename number>
1399template <
typename number>
1408template <
typename number>
1409template <
typename MatrixType>
1412 const bool elide_zero_values)
1414 reinit(M.m(), M.n(), this->saved_default_row_length, this->increment);
1419 for (
size_type row = 0; row < M.m(); ++row)
1421 const typename MatrixType::const_iterator end_row = M.end(row);
1422 for (
typename MatrixType::const_iterator entry = M.begin(row);
1425 set(row, entry->column(), entry->value(), elide_zero_values);
1431template <
typename number>
1432template <
typename MatrixType>
1445 for (
size_type row = 0; row < M.m(); ++row)
1447 const typename MatrixType::const_iterator end_row = M.end(row);
1448 for (
typename MatrixType::const_iterator entry = M.begin(row);
1451 if (entry->value() != 0)
1452 add(row, entry->column(), factor * entry->value());
1458template <
typename number>
1459template <
typename MatrixTypeA,
typename MatrixTypeB>
1462 const MatrixTypeB &B,
1481 typename MatrixTypeB::const_iterator b1 = B.begin();
1482 const typename MatrixTypeB::const_iterator b_final = B.end();
1484 while (b1 != b_final)
1488 typename MatrixTypeB::const_iterator b2 = B.begin();
1489 while (b2 != b_final)
1494 const typename MatrixTypeA::value_type a = A.el(k, l);
1497 add(i, j, a * b1->value() * b2->value());
1508 std::vector<size_type> minrow(B.n(), B.m());
1509 std::vector<size_type> maxrow(B.n(), 0);
1510 while (b1 != b_final)
1513 if (r < minrow[b1->column()])
1514 minrow[b1->column()] = r;
1515 if (r > maxrow[b1->column()])
1516 maxrow[b1->column()] = r;
1520 typename MatrixTypeA::const_iterator ai = A.begin();
1521 const typename MatrixTypeA::const_iterator ae = A.end();
1525 const typename MatrixTypeA::value_type a = ai->value();
1535 b1 = B.begin(minrow[ai->row()]);
1536 const typename MatrixTypeB::const_iterator be1 =
1537 B.end(maxrow[ai->row()]);
1538 const typename MatrixTypeB::const_iterator be2 =
1539 B.end(maxrow[ai->column()]);
1543 const double b1v = b1->value();
1548 if (b1->column() == ai->row() && (b1v != 0.))
1552 typename MatrixTypeB::const_iterator b2 =
1553 B.begin(minrow[ai->column()]);
1556 if (b2->column() == ai->column())
1559 add(i, j, a * b1v * b2->value());
1572template <
typename number>
1573template <
class StreamType>
1580 std::vector<size_type> used_by_line;
1584 out <<
"SparseMatrixEZ:used entries:" << used << std::endl
1585 <<
"SparseMatrixEZ:allocated entries:" << allocated << std::endl
1586 <<
"SparseMatrixEZ:reserved entries:" << reserved << std::endl;
1590 for (
size_type i = 0; i < used_by_line.size(); ++i)
1591 if (used_by_line[i] != 0)
1592 out <<
"SparseMatrixEZ:entries\t" << i <<
"\trows\t"
1593 << used_by_line[i] << std::endl;
unsigned short index() const
const SparseMatrixEZ< number > * matrix
Accessor(const SparseMatrixEZ< number > *matrix, const size_type row, const unsigned short index)
const Accessor & operator*() const
const_iterator(const SparseMatrixEZ< number > *matrix, const size_type row, const unsigned short index)
bool operator<(const const_iterator &) const
bool operator==(const const_iterator &) const
bool operator!=(const const_iterator &) const
const_iterator & operator++()
const Accessor * operator->() const
void block_read(std::istream &in)
SparseMatrixEZ< number > & copy_from(const MatrixType &source, const bool elide_zero_values=true)
std::vector< Entry > data
void reinit(const size_type n_rows, const size_type n_columns, size_type default_row_length=0, unsigned int default_increment=1, size_type reserve=0)
void print_statistics(StreamType &s, bool full=false)
void Tvmult_add(Vector< somenumber > &dst, const Vector< somenumber > &src) const
void block_write(std::ostream &out) const
void compute_statistics(size_type &used, size_type &allocated, size_type &reserved, std::vector< size_type > &used_by_line, const bool compute_by_line) const
SparseMatrixEZ(const SparseMatrixEZ &)
number operator()(const size_type i, const size_type j) const
const Entry * locate(const size_type row, const size_type col) const
~SparseMatrixEZ() override=default
void print_formatted(std::ostream &out, const unsigned int precision=3, const bool scientific=true, const unsigned int width=0, const char *zero_string=" ", const double denominator=1.) const
Entry * allocate(const size_type row, const size_type col)
void threaded_matrix_scalar_product(const Vector< somenumber > &u, const Vector< somenumber > &v, const size_type begin_row, const size_type end_row, somenumber *partial_sum) const
size_type get_row_length(const size_type row) const
void precondition_TSOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const
std::size_t memory_consumption() const
size_type n_nonzero_elements() const
void print(std::ostream &out) const
SparseMatrixEZ< number > & operator=(const SparseMatrixEZ< number > &)
SparseMatrixEZ(const size_type n_rows, const size_type n_columns, const size_type default_row_length=0, const unsigned int default_increment=1)
void conjugate_add(const MatrixTypeA &A, const MatrixTypeB &B, const bool transpose=false)
void Tvmult(Vector< somenumber > &dst, const Vector< somenumber > &src) const
SparseMatrixEZ< number > & operator=(const double d)
void threaded_vmult(Vector< somenumber > &dst, const Vector< somenumber > &src, const size_type begin_row, const size_type end_row) const
void vmult_add(Vector< somenumber > &dst, const Vector< somenumber > &src) const
void precondition_Jacobi(Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const
const_iterator end() const
void precondition_SSOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1., const std::vector< std::size_t > &pos_right_of_diagonal=std::vector< std::size_t >()) const
std::vector< RowInfo > row_info
number el(const size_type i, const size_type j) const
const_iterator begin() const
void threaded_matrix_norm_square(const Vector< somenumber > &v, const size_type begin_row, const size_type end_row, somenumber *partial_sum) const
unsigned int saved_default_row_length
void set(const size_type i, const size_type j, const number value, const bool elide_zero_values=true)
void add(const size_type i, const size_type j, const number value)
void precondition_SOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const
void vmult(Vector< somenumber > &dst, const Vector< somenumber > &src) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DeclException0(Exception0)
static ::ExceptionBase & ExcInvalidEntry(int arg1, int arg2)
static ::ExceptionBase & ExcNoDiagonal()
#define Assert(cond, exc)
static ::ExceptionBase & ExcIteratorPastEnd()
#define AssertIsFinite(number)
#define DeclException2(Exception2, type1, type2, outsequence)
static ::ExceptionBase & ExcEntryAllocationFailure(int arg1, int arg2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
const types::global_dof_index invalid_size_type
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index
static const size_type invalid
static const unsigned short invalid_diagonal
RowInfo(const size_type start=Entry::invalid)