Reference documentation for deal.II version 8.5.1
constraint_matrix.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 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 
17 #ifndef dealii__constraint_matrix_h
18 #define dealii__constraint_matrix_h
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/exceptions.h>
22 #include <deal.II/base/index_set.h>
23 #include <deal.II/base/subscriptor.h>
24 #include <deal.II/base/template_constraints.h>
25 
26 #include <deal.II/lac/vector.h>
27 
28 #include <vector>
29 #include <set>
30 #include <utility>
31 
32 
33 DEAL_II_NAMESPACE_OPEN
34 
35 template<int dim, class T> class Table;
36 template <typename> class FullMatrix;
37 class SparsityPattern;
41 template <typename number> class SparseMatrix;
42 template <typename number> class BlockSparseMatrix;
43 
44 namespace internals
45 {
46  class GlobalRowsFromLocal;
47 }
48 
49 
50 //TODO[WB]: We should have a function of the kind
51 // ConstraintMatrix::add_constraint (const size_type constrained_dof,
52 // const std::vector<std::pair<size_type, double> > &entries,
53 // const double inhomogeneity = 0);
54 // rather than building up constraints piecemeal through add_line/add_entry
55 // etc. This would also eliminate the possibility of accidentally changing
56 // existing constraints into something pointless, see the discussion on the
57 // mailing list on "Tiny bug in interpolate_boundary_values" in Sept. 2010.
58 
135 {
136 public:
141 
148  {
154 
161 
168  };
169 
187  explicit ConstraintMatrix (const IndexSet &local_constraints = IndexSet());
188 
192  explicit ConstraintMatrix (const ConstraintMatrix &constraint_matrix);
193 
200  void reinit (const IndexSet &local_constraints = IndexSet());
201 
208  bool can_store_line (const size_type line_index) const;
209 
216  const IndexSet &get_local_lines() const;
217 
237  void add_selected_constraints (const ConstraintMatrix &constraints_in,
238  const IndexSet &filter);
239 
249  void add_line (const size_type line);
250 
263  void add_lines (const std::vector<bool> &lines);
264 
277  void add_lines (const std::set<size_type> &lines);
278 
291  void add_lines (const IndexSet &lines);
292 
304  void add_entry (const size_type line,
305  const size_type column,
306  const double value);
307 
313  void add_entries (const size_type line,
314  const std::vector<std::pair<size_type,double> > &col_val_pairs);
315 
322  void set_inhomogeneity (const size_type line,
323  const double value);
324 
345  void close ();
346 
371  void merge (const ConstraintMatrix &other_constraints,
372  const MergeConflictBehavior merge_conflict_behavior = no_conflicts_allowed,
373  const bool allow_different_local_lines = false);
374 
387  void shift (const size_type offset);
388 
396  void clear ();
397 
411  size_type n_constraints () const;
412 
422  bool is_constrained (const size_type index) const;
423 
435  bool is_identity_constrained (const size_type index) const;
436 
443  bool are_identity_constrained (const size_type index1,
444  const size_type index2) const;
445 
457 
462  bool is_inhomogeneously_constrained (const size_type index) const;
463 
469  bool has_inhomogeneities () const;
470 
475  const std::vector<std::pair<size_type,double> > *
476  get_constraint_entries (const size_type line) const;
477 
482  double get_inhomogeneity (const size_type line) const;
483 
504  void print (std::ostream &out) const;
505 
518  void write_dot (std::ostream &) const;
519 
524  std::size_t memory_consumption () const;
525 
532  void resolve_indices(std::vector<types::global_dof_index> &indices) const;
533 
561  void condense (SparsityPattern &sparsity) const;
562 
566  void condense (BlockSparsityPattern &sparsity) const;
567 
572  void condense (DynamicSparsityPattern &sparsity) const;
573 
578  void condense (BlockDynamicSparsityPattern &sparsity) const;
579 
587  template<typename number>
588  void condense (SparseMatrix<number> &matrix) const;
589 
593  template <typename number>
594  void condense (BlockSparseMatrix<number> &matrix) const;
595 
607  template <class VectorType>
608  void condense (VectorType &vec) const;
609 
616  template <class VectorType>
617  void condense (const VectorType &vec_ghosted,
618  VectorType &output) const;
619 
632  template<typename number, class VectorType>
633  void condense (SparseMatrix<number> &matrix,
634  VectorType &vector) const;
635 
640  template <typename number, class BlockVectorType>
641  void condense (BlockSparseMatrix<number> &matrix,
642  BlockVectorType &vector) const;
643 
650  template <class VectorType>
651  void set_zero (VectorType &vec) const;
652 
704  template <class InVector, class OutVector>
705  void
706  distribute_local_to_global (const InVector &local_vector,
707  const std::vector<size_type> &local_dof_indices,
708  OutVector &global_vector) const;
709 
753  template <typename VectorType, typename LocalType>
754  void
755  distribute_local_to_global (const Vector<LocalType> &local_vector,
756  const std::vector<size_type> &local_dof_indices,
757  VectorType &global_vector,
758  const FullMatrix<LocalType> &local_matrix) const;
759 
777  template <typename VectorType, typename LocalType>
778  void
779  distribute_local_to_global (const Vector<LocalType> &local_vector,
780  const std::vector<size_type> &local_dof_indices_row,
781  const std::vector<size_type> &local_dof_indices_col,
782  VectorType &global_vector,
783  const FullMatrix<LocalType> &local_matrix,
784  bool diagonal = false) const;
785 
789  template <class VectorType>
790  void
792  const double value,
793  VectorType &global_vector) const;
794 
823  template <typename ForwardIteratorVec, typename ForwardIteratorInd,
824  class VectorType>
825  void
826  distribute_local_to_global (ForwardIteratorVec local_vector_begin,
827  ForwardIteratorVec local_vector_end,
828  ForwardIteratorInd local_indices_begin,
829  VectorType &global_vector) const;
830 
878  template <typename MatrixType>
879  void
881  const std::vector<size_type> &local_dof_indices,
882  MatrixType &global_matrix) const;
883 
911  template <typename MatrixType>
912  void
914  const std::vector<size_type> &row_indices,
915  const std::vector<size_type> &col_indices,
916  MatrixType &global_matrix) const;
917 
934  template <typename MatrixType, typename VectorType>
935  void
937  const Vector<typename VectorType::value_type> &local_vector,
938  const std::vector<size_type> &local_dof_indices,
939  MatrixType &global_matrix,
940  VectorType &global_vector,
941  bool use_inhomogeneities_for_rhs = false) const;
942 
996  template <typename SparsityPatternType>
997  void
998  add_entries_local_to_global (const std::vector<size_type> &local_dof_indices,
999  SparsityPatternType &sparsity_pattern,
1000  const bool keep_constrained_entries = true,
1001  const Table<2,bool> &dof_mask = default_empty_table) const;
1002 
1006  template <typename SparsityPatternType>
1007  void
1008  add_entries_local_to_global (const std::vector<size_type> &row_indices,
1009  const std::vector<size_type> &col_indices,
1010  SparsityPatternType &sparsity_pattern,
1011  const bool keep_constrained_entries = true,
1012  const Table<2,bool> &dof_mask = default_empty_table) const;
1013 
1033  template <typename ForwardIteratorVec, typename ForwardIteratorInd,
1034  class VectorType>
1035  void
1036  get_dof_values (const VectorType &global_vector,
1037  ForwardIteratorInd local_indices_begin,
1038  ForwardIteratorVec local_vector_begin,
1039  ForwardIteratorVec local_vector_end) const;
1040 
1062  template <class VectorType>
1063  void distribute (VectorType &vec) const;
1064 
1087  size_type,
1088  << "The specified line " << arg1
1089  << " does not exist.");
1096  size_type, size_type, double, double,
1097  << "The entry for the indices " << arg1 << " and "
1098  << arg2 << " already exists, but the values "
1099  << arg3 << " (old) and " << arg4 << " (new) differ "
1100  << "by " << (arg4-arg3) << ".");
1107  int, int,
1108  << "You tried to constrain DoF " << arg1
1109  << " to DoF " << arg2
1110  << ", but that one is also constrained. This is not allowed!");
1117  size_type,
1118  << "Degree of freedom " << arg1
1119  << " is constrained from both object in a merge operation.");
1126  size_type,
1127  << "In the given argument a degree of freedom is constrained "
1128  << "to another DoF with number " << arg1
1129  << ", which however is constrained by this object. This is not"
1130  << " allowed.");
1137  size_type,
1138  << "The index set given to this constraint matrix indicates "
1139  << "constraints for degree of freedom " << arg1
1140  << " should not be stored by this object, but a constraint "
1141  << "is being added.");
1142 
1149  size_type,
1150  size_type,
1151  << "The index set given to this constraint matrix indicates "
1152  << "constraints using degree of freedom " << arg2
1153  << " should not be stored by this object, but a constraint "
1154  << "for degree of freedom " << arg1 <<" uses it.");
1155 
1162  int, int,
1163  << "While distributing the constraint for DoF "
1164  << arg1 << ", it turns out that one of the processors "
1165  << "who own the " << arg2
1166  << " degrees of freedom that x_" << arg1
1167  << " is constrained against does not know about "
1168  << "the constraint on x_" << arg1
1169  << ". Did you not initialize the ConstraintMatrix "
1170  << "with the appropriate locally_relevant set so "
1171  << "that every processor who owns a DoF that constrains "
1172  << "another DoF also knows about this constraint?");
1173 
1174 private:
1175 
1180  {
1185  typedef std::vector<std::pair<size_type,double> > Entries;
1186 
1192 
1201 
1206 
1214  bool operator < (const ConstraintLine &) const;
1215 
1221  bool operator == (const ConstraintLine &) const;
1222 
1227  std::size_t memory_consumption () const;
1228  };
1229 
1241  std::vector<ConstraintLine> lines;
1242 
1275  std::vector<size_type> lines_cache;
1276 
1283 
1287  bool sorted;
1288 
1293  size_type calculate_line_index (const size_type line) const;
1294 
1299  static bool check_zero_weight (const std::pair<size_type, double> &p);
1300 
1306 
1311  template <typename MatrixType, typename VectorType>
1312  void
1314  const Vector<typename VectorType::value_type> &local_vector,
1315  const std::vector<size_type> &local_dof_indices,
1316  MatrixType &global_matrix,
1317  VectorType &global_vector,
1318  bool use_inhomogeneities_for_rhs,
1320 
1325  template <typename MatrixType, typename VectorType>
1326  void
1328  const Vector<typename VectorType::value_type> &local_vector,
1329  const std::vector<size_type> &local_dof_indices,
1330  MatrixType &global_matrix,
1331  VectorType &global_vector,
1332  bool use_inhomogeneities_for_rhs,
1334 
1339  template <typename SparsityPatternType>
1340  void
1341  add_entries_local_to_global (const std::vector<size_type> &local_dof_indices,
1342  SparsityPatternType &sparsity_pattern,
1343  const bool keep_constrained_entries,
1344  const Table<2,bool> &dof_mask,
1346 
1351  template <typename SparsityPatternType>
1352  void
1353  add_entries_local_to_global (const std::vector<size_type> &local_dof_indices,
1354  SparsityPatternType &sparsity_pattern,
1355  const bool keep_constrained_entries,
1356  const Table<2,bool> &dof_mask,
1358 
1366  void
1367  make_sorted_row_list (const std::vector<size_type> &local_dof_indices,
1368  internals::GlobalRowsFromLocal &global_rows) const;
1369 
1377  void
1378  make_sorted_row_list (const std::vector<size_type> &local_dof_indices,
1379  std::vector<size_type> &active_dofs) const;
1380 
1384  template <typename LocalType>
1385  LocalType
1387  const internals::GlobalRowsFromLocal &global_rows,
1388  const Vector<LocalType> &local_vector,
1389  const std::vector<size_type> &local_dof_indices,
1390  const FullMatrix<LocalType> &local_matrix) const;
1391 
1397 };
1398 
1399 
1400 
1401 /* ---------------- template and inline functions ----------------- */
1402 
1403 inline
1405  :
1406  lines (),
1407  local_lines (local_constraints),
1408  sorted (false)
1409 {
1410  // make sure the IndexSet is compressed. Otherwise this can lead to crashes
1411  // that are hard to find (only happen in release mode).
1412  // see tests/mpi/constraint_matrix_crash_01
1414 }
1415 
1416 
1417 
1418 inline
1420  :
1421  Subscriptor (),
1422  lines (constraint_matrix.lines),
1423  lines_cache (constraint_matrix.lines_cache),
1424  local_lines (constraint_matrix.local_lines),
1425  sorted (constraint_matrix.sorted)
1426 {}
1427 
1428 
1429 inline
1430 void
1432 {
1433  Assert (sorted==false, ExcMatrixIsClosed());
1434 
1435  // the following can happen when we compute with distributed meshes and dof
1436  // handlers and we constrain a degree of freedom whose number we don't have
1437  // locally. if we don't abort here the program will try to allocate several
1438  // terabytes of memory to resize the various arrays below :-)
1440  ExcInternalError());
1441  const size_type line_index = calculate_line_index (line);
1442 
1443  // check whether line already exists; it may, in which case we can just quit
1444  if (is_constrained(line))
1445  return;
1446 
1447  // if necessary enlarge vector of existing entries for cache
1448  if (line_index >= lines_cache.size())
1449  lines_cache.resize (std::max(2*static_cast<size_type>(lines_cache.size()),
1450  line_index+1),
1452 
1453  // push a new line to the end of the list
1454  lines.push_back (ConstraintLine());
1455  lines.back().line = line;
1456  lines.back().inhomogeneity = 0.;
1457  lines_cache[line_index] = lines.size()-1;
1458 }
1459 
1460 
1461 
1462 inline
1463 void
1465  const size_type column,
1466  const double value)
1467 {
1468  Assert (sorted==false, ExcMatrixIsClosed());
1469  Assert (line != column,
1470  ExcMessage ("Can't constrain a degree of freedom to itself"));
1471 
1472  // if in debug mode, check whether an entry for this column already exists
1473  // and if it's the same as the one entered at present
1474  //
1475  // in any case: exit the function if an entry for this column already
1476  // exists, since we don't want to enter it twice
1478  ExcInternalError());
1479  Assert (!local_lines.size() || local_lines.is_element(column),
1480  ExcColumnNotStoredHere(line, column));
1482  Assert (line_ptr->line == line, ExcInternalError());
1483  for (ConstraintLine::Entries::const_iterator
1484  p=line_ptr->entries.begin();
1485  p != line_ptr->entries.end(); ++p)
1486  if (p->first == column)
1487  {
1488  Assert (std::fabs(p->second - value) < 1.e-14,
1489  ExcEntryAlreadyExists(line, column, p->second, value));
1490  return;
1491  }
1492 
1493  line_ptr->entries.push_back (std::make_pair(column,value));
1494 }
1495 
1496 
1497 
1498 inline
1499 void
1501  const double value)
1502 {
1503  const size_type line_index = calculate_line_index(line);
1504  Assert( line_index < lines_cache.size() &&
1505  lines_cache[line_index] != numbers::invalid_size_type,
1506  ExcMessage("call add_line() before calling set_inhomogeneity()"));
1507  Assert(lines_cache[line_index] < lines.size(), ExcInternalError());
1508  ConstraintLine *line_ptr = &lines[lines_cache[line_index]];
1509  line_ptr->inhomogeneity = value;
1510 }
1511 
1512 
1513 
1514 inline
1517 {
1518  return lines.size();
1519 }
1520 
1521 
1522 
1523 inline
1524 bool
1526 {
1527  const size_type line_index = calculate_line_index(index);
1528  return ((line_index < lines_cache.size())
1529  &&
1530  (lines_cache[line_index] != numbers::invalid_size_type));
1531 }
1532 
1533 
1534 
1535 inline
1536 bool
1538 {
1539  // check whether the entry is constrained. could use is_constrained, but
1540  // that means computing the line index twice
1541  const size_type line_index = calculate_line_index(index);
1542  if (line_index >= lines_cache.size() ||
1543  lines_cache[line_index] == numbers::invalid_size_type)
1544  return false;
1545  else
1546  {
1547  Assert(lines_cache[line_index] < lines.size(), ExcInternalError());
1548  return !(lines[lines_cache[line_index]].inhomogeneity == 0);
1549  }
1550 }
1551 
1552 
1553 
1554 inline
1555 const std::vector<std::pair<types::global_dof_index,double> > *
1557 {
1558  // check whether the entry is constrained. could use is_constrained, but
1559  // that means computing the line index twice
1560  const size_type line_index = calculate_line_index(line);
1561  if (line_index >= lines_cache.size() ||
1562  lines_cache[line_index] == numbers::invalid_size_type)
1563  return 0;
1564  else
1565  return &lines[lines_cache[line_index]].entries;
1566 }
1567 
1568 
1569 
1570 inline
1571 double
1573 {
1574  // check whether the entry is constrained. could use is_constrained, but
1575  // that means computing the line index twice
1576  const size_type line_index = calculate_line_index(line);
1577  if (line_index >= lines_cache.size() ||
1578  lines_cache[line_index] == numbers::invalid_size_type)
1579  return 0;
1580  else
1581  return lines[lines_cache[line_index]].inhomogeneity;
1582 }
1583 
1584 
1585 
1588 {
1589  //IndexSet is unused (serial case)
1590  if (!local_lines.size())
1591  return line;
1592 
1594  ExcRowNotStoredHere(line));
1595 
1596  return local_lines.index_within_set(line);
1597 }
1598 
1599 
1600 
1601 inline bool
1603 {
1604  return !local_lines.size() || local_lines.is_element(line_index);
1605 }
1606 
1607 
1608 
1609 inline
1610 const IndexSet &
1612 {
1613  return local_lines;
1614 }
1615 
1616 
1617 
1618 template <class VectorType>
1619 inline
1621  const size_type index,
1622  const double value,
1623  VectorType &global_vector) const
1624 {
1625  Assert (lines.empty() || sorted == true, ExcMatrixNotClosed());
1626 
1627  if (is_constrained(index) == false)
1628  global_vector(index) += value;
1629  else
1630  {
1631  const ConstraintLine &position =
1633  for (size_type j=0; j<position.entries.size(); ++j)
1634  global_vector(position.entries[j].first)
1635  += value * position.entries[j].second;
1636  }
1637 }
1638 
1639 
1640 template <typename ForwardIteratorVec, typename ForwardIteratorInd,
1641  class VectorType>
1642 inline
1644  ForwardIteratorVec local_vector_begin,
1645  ForwardIteratorVec local_vector_end,
1646  ForwardIteratorInd local_indices_begin,
1647  VectorType &global_vector) const
1648 {
1649  Assert (lines.empty() || sorted == true, ExcMatrixNotClosed());
1650  for ( ; local_vector_begin != local_vector_end;
1651  ++local_vector_begin, ++local_indices_begin)
1652  {
1653  if (is_constrained(*local_indices_begin) == false)
1654  global_vector(*local_indices_begin) += *local_vector_begin;
1655  else
1656  {
1657  const ConstraintLine &position =
1658  lines[lines_cache[calculate_line_index(*local_indices_begin)]];
1659  for (size_type j=0; j<position.entries.size(); ++j)
1660  global_vector(position.entries[j].first)
1661  += *local_vector_begin * position.entries[j].second;
1662  }
1663  }
1664 }
1665 
1666 
1667 template <class InVector, class OutVector>
1668 inline
1669 void
1671  const InVector &local_vector,
1672  const std::vector<size_type> &local_dof_indices,
1673  OutVector &global_vector) const
1674 {
1675  Assert (local_vector.size() == local_dof_indices.size(),
1676  ExcDimensionMismatch(local_vector.size(), local_dof_indices.size()));
1677  distribute_local_to_global (local_vector.begin(), local_vector.end(),
1678  local_dof_indices.begin(), global_vector);
1679 }
1680 
1681 
1682 
1683 template <typename ForwardIteratorVec, typename ForwardIteratorInd,
1684  class VectorType>
1685 inline
1686 void ConstraintMatrix::get_dof_values (const VectorType &global_vector,
1687  ForwardIteratorInd local_indices_begin,
1688  ForwardIteratorVec local_vector_begin,
1689  ForwardIteratorVec local_vector_end) const
1690 {
1691  Assert (lines.empty() || sorted == true, ExcMatrixNotClosed());
1692  for ( ; local_vector_begin != local_vector_end;
1693  ++local_vector_begin, ++local_indices_begin)
1694  {
1695  if (is_constrained(*local_indices_begin) == false)
1696  *local_vector_begin = global_vector(*local_indices_begin);
1697  else
1698  {
1699  const ConstraintLine &position =
1700  lines[lines_cache[calculate_line_index(*local_indices_begin)]];
1701  typename VectorType::value_type value = position.inhomogeneity;
1702  for (size_type j=0; j<position.entries.size(); ++j)
1703  value += (global_vector(position.entries[j].first) *
1704  position.entries[j].second);
1705  *local_vector_begin = value;
1706  }
1707  }
1708 }
1709 
1710 
1711 template <typename MatrixType> class BlockMatrixBase;
1712 template <typename SparsityPatternType> class BlockSparsityPatternBase;
1713 template <typename number> class BlockSparseMatrixEZ;
1714 
1733 template <typename MatrixType>
1735 {
1736 private:
1737  struct yes_type
1738  {
1739  char c[1];
1740  };
1741  struct no_type
1742  {
1743  char c[2];
1744  };
1745 
1751  template <typename T>
1752  static yes_type check_for_block_matrix (const BlockMatrixBase<T> *);
1753 
1758  template <typename T>
1759  static yes_type check_for_block_matrix (const BlockSparsityPatternBase<T> *);
1760 
1765  template <typename T>
1766  static yes_type check_for_block_matrix (const BlockSparseMatrixEZ<T> *);
1767 
1772  static no_type check_for_block_matrix (...);
1773 
1774 public:
1780  static const bool value = (sizeof(check_for_block_matrix
1781  ((MatrixType *)0))
1782  ==
1783  sizeof(yes_type));
1784 };
1785 
1786 
1787 // instantiation of the static member
1788 template <typename MatrixType>
1790 
1791 
1792 template <typename MatrixType>
1793 inline
1794 void
1797  const std::vector<size_type> &local_dof_indices,
1798  MatrixType &global_matrix) const
1799 {
1800  // create a dummy and hand on to the function actually implementing this
1801  // feature in the cm.templates.h file.
1803  distribute_local_to_global (local_matrix, dummy, local_dof_indices,
1804  global_matrix, dummy, false,
1806 }
1807 
1808 
1809 
1810 
1811 template <typename MatrixType, typename VectorType>
1812 inline
1813 void
1816  const Vector<typename VectorType::value_type> &local_vector,
1817  const std::vector<size_type> &local_dof_indices,
1818  MatrixType &global_matrix,
1819  VectorType &global_vector,
1820  bool use_inhomogeneities_for_rhs) const
1821 {
1822  // enter the internal function with the respective block information set,
1823  // the actual implementation follows in the cm.templates.h file.
1824  distribute_local_to_global (local_matrix, local_vector, local_dof_indices,
1825  global_matrix, global_vector, use_inhomogeneities_for_rhs,
1827 }
1828 
1829 
1830 
1831 
1832 template <typename SparsityPatternType>
1833 inline
1834 void
1836 add_entries_local_to_global (const std::vector<size_type> &local_dof_indices,
1837  SparsityPatternType &sparsity_pattern,
1838  const bool keep_constrained_entries,
1839  const Table<2,bool> &dof_mask) const
1840 {
1841  // enter the internal function with the respective block information set,
1842  // the actual implementation follows in the cm.templates.h file.
1843  add_entries_local_to_global (local_dof_indices, sparsity_pattern,
1844  keep_constrained_entries, dof_mask,
1846 }
1847 
1848 
1849 DEAL_II_NAMESPACE_CLOSE
1850 
1851 #endif
const types::global_dof_index invalid_size_type
Definition: types.h:179
bool operator<(const ConstraintLine &) const
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:576
std::vector< size_type > lines_cache
void get_dof_values(const VectorType &global_vector, ForwardIteratorInd local_indices_begin, ForwardIteratorVec local_vector_begin, ForwardIteratorVec local_vector_end) const
static ::ExceptionBase & ExcDoFIsConstrainedFromBothObjects(size_type arg1)
static ::ExceptionBase & ExcEntryAlreadyExists(size_type arg1, size_type arg2, double arg3, double arg4)
void merge(const ConstraintMatrix &other_constraints, const MergeConflictBehavior merge_conflict_behavior=no_conflicts_allowed, const bool allow_different_local_lines=false)
static bool check_zero_weight(const std::pair< size_type, double > &p)
ConstraintMatrix & operator=(const ConstraintMatrix &other)
std::size_t memory_consumption() const
void set_zero(VectorType &vec) const
void add_entries(const size_type line, const std::vector< std::pair< size_type, double > > &col_val_pairs)
const IndexSet & get_local_lines() const
static ::ExceptionBase & ExcDoFIsConstrainedToConstrainedDoF(size_type arg1)
static const bool value
types::global_dof_index size_type
static yes_type check_for_block_matrix(const BlockMatrixBase< T > *)
std::vector< std::pair< size_type, double > > Entries
static ::ExceptionBase & ExcLineInexistant(size_type arg1)
size_type size() const
Definition: index_set.h:1419
void add_line(const size_type line)
bool operator==(const ConstraintLine &) const
static ::ExceptionBase & ExcMessage(std::string arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:564
void add_entry(const size_type line, const size_type column, const double value)
void make_sorted_row_list(const std::vector< size_type > &local_dof_indices, internals::GlobalRowsFromLocal &global_rows) const
unsigned int global_dof_index
Definition: types.h:88
bool is_inhomogeneously_constrained(const size_type index) const
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
void print(std::ostream &out) const
size_type index_within_set(const size_type global_index) const
Definition: index_set.h:1653
static ::ExceptionBase & ExcColumnNotStoredHere(size_type arg1, size_type arg2)
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
ConstraintMatrix(const IndexSet &local_constraints=IndexSet())
void condense(SparsityPattern &sparsity) const
size_type n_constraints() const
#define DeclException0(Exception0)
Definition: exceptions.h:541
void resolve_indices(std::vector< types::global_dof_index > &indices) const
bool has_inhomogeneities() const
static const Table< 2, bool > default_empty_table
std::size_t memory_consumption() const
void add_lines(const std::vector< bool > &lines)
size_type max_constraint_indirections() const
size_type calculate_line_index(const size_type line) const
static ::ExceptionBase & ExcIncorrectConstraint(int arg1, int arg2)
bool can_store_line(const size_type line_index) const
static ::ExceptionBase & ExcMatrixNotClosed()
LocalType resolve_vector_entry(const size_type i, const internals::GlobalRowsFromLocal &global_rows, const Vector< LocalType > &local_vector, const std::vector< size_type > &local_dof_indices, const FullMatrix< LocalType > &local_matrix) const
static ::ExceptionBase & ExcDoFConstrainedToConstrainedDoF(int arg1, int arg2)
void distribute(VectorType &vec) const
void compress() const
Definition: index_set.h:1428
void write_dot(std::ostream &) const
bool is_identity_constrained(const size_type index) const
void shift(const size_type offset)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
Definition: exceptions.h:600
void add_entries_local_to_global(const std::vector< size_type > &local_dof_indices, SparsityPatternType &sparsity_pattern, const bool keep_constrained_entries=true, const Table< 2, bool > &dof_mask=default_empty_table) const
bool are_identity_constrained(const size_type index1, const size_type index2) const
static ::ExceptionBase & ExcRowNotStoredHere(size_type arg1)
Definition: table.h:33
bool is_element(const size_type index) const
Definition: index_set.h:1489
void reinit(const IndexSet &local_constraints=IndexSet())
void add_selected_constraints(const ConstraintMatrix &constraints_in, const IndexSet &filter)
const std::vector< std::pair< size_type, double > > * get_constraint_entries(const size_type line) const
double get_inhomogeneity(const size_type line) const
void set_inhomogeneity(const size_type line, const double value)
std::vector< ConstraintLine > lines
bool is_constrained(const size_type index) const
static ::ExceptionBase & ExcMatrixIsClosed()
static ::ExceptionBase & ExcInternalError()