16 #ifndef dealii_sparse_matrix_ez_h 17 #define dealii_sparse_matrix_ez_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/subscriptor.h> 22 #include <deal.II/base/smartpointer.h> 23 #include <deal.II/lac/exceptions.h> 27 DEAL_II_NAMESPACE_OPEN
29 template <
typename number>
class Vector;
99 template <
typename number>
123 const number &
value);
167 static const unsigned short 191 const unsigned short index);
201 unsigned short index()
const;
211 number
value()
const;
241 const unsigned short index);
321 const unsigned int default_increment = 1);
353 unsigned int default_increment = 1,
405 template <
class StreamType>
420 std::vector<size_type> &used_by_line,
421 const bool compute_by_line)
const;
444 const number value,
const bool elide_zero_values =
true);
469 template <
typename number2>
470 void add (
const std::vector<size_type> &indices,
472 const bool elide_zero_values =
true);
479 template <
typename number2>
480 void add (
const std::vector<size_type> &row_indices,
481 const std::vector<size_type> &col_indices,
483 const bool elide_zero_values =
true);
494 template <
typename number2>
496 const std::vector<size_type> &col_indices,
497 const std::vector<number2> &values,
498 const bool elide_zero_values =
true);
509 template <
typename number2>
513 const number2 *values,
514 const bool elide_zero_values =
true,
515 const bool col_indices_are_sorted =
false);
538 template <
typename MatrixType>
541 const bool elide_zero_values =
true);
550 template <
typename MatrixType>
551 void add (
const number factor,
552 const MatrixType &matrix);
585 template <
typename somenumber>
586 void vmult (Vector<somenumber> &dst,
587 const Vector<somenumber> &src)
const;
594 template <
typename somenumber>
595 void Tvmult (Vector<somenumber> &dst,
596 const Vector<somenumber> &src)
const;
602 template <
typename somenumber>
604 const Vector<somenumber> &src)
const;
611 template <
typename somenumber>
613 const Vector<somenumber> &src)
const;
633 template <
typename somenumber>
635 const Vector<somenumber> &src,
636 const number omega = 1.)
const;
641 template <
typename somenumber>
643 const Vector<somenumber> &src,
644 const number om = 1.,
645 const std::vector<std::size_t> &pos_right_of_diagonal = std::vector<std::size_t>())
const;
651 template <
typename somenumber>
653 const Vector<somenumber> &src,
654 const number om = 1.)
const;
660 template <
typename somenumber>
662 const Vector<somenumber> &src,
663 const number om = 1.)
const;
673 template <
typename MatrixTypeA,
typename MatrixTypeB>
675 const MatrixTypeB &B,
712 void print (std::ostream &out)
const;
735 const unsigned int precision = 3,
736 const bool scientific =
true,
737 const unsigned int width = 0,
738 const char *zero_string =
" ",
739 const double denominator = 1.)
const;
776 <<
"The entry with index (" << arg1 <<
',' << arg2
777 <<
") does not exist.");
781 <<
"An entry with index (" << arg1 <<
',' << arg2
782 <<
") cannot be allocated.");
810 template <
typename somenumber>
812 const Vector<somenumber> &src,
821 template <
typename somenumber>
825 somenumber *partial_sum)
const;
832 template <
typename somenumber>
834 const Vector<somenumber> &v,
837 somenumber *partial_sum)
const;
870 template <
typename number>
881 template <
typename number>
890 template <
typename number>
896 diagonal(invalid_diagonal)
901 template <
typename number>
906 const unsigned short i)
914 template <
typename number>
923 template <
typename number>
928 return matrix->data[matrix->row_info[a_row].start+a_index].column;
932 template <
typename number>
942 template <
typename number>
947 return matrix->data[matrix->row_info[a_row].start+a_index].value;
951 template <
typename number>
956 const unsigned short i)
986 template <
typename number>
994 ++(accessor.a_index);
998 if (accessor.a_index >= accessor.matrix->row_info[accessor.a_row].length)
1000 accessor.a_index = 0;
1007 while (accessor.a_row < accessor.matrix->m()
1008 && accessor.matrix->row_info[accessor.a_row].length == 0);
1014 template <
typename number>
1023 template <
typename number>
1032 template <
typename number>
1043 template <
typename number>
1049 return ! (*
this == other);
1053 template <
typename number>
1066 template <
typename number>
1074 template <
typename number>
1082 template <
typename number>
1096 if (entry->
column == col)
1106 template <
typename number>
1113 return t->
locate(row,col);
1117 template <
typename number>
1136 while (i<
end &&
data[i].column < col) ++i;
1139 if (i !=
end &&
data[i].column == col)
1178 Entry temp = *entry;
1215 template <
typename number>
1220 const bool elide_zero_values)
1227 if (elide_zero_values && value == 0.)
1230 if (entry !=
nullptr)
1236 entry->
value = value;
1242 template <
typename number>
1259 entry->
value += value;
1263 template <
typename number>
1264 template <
typename number2>
1267 const bool elide_zero_values)
1270 for (
size_type i=0; i<indices.size(); ++i)
1271 for (
size_type j=0; j<indices.size(); ++j)
1272 if ((full_matrix(i,j) != 0) || (elide_zero_values ==
false))
1273 add (indices[i], indices[j], full_matrix(i,j));
1278 template <
typename number>
1279 template <
typename number2>
1281 const std::vector<size_type> &col_indices,
1283 const bool elide_zero_values)
1286 for (
size_type i=0; i<row_indices.size(); ++i)
1287 for (
size_type j=0; j<col_indices.size(); ++j)
1288 if ((full_matrix(i,j) != 0) || (elide_zero_values ==
false))
1289 add (row_indices[i], col_indices[j], full_matrix(i,j));
1295 template <
typename number>
1296 template <
typename number2>
1298 const std::vector<size_type> &col_indices,
1299 const std::vector<number2> &values,
1300 const bool elide_zero_values)
1303 for (
size_type j=0; j<col_indices.size(); ++j)
1304 if ((values[j] != 0) || (elide_zero_values ==
false))
1305 add (row, col_indices[j], values[j]);
1310 template <
typename number>
1311 template <
typename number2>
1315 const number2 *values,
1316 const bool elide_zero_values,
1321 if ((values[j] != 0) || (elide_zero_values ==
false))
1322 add (row, col_indices[j], values[j]);
1328 template <
typename number>
1335 return entry->
value;
1341 template <
typename number>
1348 return entry->
value;
1354 template <
typename number>
1363 template <
typename number>
1371 template <
typename number>
1381 template <
typename number>
1391 template <
typename number>
1392 template <
typename MatrixType>
1405 for (size_type row = 0; row < M.m(); ++row)
1407 const typename MatrixType::const_iterator end_row = M.end(row);
1408 for (
typename MatrixType::const_iterator entry = M.begin(row);
1409 entry != end_row; ++entry)
1410 set(row, entry->column(), entry->value(), elide_zero_values);
1416 template <
typename number>
1417 template <
typename MatrixType>
1421 const MatrixType &M)
1432 for (
size_type row = 0; row < M.m(); ++row)
1434 const typename MatrixType::const_iterator end_row = M.end(row);
1435 for (
typename MatrixType::const_iterator entry = M.begin(row);
1436 entry != end_row; ++entry)
1437 if (entry->value() != 0)
1438 add(row, entry->column(), factor * entry->value());
1444 template <
typename number>
1445 template <
typename MatrixTypeA,
typename MatrixTypeB>
1448 const MatrixTypeB &B,
1467 typename MatrixTypeB::const_iterator b1 = B.begin();
1468 const typename MatrixTypeB::const_iterator b_final = B.end();
1470 while (b1 != b_final)
1474 typename MatrixTypeB::const_iterator b2 = B.begin();
1475 while (b2 != b_final)
1480 const typename MatrixTypeA::value_type a = A.el(k,l);
1483 add (i, j, a * b1->value() * b2->value());
1494 std::vector<size_type> minrow(B.n(), B.m());
1495 std::vector<size_type> maxrow(B.n(), 0);
1496 while (b1 != b_final)
1499 if (r < minrow[b1->column()])
1500 minrow[b1->column()] = r;
1501 if (r > maxrow[b1->column()])
1502 maxrow[b1->column()] = r;
1506 typename MatrixTypeA::const_iterator ai = A.begin();
1507 const typename MatrixTypeA::const_iterator ae = A.end();
1511 const typename MatrixTypeA::value_type a = ai->value();
1514 if (a == 0.)
continue;
1520 b1 = B.begin(minrow[ai->row()]);
1521 const typename MatrixTypeB::const_iterator
1522 be1 = B.end(maxrow[ai->row()]);
1523 const typename MatrixTypeB::const_iterator
1524 be2 = B.end(maxrow[ai->column()]);
1528 const double b1v = b1->value();
1533 if (b1->column() == ai->row() && (b1v != 0.))
1537 typename MatrixTypeB::const_iterator
1538 b2 = B.begin(minrow[ai->column()]);
1541 if (b2->column() == ai->column())
1544 add (i, j, a * b1v * b2->value());
1557 template <
typename number>
1558 template <
class StreamType>
1566 std::vector<size_type> used_by_line;
1570 out <<
"SparseMatrixEZ:used entries:" << used << std::endl
1571 <<
"SparseMatrixEZ:allocated entries:" << allocated << std::endl
1572 <<
"SparseMatrixEZ:reserved entries:" << reserved << std::endl;
1576 for (
size_type i=0; i< used_by_line.size(); ++i)
1577 if (used_by_line[i] != 0)
1578 out <<
"SparseMatrixEZ:entries\t" << i
1579 <<
"\trows\t" << used_by_line[i]
1586 DEAL_II_NAMESPACE_CLOSE
size_type get_row_length(const size_type row) const
const types::global_dof_index invalid_size_type
void Tvmult_add(Vector< somenumber > &dst, const Vector< somenumber > &src) const
const_iterator & operator++()
static ::ExceptionBase & ExcEntryAllocationFailure(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
void add(const size_type i, const size_type j, const number value)
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
number operator()(const size_type i, const size_type j) const
SparseMatrixEZ< number > & copy_from(const MatrixType &source, const bool elide_zero_values=true)
void set(const size_type i, const size_type j, const number value, const bool elide_zero_values=true)
void conjugate_add(const MatrixTypeA &A, const MatrixTypeB &B, const bool transpose=false)
void block_read(std::istream &in)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
bool operator<(const const_iterator &) const
void vmult_add(Vector< somenumber > &dst, const Vector< somenumber > &src) const
void vmult(Vector< somenumber > &dst, const Vector< somenumber > &src) const
static const size_type invalid
const Accessor & operator*() const
std::size_t memory_consumption() const
void threaded_matrix_norm_square(const Vector< somenumber > &v, const size_type begin_row, const size_type end_row, somenumber *partial_sum) const
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
bool operator==(const const_iterator &) const
void print(std::ostream &out) const
static ::ExceptionBase & ExcNoDiagonal()
void precondition_SOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const
SparseMatrixEZ< number > & operator=(const SparseMatrixEZ< number > &)
Accessor(const SparseMatrixEZ< number > *matrix, const size_type row, const unsigned short index)
unsigned int global_dof_index
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
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)
bool operator!=(const const_iterator &) const
const_iterator(const SparseMatrixEZ< number > *matrix, const size_type row, const unsigned short index)
#define DeclException0(Exception0)
const SparseMatrixEZ< number > * matrix
std::vector< Entry > data
Entry * allocate(const size_type row, const size_type col)
SymmetricTensor< rank_, dim, Number > transpose(const SymmetricTensor< rank_, dim, Number > &t)
RowInfo(const size_type start=Entry::invalid)
void print_statistics(StreamType &s, bool full=false)
static ::ExceptionBase & ExcIteratorPastEnd()
void swap(Vector< Number > &u, Vector< Number > &v)
void block_write(std::ostream &out) const
types::global_dof_index size_type
const Accessor * operator->() const
void precondition_TSOR(Vector< somenumber > &dst, const Vector< somenumber > &src, const number om=1.) const
unsigned short index() const
number el(const size_type i, const size_type j) const
~SparseMatrixEZ()=default
const Entry * locate(const size_type row, const size_type col) const
const_iterator begin() const
static ::ExceptionBase & ExcInvalidEntry(int arg1, int arg2)
void precondition_Jacobi(Vector< somenumber > &dst, const Vector< somenumber > &src, const number omega=1.) const
void threaded_vmult(Vector< somenumber > &dst, const Vector< somenumber > &src, const size_type begin_row, const size_type end_row) const
void Tvmult(Vector< somenumber > &dst, const Vector< somenumber > &src) const
static const unsigned short invalid_diagonal
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
const_iterator end() const
#define AssertIsFinite(number)
unsigned int saved_default_row_length
std::vector< RowInfo > row_info
static ::ExceptionBase & ExcInternalError()
size_type n_nonzero_elements() 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