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;
97 template <
typename number>
121 const number &
value);
165 static const unsigned short 189 const unsigned short index);
199 unsigned short index()
const;
209 number
value()
const;
239 const unsigned short index);
324 const unsigned int default_increment = 1);
357 unsigned int default_increment = 1,
409 template <
class StreamType>
424 std::vector<size_type> &used_by_line,
425 const bool compute_by_line)
const;
448 const number value,
const bool elide_zero_values =
true);
473 template <
typename number2>
474 void add (
const std::vector<size_type> &indices,
476 const bool elide_zero_values =
true);
483 template <
typename number2>
484 void add (
const std::vector<size_type> &row_indices,
485 const std::vector<size_type> &col_indices,
487 const bool elide_zero_values =
true);
498 template <
typename number2>
500 const std::vector<size_type> &col_indices,
501 const std::vector<number2> &values,
502 const bool elide_zero_values =
true);
513 template <
typename number2>
517 const number2 *values,
518 const bool elide_zero_values =
true,
519 const bool col_indices_are_sorted =
false);
542 template <
typename MatrixType>
545 const bool elide_zero_values =
true);
554 template <
typename MatrixType>
555 void add (
const number factor,
556 const MatrixType &matrix);
589 template <
typename somenumber>
598 template <
typename somenumber>
606 template <
typename somenumber>
615 template <
typename somenumber>
637 template <
typename somenumber>
640 const number omega = 1.)
const;
645 template <
typename somenumber>
648 const number om = 1.,
649 const std::vector<std::size_t> &pos_right_of_diagonal = std::vector<std::size_t>())
const;
655 template <
typename somenumber>
658 const number om = 1.)
const;
664 template <
typename somenumber>
667 const number om = 1.)
const;
677 template <
typename MatrixTypeA,
typename MatrixTypeB>
679 const MatrixTypeB &B,
680 const bool transpose =
false);
716 void print (std::ostream &out)
const;
739 const unsigned int precision = 3,
740 const bool scientific =
true,
741 const unsigned int width = 0,
742 const char *zero_string =
" ",
743 const double denominator = 1.)
const;
780 <<
"The entry with index (" << arg1 <<
',' << arg2
781 <<
") does not exist.");
785 <<
"An entry with index (" << arg1 <<
',' << arg2
786 <<
") cannot be allocated.");
814 template <
typename somenumber>
825 template <
typename somenumber>
829 somenumber *partial_sum)
const;
836 template <
typename somenumber>
841 somenumber *partial_sum)
const;
874 template <
typename number>
885 template <
typename number>
894 template <
typename number>
900 diagonal(invalid_diagonal)
905 template <
typename number>
910 const unsigned short i)
918 template <
typename number>
927 template <
typename number>
932 return matrix->data[matrix->row_info[a_row].start+a_index].column;
936 template <
typename number>
946 template <
typename number>
951 return matrix->data[matrix->row_info[a_row].start+a_index].value;
955 template <
typename number>
960 const unsigned short i)
990 template <
typename number>
998 ++(accessor.a_index);
1002 if (accessor.a_index >= accessor.matrix->row_info[accessor.a_row].length)
1004 accessor.a_index = 0;
1011 while (accessor.a_row < accessor.matrix->m()
1012 && accessor.matrix->row_info[accessor.a_row].length == 0);
1018 template <
typename number>
1027 template <
typename number>
1036 template <
typename number>
1047 template <
typename number>
1053 return ! (*
this == other);
1057 template <
typename number>
1070 template <
typename number>
1078 template <
typename number>
1086 template <
typename number>
1100 if (entry->
column == col)
1110 template <
typename number>
1117 return t->
locate(row,col);
1121 template <
typename number>
1140 while (i<
end &&
data[i].column < col) ++i;
1143 if (i !=
end &&
data[i].column == col)
1182 Entry temp = *entry;
1208 std::swap (
data[j], temp);
1219 template <
typename number>
1224 const bool elide_zero_values)
1231 if (elide_zero_values && value == 0.)
1240 entry->
value = value;
1246 template <
typename number>
1263 entry->
value += value;
1267 template <
typename number>
1268 template <
typename number2>
1271 const bool elide_zero_values)
1274 for (
size_type i=0; i<indices.size(); ++i)
1275 for (
size_type j=0; j<indices.size(); ++j)
1276 if ((full_matrix(i,j) != 0) || (elide_zero_values ==
false))
1277 add (indices[i], indices[j], full_matrix(i,j));
1282 template <
typename number>
1283 template <
typename number2>
1285 const std::vector<size_type> &col_indices,
1287 const bool elide_zero_values)
1290 for (
size_type i=0; i<row_indices.size(); ++i)
1291 for (
size_type j=0; j<col_indices.size(); ++j)
1292 if ((full_matrix(i,j) != 0) || (elide_zero_values ==
false))
1293 add (row_indices[i], col_indices[j], full_matrix(i,j));
1299 template <
typename number>
1300 template <
typename number2>
1302 const std::vector<size_type> &col_indices,
1303 const std::vector<number2> &values,
1304 const bool elide_zero_values)
1307 for (
size_type j=0; j<col_indices.size(); ++j)
1308 if ((values[j] != 0) || (elide_zero_values ==
false))
1309 add (row, col_indices[j], values[j]);
1314 template <
typename number>
1315 template <
typename number2>
1319 const number2 *values,
1320 const bool elide_zero_values,
1325 if ((values[j] != 0) || (elide_zero_values ==
false))
1326 add (row, col_indices[j], values[j]);
1332 template <
typename number>
1339 return entry->
value;
1345 template <
typename number>
1352 return entry->
value;
1358 template <
typename number>
1367 template <
typename number>
1375 template <
typename number>
1385 template <
typename number>
1395 template<
typename number>
1396 template <
typename MatrixType>
1409 for (size_type row = 0; row < M.m(); ++row)
1411 const typename MatrixType::const_iterator end_row = M.end(row);
1412 for (
typename MatrixType::const_iterator entry = M.begin(row);
1413 entry != end_row; ++entry)
1414 set(row, entry->column(), entry->value(), elide_zero_values);
1420 template<
typename number>
1421 template <
typename MatrixType>
1425 const MatrixType &M)
1436 for (
size_type row = 0; row < M.m(); ++row)
1438 const typename MatrixType::const_iterator end_row = M.end(row);
1439 for (
typename MatrixType::const_iterator entry = M.begin(row);
1440 entry != end_row; ++entry)
1441 if (entry->value() != 0)
1442 add(row, entry->column(), factor * entry->value());
1448 template<
typename number>
1449 template <
typename MatrixTypeA,
typename MatrixTypeB>
1452 const MatrixTypeB &B,
1453 const bool transpose)
1471 typename MatrixTypeB::const_iterator b1 = B.begin();
1472 const typename MatrixTypeB::const_iterator b_final = B.end();
1474 while (b1 != b_final)
1478 typename MatrixTypeB::const_iterator b2 = B.begin();
1479 while (b2 != b_final)
1484 const typename MatrixTypeA::value_type a = A.el(k,l);
1487 add (i, j, a * b1->value() * b2->value());
1498 std::vector<size_type> minrow(B.n(), B.m());
1499 std::vector<size_type> maxrow(B.n(), 0);
1500 while (b1 != b_final)
1503 if (r < minrow[b1->column()])
1504 minrow[b1->column()] = r;
1505 if (r > maxrow[b1->column()])
1506 maxrow[b1->column()] = r;
1510 typename MatrixTypeA::const_iterator ai = A.begin();
1511 const typename MatrixTypeA::const_iterator ae = A.end();
1515 const typename MatrixTypeA::value_type a = ai->value();
1518 if (a == 0.)
continue;
1524 b1 = B.begin(minrow[ai->row()]);
1525 const typename MatrixTypeB::const_iterator
1526 be1 = B.end(maxrow[ai->row()]);
1527 const typename MatrixTypeB::const_iterator
1528 be2 = B.end(maxrow[ai->column()]);
1532 const double b1v = b1->value();
1537 if (b1->column() == ai->row() && (b1v != 0.))
1541 typename MatrixTypeB::const_iterator
1542 b2 = B.begin(minrow[ai->column()]);
1545 if (b2->column() == ai->column())
1548 add (i, j, a * b1v * b2->value());
1561 template <
typename number>
1562 template <
class StreamType>
1570 std::vector<size_type> used_by_line;
1574 out <<
"SparseMatrixEZ:used entries:" << used << std::endl
1575 <<
"SparseMatrixEZ:allocated entries:" << allocated << std::endl
1576 <<
"SparseMatrixEZ:reserved entries:" << reserved << std::endl;
1580 for (
size_type i=0; i< used_by_line.size(); ++i)
1581 if (used_by_line[i] != 0)
1582 out <<
"SparseMatrixEZ:entries\t" << i
1583 <<
"\trows\t" << used_by_line[i]
1590 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)
void print_statistics(StreamType &s, bool full=false)
static ::ExceptionBase & ExcIteratorPastEnd()
RowInfo(size_type start=Entry::invalid)
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
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