17 #ifndef dealii_constraint_matrix_h 18 #define dealii_constraint_matrix_h 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> 26 #include <deal.II/lac/vector.h> 27 #include <deal.II/lac/vector_element_access.h> 29 #include <boost/range/iterator_range.hpp> 36 DEAL_II_NAMESPACE_OPEN
38 template <
int dim,
class T>
class Table;
49 class GlobalRowsFromLocal;
346 const std::vector<std::pair<size_type,double> > &col_val_pairs);
405 const bool allow_different_local_lines =
false);
507 const std::vector<std::pair<size_type,double> > *
536 void print (std::ostream &out)
const;
564 void resolve_indices(std::vector<types::global_dof_index> &indices)
const;
619 template <
typename number>
625 template <
typename number>
639 template <
class VectorType>
640 void condense (VectorType &vec)
const;
648 template <
class VectorType>
649 void condense (
const VectorType &vec_ghosted,
650 VectorType &output)
const;
664 template <
typename number,
class VectorType>
666 VectorType &vector)
const;
672 template <
typename number,
class BlockVectorType>
674 BlockVectorType &vector)
const;
682 template <
class VectorType>
683 void set_zero (VectorType &vec)
const;
736 template <
class InVector,
class OutVector>
739 const std::vector<size_type> &local_dof_indices,
740 OutVector &global_vector)
const;
785 template <
typename VectorType,
typename LocalType>
788 const std::vector<size_type> &local_dof_indices,
789 VectorType &global_vector,
809 template <
typename VectorType,
typename LocalType>
812 const std::vector<size_type> &local_dof_indices_row,
813 const std::vector<size_type> &local_dof_indices_col,
814 VectorType &global_vector,
816 bool diagonal =
false)
const;
821 template <
class VectorType>
825 VectorType &global_vector)
const;
855 template <
typename ForwardIteratorVec,
typename ForwardIteratorInd,
859 ForwardIteratorVec local_vector_end,
860 ForwardIteratorInd local_indices_begin,
861 VectorType &global_vector)
const;
910 template <
typename MatrixType>
913 const std::vector<size_type> &local_dof_indices,
914 MatrixType &global_matrix)
const;
943 template <
typename MatrixType>
946 const std::vector<size_type> &row_indices,
947 const std::vector<size_type> &col_indices,
948 MatrixType &global_matrix)
const;
966 template <
typename MatrixType>
969 const std::vector<size_type> &row_indices,
971 const std::vector<size_type> &column_indices,
972 MatrixType &global_matrix)
const;
990 template <
typename MatrixType,
typename VectorType>
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;
1052 template <
typename SparsityPatternType>
1055 SparsityPatternType &sparsity_pattern,
1056 const bool keep_constrained_entries =
true,
1062 template <
typename SparsityPatternType>
1065 const std::vector<size_type> &col_indices,
1066 SparsityPatternType &sparsity_pattern,
1067 const bool keep_constrained_entries =
true,
1089 template <
typename ForwardIteratorVec,
typename ForwardIteratorInd,
1093 ForwardIteratorInd local_indices_begin,
1094 ForwardIteratorVec local_vector_begin,
1095 ForwardIteratorVec local_vector_end)
const;
1118 template <
class VectorType>
1136 typedef std::vector<std::pair<size_type,double> >
Entries;
1184 template <
class Archive>
1244 const IndexSet &locally_active_dofs,
1245 const MPI_Comm mpi_communicator,
1246 const bool verbose=
false)
const;
1268 <<
"The specified line " << arg1
1269 <<
" does not exist.");
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) <<
".");
1288 <<
"You tried to constrain DoF " << arg1
1289 <<
" to DoF " << arg2
1290 <<
", but that one is also constrained. This is not allowed!");
1298 <<
"Degree of freedom " << arg1
1299 <<
" is constrained from both object in a merge operation.");
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" 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.");
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.");
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?");
1437 template <
typename MatrixType,
typename VectorType>
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;
1451 template <
typename MatrixType,
typename VectorType>
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;
1465 template <
typename SparsityPatternType>
1468 SparsityPatternType &sparsity_pattern,
1469 const bool keep_constrained_entries,
1471 std::integral_constant<bool, false>)
const;
1477 template <
typename SparsityPatternType>
1480 SparsityPatternType &sparsity_pattern,
1481 const bool keep_constrained_entries,
1483 std::integral_constant<bool, true>)
const;
1494 internals::GlobalRowsFromLocal &global_rows)
const;
1505 std::vector<size_type> &active_dofs)
const;
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,
1527 local_lines (local_constraints),
1542 lines (constraint_matrix.lines),
1543 lines_cache (constraint_matrix.lines_cache),
1544 local_lines (constraint_matrix.local_lines),
1545 sorted (constraint_matrix.sorted)
1574 lines.emplace_back ();
1575 lines.back().index = line;
1576 lines.back().inhomogeneity = 0.;
1590 ExcMessage (
"Can't constrain a degree of freedom to itself"));
1595 ExcMessage(
"The current ConstraintMatrix does not contain the line " 1596 "for the current entry. Call ConstraintMatrix::add_line " 1597 "before calling this function."));
1610 for (ConstraintLine::Entries::const_iterator
1611 p=line_ptr->entries.begin();
1612 p != line_ptr->entries.end(); ++p)
1613 if (p->first == column)
1615 Assert (std::fabs(p->second - value) < 1.e-14,
1620 line_ptr->entries.emplace_back (column, value);
1633 ExcMessage(
"call add_line() before calling set_inhomogeneity()"));
1636 line_ptr->inhomogeneity = value;
1645 return lines.size();
1682 const std::vector<std::pair<types::global_dof_index,double> > *
1745 template <
class VectorType>
1750 VectorType &global_vector)
const 1755 global_vector(index) += value;
1761 global_vector(position.
entries[j].first)
1762 += value * position.
entries[j].second;
1767 template <
typename ForwardIteratorVec,
typename ForwardIteratorInd,
1771 ForwardIteratorVec local_vector_begin,
1772 ForwardIteratorVec local_vector_end,
1773 ForwardIteratorInd local_indices_begin,
1774 VectorType &global_vector)
const 1777 for ( ; local_vector_begin != local_vector_end;
1778 ++local_vector_begin, ++local_indices_begin)
1781 internal::ElementAccess<VectorType>::add(*local_vector_begin,
1782 *local_indices_begin, global_vector);
1788 internal::ElementAccess<VectorType>::add(*local_vector_begin *
1796 template <
class InVector,
class OutVector>
1800 const InVector &local_vector,
1801 const std::vector<size_type> &local_dof_indices,
1802 OutVector &global_vector)
const 1804 Assert (local_vector.size() == local_dof_indices.size(),
1807 local_dof_indices.begin(), global_vector);
1812 template <
typename ForwardIteratorVec,
typename ForwardIteratorInd,
1816 ForwardIteratorInd local_indices_begin,
1817 ForwardIteratorVec local_vector_begin,
1818 ForwardIteratorVec local_vector_end)
const 1821 for ( ; local_vector_begin != local_vector_end;
1822 ++local_vector_begin, ++local_indices_begin)
1825 *local_vector_begin = global_vector(*local_indices_begin);
1830 typename VectorType::value_type value = position.
inhomogeneity;
1832 value += (global_vector(position.
entries[j].first) *
1834 *local_vector_begin = value;
1862 template <
typename MatrixType>
1880 template <
typename T>
1887 template <
typename T>
1894 template <
typename T>
1910 ((MatrixType *)
nullptr))
1917 template <
typename MatrixType>
1921 template <
typename MatrixType>
1926 const std::vector<size_type> &local_dof_indices,
1927 MatrixType &global_matrix)
const 1931 Vector<typename MatrixType::value_type> dummy(0);
1933 global_matrix, dummy,
false,
1940 template <
typename MatrixType,
typename VectorType>
1946 const std::vector<size_type> &local_dof_indices,
1947 MatrixType &global_matrix,
1948 VectorType &global_vector,
1949 bool use_inhomogeneities_for_rhs)
const 1954 global_matrix, global_vector, use_inhomogeneities_for_rhs,
1961 template <
typename SparsityPatternType>
1966 SparsityPatternType &sparsity_pattern,
1967 const bool keep_constrained_entries,
1973 keep_constrained_entries, dof_mask,
1978 DEAL_II_NAMESPACE_CLOSE
const types::global_dof_index invalid_size_type
bool operator<(const ConstraintLine &) const
#define DeclException2(Exception2, type1, type2, outsequence)
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)
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)
void add_line(const size_type line)
bool operator==(const ConstraintLine &) const
static ::ExceptionBase & ExcMessage(std::string arg1)
#define DeclException1(Exception1, type1, outsequence)
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
ConstraintMatrix & operator=(const ConstraintMatrix &)=delete
bool is_inhomogeneously_constrained(const size_type index) const
#define Assert(cond, exc)
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
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)
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 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)
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)
bool is_element(const size_type index) const
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)