16#ifndef dealii_affine_constraints_h
17#define dealii_affine_constraints_h
31#include <boost/range/iterator_range.hpp>
49template <
typename number>
51template <
typename number>
82 Distributing(
const Distributing &in);
85 operator=(
const Distributing &in);
90 return global_row < in.global_row;
111 template <
typename number>
120 insert_new_index(
const std::pair<size_type, number> &pair);
123 append_index(
const size_type index,
124 const std::pair<size_type, number> &pair);
127 get_size(
const size_type index)
const;
129 const std::pair<size_type, number> *
130 get_entry(
const size_type index)
const;
134 std::vector<std::pair<size_type, number>> data;
136 std::vector<size_type> individual_size;
169 template <
typename number>
170 class GlobalRowsFromLocal
176 GlobalRowsFromLocal();
179 reinit(
const size_type n_local_rows);
182 insert_index(
const size_type global_row,
183 const size_type local_row,
184 const number constraint_value);
189 print(std::ostream &os);
202 size(
const size_type counter_index)
const;
208 global_row(
const size_type counter_index)
const;
214 global_row(
const size_type counter_index);
222 local_row(
const size_type counter_index)
const;
228 local_row(
const size_type counter_index);
236 local_row(
const size_type counter_index,
237 const size_type index_in_constraint)
const;
244 constraint_value(
const size_type counter_index,
245 const size_type index_in_constraint)
const;
253 have_indirect_rows()
const;
260 insert_constraint(
const size_type constrained_local_dof);
268 n_constraints()
const;
275 n_inhomogeneities()
const;
284 set_ith_constraint_inhomogeneous(
const size_type i);
291 constraint_origin(size_type i)
const;
298 std::vector<Distributing> total_row_indices;
304 DataCache<number> data_cache;
338 template <
typename number>
351 ScratchData(
const ScratchData &)
363 std::vector<size_type> columns;
368 std::vector<number>
values;
373 std::vector<size_type> block_starts;
378 std::vector<size_type> vector_indices;
383 std::vector<number> vector_values;
388 GlobalRowsFromLocal<number> global_rows;
393 GlobalRowsFromLocal<number> global_columns;
400 namespace AffineConstraintsImplementation
402 template <
class VectorType>
404 set_zero_all(
const std::vector<types::global_dof_index> &cm,
409 set_zero_all(
const std::vector<types::global_dof_index> &cm,
414 set_zero_all(
const std::vector<types::global_dof_index> &cm,
502template <
typename number =
double>
594 template <
typename other_number>
723 const number weight);
733 const std::vector<std::pair<size_type, number>> &col_weight_pairs);
817 const bool allow_different_local_lines =
false);
942 const std::vector<std::pair<size_type, number>> *
1081 template <
class VectorType>
1091 template <
class VectorType>
1093 condense(
const VectorType &vec_ghosted, VectorType &output)
const;
1107 template <
class VectorType>
1115 template <
class BlockVectorType>
1125 template <
class VectorType>
1184 template <
class InVector,
class OutVector>
1187 const std::vector<size_type> &local_dof_indices,
1188 OutVector & global_vector)
const;
1237 template <
typename VectorType>
1240 const std::vector<size_type> &local_dof_indices,
1241 VectorType & global_vector,
1263 template <
typename VectorType>
1267 const std::vector<size_type> &local_dof_indices_row,
1268 const std::vector<size_type> &local_dof_indices_col,
1269 VectorType & global_vector,
1271 bool diagonal =
false)
const;
1276 template <
class VectorType>
1280 VectorType & global_vector)
const;
1314 template <
typename ForwardIteratorVec,
1315 typename ForwardIteratorInd,
1319 ForwardIteratorVec local_vector_end,
1320 ForwardIteratorInd local_indices_begin,
1321 VectorType & global_vector)
const;
1374 template <
typename MatrixType>
1377 const std::vector<size_type> &local_dof_indices,
1378 MatrixType & global_matrix)
const;
1407 template <
typename MatrixType>
1410 const std::vector<size_type> &row_indices,
1411 const std::vector<size_type> &col_indices,
1412 MatrixType & global_matrix)
const;
1430 template <
typename MatrixType>
1433 const std::vector<size_type> &row_indices,
1435 const std::vector<size_type> &column_indices,
1436 MatrixType & global_matrix)
const;
1458 template <
typename MatrixType,
typename VectorType>
1462 const std::vector<size_type> &local_dof_indices,
1463 MatrixType & global_matrix,
1464 VectorType & global_vector,
1465 bool use_inhomogeneities_for_rhs =
false)
const;
1520 template <
typename SparsityPatternType>
1523 const std::vector<size_type> &local_dof_indices,
1524 SparsityPatternType & sparsity_pattern,
1525 const bool keep_constrained_entries =
true,
1531 template <
typename SparsityPatternType>
1534 const std::vector<size_type> &row_indices,
1535 const std::vector<size_type> &col_indices,
1536 SparsityPatternType & sparsity_pattern,
1537 const bool keep_constrained_entries =
true,
1544 template <
typename SparsityPatternType>
1547 const std::vector<size_type> & row_indices,
1549 const std::vector<size_type> & col_indices,
1550 SparsityPatternType & sparsity_pattern,
1551 const bool keep_constrained_entries =
true,
1573 template <
typename ForwardIteratorVec,
1574 typename ForwardIteratorInd,
1578 ForwardIteratorInd local_indices_begin,
1579 ForwardIteratorVec local_vector_begin,
1580 ForwardIteratorVec local_vector_end)
const;
1603 template <
class VectorType>
1620 using Entries = std::vector<std::pair<size_type, number>>;
1653 template <
typename Constra
intLineType>
1654 ConstraintLine(
const ConstraintLineType &other);
1659 template <
typename Constra
intLineType>
1693 template <
class Archive>
1752 const IndexSet & locally_active_dofs,
1754 const bool verbose =
false)
const;
1763 const IndexSet &locally_relevant_dofs,
1785 <<
"The specified line " << arg1 <<
" does not exist.");
1796 <<
"The entry for the indices " << arg1 <<
" and " << arg2
1797 <<
" already exists, but the values " << arg3 <<
" (old) and "
1798 << arg4 <<
" (new) differ "
1799 <<
"by " << (arg4 - arg3) <<
'.');
1808 <<
"You tried to constrain DoF " << arg1 <<
" to DoF " << arg2
1809 <<
", but that one is also constrained. This is not allowed!");
1817 <<
"Degree of freedom " << arg1
1818 <<
" is constrained from both object in a merge operation.");
1826 <<
"In the given argument a degree of freedom is constrained "
1827 <<
"to another DoF with number " << arg1
1828 <<
", which however is constrained by this object. This is not"
1837 <<
"The index set given to this constraints object indicates "
1838 <<
"constraints for degree of freedom " << arg1
1839 <<
" should not be stored by this object, but a constraint "
1840 <<
"is being added.");
1850 <<
"The index set given to this constraints object indicates "
1851 <<
"constraints using degree of freedom " << arg2
1852 <<
" should not be stored by this object, but a constraint "
1853 <<
"for degree of freedom " << arg1 <<
" uses it.");
1863 <<
"While distributing the constraint for DoF " << arg1
1864 <<
", it turns out that one of the processors "
1865 <<
"who own the " << arg2 <<
" degrees of freedom that x_"
1866 << arg1 <<
" is constrained against does not know about "
1867 <<
"the constraint on x_" << arg1
1868 <<
". Did you not initialize the AffineConstraints container "
1869 <<
"with the appropriate locally_relevant set so "
1870 <<
"that every processor who owns a DoF that constrains "
1871 <<
"another DoF also knows about this constraint?");
1937 internal::AffineConstraints::ScratchData<number>>
1951 template <
typename MatrixType,
typename VectorType>
1955 const std::vector<size_type> &local_dof_indices,
1956 MatrixType & global_matrix,
1957 VectorType & global_vector,
1958 const bool use_inhomogeneities_for_rhs,
1959 const std::integral_constant<bool, false>)
const;
1965 template <
typename MatrixType,
typename VectorType>
1969 const std::vector<size_type> &local_dof_indices,
1970 MatrixType & global_matrix,
1971 VectorType & global_vector,
1972 const bool use_inhomogeneities_for_rhs,
1973 const std::integral_constant<bool, true>)
const;
1979 template <
typename SparsityPatternType>
1982 SparsityPatternType & sparsity_pattern,
1983 const bool keep_constrained_entries,
1985 const std::integral_constant<bool, false>)
const;
1991 template <
typename SparsityPatternType>
1994 SparsityPatternType & sparsity_pattern,
1995 const bool keep_constrained_entries,
1997 const std::integral_constant<bool, true>)
const;
2008 internal::AffineConstraints::GlobalRowsFromLocal<number>
2009 &global_rows)
const;
2020 std::vector<size_type> & active_dofs)
const;
2025 template <
typename MatrixScalar,
typename VectorScalar>
2029 const internal::AffineConstraints::GlobalRowsFromLocal<number> &global_rows,
2031 const std::vector<size_type> & local_dof_indices,
2037template <
typename number>
2050template <
typename number>
2060template <
typename number>
2084 lines.emplace_back();
2085 lines.back().index = line_n;
2086 lines.back().inhomogeneity = 0.;
2092template <
typename number>
2096 const number weight)
2099 Assert(constrained_dof_index != column,
2100 ExcMessage(
"Can't constrain a degree of freedom to itself"));
2105 ExcMessage(
"The current AffineConstraints does not contain the line "
2106 "for the current entry. Call AffineConstraints::add_line "
2107 "before calling this function."));
2120 for (
const auto &p : line_ptr->
entries)
2121 if (p.first == column)
2125 constrained_dof_index, column, p.second, weight));
2129 line_ptr->
entries.emplace_back(column, weight);
2134template <
typename number>
2143 ExcMessage(
"call add_line() before calling set_inhomogeneity()"));
2151template <
typename number>
2152template <
class VectorType>
2158 std::vector<size_type> constrained_lines(
lines.size());
2159 for (
unsigned int i = 0; i <
lines.size(); ++i)
2160 constrained_lines[i] =
lines[i].index;
2161 internal::AffineConstraintsImplementation::set_zero_all(constrained_lines,
2165template <
typename number>
2169 return lines.size();
2172template <
typename number>
2176 return std::count_if(
lines.begin(),
2179 return (line.entries.size() == 1) &&
2180 (line.entries[0].second == number(1.));
2184template <
typename number>
2188 return std::count_if(
lines.begin(),
2191 return (line.inhomogeneity != number(0.));
2195template <
typename number>
2204template <
typename number>
2222template <
typename number>
2223inline const std::vector<std::pair<types::global_dof_index, number>> *
2236template <
typename number>
2250template <
typename number>
2263template <
typename number>
2270template <
typename number>
2277template <
typename number>
2278template <
class VectorType>
2283 VectorType & global_vector)
const
2288 global_vector(index) += value;
2294 global_vector(position.
entries[j].first) +=
2295 value * position.
entries[j].second;
2299template <
typename number>
2300template <
typename ForwardIteratorVec,
2301 typename ForwardIteratorInd,
2305 ForwardIteratorVec local_vector_begin,
2306 ForwardIteratorVec local_vector_end,
2307 ForwardIteratorInd local_indices_begin,
2308 VectorType & global_vector)
const
2311 for (; local_vector_begin != local_vector_end;
2312 ++local_vector_begin, ++local_indices_begin)
2316 *local_indices_begin,
2324 (*local_vector_begin) * position.
entries[j].second,
2331template <
typename number>
2332template <
class InVector,
class OutVector>
2335 const InVector & local_vector,
2336 const std::vector<size_type> &local_dof_indices,
2337 OutVector & global_vector)
const
2339 Assert(local_vector.size() == local_dof_indices.size(),
2343 local_dof_indices.begin(),
2347template <
typename number>
2348template <
typename ForwardIteratorVec,
2349 typename ForwardIteratorInd,
2353 const VectorType & global_vector,
2354 ForwardIteratorInd local_indices_begin,
2355 ForwardIteratorVec local_vector_begin,
2356 ForwardIteratorVec local_vector_end)
const
2359 for (; local_vector_begin != local_vector_end;
2360 ++local_vector_begin, ++local_indices_begin)
2363 *local_vector_begin = global_vector(*local_indices_begin);
2368 typename VectorType::value_type value = position.
inhomogeneity;
2370 value += (global_vector(position.
entries[j].first) *
2372 *local_vector_begin = value;
2379template <
typename MatrixType>
2381template <
typename SparsityPatternType>
2383template <
typename number>
2407 template <
typename MatrixType>
2408 struct IsBlockMatrix
2416 template <
typename T>
2417 static std::true_type
2424 template <
typename T>
2425 static std::true_type
2432 static std::false_type
2442 static const bool value =
2443 std::is_same<decltype(check(std::declval<MatrixType *>())),
2444 std::true_type>::value;
2448 template <
typename MatrixType>
2449 const bool IsBlockMatrix<MatrixType>::value;
2460 template <
typename MatrixType>
2461 struct IsBlockSparsityPattern
2468 template <
typename T>
2469 static std::true_type
2476 static std::false_type
2485 static const bool value =
2486 std::is_same<decltype(check(std::declval<MatrixType *>())),
2487 std::true_type>::value;
2491 template <
typename MatrixType>
2492 const bool IsBlockSparsityPattern<MatrixType>::value;
2500template <
typename number>
2501template <
typename other_number>
2507 lines.insert(lines.begin(), other.
lines.begin(), other.
lines.end());
2515template <
typename number>
2516template <
typename MatrixType>
2520 const std::vector<size_type> &local_dof_indices,
2521 MatrixType & global_matrix)
const
2526 distribute_local_to_global(
2533 std::integral_constant<
2535 internal::AffineConstraints::IsBlockMatrix<MatrixType>::value>());
2540template <
typename number>
2541template <
typename MatrixType,
typename VectorType>
2546 const std::vector<size_type> &local_dof_indices,
2547 MatrixType & global_matrix,
2548 VectorType & global_vector,
2549 bool use_inhomogeneities_for_rhs)
const
2553 distribute_local_to_global(
2559 use_inhomogeneities_for_rhs,
2560 std::integral_constant<
2562 internal::AffineConstraints::IsBlockMatrix<MatrixType>::value>());
2567template <
typename number>
2568template <
typename SparsityPatternType>
2571 const std::vector<size_type> &local_dof_indices,
2572 SparsityPatternType & sparsity_pattern,
2573 const bool keep_constrained_entries,
2578 add_entries_local_to_global(
2581 keep_constrained_entries,
2583 std::integral_constant<
bool,
2584 internal::AffineConstraints::IsBlockSparsityPattern<
2585 SparsityPatternType>::value>());
2590template <
typename number>
2594 const number &inhomogeneity)
2597 , inhomogeneity(inhomogeneity)
2602template <
typename number>
2603template <
typename Constra
intLineType>
2605 const ConstraintLineType &other)
2607 this->index = other.index;
2610 entries.insert(entries.begin(), other.entries.begin(), other.entries.end());
2612 this->inhomogeneity = other.inhomogeneity;
2617template <
typename number>
2618template <
typename Constra
intLineType>
2621 const ConstraintLineType &other)
2623 this->index = other.index;
2626 entries.insert(entries.begin(), other.entries.begin(), other.entries.end());
2628 this->inhomogeneity = other.inhomogeneity;
bool is_closed(const MPI_Comm &comm) const
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)
size_type n_inhomogeneities() const
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())
void make_consistent_in_parallel(const IndexSet &locally_owned_dofs, const IndexSet &locally_relevant_dofs, const MPI_Comm mpi_communicator)
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)
size_type n_identities() const
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 add_entries_local_to_global(const std::vector< size_type > &row_indices, const AffineConstraints< number > &col_constraints, 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)
types::global_dof_index size_type
bool check(const ConstraintKinds kind_in, const unsigned int dim)
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 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)