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> 33 DEAL_II_NAMESPACE_OPEN
35 template<
int dim,
class T>
class Table;
46 class GlobalRowsFromLocal;
314 const std::vector<std::pair<size_type,double> > &col_val_pairs);
373 const bool allow_different_local_lines =
false);
475 const std::vector<std::pair<size_type,double> > *
504 void print (std::ostream &out)
const;
532 void resolve_indices(std::vector<types::global_dof_index> &indices)
const;
587 template<
typename number>
593 template <
typename number>
607 template <
class VectorType>
608 void condense (VectorType &vec)
const;
616 template <
class VectorType>
617 void condense (
const VectorType &vec_ghosted,
618 VectorType &output)
const;
632 template<
typename number,
class VectorType>
634 VectorType &vector)
const;
640 template <
typename number,
class BlockVectorType>
642 BlockVectorType &vector)
const;
650 template <
class VectorType>
651 void set_zero (VectorType &vec)
const;
704 template <
class InVector,
class OutVector>
707 const std::vector<size_type> &local_dof_indices,
708 OutVector &global_vector)
const;
753 template <
typename VectorType,
typename LocalType>
756 const std::vector<size_type> &local_dof_indices,
757 VectorType &global_vector,
777 template <
typename VectorType,
typename LocalType>
780 const std::vector<size_type> &local_dof_indices_row,
781 const std::vector<size_type> &local_dof_indices_col,
782 VectorType &global_vector,
784 bool diagonal =
false)
const;
789 template <
class VectorType>
793 VectorType &global_vector)
const;
823 template <
typename ForwardIteratorVec,
typename ForwardIteratorInd,
827 ForwardIteratorVec local_vector_end,
828 ForwardIteratorInd local_indices_begin,
829 VectorType &global_vector)
const;
878 template <
typename MatrixType>
881 const std::vector<size_type> &local_dof_indices,
882 MatrixType &global_matrix)
const;
911 template <
typename MatrixType>
914 const std::vector<size_type> &row_indices,
915 const std::vector<size_type> &col_indices,
916 MatrixType &global_matrix)
const;
934 template <
typename MatrixType,
typename VectorType>
938 const std::vector<size_type> &local_dof_indices,
939 MatrixType &global_matrix,
940 VectorType &global_vector,
941 bool use_inhomogeneities_for_rhs =
false)
const;
996 template <
typename SparsityPatternType>
999 SparsityPatternType &sparsity_pattern,
1000 const bool keep_constrained_entries =
true,
1006 template <
typename SparsityPatternType>
1009 const std::vector<size_type> &col_indices,
1010 SparsityPatternType &sparsity_pattern,
1011 const bool keep_constrained_entries =
true,
1033 template <
typename ForwardIteratorVec,
typename ForwardIteratorInd,
1037 ForwardIteratorInd local_indices_begin,
1038 ForwardIteratorVec local_vector_begin,
1039 ForwardIteratorVec local_vector_end)
const;
1062 template <
class VectorType>
1088 <<
"The specified line " << arg1
1089 <<
" does not exist.");
1097 <<
"The entry for the indices " << arg1 <<
" and " 1098 << arg2 <<
" already exists, but the values " 1099 << arg3 <<
" (old) and " << arg4 <<
" (new) differ " 1100 <<
"by " << (arg4-arg3) <<
".");
1108 <<
"You tried to constrain DoF " << arg1
1109 <<
" to DoF " << arg2
1110 <<
", but that one is also constrained. This is not allowed!");
1118 <<
"Degree of freedom " << arg1
1119 <<
" is constrained from both object in a merge operation.");
1127 <<
"In the given argument a degree of freedom is constrained " 1128 <<
"to another DoF with number " << arg1
1129 <<
", which however is constrained by this object. This is not" 1138 <<
"The index set given to this constraint matrix indicates " 1139 <<
"constraints for degree of freedom " << arg1
1140 <<
" should not be stored by this object, but a constraint " 1141 <<
"is being added.");
1151 <<
"The index set given to this constraint matrix indicates " 1152 <<
"constraints using degree of freedom " << arg2
1153 <<
" should not be stored by this object, but a constraint " 1154 <<
"for degree of freedom " << arg1 <<
" uses it.");
1163 <<
"While distributing the constraint for DoF " 1164 << arg1 <<
", it turns out that one of the processors " 1165 <<
"who own the " << arg2
1166 <<
" degrees of freedom that x_" << arg1
1167 <<
" is constrained against does not know about " 1168 <<
"the constraint on x_" << arg1
1169 <<
". Did you not initialize the ConstraintMatrix " 1170 <<
"with the appropriate locally_relevant set so " 1171 <<
"that every processor who owns a DoF that constrains " 1172 <<
"another DoF also knows about this constraint?");
1185 typedef std::vector<std::pair<size_type,double> >
Entries;
1311 template <
typename MatrixType,
typename VectorType>
1315 const std::vector<size_type> &local_dof_indices,
1316 MatrixType &global_matrix,
1317 VectorType &global_vector,
1318 bool use_inhomogeneities_for_rhs,
1325 template <
typename MatrixType,
typename VectorType>
1329 const std::vector<size_type> &local_dof_indices,
1330 MatrixType &global_matrix,
1331 VectorType &global_vector,
1332 bool use_inhomogeneities_for_rhs,
1339 template <
typename SparsityPatternType>
1342 SparsityPatternType &sparsity_pattern,
1343 const bool keep_constrained_entries,
1351 template <
typename SparsityPatternType>
1354 SparsityPatternType &sparsity_pattern,
1355 const bool keep_constrained_entries,
1368 internals::GlobalRowsFromLocal &global_rows)
const;
1379 std::vector<size_type> &active_dofs)
const;
1384 template <
typename LocalType>
1387 const internals::GlobalRowsFromLocal &global_rows,
1389 const std::vector<size_type> &local_dof_indices,
1407 local_lines (local_constraints),
1422 lines (constraint_matrix.lines),
1423 lines_cache (constraint_matrix.lines_cache),
1424 local_lines (constraint_matrix.local_lines),
1425 sorted (constraint_matrix.sorted)
1455 lines.back().line = line;
1456 lines.back().inhomogeneity = 0.;
1470 ExcMessage (
"Can't constrain a degree of freedom to itself"));
1483 for (ConstraintLine::Entries::const_iterator
1484 p=line_ptr->entries.begin();
1485 p != line_ptr->entries.end(); ++p)
1486 if (p->first == column)
1488 Assert (std::fabs(p->second - value) < 1.e-14,
1493 line_ptr->entries.push_back (std::make_pair(column,value));
1506 ExcMessage(
"call add_line() before calling set_inhomogeneity()"));
1509 line_ptr->inhomogeneity = value;
1518 return lines.size();
1555 const std::vector<std::pair<types::global_dof_index,double> > *
1618 template <
class VectorType>
1623 VectorType &global_vector)
const 1628 global_vector(index) += value;
1634 global_vector(position.
entries[j].first)
1635 += value * position.
entries[j].second;
1640 template <
typename ForwardIteratorVec,
typename ForwardIteratorInd,
1644 ForwardIteratorVec local_vector_begin,
1645 ForwardIteratorVec local_vector_end,
1646 ForwardIteratorInd local_indices_begin,
1647 VectorType &global_vector)
const 1650 for ( ; local_vector_begin != local_vector_end;
1651 ++local_vector_begin, ++local_indices_begin)
1654 global_vector(*local_indices_begin) += *local_vector_begin;
1660 global_vector(position.
entries[j].first)
1661 += *local_vector_begin * position.
entries[j].second;
1667 template <
class InVector,
class OutVector>
1671 const InVector &local_vector,
1672 const std::vector<size_type> &local_dof_indices,
1673 OutVector &global_vector)
const 1675 Assert (local_vector.size() == local_dof_indices.size(),
1678 local_dof_indices.begin(), global_vector);
1683 template <
typename ForwardIteratorVec,
typename ForwardIteratorInd,
1687 ForwardIteratorInd local_indices_begin,
1688 ForwardIteratorVec local_vector_begin,
1689 ForwardIteratorVec local_vector_end)
const 1692 for ( ; local_vector_begin != local_vector_end;
1693 ++local_vector_begin, ++local_indices_begin)
1696 *local_vector_begin = global_vector(*local_indices_begin);
1701 typename VectorType::value_type value = position.
inhomogeneity;
1703 value += (global_vector(position.
entries[j].first) *
1705 *local_vector_begin = value;
1733 template <
typename MatrixType>
1751 template <
typename T>
1758 template <
typename T>
1765 template <
typename T>
1788 template <
typename MatrixType>
1792 template <
typename MatrixType>
1797 const std::vector<size_type> &local_dof_indices,
1798 MatrixType &global_matrix)
const 1804 global_matrix, dummy,
false,
1811 template <
typename MatrixType,
typename VectorType>
1817 const std::vector<size_type> &local_dof_indices,
1818 MatrixType &global_matrix,
1819 VectorType &global_vector,
1820 bool use_inhomogeneities_for_rhs)
const 1825 global_matrix, global_vector, use_inhomogeneities_for_rhs,
1832 template <
typename SparsityPatternType>
1837 SparsityPatternType &sparsity_pattern,
1838 const bool keep_constrained_entries,
1844 keep_constrained_entries, dof_mask,
1849 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)
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)
ConstraintMatrix & operator=(const ConstraintMatrix &other)
std::size_t memory_consumption() const
void set_zero(VectorType &vec) const
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
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)
void resolve_indices(std::vector< types::global_dof_index > &indices) const
bool has_inhomogeneities() 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()
LocalType resolve_vector_entry(const size_type i, const internals::GlobalRowsFromLocal &global_rows, const Vector< LocalType > &local_vector, const std::vector< size_type > &local_dof_indices, const FullMatrix< LocalType > &local_matrix) const
static ::ExceptionBase & ExcDoFConstrainedToConstrainedDoF(int arg1, int arg2)
void distribute(VectorType &vec) const
void write_dot(std::ostream &) const
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)
std::vector< ConstraintLine > lines
bool is_constrained(const size_type index) const
static ::ExceptionBase & ExcMatrixIsClosed()
static ::ExceptionBase & ExcInternalError()