Reference documentation for deal.II version 9.0.0
sparse_matrix_ez.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2002 - 2017 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 
99 template <typename number>
100 class SparseMatrixEZ : public Subscriptor
101 {
102 public:
107 
112  struct Entry
113  {
117  Entry();
118 
122  Entry (const size_type column,
123  const number &value);
124 
129 
133  number value;
134 
139  };
140 
145  struct RowInfo
146  {
151 
159  unsigned short length;
163  unsigned short diagonal;
167  static const unsigned short
168  invalid_diagonal = static_cast<unsigned short>(-1);
169  };
170 
171 public:
172 
177  {
178  private:
182  class Accessor
183  {
184  public:
190  const size_type row,
191  const unsigned short index);
192 
196  size_type row() const;
197 
201  unsigned short index() const;
202 
206  size_type column() const;
207 
211  number value() const;
212 
213  protected:
218 
223 
227  unsigned short a_index;
228 
232  friend class const_iterator;
233  };
234 
235  public:
240  const size_type row,
241  const unsigned short index);
242 
247 
251  const Accessor &operator* () const;
252 
256  const Accessor *operator-> () const;
257 
261  bool operator == (const const_iterator &) const;
265  bool operator != (const const_iterator &) const;
266 
271  bool operator < (const const_iterator &) const;
272 
273  private:
278 
286  };
287 
292  typedef number value_type;
293 
301  SparseMatrixEZ ();
302 
310  SparseMatrixEZ (const SparseMatrixEZ &);
311 
318  explicit SparseMatrixEZ (const size_type n_rows,
319  const size_type n_columns,
320  const size_type default_row_length = 0,
321  const unsigned int default_increment = 1);
322 
326  ~SparseMatrixEZ () = default;
327 
332 
341  SparseMatrixEZ<number> &operator = (const double d);
342 
350  void reinit (const size_type n_rows,
351  const size_type n_columns,
352  size_type default_row_length = 0,
353  unsigned int default_increment = 1,
354  size_type reserve = 0);
355 
360  void clear ();
362 
370  bool empty () const;
371 
376  size_type m () const;
377 
382  size_type n () const;
383 
387  size_type get_row_length (const size_type row) const;
388 
392  size_type n_nonzero_elements () const;
393 
398  std::size_t memory_consumption () const;
399 
405  template <class StreamType>
406  void print_statistics (StreamType &s, bool full = false);
407 
417  void compute_statistics (size_type &used,
418  size_type &allocated,
419  size_type &reserved,
420  std::vector<size_type> &used_by_line,
421  const bool compute_by_line) const;
423 
443  void set (const size_type i, const size_type j,
444  const number value, const bool elide_zero_values = true);
445 
451  void add (const size_type i,
452  const size_type j,
453  const number value);
454 
469  template <typename number2>
470  void add (const std::vector<size_type> &indices,
471  const FullMatrix<number2> &full_matrix,
472  const bool elide_zero_values = true);
473 
479  template <typename number2>
480  void add (const std::vector<size_type> &row_indices,
481  const std::vector<size_type> &col_indices,
482  const FullMatrix<number2> &full_matrix,
483  const bool elide_zero_values = true);
484 
494  template <typename number2>
495  void add (const size_type row,
496  const std::vector<size_type> &col_indices,
497  const std::vector<number2> &values,
498  const bool elide_zero_values = true);
499 
509  template <typename number2>
510  void add (const size_type row,
511  const size_type n_cols,
512  const size_type *col_indices,
513  const number2 *values,
514  const bool elide_zero_values = true,
515  const bool col_indices_are_sorted = false);
516 
538  template <typename MatrixType>
540  copy_from (const MatrixType &source,
541  const bool elide_zero_values = true);
542 
550  template <typename MatrixType>
551  void add (const number factor,
552  const MatrixType &matrix);
554 
567  number operator () (const size_type i,
568  const size_type j) const;
569 
574  number el (const size_type i,
575  const size_type j) const;
577 
585  template <typename somenumber>
586  void vmult (Vector<somenumber> &dst,
587  const Vector<somenumber> &src) const;
588 
594  template <typename somenumber>
595  void Tvmult (Vector<somenumber> &dst,
596  const Vector<somenumber> &src) const;
597 
602  template <typename somenumber>
603  void vmult_add (Vector<somenumber> &dst,
604  const Vector<somenumber> &src) const;
605 
611  template <typename somenumber>
612  void Tvmult_add (Vector<somenumber> &dst,
613  const Vector<somenumber> &src) const;
615 
622  number l2_norm () const;
624 
633  template <typename somenumber>
634  void precondition_Jacobi (Vector<somenumber> &dst,
635  const Vector<somenumber> &src,
636  const number omega = 1.) const;
637 
641  template <typename somenumber>
642  void precondition_SSOR (Vector<somenumber> &dst,
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;
646 
651  template <typename somenumber>
652  void precondition_SOR (Vector<somenumber> &dst,
653  const Vector<somenumber> &src,
654  const number om = 1.) const;
655 
660  template <typename somenumber>
661  void precondition_TSOR (Vector<somenumber> &dst,
662  const Vector<somenumber> &src,
663  const number om = 1.) const;
664 
673  template <typename MatrixTypeA, typename MatrixTypeB>
674  void conjugate_add (const MatrixTypeA &A,
675  const MatrixTypeB &B,
676  const bool transpose = false);
678 
685  const_iterator begin () const;
686 
690  const_iterator end () const;
691 
696  const_iterator begin (const size_type r) const;
697 
702  const_iterator end (const size_type r) const;
704 
712  void print (std::ostream &out) const;
713 
734  void print_formatted (std::ostream &out,
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;
740 
746  void block_write (std::ostream &out) const;
747 
758  void block_read (std::istream &in);
760 
770 
775  int, int,
776  << "The entry with index (" << arg1 << ',' << arg2
777  << ") does not exist.");
778 
780  int, int,
781  << "An entry with index (" << arg1 << ',' << arg2
782  << ") cannot be allocated.");
784 private:
789  const Entry *locate (const size_type row,
790  const size_type col) const;
791 
796  Entry *locate (const size_type row,
797  const size_type col);
798 
802  Entry *allocate (const size_type row,
803  const size_type col);
804 
810  template <typename somenumber>
811  void threaded_vmult (Vector<somenumber> &dst,
812  const Vector<somenumber> &src,
813  const size_type begin_row,
814  const size_type end_row) const;
815 
821  template <typename somenumber>
822  void threaded_matrix_norm_square (const Vector<somenumber> &v,
823  const size_type begin_row,
824  const size_type end_row,
825  somenumber *partial_sum) const;
826 
832  template <typename somenumber>
833  void threaded_matrix_scalar_product (const Vector<somenumber> &u,
834  const Vector<somenumber> &v,
835  const size_type begin_row,
836  const size_type end_row,
837  somenumber *partial_sum) const;
838 
843 
847  std::vector<RowInfo> row_info;
848 
852  std::vector<Entry> data;
853 
857  unsigned int increment;
858 
863 };
864 
868 /*---------------------- Inline functions -----------------------------------*/
869 
870 template <typename number>
871 inline
873  const number &value)
874  :
875  column(column),
876  value(value)
877 {}
878 
879 
880 
881 template <typename number>
882 inline
884  :
885  column(invalid),
886  value(0)
887 {}
888 
889 
890 template <typename number>
891 inline
893  :
894  start(start),
895  length(0),
896  diagonal(invalid_diagonal)
897 {}
898 
899 
900 //---------------------------------------------------------------------------
901 template <typename number>
902 inline
905  const size_type r,
906  const unsigned short i)
907  :
908  matrix(matrix),
909  a_row(r),
910  a_index(i)
911 {}
912 
913 
914 template <typename number>
915 inline
918 {
919  return a_row;
920 }
921 
922 
923 template <typename number>
924 inline
927 {
928  return matrix->data[matrix->row_info[a_row].start+a_index].column;
929 }
930 
931 
932 template <typename number>
933 inline
934 unsigned short
936 {
937  return a_index;
938 }
939 
940 
941 
942 template <typename number>
943 inline
944 number
946 {
947  return matrix->data[matrix->row_info[a_row].start+a_index].value;
948 }
949 
950 
951 template <typename number>
952 inline
955  const size_type r,
956  const unsigned short i)
957  :
958  accessor(matrix, r, i)
959 {
960  // Finish if this is the end()
961  if (r==accessor.matrix->m() && i==0) return;
962 
963  // Make sure we never construct an
964  // iterator pointing to a
965  // non-existing entry
966 
967  // If the index points beyond the
968  // end of the row, try the next
969  // row.
970  if (accessor.a_index >= accessor.matrix->row_info[accessor.a_row].length)
971  {
972  do
973  {
974  ++accessor.a_row;
975  }
976  // Beware! If the next row is
977  // empty, iterate until a
978  // non-empty row is found or we
979  // hit the end of the matrix.
980  while (accessor.a_row < accessor.matrix->m()
981  && accessor.matrix->row_info[accessor.a_row].length == 0);
982  }
983 }
984 
985 
986 template <typename number>
987 inline
990 {
991  Assert (accessor.a_row < accessor.matrix->m(), ExcIteratorPastEnd());
992 
993  // Increment column index
994  ++(accessor.a_index);
995  // If index exceeds number of
996  // entries in this row, proceed
997  // with next row.
998  if (accessor.a_index >= accessor.matrix->row_info[accessor.a_row].length)
999  {
1000  accessor.a_index = 0;
1001  // Do this loop to avoid
1002  // elements in empty rows
1003  do
1004  {
1005  ++accessor.a_row;
1006  }
1007  while (accessor.a_row < accessor.matrix->m()
1008  && accessor.matrix->row_info[accessor.a_row].length == 0);
1009  }
1010  return *this;
1011 }
1012 
1013 
1014 template <typename number>
1015 inline
1018 {
1019  return accessor;
1020 }
1021 
1022 
1023 template <typename number>
1024 inline
1027 {
1028  return &accessor;
1029 }
1030 
1031 
1032 template <typename number>
1033 inline
1034 bool
1036  const const_iterator &other) const
1037 {
1038  return (accessor.row() == other.accessor.row() &&
1039  accessor.index() == other.accessor.index());
1040 }
1041 
1042 
1043 template <typename number>
1044 inline
1045 bool
1047 operator != (const const_iterator &other) const
1048 {
1049  return ! (*this == other);
1050 }
1051 
1052 
1053 template <typename number>
1054 inline
1055 bool
1057 operator < (const const_iterator &other) const
1058 {
1059  return (accessor.row() < other.accessor.row() ||
1060  (accessor.row() == other.accessor.row() &&
1061  accessor.index() < other.accessor.index()));
1062 }
1063 
1064 
1065 //---------------------------------------------------------------------------
1066 template <typename number>
1067 inline
1069 {
1070  return row_info.size();
1071 }
1072 
1073 
1074 template <typename number>
1075 inline
1077 {
1078  return n_columns;
1079 }
1080 
1081 
1082 template <typename number>
1083 inline
1086  const size_type col)
1087 {
1088  Assert (row<m(), ExcIndexRange(row,0,m()));
1089  Assert (col<n(), ExcIndexRange(col,0,n()));
1090 
1091  const RowInfo &r = row_info[row];
1092  const size_type end = r.start + r.length;
1093  for (size_type i=r.start; i<end; ++i)
1094  {
1095  Entry *const entry = &data[i];
1096  if (entry->column == col)
1097  return entry;
1098  if (entry->column == Entry::invalid)
1099  return nullptr;
1100  }
1101  return nullptr;
1102 }
1103 
1104 
1105 
1106 template <typename number>
1107 inline
1108 const typename SparseMatrixEZ<number>::Entry *
1110  const size_type col) const
1111 {
1112  SparseMatrixEZ<number> *t = const_cast<SparseMatrixEZ<number>*> (this);
1113  return t->locate(row,col);
1114 }
1115 
1116 
1117 template <typename number>
1118 inline
1121  const size_type col)
1122 {
1123  Assert (row<m(), ExcIndexRange(row,0,m()));
1124  Assert (col<n(), ExcIndexRange(col,0,n()));
1125 
1126  RowInfo &r = row_info[row];
1127  const size_type end = r.start + r.length;
1128 
1129  size_type i = r.start;
1130  // If diagonal exists and this
1131  // column is higher, start only
1132  // after diagonal.
1133  if (r.diagonal != RowInfo::invalid_diagonal && col >= row)
1134  i += r.diagonal;
1135  // Find position of entry
1136  while (i<end && data[i].column < col) ++i;
1137 
1138  // entry found
1139  if (i != end && data[i].column == col)
1140  return &data[i];
1141 
1142  // Now, we must insert the new
1143  // entry and move all successive
1144  // entries back.
1145 
1146  // If no more space is available
1147  // for this row, insert new
1148  // elements into the vector.
1149 //TODO:[GK] We should not extend this row if i<end
1150  if (row != row_info.size()-1)
1151  {
1152  if (end >= row_info[row+1].start)
1153  {
1154  // Failure if increment 0
1156 
1157  // Insert new entries
1158  data.insert(data.begin()+end, increment, Entry());
1159  // Update starts of
1160  // following rows
1161  for (size_type rn=row+1; rn<row_info.size(); ++rn)
1162  row_info[rn].start += increment;
1163  }
1164  }
1165  else
1166  {
1167  if (end >= data.size())
1168  {
1169  // Here, appending a block
1170  // does not increase
1171  // performance.
1172  data.push_back(Entry());
1173  }
1174  }
1175 
1176  Entry *entry = &data[i];
1177  // Save original entry
1178  Entry temp = *entry;
1179  // Insert new entry here to
1180  // make sure all entries
1181  // are ordered by column
1182  // index
1183  entry->column = col;
1184  entry->value = 0;
1185  // Update row_info
1186  ++r.length;
1187  if (col == row)
1188  r.diagonal = i - r.start;
1189  else if (col<row && r.diagonal!= RowInfo::invalid_diagonal)
1190  ++r.diagonal;
1191 
1192  if (i == end)
1193  return entry;
1194 
1195  // Move all entries in this
1196  // row up by one
1197  for (size_type j = i+1; j < end; ++j)
1198  {
1199  // There should be no invalid
1200  // entry below end
1201  Assert (data[j].column != Entry::invalid, ExcInternalError());
1202 
1203 //TODO[GK]: This could be done more efficiently by moving starting at the top rather than swapping starting at the bottom
1204  std::swap (data[j], temp);
1205  }
1206  Assert (data[end].column == Entry::invalid, ExcInternalError());
1207 
1208  data[end] = temp;
1209 
1210  return entry;
1211 }
1212 
1213 
1214 
1215 template <typename number>
1216 inline
1218  const size_type j,
1219  const number value,
1220  const bool elide_zero_values)
1221 {
1222  AssertIsFinite(value);
1223 
1224  Assert (i<m(), ExcIndexRange(i,0,m()));
1225  Assert (j<n(), ExcIndexRange(j,0,n()));
1226 
1227  if (elide_zero_values && value == 0.)
1228  {
1229  Entry *entry = locate(i,j);
1230  if (entry != nullptr)
1231  entry->value = 0.;
1232  }
1233  else
1234  {
1235  Entry *entry = allocate(i,j);
1236  entry->value = value;
1237  }
1238 }
1239 
1240 
1241 
1242 template <typename number>
1243 inline
1245  const size_type j,
1246  const number value)
1247 {
1248 
1249  AssertIsFinite(value);
1250 
1251  Assert (i<m(), ExcIndexRange(i,0,m()));
1252  Assert (j<n(), ExcIndexRange(j,0,n()));
1253 
1254  // ignore zero additions
1255  if (value == 0)
1256  return;
1257 
1258  Entry *entry = allocate(i,j);
1259  entry->value += value;
1260 }
1261 
1262 
1263 template <typename number>
1264 template <typename number2>
1265 void SparseMatrixEZ<number>::add (const std::vector<size_type> &indices,
1266  const FullMatrix<number2> &full_matrix,
1267  const bool elide_zero_values)
1268 {
1269 //TODO: This function can surely be made more efficient
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));
1274 }
1275 
1276 
1277 
1278 template <typename number>
1279 template <typename number2>
1280 void SparseMatrixEZ<number>::add (const std::vector<size_type> &row_indices,
1281  const std::vector<size_type> &col_indices,
1282  const FullMatrix<number2> &full_matrix,
1283  const bool elide_zero_values)
1284 {
1285 //TODO: This function can surely be made more efficient
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));
1290 }
1291 
1292 
1293 
1294 
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)
1301 {
1302 //TODO: This function can surely be made more efficient
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]);
1306 }
1307 
1308 
1309 
1310 template <typename number>
1311 template <typename number2>
1313  const size_type n_cols,
1314  const size_type *col_indices,
1315  const number2 *values,
1316  const bool elide_zero_values,
1317  const bool /*col_indices_are_sorted*/)
1318 {
1319 //TODO: This function can surely be made more efficient
1320  for (size_type j=0; j<n_cols; ++j)
1321  if ((values[j] != 0) || (elide_zero_values == false))
1322  add (row, col_indices[j], values[j]);
1323 }
1324 
1325 
1326 
1327 
1328 template <typename number>
1329 inline
1331  const size_type j) const
1332 {
1333  const Entry *entry = locate(i,j);
1334  if (entry)
1335  return entry->value;
1336  return 0.;
1337 }
1338 
1339 
1340 
1341 template <typename number>
1342 inline
1344  const size_type j) const
1345 {
1346  const Entry *entry = locate(i,j);
1347  if (entry)
1348  return entry->value;
1349  Assert(false, ExcInvalidEntry(i,j));
1350  return 0.;
1351 }
1352 
1353 
1354 template <typename number>
1355 inline
1358 {
1359  const_iterator result(this, 0, 0);
1360  return result;
1361 }
1362 
1363 template <typename number>
1364 inline
1367 {
1368  return const_iterator(this, m(), 0);
1369 }
1370 
1371 template <typename number>
1372 inline
1375 {
1376  Assert (r<m(), ExcIndexRange(r,0,m()));
1377  const_iterator result (this, r, 0);
1378  return result;
1379 }
1380 
1381 template <typename number>
1382 inline
1385 {
1386  Assert (r<m(), ExcIndexRange(r,0,m()));
1387  const_iterator result(this, r+1, 0);
1388  return result;
1389 }
1390 
1391 template <typename number>
1392 template <typename MatrixType>
1393 inline
1395 SparseMatrixEZ<number>::copy_from (const MatrixType &M, const bool elide_zero_values)
1396 {
1397  reinit(M.m(),
1398  M.n(),
1400  this->increment);
1401 
1402  // loop over the elements of the argument matrix row by row, as suggested
1403  // in the documentation of the sparse matrix iterator class, and
1404  // copy them into the current object
1405  for (size_type row = 0; row < M.m(); ++row)
1406  {
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);
1411  }
1412 
1413  return *this;
1414 }
1415 
1416 template <typename number>
1417 template <typename MatrixType>
1418 inline
1419 void
1420 SparseMatrixEZ<number>::add (const number factor,
1421  const MatrixType &M)
1422 {
1423  Assert (M.m() == m(), ExcDimensionMismatch(M.m(), m()));
1424  Assert (M.n() == n(), ExcDimensionMismatch(M.n(), n()));
1425 
1426  if (factor == 0.)
1427  return;
1428 
1429  // loop over the elements of the argument matrix row by row, as suggested
1430  // in the documentation of the sparse matrix iterator class, and
1431  // add them into the current object
1432  for (size_type row = 0; row < M.m(); ++row)
1433  {
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());
1439  }
1440 }
1441 
1442 
1443 
1444 template <typename number>
1445 template <typename MatrixTypeA, typename MatrixTypeB>
1446 inline void
1448  const MatrixTypeB &B,
1449  const bool transpose)
1450 {
1451 // Compute the result
1452 // r_ij = \sum_kl b_ik b_jl a_kl
1453 
1454 // Assert (n() == B.m(), ExcDimensionMismatch(n(), B.m()));
1455 // Assert (m() == B.m(), ExcDimensionMismatch(m(), B.m()));
1456 // Assert (A.n() == B.n(), ExcDimensionMismatch(A.n(), B.n()));
1457 // Assert (A.m() == B.n(), ExcDimensionMismatch(A.m(), B.n()));
1458 
1459  // Somehow, we have to avoid making
1460  // this an operation of complexity
1461  // n^2. For the transpose case, we
1462  // can go through the non-zero
1463  // elements of A^-1 and use the
1464  // corresponding rows of B only.
1465  // For the non-transpose case, we
1466  // must find a trick.
1467  typename MatrixTypeB::const_iterator b1 = B.begin();
1468  const typename MatrixTypeB::const_iterator b_final = B.end();
1469  if (transpose)
1470  while (b1 != b_final)
1471  {
1472  const size_type i = b1->column();
1473  const size_type k = b1->row();
1474  typename MatrixTypeB::const_iterator b2 = B.begin();
1475  while (b2 != b_final)
1476  {
1477  const size_type j = b2->column();
1478  const size_type l = b2->row();
1479 
1480  const typename MatrixTypeA::value_type a = A.el(k,l);
1481 
1482  if (a != 0.)
1483  add (i, j, a * b1->value() * b2->value());
1484  ++b2;
1485  }
1486  ++b1;
1487  }
1488  else
1489  {
1490  // Determine minimal and
1491  // maximal row for a column in
1492  // advance.
1493 
1494  std::vector<size_type> minrow(B.n(), B.m());
1495  std::vector<size_type> maxrow(B.n(), 0);
1496  while (b1 != b_final)
1497  {
1498  const size_type r = b1->row();
1499  if (r < minrow[b1->column()])
1500  minrow[b1->column()] = r;
1501  if (r > maxrow[b1->column()])
1502  maxrow[b1->column()] = r;
1503  ++b1;
1504  }
1505 
1506  typename MatrixTypeA::const_iterator ai = A.begin();
1507  const typename MatrixTypeA::const_iterator ae = A.end();
1508 
1509  while (ai != ae)
1510  {
1511  const typename MatrixTypeA::value_type a = ai->value();
1512  // Don't do anything if
1513  // this entry is zero.
1514  if (a == 0.) continue;
1515 
1516  // Now, loop over all rows
1517  // having possibly a
1518  // nonzero entry in column
1519  // ai->row()
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()]);
1525 
1526  while (b1 != be1)
1527  {
1528  const double b1v = b1->value();
1529  // We need the product
1530  // of both. If it is
1531  // zero, we can save
1532  // the work
1533  if (b1->column() == ai->row() && (b1v != 0.))
1534  {
1535  const size_type i = b1->row();
1536 
1537  typename MatrixTypeB::const_iterator
1538  b2 = B.begin(minrow[ai->column()]);
1539  while (b2 != be2)
1540  {
1541  if (b2->column() == ai->column())
1542  {
1543  const size_type j = b2->row();
1544  add (i, j, a * b1v * b2->value());
1545  }
1546  ++b2;
1547  }
1548  }
1549  ++b1;
1550  }
1551  ++ai;
1552  }
1553  }
1554 }
1555 
1556 
1557 template <typename number>
1558 template <class StreamType>
1559 inline
1560 void
1561 SparseMatrixEZ<number>::print_statistics(StreamType &out, bool full)
1562 {
1563  size_type used;
1564  size_type allocated;
1565  size_type reserved;
1566  std::vector<size_type> used_by_line;
1567 
1568  compute_statistics (used, allocated, reserved, used_by_line, full);
1569 
1570  out << "SparseMatrixEZ:used entries:" << used << std::endl
1571  << "SparseMatrixEZ:allocated entries:" << allocated << std::endl
1572  << "SparseMatrixEZ:reserved entries:" << reserved << std::endl;
1573 
1574  if (full)
1575  {
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]
1580  << std::endl;
1581 
1582  }
1583 }
1584 
1585 
1586 DEAL_II_NAMESPACE_CLOSE
1587 
1588 #endif
1589 /*---------------------------- sparse_matrix.h ---------------------------*/
size_type get_row_length(const size_type row) const
const types::global_dof_index invalid_size_type
Definition: types.h:182
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:358
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:1142
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:323
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)
Definition: vector.h:1312
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
~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
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:1299
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