15#ifndef dealii_sparse_matrix_ez_h
16#define dealii_sparse_matrix_ez_h
32template <
typename number>
34template <
typename number>
102template <
typename number>
170 static_cast<unsigned short>(-1);
192 const unsigned short index);
244 const unsigned short index);
322 const unsigned int default_increment = 1);
357 unsigned int default_increment = 1,
416 template <
typename StreamType>
433 std::vector<size_type> &used_by_line,
434 const bool compute_by_line)
const;
461 const bool elide_zero_values =
true);
490 template <
typename number2>
492 add(
const std::vector<size_type> &indices,
494 const bool elide_zero_values =
true);
501 template <
typename number2>
503 add(
const std::vector<size_type> &row_indices,
504 const std::vector<size_type> &col_indices,
506 const bool elide_zero_values =
true);
517 template <
typename number2>
520 const std::vector<size_type> &col_indices,
521 const std::vector<number2> &values,
522 const bool elide_zero_values =
true);
533 template <
typename number2>
538 const number2 *values,
539 const bool elide_zero_values =
true,
540 const bool col_indices_are_sorted =
false);
563 template <
typename MatrixType>
565 copy_from(
const MatrixType &source,
const bool elide_zero_values =
true);
574 template <
typename MatrixType>
576 add(
const number factor,
const MatrixType &matrix);
609 template <
typename somenumber>
618 template <
typename somenumber>
626 template <
typename somenumber>
635 template <
typename somenumber>
658 template <
typename somenumber>
662 const number omega = 1.)
const;
667 template <
typename somenumber>
671 const number om = 1.,
672 const std::vector<std::size_t> &pos_right_of_diagonal =
673 std::vector<std::size_t>())
const;
679 template <
typename somenumber>
683 const number om = 1.)
const;
689 template <
typename somenumber>
693 const number om = 1.)
const;
703 template <
typename MatrixTypeA,
typename MatrixTypeB>
706 const MatrixTypeB &B,
774 const unsigned int precision = 3,
775 const bool scientific =
true,
776 const unsigned int width = 0,
777 const char *zero_string =
" ",
778 const double denominator = 1.,
779 const char *separator =
" ")
const;
819 <<
"The entry with index (" << arg1 <<
',' << arg2
820 <<
") does not exist.");
825 <<
"An entry with index (" << arg1 <<
',' << arg2
826 <<
") cannot be allocated.");
854 template <
typename somenumber>
866 template <
typename somenumber>
871 somenumber *partial_sum)
const;
878 template <
typename somenumber>
884 somenumber *partial_sum)
const;
917template <
typename number>
926template <
typename number>
933template <
typename number>
937 , diagonal(invalid_diagonal)
942template <
typename number>
946 const unsigned short i)
953template <
typename number>
961template <
typename number>
965 return matrix->data[matrix->row_info[a_row].start + a_index].column;
969template <
typename number>
978template <
typename number>
982 return matrix->data[matrix->row_info[a_row].start + a_index].value;
986template <
typename number>
990 const unsigned short i)
1020template <
typename number>
1027 ++(accessor.a_index);
1031 if (accessor.a_index >= accessor.matrix->row_info[accessor.a_row].length)
1033 accessor.a_index = 0;
1040 while (accessor.a_row < accessor.matrix->m() &&
1041 accessor.matrix->row_info[accessor.a_row].length == 0);
1047template <
typename number>
1055template <
typename number>
1063template <
typename number>
1073template <
typename number>
1078 return !(*
this == other);
1082template <
typename number>
1094template <
typename number>
1102template <
typename number>
1110template <
typename number>
1122 if (entry->
column == col)
1132template <
typename number>
1137 return t->
locate(row, col);
1141template <
typename number>
1158 while (i <
end &&
data[i].column < col)
1162 if (i !=
end &&
data[i].column == col)
1201 Entry temp = *entry;
1228 std::swap(
data[j], temp);
1239template <
typename number>
1244 const bool elide_zero_values)
1251 if (elide_zero_values && value == 0.)
1254 if (entry !=
nullptr)
1260 entry->
value = value;
1266template <
typename number>
1282 entry->
value += value;
1286template <
typename number>
1287template <
typename number2>
1291 const bool elide_zero_values)
1294 for (
size_type i = 0; i < indices.size(); ++i)
1295 for (
size_type j = 0; j < indices.size(); ++j)
1296 if ((full_matrix(i, j) != 0) || (elide_zero_values ==
false))
1297 add(indices[i], indices[j], full_matrix(i, j));
1302template <
typename number>
1303template <
typename number2>
1306 const std::vector<size_type> &col_indices,
1308 const bool elide_zero_values)
1311 for (
size_type i = 0; i < row_indices.size(); ++i)
1312 for (
size_type j = 0; j < col_indices.size(); ++j)
1313 if ((full_matrix(i, j) != 0) || (elide_zero_values ==
false))
1314 add(row_indices[i], col_indices[j], full_matrix(i, j));
1319template <
typename number>
1320template <
typename number2>
1323 const std::vector<size_type> &col_indices,
1324 const std::vector<number2> &values,
1325 const bool elide_zero_values)
1328 for (
size_type j = 0; j < col_indices.size(); ++j)
1329 if ((values[j] != 0) || (elide_zero_values ==
false))
1330 add(row, col_indices[j], values[j]);
1335template <
typename number>
1336template <
typename number2>
1341 const number2 *values,
1342 const bool elide_zero_values,
1347 if ((
std::abs(values[j]) != 0) || (elide_zero_values ==
false))
1348 add(row, col_indices[j], values[j]);
1353template <
typename number>
1359 return entry->
value;
1365template <
typename number>
1371 return entry->
value;
1377template <
typename number>
1385template <
typename number>
1392template <
typename number>
1401template <
typename number>
1410template <
typename number>
1411template <
typename MatrixType>
1414 const bool elide_zero_values)
1416 reinit(M.m(), M.n(), this->saved_default_row_length, this->increment);
1421 for (
size_type row = 0; row < M.m(); ++row)
1423 const typename MatrixType::const_iterator end_row = M.end(row);
1424 for (
typename MatrixType::const_iterator entry = M.begin(row);
1427 set(row, entry->column(), entry->value(), elide_zero_values);
1433template <
typename number>
1434template <
typename MatrixType>
1447 for (
size_type row = 0; row < M.m(); ++row)
1449 const typename MatrixType::const_iterator end_row = M.end(row);
1450 for (
typename MatrixType::const_iterator entry = M.begin(row);
1453 if (entry->value() != 0)
1454 add(row, entry->column(), factor * entry->value());
1460template <
typename number>
1461template <
typename MatrixTypeA,
typename MatrixTypeB>
1464 const MatrixTypeB &B,
1483 typename MatrixTypeB::const_iterator b1 = B.begin();
1484 const typename MatrixTypeB::const_iterator b_final = B.end();
1486 while (b1 != b_final)
1490 typename MatrixTypeB::const_iterator b2 = B.begin();
1491 while (b2 != b_final)
1496 const typename MatrixTypeA::value_type a = A.el(k, l);
1499 add(i, j, a * b1->value() * b2->value());
1510 std::vector<size_type> minrow(B.n(), B.m());
1511 std::vector<size_type> maxrow(B.n(), 0);
1512 while (b1 != b_final)
1515 if (r < minrow[b1->column()])
1516 minrow[b1->column()] = r;
1517 if (r > maxrow[b1->column()])
1518 maxrow[b1->column()] = r;
1522 typename MatrixTypeA::const_iterator ai = A.begin();
1523 const typename MatrixTypeA::const_iterator ae = A.end();
1527 const typename MatrixTypeA::value_type a = ai->value();
1537 b1 = B.begin(minrow[ai->row()]);
1538 const typename MatrixTypeB::const_iterator be1 =
1539 B.end(maxrow[ai->row()]);
1540 const typename MatrixTypeB::const_iterator be2 =
1541 B.end(maxrow[ai->column()]);
1545 const double b1v = b1->value();
1550 if (b1->column() == ai->row() && (b1v != 0.))
1554 typename MatrixTypeB::const_iterator b2 =
1555 B.begin(minrow[ai->column()]);
1558 if (b2->column() == ai->column())
1561 add(i, j, a * b1v * b2->value());
1574template <
typename number>
1575template <
typename StreamType>
1582 std::vector<size_type> used_by_line;
1586 out <<
"SparseMatrixEZ:used entries:" << used << std::endl
1587 <<
"SparseMatrixEZ:allocated entries:" << allocated << std::endl
1588 <<
"SparseMatrixEZ:reserved entries:" << reserved << std::endl;
1592 for (
size_type i = 0; i < used_by_line.size(); ++i)
1593 if (used_by_line[i] != 0)
1594 out <<
"SparseMatrixEZ:entries\t" << i <<
"\trows\t"
1595 << 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 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 char *separator=" ") 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
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)