Reference documentation for deal.II version 8.5.1
sparse_matrix_ez.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2002 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__sparse_matrix_ez_h
17 #define dealii__sparse_matrix_ez_h
18 
19 
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>
24 
25 #include <vector>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 template<typename number> class Vector;
30 template<typename number> class FullMatrix;
31 
97 template <typename number>
98 class SparseMatrixEZ : public Subscriptor
99 {
100 public:
105 
110  struct Entry
111  {
115  Entry();
116 
120  Entry (const size_type column,
121  const number &value);
122 
127 
131  number value;
132 
137  };
138 
143  struct RowInfo
144  {
149 
157  unsigned short length;
161  unsigned short diagonal;
165  static const unsigned short
166  invalid_diagonal = static_cast<unsigned short>(-1);
167  };
168 
169 public:
170 
175  {
176  private:
180  class Accessor
181  {
182  public:
188  const size_type row,
189  const unsigned short index);
190 
194  size_type row() const;
195 
199  unsigned short index() const;
200 
204  size_type column() const;
205 
209  number value() const;
210 
211  protected:
216 
221 
225  unsigned short a_index;
226 
230  friend class const_iterator;
231  };
232 
233  public:
238  const size_type row,
239  const unsigned short index);
240 
245 
250 
254  const Accessor &operator* () const;
255 
259  const Accessor *operator-> () const;
260 
264  bool operator == (const const_iterator &) const;
268  bool operator != (const const_iterator &) const;
269 
274  bool operator < (const const_iterator &) const;
275 
276  private:
281 
289  };
290 
295  typedef number value_type;
296 
304  SparseMatrixEZ ();
305 
313  SparseMatrixEZ (const SparseMatrixEZ &);
314 
321  explicit SparseMatrixEZ (const size_type n_rows,
322  const size_type n_columns,
323  const size_type default_row_length = 0,
324  const unsigned int default_increment = 1);
325 
330  ~SparseMatrixEZ ();
331 
336 
345  SparseMatrixEZ<number> &operator = (const double d);
346 
354  void reinit (const size_type n_rows,
355  const size_type n_columns,
356  size_type default_row_length = 0,
357  unsigned int default_increment = 1,
358  size_type reserve = 0);
359 
364  void clear ();
366 
374  bool empty () const;
375 
380  size_type m () const;
381 
386  size_type n () const;
387 
391  size_type get_row_length (const size_type row) const;
392 
396  size_type n_nonzero_elements () const;
397 
402  std::size_t memory_consumption () const;
403 
409  template <class StreamType>
410  void print_statistics (StreamType &s, bool full = false);
411 
421  void compute_statistics (size_type &used,
422  size_type &allocated,
423  size_type &reserved,
424  std::vector<size_type> &used_by_line,
425  const bool compute_by_line) const;
427 
447  void set (const size_type i, const size_type j,
448  const number value, const bool elide_zero_values = true);
449 
455  void add (const size_type i,
456  const size_type j,
457  const number value);
458 
473  template <typename number2>
474  void add (const std::vector<size_type> &indices,
475  const FullMatrix<number2> &full_matrix,
476  const bool elide_zero_values = true);
477 
483  template <typename number2>
484  void add (const std::vector<size_type> &row_indices,
485  const std::vector<size_type> &col_indices,
486  const FullMatrix<number2> &full_matrix,
487  const bool elide_zero_values = true);
488 
498  template <typename number2>
499  void add (const size_type row,
500  const std::vector<size_type> &col_indices,
501  const std::vector<number2> &values,
502  const bool elide_zero_values = true);
503 
513  template <typename number2>
514  void add (const size_type row,
515  const size_type n_cols,
516  const size_type *col_indices,
517  const number2 *values,
518  const bool elide_zero_values = true,
519  const bool col_indices_are_sorted = false);
520 
542  template <typename MatrixType>
544  copy_from (const MatrixType &source,
545  const bool elide_zero_values = true);
546 
554  template <typename MatrixType>
555  void add (const number factor,
556  const MatrixType &matrix);
558 
571  number operator () (const size_type i,
572  const size_type j) const;
573 
578  number el (const size_type i,
579  const size_type j) const;
581 
589  template <typename somenumber>
590  void vmult (Vector<somenumber> &dst,
591  const Vector<somenumber> &src) const;
592 
598  template <typename somenumber>
599  void Tvmult (Vector<somenumber> &dst,
600  const Vector<somenumber> &src) const;
601 
606  template <typename somenumber>
607  void vmult_add (Vector<somenumber> &dst,
608  const Vector<somenumber> &src) const;
609 
615  template <typename somenumber>
616  void Tvmult_add (Vector<somenumber> &dst,
617  const Vector<somenumber> &src) const;
619 
626  number l2_norm () const;
628 
637  template <typename somenumber>
639  const Vector<somenumber> &src,
640  const number omega = 1.) const;
641 
645  template <typename somenumber>
647  const Vector<somenumber> &src,
648  const number om = 1.,
649  const std::vector<std::size_t> &pos_right_of_diagonal = std::vector<std::size_t>()) const;
650 
655  template <typename somenumber>
657  const Vector<somenumber> &src,
658  const number om = 1.) const;
659 
664  template <typename somenumber>
666  const Vector<somenumber> &src,
667  const number om = 1.) const;
668 
677  template <typename MatrixTypeA, typename MatrixTypeB>
678  void conjugate_add (const MatrixTypeA &A,
679  const MatrixTypeB &B,
680  const bool transpose = false);
682 
689  const_iterator begin () const;
690 
694  const_iterator end () const;
695 
700  const_iterator begin (const size_type r) const;
701 
706  const_iterator end (const size_type r) const;
708 
716  void print (std::ostream &out) const;
717 
738  void print_formatted (std::ostream &out,
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;
744 
750  void block_write (std::ostream &out) const;
751 
762  void block_read (std::istream &in);
764 
774 
779  int, int,
780  << "The entry with index (" << arg1 << ',' << arg2
781  << ") does not exist.");
782 
784  int, int,
785  << "An entry with index (" << arg1 << ',' << arg2
786  << ") cannot be allocated.");
788 private:
793  const Entry *locate (const size_type row,
794  const size_type col) const;
795 
800  Entry *locate (const size_type row,
801  const size_type col);
802 
806  Entry *allocate (const size_type row,
807  const size_type col);
808 
814  template <typename somenumber>
816  const Vector<somenumber> &src,
817  const size_type begin_row,
818  const size_type end_row) const;
819 
825  template <typename somenumber>
827  const size_type begin_row,
828  const size_type end_row,
829  somenumber *partial_sum) const;
830 
836  template <typename somenumber>
838  const Vector<somenumber> &v,
839  const size_type begin_row,
840  const size_type end_row,
841  somenumber *partial_sum) const;
842 
847 
851  std::vector<RowInfo> row_info;
852 
856  std::vector<Entry> data;
857 
861  unsigned int increment;
862 
867 };
868 
872 /*---------------------- Inline functions -----------------------------------*/
873 
874 template <typename number>
875 inline
877  const number &value)
878  :
879  column(column),
880  value(value)
881 {}
882 
883 
884 
885 template <typename number>
886 inline
888  :
889  column(invalid),
890  value(0)
891 {}
892 
893 
894 template <typename number>
895 inline
897  :
898  start(start),
899  length(0),
900  diagonal(invalid_diagonal)
901 {}
902 
903 
904 //---------------------------------------------------------------------------
905 template <typename number>
906 inline
909  const size_type r,
910  const unsigned short i)
911  :
912  matrix(matrix),
913  a_row(r),
914  a_index(i)
915 {}
916 
917 
918 template <typename number>
919 inline
922 {
923  return a_row;
924 }
925 
926 
927 template <typename number>
928 inline
931 {
932  return matrix->data[matrix->row_info[a_row].start+a_index].column;
933 }
934 
935 
936 template <typename number>
937 inline
938 unsigned short
940 {
941  return a_index;
942 }
943 
944 
945 
946 template <typename number>
947 inline
948 number
950 {
951  return matrix->data[matrix->row_info[a_row].start+a_index].value;
952 }
953 
954 
955 template <typename number>
956 inline
959  const size_type r,
960  const unsigned short i)
961  :
962  accessor(matrix, r, i)
963 {
964  // Finish if this is the end()
965  if (r==accessor.matrix->m() && i==0) return;
966 
967  // Make sure we never construct an
968  // iterator pointing to a
969  // non-existing entry
970 
971  // If the index points beyond the
972  // end of the row, try the next
973  // row.
974  if (accessor.a_index >= accessor.matrix->row_info[accessor.a_row].length)
975  {
976  do
977  {
978  ++accessor.a_row;
979  }
980  // Beware! If the next row is
981  // empty, iterate until a
982  // non-empty row is found or we
983  // hit the end of the matrix.
984  while (accessor.a_row < accessor.matrix->m()
985  && accessor.matrix->row_info[accessor.a_row].length == 0);
986  }
987 }
988 
989 
990 template <typename number>
991 inline
994 {
995  Assert (accessor.a_row < accessor.matrix->m(), ExcIteratorPastEnd());
996 
997  // Increment column index
998  ++(accessor.a_index);
999  // If index exceeds number of
1000  // entries in this row, proceed
1001  // with next row.
1002  if (accessor.a_index >= accessor.matrix->row_info[accessor.a_row].length)
1003  {
1004  accessor.a_index = 0;
1005  // Do this loop to avoid
1006  // elements in empty rows
1007  do
1008  {
1009  ++accessor.a_row;
1010  }
1011  while (accessor.a_row < accessor.matrix->m()
1012  && accessor.matrix->row_info[accessor.a_row].length == 0);
1013  }
1014  return *this;
1015 }
1016 
1017 
1018 template <typename number>
1019 inline
1022 {
1023  return accessor;
1024 }
1025 
1026 
1027 template <typename number>
1028 inline
1031 {
1032  return &accessor;
1033 }
1034 
1035 
1036 template <typename number>
1037 inline
1038 bool
1040  const const_iterator &other) const
1041 {
1042  return (accessor.row() == other.accessor.row() &&
1043  accessor.index() == other.accessor.index());
1044 }
1045 
1046 
1047 template <typename number>
1048 inline
1049 bool
1051 operator != (const const_iterator &other) const
1052 {
1053  return ! (*this == other);
1054 }
1055 
1056 
1057 template <typename number>
1058 inline
1059 bool
1061 operator < (const const_iterator &other) const
1062 {
1063  return (accessor.row() < other.accessor.row() ||
1064  (accessor.row() == other.accessor.row() &&
1065  accessor.index() < other.accessor.index()));
1066 }
1067 
1068 
1069 //---------------------------------------------------------------------------
1070 template <typename number>
1071 inline
1073 {
1074  return row_info.size();
1075 }
1076 
1077 
1078 template <typename number>
1079 inline
1081 {
1082  return n_columns;
1083 }
1084 
1085 
1086 template <typename number>
1087 inline
1090  const size_type col)
1091 {
1092  Assert (row<m(), ExcIndexRange(row,0,m()));
1093  Assert (col<n(), ExcIndexRange(col,0,n()));
1094 
1095  const RowInfo &r = row_info[row];
1096  const size_type end = r.start + r.length;
1097  for (size_type i=r.start; i<end; ++i)
1098  {
1099  Entry *const entry = &data[i];
1100  if (entry->column == col)
1101  return entry;
1102  if (entry->column == Entry::invalid)
1103  return 0;
1104  }
1105  return 0;
1106 }
1107 
1108 
1109 
1110 template <typename number>
1111 inline
1112 const typename SparseMatrixEZ<number>::Entry *
1114  const size_type col) const
1115 {
1116  SparseMatrixEZ<number> *t = const_cast<SparseMatrixEZ<number>*> (this);
1117  return t->locate(row,col);
1118 }
1119 
1120 
1121 template <typename number>
1122 inline
1125  const size_type col)
1126 {
1127  Assert (row<m(), ExcIndexRange(row,0,m()));
1128  Assert (col<n(), ExcIndexRange(col,0,n()));
1129 
1130  RowInfo &r = row_info[row];
1131  const size_type end = r.start + r.length;
1132 
1133  size_type i = r.start;
1134  // If diagonal exists and this
1135  // column is higher, start only
1136  // after diagonal.
1137  if (r.diagonal != RowInfo::invalid_diagonal && col >= row)
1138  i += r.diagonal;
1139  // Find position of entry
1140  while (i<end && data[i].column < col) ++i;
1141 
1142  // entry found
1143  if (i != end && data[i].column == col)
1144  return &data[i];
1145 
1146  // Now, we must insert the new
1147  // entry and move all successive
1148  // entries back.
1149 
1150  // If no more space is available
1151  // for this row, insert new
1152  // elements into the vector.
1153 //TODO:[GK] We should not extend this row if i<end
1154  if (row != row_info.size()-1)
1155  {
1156  if (end >= row_info[row+1].start)
1157  {
1158  // Failure if increment 0
1160 
1161  // Insert new entries
1162  data.insert(data.begin()+end, increment, Entry());
1163  // Update starts of
1164  // following rows
1165  for (size_type rn=row+1; rn<row_info.size(); ++rn)
1166  row_info[rn].start += increment;
1167  }
1168  }
1169  else
1170  {
1171  if (end >= data.size())
1172  {
1173  // Here, appending a block
1174  // does not increase
1175  // performance.
1176  data.push_back(Entry());
1177  }
1178  }
1179 
1180  Entry *entry = &data[i];
1181  // Save original entry
1182  Entry temp = *entry;
1183  // Insert new entry here to
1184  // make sure all entries
1185  // are ordered by column
1186  // index
1187  entry->column = col;
1188  entry->value = 0;
1189  // Update row_info
1190  ++r.length;
1191  if (col == row)
1192  r.diagonal = i - r.start;
1193  else if (col<row && r.diagonal!= RowInfo::invalid_diagonal)
1194  ++r.diagonal;
1195 
1196  if (i == end)
1197  return entry;
1198 
1199  // Move all entries in this
1200  // row up by one
1201  for (size_type j = i+1; j < end; ++j)
1202  {
1203  // There should be no invalid
1204  // entry below end
1205  Assert (data[j].column != Entry::invalid, ExcInternalError());
1206 
1207 //TODO[GK]: This could be done more efficiently by moving starting at the top rather than swapping starting at the bottom
1208  std::swap (data[j], temp);
1209  }
1210  Assert (data[end].column == Entry::invalid, ExcInternalError());
1211 
1212  data[end] = temp;
1213 
1214  return entry;
1215 }
1216 
1217 
1218 
1219 template <typename number>
1220 inline
1222  const size_type j,
1223  const number value,
1224  const bool elide_zero_values)
1225 {
1226  AssertIsFinite(value);
1227 
1228  Assert (i<m(), ExcIndexRange(i,0,m()));
1229  Assert (j<n(), ExcIndexRange(j,0,n()));
1230 
1231  if (elide_zero_values && value == 0.)
1232  {
1233  Entry *entry = locate(i,j);
1234  if (entry != 0)
1235  entry->value = 0.;
1236  }
1237  else
1238  {
1239  Entry *entry = allocate(i,j);
1240  entry->value = value;
1241  }
1242 }
1243 
1244 
1245 
1246 template <typename number>
1247 inline
1249  const size_type j,
1250  const number value)
1251 {
1252 
1253  AssertIsFinite(value);
1254 
1255  Assert (i<m(), ExcIndexRange(i,0,m()));
1256  Assert (j<n(), ExcIndexRange(j,0,n()));
1257 
1258  // ignore zero additions
1259  if (value == 0)
1260  return;
1261 
1262  Entry *entry = allocate(i,j);
1263  entry->value += value;
1264 }
1265 
1266 
1267 template <typename number>
1268 template <typename number2>
1269 void SparseMatrixEZ<number>::add (const std::vector<size_type> &indices,
1270  const FullMatrix<number2> &full_matrix,
1271  const bool elide_zero_values)
1272 {
1273 //TODO: This function can surely be made more efficient
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));
1278 }
1279 
1280 
1281 
1282 template <typename number>
1283 template <typename number2>
1284 void SparseMatrixEZ<number>::add (const std::vector<size_type> &row_indices,
1285  const std::vector<size_type> &col_indices,
1286  const FullMatrix<number2> &full_matrix,
1287  const bool elide_zero_values)
1288 {
1289 //TODO: This function can surely be made more efficient
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));
1294 }
1295 
1296 
1297 
1298 
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)
1305 {
1306 //TODO: This function can surely be made more efficient
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]);
1310 }
1311 
1312 
1313 
1314 template <typename number>
1315 template <typename number2>
1317  const size_type n_cols,
1318  const size_type *col_indices,
1319  const number2 *values,
1320  const bool elide_zero_values,
1321  const bool /*col_indices_are_sorted*/)
1322 {
1323 //TODO: This function can surely be made more efficient
1324  for (size_type j=0; j<n_cols; ++j)
1325  if ((values[j] != 0) || (elide_zero_values == false))
1326  add (row, col_indices[j], values[j]);
1327 }
1328 
1329 
1330 
1331 
1332 template <typename number>
1333 inline
1335  const size_type j) const
1336 {
1337  const Entry *entry = locate(i,j);
1338  if (entry)
1339  return entry->value;
1340  return 0.;
1341 }
1342 
1343 
1344 
1345 template <typename number>
1346 inline
1348  const size_type j) const
1349 {
1350  const Entry *entry = locate(i,j);
1351  if (entry)
1352  return entry->value;
1353  Assert(false, ExcInvalidEntry(i,j));
1354  return 0.;
1355 }
1356 
1357 
1358 template <typename number>
1359 inline
1362 {
1363  const_iterator result(this, 0, 0);
1364  return result;
1365 }
1366 
1367 template <typename number>
1368 inline
1371 {
1372  return const_iterator(this, m(), 0);
1373 }
1374 
1375 template <typename number>
1376 inline
1379 {
1380  Assert (r<m(), ExcIndexRange(r,0,m()));
1381  const_iterator result (this, r, 0);
1382  return result;
1383 }
1384 
1385 template <typename number>
1386 inline
1389 {
1390  Assert (r<m(), ExcIndexRange(r,0,m()));
1391  const_iterator result(this, r+1, 0);
1392  return result;
1393 }
1394 
1395 template<typename number>
1396 template <typename MatrixType>
1397 inline
1399 SparseMatrixEZ<number>::copy_from (const MatrixType &M, const bool elide_zero_values)
1400 {
1401  reinit(M.m(),
1402  M.n(),
1404  this->increment);
1405 
1406  // loop over the elements of the argument matrix row by row, as suggested
1407  // in the documentation of the sparse matrix iterator class, and
1408  // copy them into the current object
1409  for (size_type row = 0; row < M.m(); ++row)
1410  {
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);
1415  }
1416 
1417  return *this;
1418 }
1419 
1420 template<typename number>
1421 template <typename MatrixType>
1422 inline
1423 void
1424 SparseMatrixEZ<number>::add (const number factor,
1425  const MatrixType &M)
1426 {
1427  Assert (M.m() == m(), ExcDimensionMismatch(M.m(), m()));
1428  Assert (M.n() == n(), ExcDimensionMismatch(M.n(), n()));
1429 
1430  if (factor == 0.)
1431  return;
1432 
1433  // loop over the elements of the argument matrix row by row, as suggested
1434  // in the documentation of the sparse matrix iterator class, and
1435  // add them into the current object
1436  for (size_type row = 0; row < M.m(); ++row)
1437  {
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());
1443  }
1444 }
1445 
1446 
1447 
1448 template<typename number>
1449 template <typename MatrixTypeA, typename MatrixTypeB>
1450 inline void
1452  const MatrixTypeB &B,
1453  const bool transpose)
1454 {
1455 // Compute the result
1456 // r_ij = \sum_kl b_ik b_jl a_kl
1457 
1458 // Assert (n() == B.m(), ExcDimensionMismatch(n(), B.m()));
1459 // Assert (m() == B.m(), ExcDimensionMismatch(m(), B.m()));
1460 // Assert (A.n() == B.n(), ExcDimensionMismatch(A.n(), B.n()));
1461 // Assert (A.m() == B.n(), ExcDimensionMismatch(A.m(), B.n()));
1462 
1463  // Somehow, we have to avoid making
1464  // this an operation of complexity
1465  // n^2. For the transpose case, we
1466  // can go through the non-zero
1467  // elements of A^-1 and use the
1468  // corresponding rows of B only.
1469  // For the non-transpose case, we
1470  // must find a trick.
1471  typename MatrixTypeB::const_iterator b1 = B.begin();
1472  const typename MatrixTypeB::const_iterator b_final = B.end();
1473  if (transpose)
1474  while (b1 != b_final)
1475  {
1476  const size_type i = b1->column();
1477  const size_type k = b1->row();
1478  typename MatrixTypeB::const_iterator b2 = B.begin();
1479  while (b2 != b_final)
1480  {
1481  const size_type j = b2->column();
1482  const size_type l = b2->row();
1483 
1484  const typename MatrixTypeA::value_type a = A.el(k,l);
1485 
1486  if (a != 0.)
1487  add (i, j, a * b1->value() * b2->value());
1488  ++b2;
1489  }
1490  ++b1;
1491  }
1492  else
1493  {
1494  // Determine minimal and
1495  // maximal row for a column in
1496  // advance.
1497 
1498  std::vector<size_type> minrow(B.n(), B.m());
1499  std::vector<size_type> maxrow(B.n(), 0);
1500  while (b1 != b_final)
1501  {
1502  const size_type r = b1->row();
1503  if (r < minrow[b1->column()])
1504  minrow[b1->column()] = r;
1505  if (r > maxrow[b1->column()])
1506  maxrow[b1->column()] = r;
1507  ++b1;
1508  }
1509 
1510  typename MatrixTypeA::const_iterator ai = A.begin();
1511  const typename MatrixTypeA::const_iterator ae = A.end();
1512 
1513  while (ai != ae)
1514  {
1515  const typename MatrixTypeA::value_type a = ai->value();
1516  // Don't do anything if
1517  // this entry is zero.
1518  if (a == 0.) continue;
1519 
1520  // Now, loop over all rows
1521  // having possibly a
1522  // nonzero entry in column
1523  // ai->row()
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()]);
1529 
1530  while (b1 != be1)
1531  {
1532  const double b1v = b1->value();
1533  // We need the product
1534  // of both. If it is
1535  // zero, we can save
1536  // the work
1537  if (b1->column() == ai->row() && (b1v != 0.))
1538  {
1539  const size_type i = b1->row();
1540 
1541  typename MatrixTypeB::const_iterator
1542  b2 = B.begin(minrow[ai->column()]);
1543  while (b2 != be2)
1544  {
1545  if (b2->column() == ai->column())
1546  {
1547  const size_type j = b2->row();
1548  add (i, j, a * b1v * b2->value());
1549  }
1550  ++b2;
1551  }
1552  }
1553  ++b1;
1554  }
1555  ++ai;
1556  }
1557  }
1558 }
1559 
1560 
1561 template <typename number>
1562 template <class StreamType>
1563 inline
1564 void
1565 SparseMatrixEZ<number>::print_statistics(StreamType &out, bool full)
1566 {
1567  size_type used;
1568  size_type allocated;
1569  size_type reserved;
1570  std::vector<size_type> used_by_line;
1571 
1572  compute_statistics (used, allocated, reserved, used_by_line, full);
1573 
1574  out << "SparseMatrixEZ:used entries:" << used << std::endl
1575  << "SparseMatrixEZ:allocated entries:" << allocated << std::endl
1576  << "SparseMatrixEZ:reserved entries:" << reserved << std::endl;
1577 
1578  if (full)
1579  {
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]
1584  << std::endl;
1585 
1586  }
1587 }
1588 
1589 
1590 DEAL_II_NAMESPACE_CLOSE
1591 
1592 #endif
1593 /*---------------------------- sparse_matrix.h ---------------------------*/
size_type get_row_length(const size_type row) const
const types::global_dof_index invalid_size_type
Definition: types.h:179
void Tvmult_add(Vector< somenumber > &dst, const Vector< somenumber > &src) const
static ::ExceptionBase & ExcEntryAllocationFailure(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
void add(const size_type i, const size_type j, const number value)
unsigned int increment
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)
size_type n() const
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
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:313
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)
Definition: exceptions.h:541
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
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
bool empty() const
static const unsigned short invalid_diagonal
number l2_norm() 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
const_iterator end() const
#define AssertIsFinite(number)
Definition: exceptions.h:1197
unsigned int saved_default_row_length
size_type m() const
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