16#ifndef dealii_affine_constraints_h
17#define dealii_affine_constraints_h
31#include <boost/range/iterator_range.hpp>
48template <
typename number>
50template <
typename number>
81 Distributing(
const Distributing &in);
84 operator=(
const Distributing &in);
89 return global_row < in.global_row;
110 template <
typename number>
119 insert_new_index(
const std::pair<size_type, number> &pair);
123 const std::pair<size_type, number> &pair);
128 const std::pair<size_type, number> *
133 std::vector<std::pair<size_type, number>> data;
135 std::vector<size_type> individual_size;
168 template <
typename number>
169 class GlobalRowsFromLocal
175 GlobalRowsFromLocal();
183 const number constraint_value);
188 print(std::ostream &os);
201 size(
const size_type counter_index)
const;
207 global_row(
const size_type counter_index)
const;
213 global_row(
const size_type counter_index);
221 local_row(
const size_type counter_index)
const;
227 local_row(
const size_type counter_index);
236 const size_type index_in_constraint)
const;
243 constraint_value(
const size_type counter_index,
244 const size_type index_in_constraint)
const;
252 have_indirect_rows()
const;
259 insert_constraint(
const size_type constrained_local_dof);
267 n_constraints()
const;
274 n_inhomogeneities()
const;
283 set_ith_constraint_inhomogeneous(
const size_type i);
297 std::vector<Distributing> total_row_indices;
303 DataCache<number> data_cache;
337 template <
typename number>
350 ScratchData(
const ScratchData &)
362 std::vector<size_type> columns;
367 std::vector<number>
values;
372 std::vector<size_type> block_starts;
377 std::vector<size_type> vector_indices;
382 std::vector<number> vector_values;
387 GlobalRowsFromLocal<number> global_rows;
392 GlobalRowsFromLocal<number> global_columns;
399 namespace AffineConstraintsImplementation
401 template <
class VectorType>
403 set_zero_all(
const std::vector<types::global_dof_index> &cm,
408 set_zero_all(
const std::vector<types::global_dof_index> &cm,
413 set_zero_all(
const std::vector<types::global_dof_index> &cm,
419template <
typename number>
505template <
typename number =
double>
597 template <
typename other_number>
726 const number weight);
736 const std::vector<std::pair<size_type, number>> &col_weight_pairs);
804 const bool allow_different_local_lines =
false);
914 const std::vector<std::pair<size_type, number>> *
1053 template <
class VectorType>
1063 template <
class VectorType>
1065 condense(
const VectorType &vec_ghosted, VectorType &output)
const;
1079 template <
class VectorType>
1087 template <
class BlockVectorType>
1097 template <
class VectorType>
1156 template <
class InVector,
class OutVector>
1159 const std::vector<size_type> &local_dof_indices,
1160 OutVector & global_vector)
const;
1209 template <
typename VectorType>
1212 const std::vector<size_type> &local_dof_indices,
1213 VectorType & global_vector,
1235 template <
typename VectorType>
1239 const std::vector<size_type> &local_dof_indices_row,
1240 const std::vector<size_type> &local_dof_indices_col,
1241 VectorType & global_vector,
1248 template <
class VectorType>
1252 VectorType & global_vector)
const;
1286 template <
typename ForwardIteratorVec,
1287 typename ForwardIteratorInd,
1291 ForwardIteratorVec local_vector_end,
1292 ForwardIteratorInd local_indices_begin,
1293 VectorType & global_vector)
const;
1346 template <
typename MatrixType>
1349 const std::vector<size_type> &local_dof_indices,
1350 MatrixType & global_matrix)
const;
1379 template <
typename MatrixType>
1382 const std::vector<size_type> &row_indices,
1383 const std::vector<size_type> &col_indices,
1384 MatrixType & global_matrix)
const;
1402 template <
typename MatrixType>
1405 const std::vector<size_type> &row_indices,
1407 const std::vector<size_type> &column_indices,
1408 MatrixType & global_matrix)
const;
1430 template <
typename MatrixType,
typename VectorType>
1434 const std::vector<size_type> &local_dof_indices,
1435 MatrixType & global_matrix,
1436 VectorType & global_vector,
1437 bool use_inhomogeneities_for_rhs =
false)
const;
1492 template <
typename SparsityPatternType>
1495 const std::vector<size_type> &local_dof_indices,
1496 SparsityPatternType & sparsity_pattern,
1497 const bool keep_constrained_entries =
true,
1503 template <
typename SparsityPatternType>
1506 const std::vector<size_type> &row_indices,
1507 const std::vector<size_type> &col_indices,
1508 SparsityPatternType & sparsity_pattern,
1509 const bool keep_constrained_entries =
true,
1531 template <
typename ForwardIteratorVec,
1532 typename ForwardIteratorInd,
1536 ForwardIteratorInd local_indices_begin,
1537 ForwardIteratorVec local_vector_begin,
1538 ForwardIteratorVec local_vector_end)
const;
1561 template <
class VectorType>
1578 using Entries = std::vector<std::pair<size_type, number>>;
1611 template <
typename Constra
intLineType>
1617 template <
typename Constra
intLineType>
1651 template <
class Archive>
1710 const IndexSet & locally_active_dofs,
1712 const bool verbose =
false)
const;
1733 <<
"The specified line " << arg1 <<
" does not exist.");
1744 <<
"The entry for the indices " << arg1 <<
" and " << arg2
1745 <<
" already exists, but the values " << arg3 <<
" (old) and "
1746 << arg4 <<
" (new) differ "
1747 <<
"by " << (arg4 - arg3) <<
".");
1756 <<
"You tried to constrain DoF " << arg1 <<
" to DoF " << arg2
1757 <<
", but that one is also constrained. This is not allowed!");
1765 <<
"Degree of freedom " << arg1
1766 <<
" is constrained from both object in a merge operation.");
1774 <<
"In the given argument a degree of freedom is constrained "
1775 <<
"to another DoF with number " << arg1
1776 <<
", which however is constrained by this object. This is not"
1785 <<
"The index set given to this constraints object indicates "
1786 <<
"constraints for degree of freedom " << arg1
1787 <<
" should not be stored by this object, but a constraint "
1788 <<
"is being added.");
1798 <<
"The index set given to this constraints object indicates "
1799 <<
"constraints using degree of freedom " << arg2
1800 <<
" should not be stored by this object, but a constraint "
1801 <<
"for degree of freedom " << arg1 <<
" uses it.");
1811 <<
"While distributing the constraint for DoF " << arg1
1812 <<
", it turns out that one of the processors "
1813 <<
"who own the " << arg2 <<
" degrees of freedom that x_"
1814 << arg1 <<
" is constrained against does not know about "
1815 <<
"the constraint on x_" << arg1
1816 <<
". Did you not initialize the AffineConstraints container "
1817 <<
"with the appropriate locally_relevant set so "
1818 <<
"that every processor who owns a DoF that constrains "
1819 <<
"another DoF also knows about this constraint?");
1885 internal::AffineConstraints::ScratchData<number>>
1899 template <
typename MatrixType,
typename VectorType>
1903 const std::vector<size_type> &local_dof_indices,
1904 MatrixType & global_matrix,
1905 VectorType & global_vector,
1906 const bool use_inhomogeneities_for_rhs,
1907 const std::integral_constant<bool, false>)
const;
1913 template <
typename MatrixType,
typename VectorType>
1917 const std::vector<size_type> &local_dof_indices,
1918 MatrixType & global_matrix,
1919 VectorType & global_vector,
1920 const bool use_inhomogeneities_for_rhs,
1921 const std::integral_constant<bool, true>)
const;
1927 template <
typename SparsityPatternType>
1930 SparsityPatternType & sparsity_pattern,
1931 const bool keep_constrained_entries,
1933 const std::integral_constant<bool, false>)
const;
1939 template <
typename SparsityPatternType>
1942 SparsityPatternType & sparsity_pattern,
1943 const bool keep_constrained_entries,
1945 const std::integral_constant<bool, true>)
const;
1956 internal::AffineConstraints::GlobalRowsFromLocal<number>
1957 &global_rows)
const;
1968 std::vector<size_type> & active_dofs)
const;
1973 template <
typename MatrixScalar,
typename VectorScalar>
1977 const internal::AffineConstraints::GlobalRowsFromLocal<number> &global_rows,
1979 const std::vector<size_type> & local_dof_indices,
1985template <
typename number>
1998template <
typename number>
2008template <
typename number>
2032 lines.emplace_back();
2033 lines.back().index = line_n;
2034 lines.back().inhomogeneity = 0.;
2040template <
typename number>
2044 const number weight)
2047 Assert(constrained_dof_index != column,
2048 ExcMessage(
"Can't constrain a degree of freedom to itself"));
2053 ExcMessage(
"The current AffineConstraints does not contain the line "
2054 "for the current entry. Call AffineConstraints::add_line "
2055 "before calling this function."));
2068 for (
const auto &p : line_ptr->
entries)
2069 if (p.first == column)
2073 constrained_dof_index, column, p.second, weight));
2077 line_ptr->
entries.emplace_back(column, weight);
2082template <
typename number>
2091 ExcMessage(
"call add_line() before calling set_inhomogeneity()"));
2099template <
typename number>
2100template <
class VectorType>
2106 std::vector<size_type> constrained_lines(
lines.size());
2107 for (
unsigned int i = 0; i <
lines.size(); ++i)
2108 constrained_lines[i] =
lines[i].index;
2109 internal::AffineConstraintsImplementation::set_zero_all(constrained_lines,
2113template <
typename number>
2117 return lines.size();
2120template <
typename number>
2129template <
typename number>
2147template <
typename number>
2148inline const std::vector<std::pair<types::global_dof_index, number>> *
2161template <
typename number>
2175template <
typename number>
2188template <
typename number>
2195template <
typename number>
2202template <
typename number>
2203template <
class VectorType>
2208 VectorType & global_vector)
const
2213 global_vector(index) += value;
2219 global_vector(position.
entries[j].first) +=
2220 value * position.
entries[j].second;
2224template <
typename number>
2225template <
typename ForwardIteratorVec,
2226 typename ForwardIteratorInd,
2230 ForwardIteratorVec local_vector_begin,
2231 ForwardIteratorVec local_vector_end,
2232 ForwardIteratorInd local_indices_begin,
2233 VectorType & global_vector)
const
2236 for (; local_vector_begin != local_vector_end;
2237 ++local_vector_begin, ++local_indices_begin)
2241 *local_indices_begin,
2249 (*local_vector_begin) * position.
entries[j].second,
2256template <
typename number>
2257template <
class InVector,
class OutVector>
2260 const InVector & local_vector,
2261 const std::vector<size_type> &local_dof_indices,
2262 OutVector & global_vector)
const
2264 Assert(local_vector.size() == local_dof_indices.size(),
2268 local_dof_indices.begin(),
2272template <
typename number>
2273template <
typename ForwardIteratorVec,
2274 typename ForwardIteratorInd,
2278 const VectorType & global_vector,
2279 ForwardIteratorInd local_indices_begin,
2280 ForwardIteratorVec local_vector_begin,
2281 ForwardIteratorVec local_vector_end)
const
2284 for (; local_vector_begin != local_vector_end;
2285 ++local_vector_begin, ++local_indices_begin)
2288 *local_vector_begin = global_vector(*local_indices_begin);
2293 typename VectorType::value_type value = position.
inhomogeneity;
2295 value += (global_vector(position.
entries[j].first) *
2297 *local_vector_begin = value;
2302template <
typename MatrixType>
2304template <
typename SparsityPatternType>
2306template <
typename number>
2330 template <
typename MatrixType>
2339 template <
typename T>
2340 static std::true_type
2347 template <
typename T>
2348 static std::true_type
2355 static std::false_type
2366 std::is_same<decltype(check(std::declval<MatrixType *>())),
2371 template <
typename MatrixType>
2383 template <
typename MatrixType>
2391 template <
typename T>
2392 static std::true_type
2399 static std::false_type
2409 std::is_same<decltype(check(std::declval<MatrixType *>())),
2414 template <
typename MatrixType>
2422template <
typename number>
2423template <
typename other_number>
2429 lines.insert(lines.begin(), other.
lines.begin(), other.
lines.end());
2437template <
typename number>
2438template <
typename MatrixType>
2442 const std::vector<size_type> &local_dof_indices,
2443 MatrixType & global_matrix)
const
2448 distribute_local_to_global(
2455 std::integral_constant<
2462template <
typename number>
2463template <
typename MatrixType,
typename VectorType>
2468 const std::vector<size_type> &local_dof_indices,
2469 MatrixType & global_matrix,
2470 VectorType & global_vector,
2471 bool use_inhomogeneities_for_rhs)
const
2475 distribute_local_to_global(
2481 use_inhomogeneities_for_rhs,
2482 std::integral_constant<
2489template <
typename number>
2490template <
typename SparsityPatternType>
2493 const std::vector<size_type> &local_dof_indices,
2494 SparsityPatternType & sparsity_pattern,
2495 const bool keep_constrained_entries,
2500 add_entries_local_to_global(
2503 keep_constrained_entries,
2505 std::integral_constant<
bool,
2507 SparsityPatternType>::value>());
2512template <
typename number>
2516 const number &inhomogeneity)
2519 , inhomogeneity(inhomogeneity)
2524template <
typename number>
2525template <
typename Constra
intLineType>
2527 const ConstraintLineType &other)
2529 this->index = other.index;
2532 entries.insert(entries.begin(), other.entries.begin(), other.entries.end());
2534 this->inhomogeneity = other.inhomogeneity;
2539template <
typename number>
2540template <
typename Constra
intLineType>
2543operator=(
const ConstraintLineType &other)
2545 this->index = other.index;
2548 entries.insert(entries.begin(), other.entries.begin(), other.entries.end());
2550 this->inhomogeneity = other.inhomogeneity;
const LineRange get_lines() const
void add_entries(const size_type constrained_dof_index, const std::vector< std::pair< size_type, number > > &col_weight_pairs)
void condense(DynamicSparsityPattern &sparsity) const
void distribute_local_to_global(ForwardIteratorVec local_vector_begin, ForwardIteratorVec local_vector_end, ForwardIteratorInd local_indices_begin, VectorType &global_vector) const
void merge(const AffineConstraints &other_constraints, const MergeConflictBehavior merge_conflict_behavior=no_conflicts_allowed, const bool allow_different_local_lines=false)
bool has_inhomogeneities() const
void add_line(const size_type line_n)
void distribute_local_to_global(const FullMatrix< number > &local_matrix, const Vector< number > &local_vector, const std::vector< size_type > &local_dof_indices, MatrixType &global_matrix, VectorType &global_vector, const bool use_inhomogeneities_for_rhs, const std::integral_constant< bool, true >) const
void add_selected_constraints(const AffineConstraints &constraints_in, const IndexSet &filter)
void condense(VectorType &vec) const
bool can_store_line(const size_type line_n) const
AffineConstraints(AffineConstraints &&affine_constraints) noexcept=default
Threads::ThreadLocalStorage< internal::AffineConstraints::ScratchData< number > > scratch_data
friend class AffineConstraints
void add_entry(const size_type constrained_dof_index, const size_type column, const number weight)
void reinit(const IndexSet &local_constraints=IndexSet())
number get_inhomogeneity(const size_type line_n) const
void make_sorted_row_list(const std::vector< size_type > &local_dof_indices, std::vector< size_type > &active_dofs) const
void get_dof_values(const VectorType &global_vector, ForwardIteratorInd local_indices_begin, ForwardIteratorVec local_vector_begin, ForwardIteratorVec local_vector_end) const
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
void make_sorted_row_list(const std::vector< size_type > &local_dof_indices, internal::AffineConstraints::GlobalRowsFromLocal< number > &global_rows) const
const IndexSet & get_local_lines() const
void shift(const size_type offset)
std::size_t memory_consumption() const
void distribute_local_to_global(const Vector< number > &local_vector, const std::vector< size_type > &local_dof_indices_row, const std::vector< size_type > &local_dof_indices_col, VectorType &global_vector, const FullMatrix< number > &local_matrix, bool diagonal=false) const
typename std::vector< ConstraintLine >::const_iterator const_iterator
void condense(SparseMatrix< number > &matrix, VectorType &vector) const
void set_inhomogeneity(const size_type constrained_dof_index, const number value)
void distribute_local_to_global(const FullMatrix< number > &local_matrix, const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, MatrixType &global_matrix) const
void condense(BlockSparsityPattern &sparsity) const
void distribute_local_to_global(const FullMatrix< number > &local_matrix, const Vector< number > &local_vector, const std::vector< size_type > &local_dof_indices, MatrixType &global_matrix, VectorType &global_vector, const bool use_inhomogeneities_for_rhs, const std::integral_constant< bool, false >) const
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 condense(SparsityPattern &sparsity) const
void add_entries_local_to_global(const std::vector< size_type > &row_indices, const std::vector< size_type > &col_indices, SparsityPatternType &sparsity_pattern, const bool keep_constrained_entries=true, const Table< 2, bool > &dof_mask=Table< 2, bool >()) const
void distribute_local_to_global(const size_type index, const number value, VectorType &global_vector) const
size_type max_constraint_indirections() const
void add_entries_local_to_global(const std::vector< size_type > &local_dof_indices, SparsityPatternType &sparsity_pattern, const bool keep_constrained_entries, const Table< 2, bool > &dof_mask, const std::integral_constant< bool, false >) const
ProductType< VectorScalar, MatrixScalar >::type resolve_vector_entry(const size_type i, const internal::AffineConstraints::GlobalRowsFromLocal< number > &global_rows, const Vector< VectorScalar > &local_vector, const std::vector< size_type > &local_dof_indices, const FullMatrix< MatrixScalar > &local_matrix) const
void add_lines(const std::set< size_type > &lines)
AffineConstraints(const AffineConstraints &affine_constraints)
void distribute(VectorType &vec) const
void add_lines(const IndexSet &lines)
void condense(BlockSparseMatrix< number > &matrix) const
void resolve_indices(std::vector< types::global_dof_index > &indices) const
std::vector< ConstraintLine > lines
bool is_constrained(const size_type line_n) const
void condense(BlockDynamicSparsityPattern &sparsity) const
void distribute_local_to_global(const FullMatrix< number > &local_matrix, const std::vector< size_type > &local_dof_indices, MatrixType &global_matrix) const
const std::vector< std::pair< size_type, number > > * get_constraint_entries(const size_type line_n) const
size_type calculate_line_index(const size_type line_n) const
AffineConstraints & operator=(AffineConstraints &&affine_constraints) noexcept=default
void condense(SparseMatrix< number > &matrix) const
void add_entries_local_to_global(const std::vector< size_type > &local_dof_indices, SparsityPatternType &sparsity_pattern, const bool keep_constrained_entries, const Table< 2, bool > &dof_mask, const std::integral_constant< bool, true >) const
bool is_inhomogeneously_constrained(const size_type index) const
void distribute_local_to_global(const FullMatrix< number > &local_matrix, const Vector< number > &local_vector, const std::vector< size_type > &local_dof_indices, MatrixType &global_matrix, VectorType &global_vector, bool use_inhomogeneities_for_rhs=false) const
void condense(const VectorType &vec_ghosted, VectorType &output) const
void copy_from(const AffineConstraints< other_number > &other)
std::vector< size_type > lines_cache
AffineConstraints(const IndexSet &local_constraints=IndexSet())
void condense(BlockSparseMatrix< number > &matrix, BlockVectorType &vector) 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
bool is_identity_constrained(const size_type line_n) const
size_type n_constraints() const
void print(std::ostream &out) const
void write_dot(std::ostream &) const
void distribute_local_to_global(const Vector< number > &local_vector, const std::vector< size_type > &local_dof_indices, VectorType &global_vector, const FullMatrix< number > &local_matrix) const
void set_zero(VectorType &vec) const
void distribute_local_to_global(const FullMatrix< number > &local_matrix, const std::vector< size_type > &row_indices, const AffineConstraints &column_affine_constraints, const std::vector< size_type > &column_indices, MatrixType &global_matrix) const
void add_lines(const std::vector< bool > &lines)
boost::iterator_range< const_iterator > LineRange
bool are_identity_constrained(const size_type line_n_1, const size_type line_n_2) const
AffineConstraints & operator=(const AffineConstraints &)=delete
size_type index_within_set(const size_type global_index) const
bool is_element(const size_type index) const
A class that provides a separate storage location on each thread that accesses the object.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DeclException0(Exception0)
static ::ExceptionBase & ExcMatrixNotClosed()
static ::ExceptionBase & ExcLineInexistant(size_type arg1)
#define DeclException4(Exception4, type1, type2, type3, type4, outsequence)
static ::ExceptionBase & ExcDoFIsConstrainedFromBothObjects(size_type arg1)
static ::ExceptionBase & ExcDoFConstrainedToConstrainedDoF(int arg1, int arg2)
#define Assert(cond, exc)
static ::ExceptionBase & ExcMatrixIsClosed()
#define DeclException2(Exception2, type1, type2, outsequence)
static ::ExceptionBase & ExcDoFIsConstrainedToConstrainedDoF(size_type arg1)
static ::ExceptionBase & ExcColumnNotStoredHere(size_type arg1, size_type arg2)
static ::ExceptionBase & ExcRowNotStoredHere(size_type arg1)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcIncorrectConstraint(int arg1, int arg2)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcEntryAlreadyExists(size_type arg1, size_type arg2, number arg3, number arg4)
static ::ExceptionBase & ExcMessage(std::string arg1)
@ matrix
Contents is actually a matrix.
@ diagonal
Matrix is diagonal.
types::global_dof_index size_type
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
const types::global_dof_index invalid_dof_index
const types::global_dof_index invalid_size_type
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index
bool operator<(const ConstraintLine &) const
void serialize(Archive &ar, const unsigned int)
ConstraintLine & operator=(const ConstraintLineType &other)
ConstraintLine(const size_type &index=numbers::invalid_dof_index, const Entries &entries={}, const number &inhomogeneity=0.0)
std::vector< std::pair< size_type, number > > Entries
std::size_t memory_consumption() const
bool operator==(const ConstraintLine &) const
typename internal::ProductTypeImpl< typename std::decay< T >::type, typename std::decay< U >::type >::type type
static std::false_type check(...)
static std::true_type check(const BlockMatrixBase< T > *)
static std::true_type check(const BlockSparseMatrixEZ< T > *)
static std::true_type check(const BlockSparsityPatternBase< T > *)
static std::false_type check(...)
static void add(const typename VectorType::value_type value, const types::global_dof_index i, VectorType &V)
bool operator<(const SynchronousIterators< Iterators > &a, const SynchronousIterators< Iterators > &b)