Reference documentation for deal.II version 9.0.0
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 #include <deal.II/lac/vector_element_access.h>
28 
29 #include <boost/range/iterator_range.hpp>
30 
31 #include <vector>
32 #include <set>
33 #include <utility>
34 
35 
36 DEAL_II_NAMESPACE_OPEN
37 
38 template <int dim, class T> class Table;
39 template <typename> class FullMatrix;
40 class SparsityPattern;
44 template <typename number> class SparseMatrix;
45 template <typename number> class BlockSparseMatrix;
46 
47 namespace internals
48 {
49  class GlobalRowsFromLocal;
50 }
51 
52 
53 //TODO[WB]: We should have a function of the kind
54 // ConstraintMatrix::add_constraint (const size_type constrained_dof,
55 // const std::vector<std::pair<size_type, double> > &entries,
56 // const double inhomogeneity = 0);
57 // rather than building up constraints piecemeal through add_line/add_entry
58 // etc. This would also eliminate the possibility of accidentally changing
59 // existing constraints into something pointless, see the discussion on the
60 // mailing list on "Tiny bug in interpolate_boundary_values" in Sept. 2010.
61 
138 {
139 public:
144 
151  {
157 
164 
171  };
172 
190  explicit ConstraintMatrix (const IndexSet &local_constraints = IndexSet());
191 
195  explicit ConstraintMatrix (const ConstraintMatrix &constraint_matrix);
196 
200  ConstraintMatrix (ConstraintMatrix &&constraint_matrix) = default;
201 
211  ConstraintMatrix &operator= (const ConstraintMatrix &) = delete;
212 
216  ConstraintMatrix &operator= (ConstraintMatrix &&constraint_matrix) = default;
217 
224  void copy_from (const ConstraintMatrix &other);
225 
232  void reinit (const IndexSet &local_constraints = IndexSet());
233 
240  bool can_store_line (const size_type line_index) const;
241 
248  const IndexSet &get_local_lines() const;
249 
269  void add_selected_constraints (const ConstraintMatrix &constraints_in,
270  const IndexSet &filter);
271 
281  void add_line (const size_type line);
282 
295  void add_lines (const std::vector<bool> &lines);
296 
309  void add_lines (const std::set<size_type> &lines);
310 
323  void add_lines (const IndexSet &lines);
324 
336  void add_entry (const size_type line,
337  const size_type column,
338  const double value);
339 
345  void add_entries (const size_type line,
346  const std::vector<std::pair<size_type,double> > &col_val_pairs);
347 
354  void set_inhomogeneity (const size_type line,
355  const double value);
356 
377  void close ();
378 
403  void merge (const ConstraintMatrix &other_constraints,
404  const MergeConflictBehavior merge_conflict_behavior = no_conflicts_allowed,
405  const bool allow_different_local_lines = false);
406 
419  void shift (const size_type offset);
420 
428  void clear ();
429 
443  size_type n_constraints () const;
444 
454  bool is_constrained (const size_type index) const;
455 
467  bool is_identity_constrained (const size_type index) const;
468 
475  bool are_identity_constrained (const size_type index1,
476  const size_type index2) const;
477 
489 
494  bool is_inhomogeneously_constrained (const size_type index) const;
495 
501  bool has_inhomogeneities () const;
502 
507  const std::vector<std::pair<size_type,double> > *
508  get_constraint_entries (const size_type line) const;
509 
514  double get_inhomogeneity (const size_type line) const;
515 
536  void print (std::ostream &out) const;
537 
550  void write_dot (std::ostream &) const;
551 
556  std::size_t memory_consumption () const;
557 
564  void resolve_indices(std::vector<types::global_dof_index> &indices) const;
565 
593  void condense (SparsityPattern &sparsity) const;
594 
598  void condense (BlockSparsityPattern &sparsity) const;
599 
604  void condense (DynamicSparsityPattern &sparsity) const;
605 
610  void condense (BlockDynamicSparsityPattern &sparsity) const;
611 
619  template <typename number>
620  void condense (SparseMatrix<number> &matrix) const;
621 
625  template <typename number>
626  void condense (BlockSparseMatrix<number> &matrix) const;
627 
639  template <class VectorType>
640  void condense (VectorType &vec) const;
641 
648  template <class VectorType>
649  void condense (const VectorType &vec_ghosted,
650  VectorType &output) const;
651 
664  template <typename number, class VectorType>
665  void condense (SparseMatrix<number> &matrix,
666  VectorType &vector) const;
667 
672  template <typename number, class BlockVectorType>
673  void condense (BlockSparseMatrix<number> &matrix,
674  BlockVectorType &vector) const;
675 
682  template <class VectorType>
683  void set_zero (VectorType &vec) const;
684 
736  template <class InVector, class OutVector>
737  void
738  distribute_local_to_global (const InVector &local_vector,
739  const std::vector<size_type> &local_dof_indices,
740  OutVector &global_vector) const;
741 
785  template <typename VectorType, typename LocalType>
786  void
787  distribute_local_to_global (const Vector<LocalType> &local_vector,
788  const std::vector<size_type> &local_dof_indices,
789  VectorType &global_vector,
790  const FullMatrix<LocalType> &local_matrix) const;
791 
809  template <typename VectorType, typename LocalType>
810  void
811  distribute_local_to_global (const Vector<LocalType> &local_vector,
812  const std::vector<size_type> &local_dof_indices_row,
813  const std::vector<size_type> &local_dof_indices_col,
814  VectorType &global_vector,
815  const FullMatrix<LocalType> &local_matrix,
816  bool diagonal = false) const;
817 
821  template <class VectorType>
822  void
824  const double value,
825  VectorType &global_vector) const;
826 
855  template <typename ForwardIteratorVec, typename ForwardIteratorInd,
856  class VectorType>
857  void
858  distribute_local_to_global (ForwardIteratorVec local_vector_begin,
859  ForwardIteratorVec local_vector_end,
860  ForwardIteratorInd local_indices_begin,
861  VectorType &global_vector) const;
862 
910  template <typename MatrixType>
911  void
913  const std::vector<size_type> &local_dof_indices,
914  MatrixType &global_matrix) const;
915 
943  template <typename MatrixType>
944  void
946  const std::vector<size_type> &row_indices,
947  const std::vector<size_type> &col_indices,
948  MatrixType &global_matrix) const;
949 
966  template <typename MatrixType>
967  void
969  const std::vector<size_type> &row_indices,
970  const ConstraintMatrix &column_constraint_matrix,
971  const std::vector<size_type> &column_indices,
972  MatrixType &global_matrix) const;
973 
990  template <typename MatrixType, typename VectorType>
991  void
993  const Vector<typename VectorType::value_type> &local_vector,
994  const std::vector<size_type> &local_dof_indices,
995  MatrixType &global_matrix,
996  VectorType &global_vector,
997  bool use_inhomogeneities_for_rhs = false) const;
998 
1052  template <typename SparsityPatternType>
1053  void
1054  add_entries_local_to_global (const std::vector<size_type> &local_dof_indices,
1055  SparsityPatternType &sparsity_pattern,
1056  const bool keep_constrained_entries = true,
1057  const Table<2,bool> &dof_mask = default_empty_table) const;
1058 
1062  template <typename SparsityPatternType>
1063  void
1064  add_entries_local_to_global (const std::vector<size_type> &row_indices,
1065  const std::vector<size_type> &col_indices,
1066  SparsityPatternType &sparsity_pattern,
1067  const bool keep_constrained_entries = true,
1068  const Table<2,bool> &dof_mask = default_empty_table) const;
1069 
1089  template <typename ForwardIteratorVec, typename ForwardIteratorInd,
1090  class VectorType>
1091  void
1092  get_dof_values (const VectorType &global_vector,
1093  ForwardIteratorInd local_indices_begin,
1094  ForwardIteratorVec local_vector_begin,
1095  ForwardIteratorVec local_vector_end) const;
1096 
1118  template <class VectorType>
1119  void distribute (VectorType &vec) const;
1120 
1131  {
1136  typedef std::vector<std::pair<size_type,double> > Entries;
1137 
1144 
1153 
1158 
1166  bool operator < (const ConstraintLine &) const;
1167 
1173  bool operator == (const ConstraintLine &) const;
1174 
1179  std::size_t memory_consumption () const;
1180 
1184  template <class Archive>
1185  void serialize(Archive &ar, const unsigned int)
1186  {
1187  ar &index &entries &inhomogeneity;
1188  }
1189 
1190  };
1191 
1192 
1196  typedef std::vector<ConstraintLine>::const_iterator const_iterator;
1197 
1198 
1202  typedef boost::iterator_range<const_iterator> LineRange;
1203 
1204 
1213  const LineRange get_lines() const;
1214 
1215 
1243  bool is_consistent_in_parallel(const std::vector<IndexSet> &locally_owned_dofs,
1244  const IndexSet &locally_active_dofs,
1245  const MPI_Comm mpi_communicator,
1246  const bool verbose=false) const;
1247 
1248 
1267  size_type,
1268  << "The specified line " << arg1
1269  << " does not exist.");
1276  size_type, size_type, double, double,
1277  << "The entry for the indices " << arg1 << " and "
1278  << arg2 << " already exists, but the values "
1279  << arg3 << " (old) and " << arg4 << " (new) differ "
1280  << "by " << (arg4-arg3) << ".");
1287  int, int,
1288  << "You tried to constrain DoF " << arg1
1289  << " to DoF " << arg2
1290  << ", but that one is also constrained. This is not allowed!");
1297  size_type,
1298  << "Degree of freedom " << arg1
1299  << " is constrained from both object in a merge operation.");
1306  size_type,
1307  << "In the given argument a degree of freedom is constrained "
1308  << "to another DoF with number " << arg1
1309  << ", which however is constrained by this object. This is not"
1310  << " allowed.");
1317  size_type,
1318  << "The index set given to this constraint matrix indicates "
1319  << "constraints for degree of freedom " << arg1
1320  << " should not be stored by this object, but a constraint "
1321  << "is being added.");
1322 
1329  size_type,
1330  size_type,
1331  << "The index set given to this constraint matrix indicates "
1332  << "constraints using degree of freedom " << arg2
1333  << " should not be stored by this object, but a constraint "
1334  << "for degree of freedom " << arg1 <<" uses it.");
1335 
1342  int, int,
1343  << "While distributing the constraint for DoF "
1344  << arg1 << ", it turns out that one of the processors "
1345  << "who own the " << arg2
1346  << " degrees of freedom that x_" << arg1
1347  << " is constrained against does not know about "
1348  << "the constraint on x_" << arg1
1349  << ". Did you not initialize the ConstraintMatrix "
1350  << "with the appropriate locally_relevant set so "
1351  << "that every processor who owns a DoF that constrains "
1352  << "another DoF also knows about this constraint?");
1353 
1354 private:
1355 
1367  std::vector<ConstraintLine> lines;
1368 
1401  std::vector<size_type> lines_cache;
1402 
1409 
1413  bool sorted;
1414 
1419  size_type calculate_line_index (const size_type line) const;
1420 
1425  static bool check_zero_weight (const std::pair<size_type, double> &p);
1426 
1432 
1437  template <typename MatrixType, typename VectorType>
1438  void
1440  const Vector<typename VectorType::value_type> &local_vector,
1441  const std::vector<size_type> &local_dof_indices,
1442  MatrixType &global_matrix,
1443  VectorType &global_vector,
1444  bool use_inhomogeneities_for_rhs,
1445  std::integral_constant<bool, false>) const;
1446 
1451  template <typename MatrixType, typename VectorType>
1452  void
1454  const Vector<typename VectorType::value_type> &local_vector,
1455  const std::vector<size_type> &local_dof_indices,
1456  MatrixType &global_matrix,
1457  VectorType &global_vector,
1458  bool use_inhomogeneities_for_rhs,
1459  std::integral_constant<bool, true>) const;
1460 
1465  template <typename SparsityPatternType>
1466  void
1467  add_entries_local_to_global (const std::vector<size_type> &local_dof_indices,
1468  SparsityPatternType &sparsity_pattern,
1469  const bool keep_constrained_entries,
1470  const Table<2,bool> &dof_mask,
1471  std::integral_constant<bool, false>) const;
1472 
1477  template <typename SparsityPatternType>
1478  void
1479  add_entries_local_to_global (const std::vector<size_type> &local_dof_indices,
1480  SparsityPatternType &sparsity_pattern,
1481  const bool keep_constrained_entries,
1482  const Table<2,bool> &dof_mask,
1483  std::integral_constant<bool, true>) const;
1484 
1492  void
1493  make_sorted_row_list (const std::vector<size_type> &local_dof_indices,
1494  internals::GlobalRowsFromLocal &global_rows) const;
1495 
1503  void
1504  make_sorted_row_list (const std::vector<size_type> &local_dof_indices,
1505  std::vector<size_type> &active_dofs) const;
1506 
1510  template <typename MatrixScalar, typename VectorScalar>
1511  typename ProductType<VectorScalar,MatrixScalar>::type
1513  const internals::GlobalRowsFromLocal &global_rows,
1514  const Vector<VectorScalar> &local_vector,
1515  const std::vector<size_type> &local_dof_indices,
1516  const FullMatrix<MatrixScalar> &local_matrix) const;
1517 };
1518 
1519 
1520 
1521 /* ---------------- template and inline functions ----------------- */
1522 
1523 inline
1525  :
1526  lines (),
1527  local_lines (local_constraints),
1528  sorted (false)
1529 {
1530  // make sure the IndexSet is compressed. Otherwise this can lead to crashes
1531  // that are hard to find (only happen in release mode).
1532  // see tests/mpi/constraint_matrix_crash_01
1534 }
1535 
1536 
1537 
1538 inline
1540  :
1541  Subscriptor (),
1542  lines (constraint_matrix.lines),
1543  lines_cache (constraint_matrix.lines_cache),
1544  local_lines (constraint_matrix.local_lines),
1545  sorted (constraint_matrix.sorted)
1546 {}
1547 
1548 
1549 inline
1550 void
1552 {
1553  Assert (sorted==false, ExcMatrixIsClosed());
1554 
1555  // the following can happen when we compute with distributed meshes and dof
1556  // handlers and we constrain a degree of freedom whose number we don't have
1557  // locally. if we don't abort here the program will try to allocate several
1558  // terabytes of memory to resize the various arrays below :-)
1560  ExcInternalError());
1561  const size_type line_index = calculate_line_index (line);
1562 
1563  // check whether line already exists; it may, in which case we can just quit
1564  if (is_constrained(line))
1565  return;
1566 
1567  // if necessary enlarge vector of existing entries for cache
1568  if (line_index >= lines_cache.size())
1569  lines_cache.resize (std::max(2*static_cast<size_type>(lines_cache.size()),
1570  line_index+1),
1572 
1573  // push a new line to the end of the list
1574  lines.emplace_back ();
1575  lines.back().index = line;
1576  lines.back().inhomogeneity = 0.;
1577  lines_cache[line_index] = lines.size()-1;
1578 }
1579 
1580 
1581 
1582 inline
1583 void
1585  const size_type column,
1586  const double value)
1587 {
1588  Assert (sorted==false, ExcMatrixIsClosed());
1589  Assert (line != column,
1590  ExcMessage ("Can't constrain a degree of freedom to itself"));
1591 
1592  // Ensure that the current line is present in the cache:
1593  const size_type line_index = calculate_line_index(line);
1594  Assert (line_index < lines_cache.size(),
1595  ExcMessage("The current ConstraintMatrix does not contain the line "
1596  "for the current entry. Call ConstraintMatrix::add_line "
1597  "before calling this function."));
1598 
1599  // if in debug mode, check whether an entry for this column already exists
1600  // and if it's the same as the one entered at present
1601  //
1602  // in any case: exit the function if an entry for this column already
1603  // exists, since we don't want to enter it twice
1605  ExcInternalError());
1606  Assert (!local_lines.size() || local_lines.is_element(column),
1607  ExcColumnNotStoredHere(line, column));
1608  ConstraintLine *line_ptr = &lines[lines_cache[line_index]];
1609  Assert (line_ptr->index == line, ExcInternalError());
1610  for (ConstraintLine::Entries::const_iterator
1611  p=line_ptr->entries.begin();
1612  p != line_ptr->entries.end(); ++p)
1613  if (p->first == column)
1614  {
1615  Assert (std::fabs(p->second - value) < 1.e-14,
1616  ExcEntryAlreadyExists(line, column, p->second, value));
1617  return;
1618  }
1619 
1620  line_ptr->entries.emplace_back (column, value);
1621 }
1622 
1623 
1624 
1625 inline
1626 void
1628  const double value)
1629 {
1630  const size_type line_index = calculate_line_index(line);
1631  Assert( line_index < lines_cache.size() &&
1632  lines_cache[line_index] != numbers::invalid_size_type,
1633  ExcMessage("call add_line() before calling set_inhomogeneity()"));
1634  Assert(lines_cache[line_index] < lines.size(), ExcInternalError());
1635  ConstraintLine *line_ptr = &lines[lines_cache[line_index]];
1636  line_ptr->inhomogeneity = value;
1637 }
1638 
1639 
1640 
1641 inline
1644 {
1645  return lines.size();
1646 }
1647 
1648 
1649 
1650 inline
1651 bool
1653 {
1654  const size_type line_index = calculate_line_index(index);
1655  return ((line_index < lines_cache.size())
1656  &&
1657  (lines_cache[line_index] != numbers::invalid_size_type));
1658 }
1659 
1660 
1661 
1662 inline
1663 bool
1665 {
1666  // check whether the entry is constrained. could use is_constrained, but
1667  // that means computing the line index twice
1668  const size_type line_index = calculate_line_index(index);
1669  if (line_index >= lines_cache.size() ||
1670  lines_cache[line_index] == numbers::invalid_size_type)
1671  return false;
1672  else
1673  {
1674  Assert(lines_cache[line_index] < lines.size(), ExcInternalError());
1675  return !(lines[lines_cache[line_index]].inhomogeneity == 0);
1676  }
1677 }
1678 
1679 
1680 
1681 inline
1682 const std::vector<std::pair<types::global_dof_index,double> > *
1684 {
1685  // check whether the entry is constrained. could use is_constrained, but
1686  // that means computing the line index twice
1687  const size_type line_index = calculate_line_index(line);
1688  if (line_index >= lines_cache.size() ||
1689  lines_cache[line_index] == numbers::invalid_size_type)
1690  return nullptr;
1691  else
1692  return &lines[lines_cache[line_index]].entries;
1693 }
1694 
1695 
1696 
1697 inline
1698 double
1700 {
1701  // check whether the entry is constrained. could use is_constrained, but
1702  // that means computing the line index twice
1703  const size_type line_index = calculate_line_index(line);
1704  if (line_index >= lines_cache.size() ||
1705  lines_cache[line_index] == numbers::invalid_size_type)
1706  return 0;
1707  else
1708  return lines[lines_cache[line_index]].inhomogeneity;
1709 }
1710 
1711 
1712 
1715 {
1716  //IndexSet is unused (serial case)
1717  if (!local_lines.size())
1718  return line;
1719 
1721  ExcRowNotStoredHere(line));
1722 
1723  return local_lines.index_within_set(line);
1724 }
1725 
1726 
1727 
1728 inline bool
1730 {
1731  return !local_lines.size() || local_lines.is_element(line_index);
1732 }
1733 
1734 
1735 
1736 inline
1737 const IndexSet &
1739 {
1740  return local_lines;
1741 }
1742 
1743 
1744 
1745 template <class VectorType>
1746 inline
1748  const size_type index,
1749  const double value,
1750  VectorType &global_vector) const
1751 {
1752  Assert (lines.empty() || sorted == true, ExcMatrixNotClosed());
1753 
1754  if (is_constrained(index) == false)
1755  global_vector(index) += value;
1756  else
1757  {
1758  const ConstraintLine &position =
1760  for (size_type j=0; j<position.entries.size(); ++j)
1761  global_vector(position.entries[j].first)
1762  += value * position.entries[j].second;
1763  }
1764 }
1765 
1766 
1767 template <typename ForwardIteratorVec, typename ForwardIteratorInd,
1768  class VectorType>
1769 inline
1771  ForwardIteratorVec local_vector_begin,
1772  ForwardIteratorVec local_vector_end,
1773  ForwardIteratorInd local_indices_begin,
1774  VectorType &global_vector) const
1775 {
1776  Assert (lines.empty() || sorted == true, ExcMatrixNotClosed());
1777  for ( ; local_vector_begin != local_vector_end;
1778  ++local_vector_begin, ++local_indices_begin)
1779  {
1780  if (is_constrained(*local_indices_begin) == false)
1781  internal::ElementAccess<VectorType>::add(*local_vector_begin,
1782  *local_indices_begin, global_vector);
1783  else
1784  {
1785  const ConstraintLine &position =
1786  lines[lines_cache[calculate_line_index(*local_indices_begin)]];
1787  for (size_type j=0; j<position.entries.size(); ++j)
1788  internal::ElementAccess<VectorType>::add(*local_vector_begin *
1789  position.entries[j].second, position.entries[j].first,
1790  global_vector);
1791  }
1792  }
1793 }
1794 
1795 
1796 template <class InVector, class OutVector>
1797 inline
1798 void
1800  const InVector &local_vector,
1801  const std::vector<size_type> &local_dof_indices,
1802  OutVector &global_vector) const
1803 {
1804  Assert (local_vector.size() == local_dof_indices.size(),
1805  ExcDimensionMismatch(local_vector.size(), local_dof_indices.size()));
1806  distribute_local_to_global (local_vector.begin(), local_vector.end(),
1807  local_dof_indices.begin(), global_vector);
1808 }
1809 
1810 
1811 
1812 template <typename ForwardIteratorVec, typename ForwardIteratorInd,
1813  class VectorType>
1814 inline
1815 void ConstraintMatrix::get_dof_values (const VectorType &global_vector,
1816  ForwardIteratorInd local_indices_begin,
1817  ForwardIteratorVec local_vector_begin,
1818  ForwardIteratorVec local_vector_end) const
1819 {
1820  Assert (lines.empty() || sorted == true, ExcMatrixNotClosed());
1821  for ( ; local_vector_begin != local_vector_end;
1822  ++local_vector_begin, ++local_indices_begin)
1823  {
1824  if (is_constrained(*local_indices_begin) == false)
1825  *local_vector_begin = global_vector(*local_indices_begin);
1826  else
1827  {
1828  const ConstraintLine &position =
1829  lines[lines_cache[calculate_line_index(*local_indices_begin)]];
1830  typename VectorType::value_type value = position.inhomogeneity;
1831  for (size_type j=0; j<position.entries.size(); ++j)
1832  value += (global_vector(position.entries[j].first) *
1833  position.entries[j].second);
1834  *local_vector_begin = value;
1835  }
1836  }
1837 }
1838 
1839 
1840 template <typename MatrixType> class BlockMatrixBase;
1841 template <typename SparsityPatternType> class BlockSparsityPatternBase;
1842 template <typename number> class BlockSparseMatrixEZ;
1843 
1862 template <typename MatrixType>
1864 {
1865 private:
1866  struct yes_type
1867  {
1868  char c[1];
1869  };
1870  struct no_type
1871  {
1872  char c[2];
1873  };
1874 
1880  template <typename T>
1881  static yes_type check_for_block_matrix (const BlockMatrixBase<T> *);
1882 
1887  template <typename T>
1888  static yes_type check_for_block_matrix (const BlockSparsityPatternBase<T> *);
1889 
1894  template <typename T>
1895  static yes_type check_for_block_matrix (const BlockSparseMatrixEZ<T> *);
1896 
1901  static no_type check_for_block_matrix (...);
1902 
1903 public:
1909  static const bool value = (sizeof(check_for_block_matrix
1910  ((MatrixType *)nullptr))
1911  ==
1912  sizeof(yes_type));
1913 };
1914 
1915 
1916 // instantiation of the static member
1917 template <typename MatrixType>
1919 
1920 
1921 template <typename MatrixType>
1922 inline
1923 void
1926  const std::vector<size_type> &local_dof_indices,
1927  MatrixType &global_matrix) const
1928 {
1929  // create a dummy and hand on to the function actually implementing this
1930  // feature in the cm.templates.h file.
1931  Vector<typename MatrixType::value_type> dummy(0);
1932  distribute_local_to_global (local_matrix, dummy, local_dof_indices,
1933  global_matrix, dummy, false,
1934  std::integral_constant<bool, IsBlockMatrix<MatrixType>::value>());
1935 }
1936 
1937 
1938 
1939 
1940 template <typename MatrixType, typename VectorType>
1941 inline
1942 void
1945  const Vector<typename VectorType::value_type> &local_vector,
1946  const std::vector<size_type> &local_dof_indices,
1947  MatrixType &global_matrix,
1948  VectorType &global_vector,
1949  bool use_inhomogeneities_for_rhs) const
1950 {
1951  // enter the internal function with the respective block information set,
1952  // the actual implementation follows in the cm.templates.h file.
1953  distribute_local_to_global (local_matrix, local_vector, local_dof_indices,
1954  global_matrix, global_vector, use_inhomogeneities_for_rhs,
1955  std::integral_constant<bool, IsBlockMatrix<MatrixType>::value>());
1956 }
1957 
1958 
1959 
1960 
1961 template <typename SparsityPatternType>
1962 inline
1963 void
1965 add_entries_local_to_global (const std::vector<size_type> &local_dof_indices,
1966  SparsityPatternType &sparsity_pattern,
1967  const bool keep_constrained_entries,
1968  const Table<2,bool> &dof_mask) const
1969 {
1970  // enter the internal function with the respective block information set,
1971  // the actual implementation follows in the cm.templates.h file.
1972  add_entries_local_to_global (local_dof_indices, sparsity_pattern,
1973  keep_constrained_entries, dof_mask,
1974  std::integral_constant<bool, IsBlockMatrix<SparsityPatternType>::value>());
1975 }
1976 
1977 
1978 DEAL_II_NAMESPACE_CLOSE
1979 
1980 #endif
const types::global_dof_index invalid_size_type
Definition: types.h:182
bool operator<(const ConstraintLine &) const
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:358
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)
std::vector< ConstraintLine >::const_iterator const_iterator
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)
std::size_t memory_consumption() const
void set_zero(VectorType &vec) const
boost::iterator_range< const_iterator > LineRange
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:1575
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:346
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
ConstraintMatrix & operator=(const ConstraintMatrix &)=delete
bool is_inhomogeneously_constrained(const size_type index) const
#define Assert(cond, exc)
Definition: exceptions.h:1142
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:1811
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:323
bool is_consistent_in_parallel(const std::vector< IndexSet > &locally_owned_dofs, const IndexSet &locally_active_dofs, const MPI_Comm mpi_communicator, const bool verbose=false) const
void resolve_indices(std::vector< types::global_dof_index > &indices) const
bool has_inhomogeneities() const
ProductType< VectorScalar, MatrixScalar >::type resolve_vector_entry(const size_type i, const internals::GlobalRowsFromLocal &global_rows, const Vector< VectorScalar > &local_vector, const std::vector< size_type > &local_dof_indices, const FullMatrix< MatrixScalar > &local_matrix) 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()
static ::ExceptionBase & ExcDoFConstrainedToConstrainedDoF(int arg1, int arg2)
void distribute(VectorType &vec) const
void compress() const
Definition: index_set.h:1584
void write_dot(std::ostream &) const
void copy_from(const ConstraintMatrix &other)
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:382
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:34
bool is_element(const size_type index) const
Definition: index_set.h:1645
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)
const LineRange get_lines() const
std::vector< ConstraintLine > lines
bool is_constrained(const size_type index) const
static ::ExceptionBase & ExcMatrixIsClosed()
static ::ExceptionBase & ExcInternalError()
void serialize(Archive &ar, const unsigned int)