16 #ifndef dealii_affine_constraints_h 17 #define dealii_affine_constraints_h 19 #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/table.h> 25 #include <deal.II/base/template_constraints.h> 27 #include <deal.II/lac/vector.h> 28 #include <deal.II/lac/vector_element_access.h> 30 #include <boost/range/iterator_range.hpp> 36 DEAL_II_NAMESPACE_OPEN
44 template <
typename number>
46 template <
typename number>
51 template <
typename number>
52 class GlobalRowsFromLocal;
56 template <
typename number>
155 template <
typename number =
double>
376 const std::vector<std::pair<size_type, number>> &col_val_pairs);
438 const bool allow_different_local_lines =
false);
548 const std::vector<std::pair<size_type, number>> *
579 print(std::ostream &out)
const;
687 template <
class VectorType>
697 template <
class VectorType>
699 condense(
const VectorType &vec_ghosted, VectorType &output)
const;
713 template <
class VectorType>
721 template <
class BlockVectorType>
731 template <
class VectorType>
786 template <
class InVector,
class OutVector>
789 const std::vector<size_type> &local_dof_indices,
790 OutVector & global_vector)
const;
835 template <
typename VectorType>
838 const std::vector<size_type> &local_dof_indices,
839 VectorType & global_vector,
861 template <
typename VectorType>
865 const std::vector<size_type> &local_dof_indices_row,
866 const std::vector<size_type> &local_dof_indices_col,
867 VectorType & global_vector,
869 bool diagonal =
false)
const;
874 template <
class VectorType>
878 VectorType & global_vector)
const;
908 template <
typename ForwardIteratorVec,
909 typename ForwardIteratorInd,
913 ForwardIteratorVec local_vector_end,
914 ForwardIteratorInd local_indices_begin,
915 VectorType & global_vector)
const;
964 template <
typename MatrixType>
967 const std::vector<size_type> &local_dof_indices,
968 MatrixType & global_matrix)
const;
997 template <
typename MatrixType>
1000 const std::vector<size_type> &row_indices,
1001 const std::vector<size_type> &col_indices,
1002 MatrixType & global_matrix)
const;
1020 template <
typename MatrixType>
1023 const std::vector<size_type> &row_indices,
1025 const std::vector<size_type> &column_indices,
1026 MatrixType & global_matrix)
const;
1044 template <
typename MatrixType,
typename VectorType>
1048 const std::vector<size_type> &local_dof_indices,
1049 MatrixType & global_matrix,
1050 VectorType & global_vector,
1051 bool use_inhomogeneities_for_rhs =
false)
const;
1106 template <
typename SparsityPatternType>
1109 const std::vector<size_type> &local_dof_indices,
1110 SparsityPatternType & sparsity_pattern,
1111 const bool keep_constrained_entries =
true,
1117 template <
typename SparsityPatternType>
1120 const std::vector<size_type> &row_indices,
1121 const std::vector<size_type> &col_indices,
1122 SparsityPatternType & sparsity_pattern,
1123 const bool keep_constrained_entries =
true,
1145 template <
typename ForwardIteratorVec,
1146 typename ForwardIteratorInd,
1150 ForwardIteratorInd local_indices_begin,
1151 ForwardIteratorVec local_vector_begin,
1152 ForwardIteratorVec local_vector_end)
const;
1175 template <
class VectorType>
1192 using Entries = std::vector<std::pair<size_type, number>>;
1243 template <
class Archive>
1301 const IndexSet & locally_active_dofs,
1302 const MPI_Comm mpi_communicator,
1303 const bool verbose =
false)
const;
1324 <<
"The specified line " << arg1 <<
" does not exist.");
1335 <<
"The entry for the indices " << arg1 <<
" and " << arg2
1336 <<
" already exists, but the values " << arg3 <<
" (old) and " 1337 << arg4 <<
" (new) differ " 1338 <<
"by " << (arg4 - arg3) <<
".");
1347 <<
"You tried to constrain DoF " << arg1 <<
" to DoF " << arg2
1348 <<
", but that one is also constrained. This is not allowed!");
1356 <<
"Degree of freedom " << arg1
1357 <<
" is constrained from both object in a merge operation.");
1365 <<
"In the given argument a degree of freedom is constrained " 1366 <<
"to another DoF with number " << arg1
1367 <<
", which however is constrained by this object. This is not" 1376 <<
"The index set given to this constraints object indicates " 1377 <<
"constraints for degree of freedom " << arg1
1378 <<
" should not be stored by this object, but a constraint " 1379 <<
"is being added.");
1389 <<
"The index set given to this constraints object indicates " 1390 <<
"constraints using degree of freedom " << arg2
1391 <<
" should not be stored by this object, but a constraint " 1392 <<
"for degree of freedom " << arg1 <<
" uses it.");
1402 <<
"While distributing the constraint for DoF " << arg1
1403 <<
", it turns out that one of the processors " 1404 <<
"who own the " << arg2 <<
" degrees of freedom that x_" 1405 << arg1 <<
" is constrained against does not know about " 1406 <<
"the constraint on x_" << arg1
1407 <<
". Did you not initialize the AffineConstraints container " 1408 <<
"with the appropriate locally_relevant set so " 1409 <<
"that every processor who owns a DoF that constrains " 1410 <<
"another DoF also knows about this constraint?");
1483 template <
typename MatrixType,
typename VectorType>
1487 const std::vector<size_type> &local_dof_indices,
1488 MatrixType & global_matrix,
1489 VectorType & global_vector,
1490 bool use_inhomogeneities_for_rhs,
1491 std::integral_constant<bool, false>)
const;
1497 template <
typename MatrixType,
typename VectorType>
1501 const std::vector<size_type> &local_dof_indices,
1502 MatrixType & global_matrix,
1503 VectorType & global_vector,
1504 bool use_inhomogeneities_for_rhs,
1505 std::integral_constant<bool, true>)
const;
1511 template <
typename SparsityPatternType>
1514 SparsityPatternType & sparsity_pattern,
1515 const bool keep_constrained_entries,
1517 std::integral_constant<bool, false>)
const;
1523 template <
typename SparsityPatternType>
1526 SparsityPatternType & sparsity_pattern,
1527 const bool keep_constrained_entries,
1529 std::integral_constant<bool, true>)
const;
1540 const std::vector<size_type> & local_dof_indices,
1541 internals::GlobalRowsFromLocal<number> &global_rows)
const;
1552 std::vector<size_type> & active_dofs)
const;
1557 template <
typename MatrixScalar,
typename VectorScalar>
1558 typename ProductType<VectorScalar, MatrixScalar>::type
1561 const internals::GlobalRowsFromLocal<number> &global_rows,
1562 const Vector<VectorScalar> & local_vector,
1563 const std::vector<size_type> & local_dof_indices,
1569 template <
typename number>
1573 , local_lines(local_constraints)
1582 template <
typename number>
1586 , lines(affine_constraints.lines)
1587 , lines_cache(affine_constraints.lines_cache)
1588 , local_lines(affine_constraints.local_lines)
1589 , sorted(affine_constraints.sorted)
1592 template <
typename number>
1596 Assert(sorted ==
false, ExcMatrixIsClosed());
1603 const size_type line_index = calculate_line_index(line_n);
1606 if (is_constrained(line_n))
1610 if (line_index >= lines_cache.size())
1611 lines_cache.resize(std::max(2 * static_cast<size_type>(lines_cache.size()),
1616 lines.emplace_back();
1617 lines.back().index = line_n;
1618 lines.back().inhomogeneity = 0.;
1619 lines_cache[line_index] = lines.size() - 1;
1622 template <
typename number>
1628 Assert(sorted ==
false, ExcMatrixIsClosed());
1630 ExcMessage(
"Can't constrain a degree of freedom to itself"));
1633 const size_type line_index = calculate_line_index(line_n);
1634 Assert(line_index < lines_cache.size(),
1635 ExcMessage(
"The current AffineConstraints does not contain the line " 1636 "for the current entry. Call AffineConstraints::add_line " 1637 "before calling this function."));
1646 Assert(!local_lines.size() || local_lines.is_element(column),
1647 ExcColumnNotStoredHere(line_n, column));
1648 ConstraintLine *line_ptr = &lines[lines_cache[line_index]];
1650 for (
const auto &p : line_ptr->entries)
1651 if (p.first == column)
1653 Assert(std::abs(p.second - value) < 1.e-14,
1654 ExcEntryAlreadyExists(line_n, column, p.second, value));
1658 line_ptr->entries.emplace_back(column, value);
1661 template <
typename number>
1666 const size_type line_index = calculate_line_index(line_n);
1667 Assert(line_index < lines_cache.size() &&
1669 ExcMessage(
"call add_line() before calling set_inhomogeneity()"));
1671 ConstraintLine *line_ptr = &lines[lines_cache[line_index]];
1672 line_ptr->inhomogeneity = value;
1675 template <
typename number>
1679 return lines.size();
1682 template <
typename number>
1686 const size_type line_index = calculate_line_index(index);
1687 return ((line_index < lines_cache.size()) &&
1691 template <
typename number>
1698 const size_type line_index = calculate_line_index(line_n);
1699 if (line_index >= lines_cache.size() ||
1705 return !(lines[lines_cache[line_index]].inhomogeneity == number(0.));
1709 template <
typename number>
1710 inline const std::vector<std::pair<types::global_dof_index, number>> *
1715 const size_type line_index = calculate_line_index(line_n);
1716 if (line_index >= lines_cache.size() ||
1720 return &lines[lines_cache[line_index]].entries;
1723 template <
typename number>
1729 const size_type line_index = calculate_line_index(line_n);
1730 if (line_index >= lines_cache.size() ||
1734 return lines[lines_cache[line_index]].inhomogeneity;
1737 template <
typename number>
1742 if (!local_lines.size())
1745 Assert(local_lines.is_element(line_n), ExcRowNotStoredHere(line_n));
1747 return local_lines.index_within_set(line_n);
1750 template <
typename number>
1754 return !local_lines.size() || local_lines.is_element(line_n);
1757 template <
typename number>
1764 template <
typename number>
1765 template <
class VectorType>
1770 VectorType & global_vector)
const 1772 Assert(lines.empty() || sorted ==
true, ExcMatrixNotClosed());
1774 if (is_constrained(index) ==
false)
1775 global_vector(index) += value;
1778 const ConstraintLine &position =
1779 lines[lines_cache[calculate_line_index(index)]];
1780 for (
size_type j = 0; j < position.entries.size(); ++j)
1781 global_vector(position.entries[j].first) +=
1782 value * position.entries[j].second;
1786 template <
typename number>
1787 template <
typename ForwardIteratorVec,
1788 typename ForwardIteratorInd,
1792 ForwardIteratorVec local_vector_begin,
1793 ForwardIteratorVec local_vector_end,
1794 ForwardIteratorInd local_indices_begin,
1795 VectorType & global_vector)
const 1797 Assert(lines.empty() || sorted ==
true, ExcMatrixNotClosed());
1798 for (; local_vector_begin != local_vector_end;
1799 ++local_vector_begin, ++local_indices_begin)
1801 if (is_constrained(*local_indices_begin) ==
false)
1802 internal::ElementAccess<VectorType>::add(*local_vector_begin,
1803 *local_indices_begin,
1807 const ConstraintLine &position =
1808 lines[lines_cache[calculate_line_index(*local_indices_begin)]];
1809 for (
size_type j = 0; j < position.entries.size(); ++j)
1810 internal::ElementAccess<VectorType>::add(
1811 (*local_vector_begin) * position.entries[j].second,
1812 position.entries[j].first,
1818 template <
typename number>
1819 template <
class InVector,
class OutVector>
1822 const InVector & local_vector,
1823 const std::vector<size_type> &local_dof_indices,
1824 OutVector & global_vector)
const 1826 Assert(local_vector.size() == local_dof_indices.size(),
1828 distribute_local_to_global(local_vector.begin(),
1830 local_dof_indices.begin(),
1834 template <
typename number>
1835 template <
typename ForwardIteratorVec,
1836 typename ForwardIteratorInd,
1840 const VectorType & global_vector,
1841 ForwardIteratorInd local_indices_begin,
1842 ForwardIteratorVec local_vector_begin,
1843 ForwardIteratorVec local_vector_end)
const 1845 Assert(lines.empty() || sorted ==
true, ExcMatrixNotClosed());
1846 for (; local_vector_begin != local_vector_end;
1847 ++local_vector_begin, ++local_indices_begin)
1849 if (is_constrained(*local_indices_begin) ==
false)
1850 *local_vector_begin = global_vector(*local_indices_begin);
1853 const ConstraintLine &position =
1854 lines[lines_cache[calculate_line_index(*local_indices_begin)]];
1855 typename VectorType::value_type value = position.inhomogeneity;
1856 for (
size_type j = 0; j < position.entries.size(); ++j)
1857 value += (global_vector(position.entries[j].first) *
1858 position.entries[j].second);
1859 *local_vector_begin = value;
1864 template <
typename MatrixType>
1866 template <
typename SparsityPatternType>
1868 template <
typename number>
1889 template <
typename MatrixType>
1907 template <
typename T>
1915 template <
typename T>
1923 template <
typename T>
1946 template <
typename MatrixType>
1949 template <
typename number>
1950 template <
typename MatrixType>
1954 const std::vector<size_type> &local_dof_indices,
1955 MatrixType & global_matrix)
const 1959 Vector<typename MatrixType::value_type> dummy(0);
1960 distribute_local_to_global(
1970 template <
typename number>
1971 template <
typename MatrixType,
typename VectorType>
1976 const std::vector<size_type> &local_dof_indices,
1977 MatrixType & global_matrix,
1978 VectorType & global_vector,
1979 bool use_inhomogeneities_for_rhs)
const 1983 distribute_local_to_global(
1989 use_inhomogeneities_for_rhs,
1993 template <
typename number>
1994 template <
typename SparsityPatternType>
1997 const std::vector<size_type> &local_dof_indices,
1998 SparsityPatternType & sparsity_pattern,
1999 const bool keep_constrained_entries,
2004 add_entries_local_to_global(
2007 keep_constrained_entries,
2012 DEAL_II_NAMESPACE_CLOSE
const types::global_dof_index invalid_size_type
std::vector< ConstraintLine > lines
void set_zero(VectorType &vec) const
bool is_constrained(const size_type line_n) const
#define DeclException2(Exception2, type1, type2, outsequence)
bool can_store_line(const size_type line_n) const
void print(std::ostream &out) const
void add_entries(const size_type line_n, const std::vector< std::pair< size_type, number >> &col_val_pairs)
void add_entry(const size_type line_n, const size_type column, const number value)
bool are_identity_constrained(const size_type line_n_1, const size_type line_n_2) const
const LineRange get_lines() const
std::size_t memory_consumption() const
static ::ExceptionBase & ExcColumnNotStoredHere(size_type arg1, size_type arg2)
types::global_dof_index size_type
void merge(const AffineConstraints &other_constraints, const MergeConflictBehavior merge_conflict_behavior=no_conflicts_allowed, const bool allow_different_local_lines=false)
static ::ExceptionBase & ExcRowNotStoredHere(size_type arg1)
void add_selected_constraints(const AffineConstraints &constraints_in, const IndexSet &filter)
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
static ::ExceptionBase & ExcMatrixNotClosed()
static yes_type check_for_block_matrix(const BlockMatrixBase< T > *)
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
bool operator<(const ConstraintLine &) const
static ::ExceptionBase & ExcMatrixIsClosed()
bool operator==(const ConstraintLine &) const
std::size_t memory_consumption() const
void serialize(Archive &ar, const unsigned int)
size_type n_constraints() const
static ::ExceptionBase & ExcMessage(std::string arg1)
void make_sorted_row_list(const std::vector< size_type > &local_dof_indices, internals::GlobalRowsFromLocal< number > &global_rows) const
#define DeclException1(Exception1, type1, outsequence)
bool has_inhomogeneities() const
static ::ExceptionBase & ExcDoFIsConstrainedToConstrainedDoF(size_type arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcLineInexistant(size_type arg1)
void shift(const size_type offset)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcIncorrectConstraint(int arg1, int arg2)
void write_dot(std::ostream &) const
#define DeclException0(Exception0)
AffineConstraints & operator=(const AffineConstraints &)=delete
static ::ExceptionBase & ExcDoFIsConstrainedFromBothObjects(size_type arg1)
static ::ExceptionBase & ExcDoFConstrainedToConstrainedDoF(int arg1, int arg2)
void resolve_indices(std::vector< types::global_dof_index > &indices) const
AffineConstraints(const IndexSet &local_constraints=IndexSet())
bool is_identity_constrained(const size_type line_n) const
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=Table< 2, bool >()) const
size_type calculate_line_index(const size_type line_n) const
std::vector< std::pair< size_type, number > > Entries
const IndexSet & get_local_lines() const
void condense(SparsityPattern &sparsity) const
std::vector< size_type > lines_cache
boost::iterator_range< const_iterator > LineRange
const std::vector< std::pair< size_type, number > > * get_constraint_entries(const size_type line_n) const
void reinit(const IndexSet &local_constraints=IndexSet())
bool is_inhomogeneously_constrained(const size_type index) const
unsigned int global_dof_index
void get_dof_values(const VectorType &global_vector, ForwardIteratorInd local_indices_begin, ForwardIteratorVec local_vector_begin, ForwardIteratorVec local_vector_end) const
number get_inhomogeneity(const size_type line_n) const
void distribute(VectorType &vec) const
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
ProductType< VectorScalar, MatrixScalar >::type resolve_vector_entry(const size_type i, const internals::GlobalRowsFromLocal< number > &global_rows, const Vector< VectorScalar > &local_vector, const std::vector< size_type > &local_dof_indices, const FullMatrix< MatrixScalar > &local_matrix) const
void copy_from(const AffineConstraints &other)
size_type max_constraint_indirections() const
void add_line(const size_type line_n)
void set_inhomogeneity(const size_type line_n, const number value)
typename std::vector< ConstraintLine >::const_iterator const_iterator
static ::ExceptionBase & ExcEntryAlreadyExists(size_type arg1, size_type arg2, number arg3, number arg4)
void add_lines(const std::vector< bool > &lines)
static ::ExceptionBase & ExcInternalError()