Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | Friends | List of all members
AffineConstraints< number > Class Template Reference

#include <deal.II/lac/affine_constraints.h>

Inheritance diagram for AffineConstraints< number >:
Inheritance graph
[legend]

Classes

struct  ConstraintLine
 

Public Types

enum  MergeConflictBehavior { no_conflicts_allowed , left_object_wins , right_object_wins }
 
using size_type = types::global_dof_index
 
using const_iterator = typename std::vector< ConstraintLine >::const_iterator
 
using LineRange = boost::iterator_range< const_iterator >
 

Public Member Functions

 AffineConstraints ()
 
 AffineConstraints (const IndexSet &locally_stored_constraints)
 
 AffineConstraints (const IndexSet &locally_owned_dofs, const IndexSet &locally_stored_constraints)
 
 AffineConstraints (const AffineConstraints &affine_constraints)
 
 AffineConstraints (AffineConstraints &&affine_constraints) noexcept=default
 
AffineConstraintsoperator= (const AffineConstraints &)=delete
 
AffineConstraintsoperator= (AffineConstraints &&affine_constraints) noexcept=default
 
template<typename other_number >
void copy_from (const AffineConstraints< other_number > &other)
 
void reinit ()
 
void reinit (const IndexSet &locally_stored_constraints)
 
void reinit (const IndexSet &locally_owned_dofs, const IndexSet &locally_stored_constraints)
 
bool can_store_line (const size_type line_n) const
 
const IndexSetget_locally_owned_indices () const
 
const IndexSetget_local_lines () const
 
void add_selected_constraints (const AffineConstraints &constraints_in, const IndexSet &filter)
 
LineRange get_lines () 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 make_consistent_in_parallel (const IndexSet &locally_owned_dofs, const IndexSet &locally_relevant_dofs, const MPI_Comm mpi_communicator)
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Adding constraints
void add_constraint (const size_type constrained_dof, const ArrayView< const std::pair< size_type, number > > &dependencies, const number inhomogeneity=0)
 
void constrain_dof_to_zero (const size_type constrained_dof)
 
void add_line (const size_type line_n)
 
void add_lines (const std::vector< bool > &lines)
 
void add_lines (const std::set< size_type > &lines)
 
void add_lines (const IndexSet &lines)
 
void add_entry (const size_type constrained_dof_index, const size_type column, const number weight)
 
void add_entries (const size_type constrained_dof_index, const std::vector< std::pair< size_type, number > > &col_weight_pairs)
 
void set_inhomogeneity (const size_type constrained_dof_index, const number value)
 
void close ()
 
bool is_closed () const
 
bool is_closed (const MPI_Comm comm) const
 
template<typename other_number >
void merge (const AffineConstraints< other_number > &other_constraints, const MergeConflictBehavior merge_conflict_behavior=no_conflicts_allowed, const bool allow_different_local_lines=false)
 
void shift (const size_type offset)
 
AffineConstraints get_view (const IndexSet &mask) const
 
void clear ()
 
Querying constraints
size_type n_constraints () const
 
size_type n_identities () const
 
size_type n_inhomogeneities () const
 
bool is_constrained (const size_type line_n) const
 
bool is_identity_constrained (const size_type line_n) const
 
bool are_identity_constrained (const size_type line_n_1, const size_type line_n_2) const
 
size_type max_constraint_indirections () const
 
bool is_inhomogeneously_constrained (const size_type index) const
 
bool has_inhomogeneities () const
 
const std::vector< std::pair< size_type, number > > * get_constraint_entries (const size_type line_n) const
 
number get_inhomogeneity (const size_type line_n) const
 
void print (std::ostream &out) const
 
void write_dot (std::ostream &) const
 
std::size_t memory_consumption () const
 
void resolve_indices (std::vector< types::global_dof_index > &indices) const
 
Eliminating constraints from linear systems after their creation
void condense (SparsityPattern &sparsity) const
 
void condense (BlockSparsityPattern &sparsity) const
 
void condense (DynamicSparsityPattern &sparsity) const
 
void condense (BlockDynamicSparsityPattern &sparsity) const
 
void condense (SparseMatrix< number > &matrix) const
 
void condense (BlockSparseMatrix< number > &matrix) const
 
template<typename VectorType >
void condense (VectorType &vec) const
 
template<typename VectorType >
void condense (const VectorType &vec_ghosted, VectorType &output) const
 
template<typename VectorType >
void condense (SparseMatrix< number > &matrix, VectorType &vector) const
 
template<typename BlockVectorType >
void condense (BlockSparseMatrix< number > &matrix, BlockVectorType &vector) const
 
template<typename VectorType >
void set_zero (VectorType &vec) const
 
Eliminating constraints from linear systems during their creation
template<class InVector , class OutVector >
void distribute_local_to_global (const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
 
template<typename VectorType >
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
 
template<typename VectorType >
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
 
template<typename VectorType >
void distribute_local_to_global (const size_type index, const number value, VectorType &global_vector) const
 
template<typename ForwardIteratorVec , typename ForwardIteratorInd , typename VectorType >
void distribute_local_to_global (ForwardIteratorVec local_vector_begin, ForwardIteratorVec local_vector_end, ForwardIteratorInd local_indices_begin, VectorType &global_vector) const
 
template<typename MatrixType >
void distribute_local_to_global (const FullMatrix< number > &local_matrix, const std::vector< size_type > &local_dof_indices, MatrixType &global_matrix) const
 
template<typename MatrixType >
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
 
template<typename MatrixType >
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
 
template<typename MatrixType , typename VectorType >
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 add_entries_local_to_global (const std::vector< size_type > &local_dof_indices, SparsityPatternBase &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 std::vector< size_type > &col_indices, SparsityPatternBase &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, SparsityPatternBase &sparsity_pattern, const bool keep_constrained_entries=true, const Table< 2, bool > &dof_mask=Table< 2, bool >()) const
 
template<typename ForwardIteratorVec , typename ForwardIteratorInd , typename VectorType >
void get_dof_values (const VectorType &global_vector, ForwardIteratorInd local_indices_begin, ForwardIteratorVec local_vector_begin, ForwardIteratorVec local_vector_end) const
 
Dealing with constraints after solving a linear system
template<typename VectorType >
void distribute (VectorType &vec) const
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Static Public Member Functions

static ::ExceptionBaseExcMatrixIsClosed ()
 
static ::ExceptionBaseExcMatrixNotClosed ()
 
static ::ExceptionBaseExcLineInexistant (size_type arg1)
 
static ::ExceptionBaseExcEntryAlreadyExists (size_type arg1, size_type arg2, number arg3, number arg4)
 
static ::ExceptionBaseExcDoFConstrainedToConstrainedDoF (int arg1, int arg2)
 
static ::ExceptionBaseExcDoFIsConstrainedFromBothObjects (size_type arg1)
 
static ::ExceptionBaseExcDoFIsConstrainedToConstrainedDoF (size_type arg1)
 
static ::ExceptionBaseExcRowNotStoredHere (size_type arg1)
 
static ::ExceptionBaseExcColumnNotStoredHere (size_type arg1, size_type arg2)
 
static ::ExceptionBaseExcIncorrectConstraint (int arg1, int arg2)
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

size_type calculate_line_index (const size_type line_n) const
 
template<typename MatrixType , typename VectorType >
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::bool_constant< false >) const
 
template<typename MatrixType , typename VectorType >
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::bool_constant< true >) const
 
void make_sorted_row_list (const std::vector< size_type > &local_dof_indices, internal::AffineConstraints::GlobalRowsFromLocal< number > &global_rows) const
 
void make_sorted_row_list (const std::vector< size_type > &local_dof_indices, std::vector< size_type > &active_dofs) const
 
template<typename MatrixScalar , typename VectorScalar >
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 check_no_subscribers () const noexcept
 

Private Attributes

std::vector< ConstraintLinelines
 
std::vector< size_typelines_cache
 
IndexSet locally_owned_dofs
 
IndexSet local_lines
 
IndexSet needed_elements_for_distribute
 
bool sorted
 
Threads::ThreadLocalStorage< internal::AffineConstraints::ScratchData< number > > scratch_data
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Friends

template<typename >
class AffineConstraints
 

Detailed Description

template<typename number = double>
class AffineConstraints< number >

This class implements dealing with linear (possibly inhomogeneous) constraints on degrees of freedom. The concept and origin of such constraints is extensively described in the Constraints on degrees of freedom module. The class is meant to deal with a limited number of constraints relative to the total number of degrees of freedom, for example a few per cent up to maybe 30 per cent; and with a linear combination of M other degrees of freedom where M is also relatively small (no larger than at most around the average number of entries per row of a linear system). It is not meant to describe full rank linear systems.

The algorithms used in the implementation of this class are described in some detail in the hp-paper. There is also a significant amount of documentation on how to use this class in the Constraints on degrees of freedom module.

Description of constraints

Each "line" in objects of this class corresponds to one constrained degree of freedom, with the number of the line being i, entered by using add_line() or add_lines(). The entries in this line are pairs of the form (j,aij), which are added by add_entry() or add_entries(). The organization is essentially a SparsityPattern, but with only a few lines containing nonzero elements, and therefore no data wasted on the others. For each line, which has been added by the mechanism above, an elimination of the constrained degree of freedom of the form

\[ x_i = \sum_j a_{ij} x_j + b_i \]

is performed, where bi is optional and set by set_inhomogeneity(). Thus, if a constraint is formulated for instance as a zero mean value of several degrees of freedom, one of the degrees has to be chosen to be eliminated.

Note that the constraints are linear in the xi, and that there might be a constant (non-homogeneous) term in the constraint. This is exactly the form we need for hanging node constraints, where we need to constrain one degree of freedom in terms of others. There are other conditions of this form possible, for example for implementing mean value conditions as is done in the step-11 tutorial program. The name of the class stems from the fact that these constraints can be represented in matrix form as X x = b, and this object then describes the matrix X and the vector b. The most frequent way to create/fill objects of this type is using the DoFTools::make_hanging_node_constraints() function. The use of these objects is first explained in step-6.

Objects of the present type are organized in lines (rows), but only those lines are stored where constraints are present. New constraints are added by adding new lines using the add_line() function, and then populating it using the add_entry() function to a given line, or add_entries() to add more than one entry at a time. The right hand side element, if nonzero, can be set using the set_inhomogeneity() function. After all constraints have been added, you need to call close(), which compresses the storage format and sorts the entries.

Note
Many of the algorithms this class implements are discussed in the hp_paper. The algorithms are also related to those shown in [175] , with the difference that the algorithms shown there completely eliminate constrained degrees of freedom, whereas we usually keep them as part of the linear system (but zero out the row and column to decouple this degree of freedom from the rest of the linear system).

Definition at line 506 of file affine_constraints.h.

Member Typedef Documentation

◆ size_type

template<typename number = double>
using AffineConstraints< number >::size_type = types::global_dof_index

Declare the type for container size.

Definition at line 512 of file affine_constraints.h.

◆ const_iterator

template<typename number = double>
using AffineConstraints< number >::const_iterator = typename std::vector<ConstraintLine>::const_iterator

Alias for the iterator type that is used in the LineRange container.

Definition at line 1940 of file affine_constraints.h.

◆ LineRange

template<typename number = double>
using AffineConstraints< number >::LineRange = boost::iterator_range<const_iterator>

Alias for the return type used by get_lines().

Definition at line 1945 of file affine_constraints.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Member Enumeration Documentation

◆ MergeConflictBehavior

template<typename number = double>
enum AffineConstraints::MergeConflictBehavior

An enum that describes what should happen if the two AffineConstraints objects involved in a call to the merge() function happen to have constraints on the same degrees of freedom.

Enumerator
no_conflicts_allowed 

Throw an exception if the two objects concerned have conflicting constraints on the same degree of freedom.

left_object_wins 

In an operation cm1.merge(cm2), if cm1 and cm2 have constraints on the same degree of freedom, take the one from cm1.

right_object_wins 

In an operation cm1.merge(cm2), if cm1 and cm2 have constraints on the same degree of freedom, take the one from cm2.

Definition at line 519 of file affine_constraints.h.

Constructor & Destructor Documentation

◆ AffineConstraints() [1/5]

template<typename number >
AffineConstraints< number >::AffineConstraints ( )
inline

Default constructor.

This constructor sets up an object that stores all constraints eventually given to it. (This is in contrast to the following constructor that restricts which degrees of freedom we should care about – for example, in a parallel context one would typically only store constraints for the locally relevant DoFs, see GlossLocallyRelevantDof .) Consequently, this constructor creates internal data structures for all possible indices, leading to memory consumption on every processor that is proportional to the overall size of the problem, not just proportional to the size of the portion of the overall problem that is handled by the current processor. Calling this constructor in a parallel program is a common bug: Call the other one, with an index set, instead.

Definition at line 2280 of file affine_constraints.h.

◆ AffineConstraints() [2/5]

template<typename number >
AffineConstraints< number >::AffineConstraints ( const IndexSet locally_stored_constraints)
inlineexplicit

Constructor. The supplied IndexSet defines for which indices this object will store constraints. In a calculation with a DoFHandler object based on parallel::distributed::Triangulation or parallel::shared::Triangulation, one should use the set of locally relevant DoFs (see GlossLocallyRelevantDof).

The given IndexSet allows the AffineConstraints container to save memory by just not caring about degrees of freedom that are not of importance to the current processor. In contrast, in parallel computations, if you do not provide such an index set (here, or using the reinit() function that takes such an argument), the current object will allocate memory proportional to the total number of degrees of freedom (accumulated over all processes), which is clearly wasteful and not efficient – and should be considered a bug.

Deprecated:
This constructor is equivalent to calling the following one with both of its arguments equal to the index set provided here. This is not wrong, but inefficient. Use the following constructor instead.

Definition at line 2287 of file affine_constraints.h.

◆ AffineConstraints() [3/5]

template<typename number >
AffineConstraints< number >::AffineConstraints ( const IndexSet locally_owned_dofs,
const IndexSet locally_stored_constraints 
)
inline

Constructor. The supplied IndexSet objects define which of the indices the object may encounter are locally owned, and for which indices this object will store constraints. In a calculation with a DoFHandler object based on parallel::distributed::Triangulation or parallel::shared::Triangulation, the first of these arguments should be the set of locally owned indices (see GlossLocallyOwnedDof) and the second argument should use the set of locally relevant DoFs (see GlossLocallyRelevantDof).

The locally_owned_dofs argument indicates which elements of vectors that will be passed to the member functions of this object are locally owned (which is generally equal to the set of locally owned degrees of freedom – thus the name of the argument). Knowledge of this set makes it possible to up front, in close(), determine which vector elements need to be imported when you later call the distribute() or related functions.

The locally_stored_constraints IndexSet allows the AffineConstraints container to save memory by just not caring about degrees of freedom that are not of importance to the current processor. In contrast, in parallel computations, if you do not provide such an index set (here, or using the reinit() function that takes such an argument), the current object will allocate memory proportional to the total number of degrees of freedom (accumulated over all processes), which is clearly wasteful and not efficient – and should be considered a bug. As mentioned above, one typically wants to provide the set of locally relevant degrees of freedom for this argument.

Definition at line 2296 of file affine_constraints.h.

◆ AffineConstraints() [4/5]

template<typename number >
AffineConstraints< number >::AffineConstraints ( const AffineConstraints< number > &  affine_constraints)
inlineexplicit

Copy constructor

Definition at line 2317 of file affine_constraints.h.

◆ AffineConstraints() [5/5]

template<typename number = double>
AffineConstraints< number >::AffineConstraints ( AffineConstraints< number > &&  affine_constraints)
defaultnoexcept

Move constructor

Member Function Documentation

◆ operator=() [1/2]

template<typename number = double>
AffineConstraints & AffineConstraints< number >::operator= ( const AffineConstraints< number > &  )
delete

Copy operator. Like for many other large objects, this operator is deleted to avoid its inadvertent use in places such as accidentally declaring a AffineConstraints object as a function argument by value, rather than by reference.

However, you can use the copy_from() function to explicitly copy AffineConstraints objects.

◆ operator=() [2/2]

template<typename number = double>
AffineConstraints & AffineConstraints< number >::operator= ( AffineConstraints< number > &&  affine_constraints)
defaultnoexcept

Move assignment operator

◆ copy_from()

template<typename number >
template<typename other_number >
void AffineConstraints< number >::copy_from ( const AffineConstraints< other_number > &  other)
inline

Copy the given object to the current one.

This function exists because operator=() is explicitly disabled.

Definition at line 2753 of file affine_constraints.h.

◆ reinit() [1/3]

template<typename number = double>
void AffineConstraints< number >::reinit ( )

Reset the AffineConstraints object. This results in an object that is as if it had been default-constructed. See the discussion in the documentation of the various constructors of this class for whether this is what you want – in particular, for parallel computations, you do not want to call this function, but rather one of the other reinit() functions that take index sets.

◆ reinit() [2/3]

template<typename number = double>
void AffineConstraints< number >::reinit ( const IndexSet locally_stored_constraints)

clear() the AffineConstraints object and supply an IndexSet that describes for which degrees of freedom this object can store constraints. See the discussion in the documentation of the constructor of this class that takes a single index set as argument.

Deprecated:
Use the reinit() function with two index set arguments instead.

◆ reinit() [3/3]

template<typename number = double>
void AffineConstraints< number >::reinit ( const IndexSet locally_owned_dofs,
const IndexSet locally_stored_constraints 
)

clear() the AffineConstraints object and supply both an IndexSet object that describes which degrees of freedom the current process owns, and an IndexSet that describes for which degrees of freedom this object can store constraints. See the discussion in the documentation of the constructor of this class that takes two index sets as argument for more information.

◆ can_store_line()

template<typename number >
bool AffineConstraints< number >::can_store_line ( const size_type  line_n) const
inline

Determines if we can store a constraint for the given line_n. This routine only matters in the distributed case and checks if the IndexSet provided to the constructor or reinit() function that determines which constraints to store, allows storage of the line provides as argument to this function. As a consequence, it will always return true if not in the distributed case.

Definition at line 2544 of file affine_constraints.h.

◆ get_locally_owned_indices()

template<typename number >
const IndexSet & AffineConstraints< number >::get_locally_owned_indices ( ) const
inline

Return the index set describing which part of the degrees of freedom to which this object stores constraints are "locally owned". Typically, these would be the locally owned degrees of freedom. This function returns the corresponding index set provided to either the constructor or the reinit() function of this class.

Definition at line 2553 of file affine_constraints.h.

◆ get_local_lines()

template<typename number >
const IndexSet & AffineConstraints< number >::get_local_lines ( ) const
inline

Return the index set describing locally relevant lines if any are present. Note that if no local lines were given, this represents an empty IndexSet, whereas otherwise it contains the global problem size and the local range.

Definition at line 2562 of file affine_constraints.h.

◆ add_selected_constraints()

template<typename number = double>
void AffineConstraints< number >::add_selected_constraints ( const AffineConstraints< number > &  constraints_in,
const IndexSet filter 
)

This function copies the content of constraints_in with DoFs that are element of the IndexSet filter. Elements that are not present in the IndexSet are ignored. All DoFs will be transformed to local index space of the filter, both the constrained DoFs and the other DoFs these entries are constrained to. The local index space of the filter is a contiguous numbering of all (global) DoFs that are elements in the filter.

If, for example, the filter represents the range [10,20), and the constraints object constraints_in includes the global indices {7,13,14}, the indices {3,4} are added to the calling constraints object (since 13 and 14 are elements in the filter and element 13 is the fourth element in the index, and 14 is the fifth).

This function provides an easy way to create a AffineConstraints for certain vector components in a vector-valued problem from a full AffineConstraints, i.e. extracting a diagonal subblock from a larger AffineConstraints. The block is specified by the IndexSet argument.

Deprecated:
This function is a combination of the get_view() function and merge() in that it selects a subset of constraints from another constraints object that is then merged into the current one. But the current function does not deal well with index sets. Furthermore, it simply discards parts of constraints that constrain one degree of freedom against ones that are not selected in the filter – something that should probably be considered a bug. Use get_view() and merge() instead.

◆ add_constraint()

template<typename number = double>
void AffineConstraints< number >::add_constraint ( const size_type  constrained_dof,
const ArrayView< const std::pair< size_type, number > > &  dependencies,
const number  inhomogeneity = 0 
)

Add a constraint to this object. This function adds a constraint of the form

\[ x_i = \sum_{j=1}^n w_j x_{k_j} + b \]

where \(i\) is the number of the degree of freedom to be constrained and is provided by the constrained_dofs argument. The weights \(w_j\) and indices \(k_j\) are provided as pairs in the dependencies argument, and the inhomogeneity \(b\) is provided by the last argument.

As an example, if you want to add the constraint

\[ x_{42} = 0.5 x_{12} + 0.5 x_{36} + 27 \]

you would call this function as follows:

constraints.add_constraint (42, {{12, 0.5}, {36, 0.5}}, 27.0);

On the other hand, if (as one often wants to) you need a constraint of the kind

\[ x_{42} = 27 \]

you would call this function as follows:

constraints.add_constraint (42, {}, 27.0);

If you want to constrain a degree of freedom to zero, i.e., require that

\[ x_{42} = 0 \]

you would call this function as follows:

constraints.add_constraint (42, {}, 0.0);

That said, this special case can be achieved in a more obvious way by calling

constraints.constrain_dof_to_zero (42);

instead.

◆ constrain_dof_to_zero()

template<typename number = double>
void AffineConstraints< number >::constrain_dof_to_zero ( const size_type  constrained_dof)

Constrain the given degree of freedom to be zero, i.e., require a constraint like

\[ x_{42} = 0. \]

Calling this function is equivalent to, but more readable than, saying

constraints.add_constraint (42, {}, 0.0);

It is not an error to call this function more than once on the same degree of freedom, but it is an error to call this function on a degree of freedom that has previously been constrained to either a different value than zero, or to a linear combination of degrees of freedom via the add_constraint() function.

◆ add_line()

template<typename number >
void AffineConstraints< number >::add_line ( const size_type  line_n)
inline

Add a new line to the matrix. If the line already exists, then the function simply returns without doing anything.

Definition at line 2333 of file affine_constraints.h.

◆ add_lines() [1/3]

template<typename number = double>
void AffineConstraints< number >::add_lines ( const std::vector< bool > &  lines)

Call the first add_line() function for every index i for which lines[i] is true.

This function essentially exists to allow adding several constraints of the form xi=0 all at once, where the set of indices i for which these constraints should be added are given by the argument of this function. On the other hand, just as if the single-argument add_line() function were called repeatedly, the constraints can later be modified to include linear dependencies using the add_entry() function as well as inhomogeneities using set_inhomogeneity().

◆ add_lines() [2/3]

template<typename number = double>
void AffineConstraints< number >::add_lines ( const std::set< size_type > &  lines)

Call the first add_line() function for every index i that appears in the argument.

This function essentially exists to allow adding several constraints of the form xi=0 all at once, where the set of indices i for which these constraints should be added are given by the argument of this function. On the other hand, just as if the single-argument add_line() function were called repeatedly, the constraints can later be modified to include linear dependencies using the add_entry() function as well as inhomogeneities using set_inhomogeneity().

◆ add_lines() [3/3]

template<typename number = double>
void AffineConstraints< number >::add_lines ( const IndexSet lines)

Call the first add_line() function for every index i that appears in the argument.

This function essentially exists to allow adding several constraints of the form xi=0 all at once, where the set of indices i for which these constraints should be added are given by the argument of this function. On the other hand, just as if the single-argument add_line() function were called repeatedly, the constraints can later be modified to include linear dependencies using the add_entry() function as well as inhomogeneities using set_inhomogeneity().

◆ add_entry()

template<typename number >
void AffineConstraints< number >::add_entry ( const size_type  constrained_dof_index,
const size_type  column,
const number  weight 
)
inline

Add an entry to a given line. In other words, this function adds a term \(a_{ij} x_j\) to the constraints for the \(i\)th degree of freedom.

If an entry with the same indices as the one this function call denotes already exists, then this function simply returns provided that the value of the entry is the same. Thus, it does no harm to enter a constraint twice.

Parameters
[in]constrained_dof_indexThe index \(i\) of the degree of freedom that is being constrained.
[in]columnThe index \(j\) of the degree of freedom being entered into the constraint for degree of freedom \(i\).
[in]weightThe factor \(a_{ij}\) that multiplies \(x_j\).

Definition at line 2365 of file affine_constraints.h.

◆ add_entries()

template<typename number = double>
void AffineConstraints< number >::add_entries ( const size_type  constrained_dof_index,
const std::vector< std::pair< size_type, number > > &  col_weight_pairs 
)

Add a whole series of entries, denoted by pairs of column indices and weight values, to a line of constraints. This function is equivalent to calling the preceding function several times, but is faster.

◆ set_inhomogeneity()

template<typename number >
void AffineConstraints< number >::set_inhomogeneity ( const size_type  constrained_dof_index,
const number  value 
)
inline

Set an inhomogeneity to the constraint for a degree of freedom. In other words, it adds a constant \(b_i\) to the constraint for degree of freedom \(i\). For this to work, you need to call add_line() first for the given degree of freedom.

Parameters
[in]constrained_dof_indexThe index \(i\) of the degree of freedom that is being constrained.
[in]valueThe right hand side value \(b_i\) for the constraint on the degree of freedom \(i\).

Definition at line 2407 of file affine_constraints.h.

◆ close()

template<typename number = double>
void AffineConstraints< number >::close ( )

Close the filling of entries. Since the lines of a matrix of this type are usually filled in an arbitrary order and since we do not want to use associative constrainers to store the lines, we need to sort the lines and within the lines the columns before usage of the matrix. This is done through this function.

Also, zero entries are discarded, since they are not needed.

After closing, no more entries are accepted. If the object was already closed, then this function returns immediately.

This function also resolves chains of constraints. For example, degree of freedom 13 may be constrained to \(u_{13} = \frac{u_3}{2} + \frac{u_7}{2}\) while degree of freedom 7 is itself constrained as \(u_{7} = \frac{u_2}{2} + \frac{u_4}{2}\). Then, the resolution will be that \(u_{13} = \frac{u_3}{2} + \frac{u_2}{4} + \frac{u_4}{4}\). Note, however, that cycles in this graph of constraints are not allowed, i.e., for example \(u_4\) may not itself be constrained, directly or indirectly, to \(u_{13}\) again.

◆ is_closed() [1/2]

template<typename number = double>
bool AffineConstraints< number >::is_closed ( ) const

Check if the function close() was called or there are no constraints locally, which is normally the case if a dummy AffineConstraints was created for the DG case.

◆ is_closed() [2/2]

template<typename number = double>
bool AffineConstraints< number >::is_closed ( const MPI_Comm  comm) const

Check if the function close() was called or there are no constraints globally, which is normally the case if a dummy AffineConstraints was created for the DG case.

◆ merge()

template<typename number >
template<typename other_number >
void AffineConstraints< number >::merge ( const AffineConstraints< other_number > &  other_constraints,
const MergeConflictBehavior  merge_conflict_behavior = no_conflicts_allowed,
const bool  allow_different_local_lines = false 
)

Merge the constraints represented by the object given as argument into the constraints represented by this object. Both objects may or may not be closed (by having their function close() called before). If this object was closed before, then it will be closed afterwards as well. Note, however, that if the other argument is closed, then merging may be significantly faster.

Using the default value of the second arguments, the constraints in each of the two objects (the old one represented by this object and the argument) may not refer to the same degree of freedom, i.e. a degree of freedom that is constrained in one object may not be constrained in the second. If this is nevertheless the case, an exception is thrown. However, this behavior can be changed by providing a different value for the second argument.

By default, merging two AffineConstraints objects that are initialized with different IndexSet objects is not allowed. This behavior can be altered by setting allow_different_local_lines appropriately.

Merging an AffineConstraints that is initialized with an IndexSet and one that is not initialized with an IndexSet is not implemented.

Definition at line 2777 of file affine_constraints.h.

◆ shift()

template<typename number = double>
void AffineConstraints< number >::shift ( const size_type  offset)

Shift all entries of this matrix down offset rows and over offset columns. If this object is initialized with an IndexSet, local_lines are shifted as well.

This function is useful if you are building block matrices, where all blocks are built by the same DoFHandler object, i.e. the matrix size is larger than the number of degrees of freedom. Since several matrix rows and columns correspond to the same degrees of freedom, you'd generate several constraint objects, then shift them, and finally merge() them together again.

◆ get_view()

template<typename number = double>
AffineConstraints AffineConstraints< number >::get_view ( const IndexSet mask) const

This function provides a "view" into a constraint object. Specifically, given a "mask" index set that describes which constraints we are interested in, it returns an AffineConstraints object that contains only those constraints that correspond to degrees of freedom that are listed in the mask, with indices shifted so that they correspond to the position within the mask. This process is the same as how IndexSet::get_view() computes the shifted indices. The function is typically used to extract from an AffineConstraints object corresponding to a DoFHandler only those constraints that correspond to a specific variable (say, to the velocity in a Stokes system) so that the resulting AffineConstraints object can be applied to a single block of a block vector of solutions; in this case, the mask would be the index set of velocity degrees of freedom, as a subset of all degrees of freedom.

This function can only work if the degrees of freedom selected by the mask are constrained only against other degrees of freedom that are listed in the mask. In the example above, this means that constraints for the selected velocity degrees of freedom are only against other velocity degrees of freedom, but not against any pressure degrees of freedom. If that is not so, an assertion will be triggered.

A typical case where this function is useful is as follows. Say, you have a block linear system in which you have blocks corresponding to variables \((u,p,T,c)\) (which you can think of as velocity, pressure, temperature, and chemical composition – or whatever other kind of problem you are currently considering in your own work). Let's assume we have developed a linear solver or preconditioner that first solves the coupled \(u\)- \(T\) system, and once that is done, solves the \(p\)- \(c\) system. In this case, it is often useful to set up block vectors with only two components corresponding to the \(u\) and \(T\) components, and later for only the \(p\)- \(c\) components of the solution. As part of this, you will want to apply constraints (using the distribute() function of this class) to only the 2-block vector, but for this you need to obtain an AffineConstraints object that represents only those constraints that correspond to the variables in question, and in the order in which they appear in the 2-block vector rather than in global 4-block vectors. This function allows you to extract such an object corresponding to a subset of constraints by applying a mask to the global constraints object that corresponds to the variables we're currently interested in. For the \(u\)- \(T\) system, we need a mask that contains all indices of \(u\) degrees of freedom as well as all indices of \(T\) degrees of freedom.

◆ clear()

template<typename number = double>
void AffineConstraints< number >::clear ( )

Clear all entries of this matrix. Reset the flag determining whether new entries are accepted or not.

This function may be called also on objects which are empty or already cleared.

◆ n_constraints()

template<typename number >
types::global_dof_index AffineConstraints< number >::n_constraints ( ) const
inline

Return number of constraints stored in this matrix.

Definition at line 2438 of file affine_constraints.h.

◆ n_identities()

template<typename number >
types::global_dof_index AffineConstraints< number >::n_identities ( ) const
inline

Return number of constraints stored in this matrix that are identities, i.e., those constraints with only one degree of freedom and weight one.

Definition at line 2445 of file affine_constraints.h.

◆ n_inhomogeneities()

template<typename number >
types::global_dof_index AffineConstraints< number >::n_inhomogeneities ( ) const
inline

Return number of constraints stored in this matrix that have an inhomogenity, i.e., those constraints with a non-trivial inhomogeneous value.

Definition at line 2457 of file affine_constraints.h.

◆ is_constrained()

template<typename number >
bool AffineConstraints< number >::is_constrained ( const size_type  line_n) const
inline

Return whether the degree of freedom with number line_n is a constrained one.

Note that if close() was called before, then this function is significantly faster, since then the constrained degrees of freedom are sorted and we can do a binary search, while before close() was called, we have to perform a linear search through all entries.

Definition at line 2468 of file affine_constraints.h.

◆ is_identity_constrained()

template<typename number = double>
bool AffineConstraints< number >::is_identity_constrained ( const size_type  line_n) const

Return whether the dof is constrained, and whether it is constrained to only one other degree of freedom with weight one. The function therefore returns whether the degree of freedom would simply be eliminated in favor of exactly one other degree of freedom.

The function returns false if either the degree of freedom is not constrained at all, or if it is constrained to more than one other degree of freedom, or if it is constrained to only one degree of freedom but with a weight different from one.

◆ are_identity_constrained()

template<typename number = double>
bool AffineConstraints< number >::are_identity_constrained ( const size_type  line_n_1,
const size_type  line_n_2 
) const

Return whether the two given degrees of freedom are linked by an equality constraint that either constrains index1 to be so that index1=index2 or constrains index2 so that index2=index1.

◆ max_constraint_indirections()

template<typename number = double>
size_type AffineConstraints< number >::max_constraint_indirections ( ) const

Return the maximum number of other dofs that one dof is constrained to. For example, in 2d a hanging node is constrained only to its two neighbors, so the returned value would be 2. However, for higher order elements and/or higher dimensions, or other types of constraints, this number is no more obvious.

The name indicates that within the system matrix, references to a constrained node are indirected to the nodes it is constrained to.

◆ is_inhomogeneously_constrained()

template<typename number >
bool AffineConstraints< number >::is_inhomogeneously_constrained ( const size_type  index) const
inline

Return true in case the dof is constrained and there is a non-trivial inhomogeneous values set to the dof.

Definition at line 2480 of file affine_constraints.h.

◆ has_inhomogeneities()

template<typename number = double>
bool AffineConstraints< number >::has_inhomogeneities ( ) const

Return false if all constraints in the AffineConstraints are homogeneous ones, and true if there is at least one inhomogeneity.

◆ get_constraint_entries()

template<typename number >
const std::vector< std::pair< types::global_dof_index, number > > * AffineConstraints< number >::get_constraint_entries ( const size_type  line_n) const
inline

Return a pointer to the vector of entries if a line is constrained, and a zero pointer in case the dof is not constrained.

Definition at line 2498 of file affine_constraints.h.

◆ get_inhomogeneity()

template<typename number >
number AffineConstraints< number >::get_inhomogeneity ( const size_type  line_n) const
inline

Return the value of the inhomogeneity stored in the constrained dof line_n. Unconstrained dofs also return a zero value.

Definition at line 2515 of file affine_constraints.h.

◆ print()

template<typename number = double>
void AffineConstraints< number >::print ( std::ostream &  out) const

Print the constraints represented by the current object to the given stream.

For each constraint of the form

\[ x_{42} = 0.5 x_2 + 0.25 x_{14} + 2.75 \]

this function will write a sequence of lines that look like this:

42 2 : 0.5
42 14 : 0.25
42 : 2.75

The last line is only shown if the inhomogeneity (here: 2.75) is nonzero.

A block of lines such as the one above is repeated for each constrained degree of freedom.

◆ write_dot()

template<typename number = double>
void AffineConstraints< number >::write_dot ( std::ostream &  ) const

Write the graph of constraints in 'dot' format. 'dot' is a program that can take a list of nodes and produce a graphical representation of the graph of constrained degrees of freedom and the degrees of freedom they are constrained to.

The output of this function can be used as input to the 'dot' program that can convert the graph into a graphical representation in postscript, png, xfig, and a number of other formats.

This function exists mostly for debugging purposes.

◆ memory_consumption()

template<typename number = double>
std::size_t AffineConstraints< number >::memory_consumption ( ) const

Determine an estimate for the memory consumption (in bytes) of this object.

◆ resolve_indices()

template<typename number = double>
void AffineConstraints< number >::resolve_indices ( std::vector< types::global_dof_index > &  indices) const

Add the constraint indices associated to the indices in the given vector. After a call to this function, the indices vector contains the initial elements and all the associated constrained indices. This function sorts the elements and suppresses duplicates.

◆ condense() [1/10]

template<typename number = double>
void AffineConstraints< number >::condense ( SparsityPattern sparsity) const

Condense a sparsity pattern. The name of the function mimics the name of the function we use to condense linear systems, but it is a bit of a misnomer for the current context. This is because in the context of linear systems, we eliminate certain rows and columns of the linear system, i.e., we "reduce" or "condense" the linear system. On the other hand, in the current context, the functions does not remove nonzero entries from the sparsity pattern. Rather, it adds those nonzero entry locations to the sparsity pattern that will later be needed for the process of condensation of constrained degrees of freedom from a linear system.

Since this function adds new nonzero entries to the sparsity pattern, the given sparsity pattern must not be compressed. The current object must be closed. The sparsity pattern is compressed at the end of the function.

◆ condense() [2/10]

template<typename number = double>
void AffineConstraints< number >::condense ( BlockSparsityPattern sparsity) const

Same function as above, but condenses square block sparsity patterns.

◆ condense() [3/10]

template<typename number = double>
void AffineConstraints< number >::condense ( DynamicSparsityPattern sparsity) const

Same function as above, but condenses square compressed sparsity patterns.

◆ condense() [4/10]

template<typename number = double>
void AffineConstraints< number >::condense ( BlockDynamicSparsityPattern sparsity) const

Same function as above, but condenses square compressed sparsity patterns.

◆ condense() [5/10]

template<typename number = double>
void AffineConstraints< number >::condense ( SparseMatrix< number > &  matrix) const

Condense a given matrix, i.e., eliminate the rows and columns of the matrix that correspond to constrained degrees of freedom.

See the general documentation of this class for more detailed information.

◆ condense() [6/10]

template<typename number = double>
void AffineConstraints< number >::condense ( BlockSparseMatrix< number > &  matrix) const

Same function as above, but condenses square block sparse matrices.

◆ condense() [7/10]

template<typename number = double>
template<typename VectorType >
void AffineConstraints< number >::condense ( VectorType &  vec) const

Condense the given vector in-place. The VectorType may be a Vector<float>, Vector<number>, BlockVector<...>, a PETSc or Trilinos vector wrapper class, or any other type having the same interface. Note that this function does not take any inhomogeneity into account and throws an exception in case there are any inhomogeneities. Use the function using both a matrix and vector for that case.

Note
This function does not work for MPI vectors. Use condense() with two vector arguments instead.

◆ condense() [8/10]

template<typename number = double>
template<typename VectorType >
void AffineConstraints< number >::condense ( const VectorType &  vec_ghosted,
VectorType &  output 
) const

The function copies and condenses values from vec_ghosted into output. In a serial code it is equivalent to calling condense (vec). If called in parallel, vec_ghosted is supposed to contain ghost elements while output should not.

◆ condense() [9/10]

template<typename number = double>
template<typename VectorType >
void AffineConstraints< number >::condense ( SparseMatrix< number > &  matrix,
VectorType &  vector 
) const

Condense a given matrix and a given vector by eliminating rows and columns of the linear system that correspond to constrained degrees of freedom. The sparsity pattern associated with the matrix needs to be condensed and compressed. This function is the appropriate choice for applying inhomogeneous constraints.

The current object must be closed to call this function.

See the general documentation of this class for more detailed information.

◆ condense() [10/10]

template<typename number = double>
template<typename BlockVectorType >
void AffineConstraints< number >::condense ( BlockSparseMatrix< number > &  matrix,
BlockVectorType &  vector 
) const

Same function as above, but condenses square block sparse matrices and vectors.

◆ set_zero()

template<typename number >
template<typename VectorType >
void AffineConstraints< number >::set_zero ( VectorType &  vec) const
inline

Set the values of all constrained DoFs in a vector to zero. The VectorType may be a Vector<float>, Vector<number>, BlockVector<...>, a PETSc or Trilinos vector wrapper class, or any other type having the same interface.

Definition at line 2425 of file affine_constraints.h.

◆ distribute_local_to_global() [1/11]

template<typename number >
template<class InVector , class OutVector >
void AffineConstraints< number >::distribute_local_to_global ( const InVector &  local_vector,
const std::vector< size_type > &  local_dof_indices,
OutVector &  global_vector 
) const
inline

This function takes a vector of local contributions (local_vector) corresponding to the degrees of freedom indices given in local_dof_indices and distributes them to the global vector. In other words, this function implements a scatter operation. In most cases, these local contributions will be the result of an integration over a cell or face of a cell. However, as long as local_vector and local_dof_indices have the same number of elements, this function is happy with whatever it is given.

In contrast to the similar function in the DoFAccessor class, this function also takes care of constraints, i.e. if one of the elements of local_dof_indices belongs to a constrained node, then rather than writing the corresponding element of local_vector into global_vector, the element is distributed to the entries in the global vector to which this particular degree of freedom is constrained.

Thus, by using this function to distribute local contributions to the global object, one saves the call to the condense function after the vectors and matrices are fully assembled. On the other hand, by consequence, the function does not only write into the entries enumerated by the local_dof_indices array, but also (possibly) others as necessary.

Note that this function will apply all constraints as if they were homogeneous. For correctly setting inhomogeneous constraints, use the similar function with a matrix argument or the function with both matrix and vector arguments.

Note
This function in itself is thread-safe, i.e., it works properly also when several threads call it simultaneously. However, the function call is only thread-safe if the underlying global vector allows for simultaneous access and the access is not to rows with the same global index at the same time. This needs to be made sure from the caller's site. There is no locking mechanism inside this method to prevent data races.
Parameters
[in]local_vectorVector of local contributions.
[in]local_dof_indicesLocal degrees of freedom indices corresponding to the vector of local contributions.
[out]global_vectorThe global vector to which all local contributions will be added.

Definition at line 2626 of file affine_constraints.h.

◆ distribute_local_to_global() [2/11]

template<typename number = double>
template<typename VectorType >
void AffineConstraints< number >::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

This function takes a vector of local contributions (local_vector) corresponding to the degrees of freedom indices given in local_dof_indices and distributes them to the global vector. In other words, this function implements a scatter operation. In most cases, these local contributions will be the result of an integration over a cell or face of a cell. However, as long as local_vector and local_dof_indices have the same number of elements, this function is happy with whatever it is given.

In contrast to the similar function in the DoFAccessor class, this function also takes care of constraints, i.e. if one of the elements of local_dof_indices belongs to a constrained node, then rather than writing the corresponding element of local_vector into global_vector, the element is distributed to the entries in the global vector to which this particular degree of freedom is constrained.

Thus, by using this function to distribute local contributions to the global object, one saves the call to the condense function after the vectors and matrices are fully assembled. On the other hand, by consequence, the function does not only write into the entries enumerated by the local_dof_indices array, but also (possibly) others as necessary. This includes writing into diagonal elements of the matrix if the corresponding degree of freedom is constrained.

The fourth argument local_matrix is intended to be used in case one wants to apply inhomogeneous constraints on the vector only. Such a situation could be where one wants to assemble of a right hand side vector on a problem with inhomogeneous constraints, but the global matrix has been assembled previously. A typical example of this is a time stepping algorithm where the stiffness matrix is assembled once, and the right hand side updated every time step. Note that, however, the entries in the columns of the local matrix have to be exactly the same as those that have been written into the global matrix. Otherwise, this function will not be able to correctly handle inhomogeneities.

Note
This function in itself is thread-safe, i.e., it works properly also when several threads call it simultaneously. However, the function call is only thread-safe if the underlying global vector allows for simultaneous access and the access is not to rows with the same global index at the same time. This needs to be made sure from the caller's site. There is no locking mechanism inside this method to prevent data races.

◆ distribute_local_to_global() [3/11]

template<typename number = double>
template<typename VectorType >
void AffineConstraints< number >::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

Same as the previous function, except that it uses two (possibly) different index sets to correctly handle inhomogeneities when the local matrix is computed from a combination of two neighboring elements, for example for an edge integral term in DG. Note that in the case that these two elements have different polynomial degree, the local matrix is rectangular.

local_dof_indices_row is the set of row indices and local_dof_indices_col is the set of column indices of the local matrix. diagonal=false says whether the two index sets are equal or not.

If both index sets are equal, diagonal must be set to true or we simply use the previous function. If both index sets are different (diagonal=false) the global_vector is modified to handle inhomogeneities but no entries from local_vector are added. Note that the edge integrals for inner edged for DG do not contribute any values to the right hand side.

◆ distribute_local_to_global() [4/11]

template<typename number >
template<typename VectorType >
void AffineConstraints< number >::distribute_local_to_global ( const size_type  index,
const number  value,
VectorType &  global_vector 
) const
inline

Enter a single value into a result vector, obeying constraints.

Definition at line 2572 of file affine_constraints.h.

◆ distribute_local_to_global() [5/11]

template<typename number >
template<typename ForwardIteratorVec , typename ForwardIteratorInd , typename VectorType >
void AffineConstraints< number >::distribute_local_to_global ( ForwardIteratorVec  local_vector_begin,
ForwardIteratorVec  local_vector_end,
ForwardIteratorInd  local_indices_begin,
VectorType &  global_vector 
) const
inline

This function takes a pointer to a vector of local contributions (local_vector) corresponding to the degrees of freedom indices given in local_dof_indices and distributes them to the global vector. In other words, this function implements a scatter operation. In most cases, these local contributions will be the result of an integration over a cell or face of a cell. However, as long as the entries in local_dof_indices indicate reasonable global vector entries, this function is happy with whatever it is given.

If one of the elements of local_dof_indices belongs to a constrained node, then rather than writing the corresponding element of local_vector into global_vector, the element is distributed to the entries in the global vector to which this particular degree of freedom is constrained.

Thus, by using this function to distribute local contributions to the global object, one saves the call to the condense function after the vectors and matrices are fully assembled. Note that this function completely ignores inhomogeneous constraints.

Note
This function in itself is thread-safe, i.e., it works properly also when several threads call it simultaneously. However, the function call is only thread-safe if the underlying global vector allows for simultaneous access and the access is not to rows with the same global index at the same time. This needs to be made sure from the caller's site. There is no locking mechanism inside this method to prevent data races.

Definition at line 2596 of file affine_constraints.h.

◆ distribute_local_to_global() [6/11]

template<typename number >
template<typename MatrixType >
void AffineConstraints< number >::distribute_local_to_global ( const FullMatrix< number > &  local_matrix,
const std::vector< size_type > &  local_dof_indices,
MatrixType &  global_matrix 
) const
inline

This function takes a matrix of local contributions (local_matrix) corresponding to the degrees of freedom indices given in local_dof_indices and distributes them to the global matrix. In other words, this function implements a scatter operation. In most cases, these local contributions will be the result of an integration over a cell or face of a cell. However, as long as local_matrix and local_dof_indices have the same number of elements, this function is happy with whatever it is given.

In contrast to the similar function in the DoFAccessor class, this function also takes care of constraints, i.e. if one of the elements of local_dof_indices belongs to a constrained node, then rather than writing the corresponding element of local_matrix into global_matrix, the element is distributed to the entries in the global matrix to which this particular degree of freedom is constrained.

With this scheme, we never write into rows or columns of constrained degrees of freedom. In order to make sure that the resulting matrix can still be inverted, we need to do something with the diagonal elements corresponding to constrained nodes. Thus, if a degree of freedom in local_dof_indices is constrained, we distribute the corresponding entries in the matrix, but also add the absolute value of the diagonal entry of the local matrix to the corresponding entry in the global matrix. Assuming the discretized operator is positive definite, this guarantees that the diagonal entry is always non-zero, positive, and of the same order of magnitude as the other entries of the matrix. On the other hand, when solving a source problem \(Au=f\) the exact value of the diagonal element is not important, since the value of the respective degree of freedom will be overwritten by the distribute() call later on anyway.

Note
The procedure described above adds an unforeseeable number of artificial eigenvalues to the spectrum of the matrix. Therefore, it is recommended to use the equivalent function with two local index vectors in such a case.

By using this function to distribute local contributions to the global object, one saves the call to the condense function after the vectors and matrices are fully assembled.

Note
This function in itself is thread-safe, i.e., it works properly also when several threads call it simultaneously. However, the function call is only thread-safe if the underlying global matrix allows for simultaneous access and the access is not to rows with the same global index at the same time. This needs to be made sure from the caller's site. There is no locking mechanism inside this method to prevent data races.

Definition at line 2931 of file affine_constraints.h.

◆ distribute_local_to_global() [7/11]

template<typename number = double>
template<typename MatrixType >
void AffineConstraints< number >::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

This function does almost the same as the function above but can treat general rectangular matrices. The main difference to achieve this is that the diagonal entries in constrained rows are left untouched instead of being filled with arbitrary values.

Since the diagonal entries corresponding to eliminated degrees of freedom are not set, the result may have a zero eigenvalue, if applied to a square matrix. This has to be considered when solving the resulting problems. For solving a source problem \(Au=f\), it is possible to set the diagonal entry after building the matrix by a piece of code of the form

for (unsigned int i=0;i<matrix.m();++i)
if (constraints.is_constrained(i))
matrix.diag_element(i) = 1.;

The value of one which is used here is arbitrary, but in the context of Krylov space methods uncritical, since it corresponds to an invariant subspace. If the other matrix entries are smaller or larger by a factor close to machine accuracy, it may be advisable to adjust it.

For solving eigenvalue problems, this will only add one spurious zero eigenvalue (with a multiplicity that is possibly greater than one). Taking this into account, nothing else has to be changed.

◆ distribute_local_to_global() [8/11]

template<typename number = double>
template<typename MatrixType >
void AffineConstraints< number >::distribute_local_to_global ( const FullMatrix< number > &  local_matrix,
const std::vector< size_type > &  row_indices,
const AffineConstraints< number > &  column_affine_constraints,
const std::vector< size_type > &  column_indices,
MatrixType &  global_matrix 
) const

This function does almost the same as the function above for general rectangular matrices but uses different AffineConstraints objects on the row and column indices. The convention is that row indices are constrained according to the calling AffineConstraints *this, whereas column indices are constrained according to the given AffineConstraints column_affine_constraints. This function allows to handle the case where rows and columns of a matrix are represented by different function spaces with their own enumeration of indices, as e.g. in mixed finite element problems with separate DoFHandler objects or for flux matrices between different levels in multigrid methods.

Like the other method with separate slots for row and column indices, this method does not add diagonal entries to eliminated degrees of freedom. See there for a more elaborate description.

◆ distribute_local_to_global() [9/11]

template<typename number >
template<typename MatrixType , typename VectorType >
void AffineConstraints< number >::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
inline

This function simultaneously writes elements into matrix and vector, according to the constraints specified by the calling AffineConstraints. In other words, it performs the scatter operation of the corresponding functions for matrices and vectors at the same time. This function can correctly handle inhomogeneous constraints as well. For the parameter use_inhomogeneities_for_rhs see the documentation in Constraints on degrees of freedom module.

Note
This function in itself is thread-safe, i.e., it works properly also when several threads call it simultaneously. However, the function call is only thread-safe if the underlying global matrix and vector allow for simultaneous access and the access is not to rows with the same global index at the same time. This needs to be made sure from the caller's site. There is no locking mechanism inside this method to prevent data races.

Definition at line 2956 of file affine_constraints.h.

◆ add_entries_local_to_global() [1/3]

template<typename number = double>
void AffineConstraints< number >::add_entries_local_to_global ( const std::vector< size_type > &  local_dof_indices,
SparsityPatternBase sparsity_pattern,
const bool  keep_constrained_entries = true,
const Table< 2, bool > &  dof_mask = Table< 2, bool >() 
) const

Do a similar operation as the distribute_local_to_global() function that distributes writing entries into a matrix for constrained degrees of freedom, except that here we don't write into a matrix but only allocate sparsity pattern entries.

As explained in the hp-paper and in step-27, first allocating a sparsity pattern and later coming back and allocating additional entries for those matrix entries that will be written to due to the elimination of constrained degrees of freedom (using AffineConstraints::condense() ), can be a very expensive procedure. It is cheaper to allocate these entries right away without having to do a second pass over the sparsity pattern object. This function does exactly that.

Because the function only allocates entries in a sparsity pattern, all it needs to know are the degrees of freedom that couple to each other. Unlike the previous function, no actual values are written, so the second input argument is not necessary here.

The third argument to this function, keep_constrained_entries determines whether the function shall allocate entries in the sparsity pattern at all for entries that will later be set to zero upon condensation of the matrix. These entries are necessary if the matrix is built unconstrained, and only later condensed. They are not necessary if the matrix is built using the distribute_local_to_global() function of this class which distributes entries right away when copying a local matrix into a global object. The default of this argument is true, meaning to allocate the few entries that may later be set to zero.

By default, the function adds entries for all pairs of indices given in the first argument to the sparsity pattern (unless keep_constrained_entries is false). However, sometimes one would like to only add a subset of all of these pairs. In that case, the last argument can be used which specifies a boolean mask which of the pairs of indices should be considered. If the mask is false for a pair of indices, then no entry will be added to the sparsity pattern for this pair, irrespective of whether one or both of the indices correspond to constrained degrees of freedom.

This function is not typically called from user code, but is used in the DoFTools::make_sparsity_pattern() function when passed an AffineConstraints object.

Note
This function in itself is thread-safe, i.e., it works properly also when several threads call it simultaneously. However, the function call is only thread-safe if the underlying global sparsity pattern allows for simultaneous access and the access is not to rows with the same global index at the same time. This needs to be made sure from the caller's site. There is no locking mechanism inside this method to prevent data races.

◆ add_entries_local_to_global() [2/3]

template<typename number = double>
void AffineConstraints< number >::add_entries_local_to_global ( const std::vector< size_type > &  row_indices,
const std::vector< size_type > &  col_indices,
SparsityPatternBase sparsity_pattern,
const bool  keep_constrained_entries = true,
const Table< 2, bool > &  dof_mask = Table< 2, bool >() 
) const

Similar to the other function, but for non-quadratic sparsity patterns.

◆ add_entries_local_to_global() [3/3]

template<typename number = double>
void AffineConstraints< number >::add_entries_local_to_global ( const std::vector< size_type > &  row_indices,
const AffineConstraints< number > &  col_constraints,
const std::vector< size_type > &  col_indices,
SparsityPatternBase sparsity_pattern,
const bool  keep_constrained_entries = true,
const Table< 2, bool > &  dof_mask = Table< 2, bool >() 
) const

Similar to the other function, but for non-quadratic sparsity patterns, and for different constraints in the column space.

◆ get_dof_values()

template<typename number >
template<typename ForwardIteratorVec , typename ForwardIteratorInd , typename VectorType >
void AffineConstraints< number >::get_dof_values ( const VectorType &  global_vector,
ForwardIteratorInd  local_indices_begin,
ForwardIteratorVec  local_vector_begin,
ForwardIteratorVec  local_vector_end 
) const
inline

This function imports values from a global vector (global_vector) by applying the constraints to a vector of local values, expressed in iterator format. In most cases, the local values will be identified by the local dof values on a cell. However, as long as the entries in local_dof_indices indicate reasonable global vector entries, this function is happy with whatever it is given.

If one of the elements of local_dof_indices belongs to a constrained node, then rather than writing the corresponding element of global_vector into local_vector, the constraints are resolved as the respective distribute function does, i.e., the local entry is constructed from the global entries to which this particular degree of freedom is constrained.

In contrast to the similar function get_dof_values in the DoFAccessor class, this function does not need the constrained values to be correctly set (i.e., distribute to be called).

Definition at line 2645 of file affine_constraints.h.

◆ distribute()

template<typename number = double>
template<typename VectorType >
void AffineConstraints< number >::distribute ( VectorType &  vec) const

Given a vector, set all constrained degrees of freedom to values so that the constraints are satisfied. For example, if the current object stores the constraint \(x_3=\frac 12 x_1 + \frac 12 x_2\), then this function will read the values of \(x_1\) and \(x_2\) from the given vector and set the element \(x_3\) according to this constraints. Similarly, if the current object stores the constraint \(x_{42}=208\), then this function will set the 42nd element of the given vector to 208.

Note
If this function is called with a parallel vector vec, then the vector must not contain ghost elements.

◆ get_lines()

template<typename number = double>
LineRange AffineConstraints< number >::get_lines ( ) const

Return a range object containing (const) iterators to all line entries stored in the AffineConstraints container. Such a range is useful to initialize range-based for loops as supported by C++11.

Returns
A range object for the half open range [this->begin(), this->end()) of line entries.

◆ is_consistent_in_parallel()

template<typename number = double>
bool AffineConstraints< number >::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

Check if the current object is consistent on all processors in a distributed computation.

This method checks if all processors agree on the constraints for their local lines as given by locally_active_dofs. This method is a collective operation and will return true only if all processors are consistent.

Please supply the owned DoFs per processor as returned by Utilities::MPI::all_gather(MPI_Comm, DoFHandler::locally_owned_dofs()) as locally_owned_dofs and the result of DoFTools::extract_locally_active_dofs() as locally_active_dofs. The former is used to determine ownership of the specific DoF, while the latter is used as the set of rows that need to be checked.

If verbose is set to true, additional debug information is written to std::cout.

Note
This method exchanges all constraint information of locally active lines and is as such slow for large computations and should probably only be used in debug mode. We do not check all lines returned by get_local_lines() but only the locally active ones, as we allow processors to not know about some locally relevant rows.
Returns
Whether all AffineConstraints objects are consistent. Returns the same value on all processors.

◆ make_consistent_in_parallel()

template<typename number = double>
void AffineConstraints< number >::make_consistent_in_parallel ( const IndexSet locally_owned_dofs,
const IndexSet locally_relevant_dofs,
const MPI_Comm  mpi_communicator 
)

Make the current object consistent on all processors in a distributed computation. One should call this function before calling close().

◆ calculate_line_index()

template<typename number >
types::global_dof_index AffineConstraints< number >::calculate_line_index ( const size_type  line_n) const
inlineprivate

Internal function to calculate the index of line line_n in the vector lines_cache using local_lines.

Definition at line 2529 of file affine_constraints.h.

◆ distribute_local_to_global() [10/11]

template<typename number = double>
template<typename MatrixType , typename VectorType >
void AffineConstraints< number >::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::bool_constant< false >   
) const
private

This function actually implements the local_to_global function for standard (non-block) matrices.

◆ distribute_local_to_global() [11/11]

template<typename number = double>
template<typename MatrixType , typename VectorType >
void AffineConstraints< number >::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::bool_constant< true >   
) const
private

This function actually implements the local_to_global function for block matrices.

◆ make_sorted_row_list() [1/2]

template<typename number = double>
void AffineConstraints< number >::make_sorted_row_list ( const std::vector< size_type > &  local_dof_indices,
internal::AffineConstraints< number >::GlobalRowsFromLocal< number > &  global_rows 
) const
private

Internal helper function for distribute_local_to_global function.

Creates a list of affected global rows for distribution, including the local rows where the entries come from. The list is sorted according to the global row indices.

◆ make_sorted_row_list() [2/2]

template<typename number = double>
void AffineConstraints< number >::make_sorted_row_list ( const std::vector< size_type > &  local_dof_indices,
std::vector< size_type > &  active_dofs 
) const
private

Internal helper function for add_entries_local_to_global function.

Creates a list of affected rows for distribution without any additional information, otherwise similar to the other make_sorted_row_list() function.

◆ resolve_vector_entry()

template<typename number = double>
template<typename MatrixScalar , typename VectorScalar >
ProductType< VectorScalar, MatrixScalar >::type AffineConstraints< number >::resolve_vector_entry ( const size_type  i,
const internal::AffineConstraints< number >::GlobalRowsFromLocal< number > &  global_rows,
const Vector< VectorScalar > &  local_vector,
const std::vector< size_type > &  local_dof_indices,
const FullMatrix< MatrixScalar > &  local_matrix 
) const
private

Internal helper function for distribute_local_to_global function.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 135 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
inlineinherited

Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.

Definition at line 309 of file subscriptor.h.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Friends And Related Symbol Documentation

◆ AffineConstraints

template<typename number = double>
template<typename >
friend class AffineConstraints
friend

Definition at line 2122 of file affine_constraints.h.

Member Data Documentation

◆ lines

template<typename number = double>
std::vector<ConstraintLine> AffineConstraints< number >::lines
private

Store the lines of the matrix. Entries are usually appended in an arbitrary order and insertion into a vector is done best at the end, so the order is unspecified after all entries are inserted. Sorting of the entries takes place when calling the close() function.

We could, instead of using a vector, use an associative array, like a map to store the lines. This, however, would mean a much more fragmented heap since it allocates many small objects, and would additionally make usage of this matrix much slower.

Definition at line 2136 of file affine_constraints.h.

◆ lines_cache

template<typename number = double>
std::vector<size_type> AffineConstraints< number >::lines_cache
private

A list of size_type that contains the position of the ConstraintLine of a constrained degree of freedom, or numbers::invalid_size_type if the degree of freedom is not constrained. The numbers::invalid_size_type return value returns thus whether there is a constraint line for a given degree of freedom index. Note that this class has no notion of how many degrees of freedom there really are, so if we check whether there is a constraint line for a given degree of freedom, then this vector may actually be shorter than the index of the DoF we check for.

This field exists since when adding a new constraint line we have to figure out whether it already exists. Previously, we would simply walk the unsorted list of constraint lines until we either hit the end or found it. This algorithm is O(N) if N is the number of constraints, which makes it O(N^2) when inserting all constraints. For large problems with many constraints, this could easily take 5-10 per cent of the total run time. With this field, we can save this time since we find any constraint in O(1) time or get to know that it a certain degree of freedom is not constrained.

To make things worse, traversing the list of existing constraints requires reads from many different places in memory. Thus, in large 3d applications, the add_line() function showed up very prominently in the overall compute time, mainly because it generated a lot of cache misses. This should also be fixed by using the O(1) algorithm to access the fields of this array.

The field is useful in a number of other contexts as well, e.g. when one needs random access to the constraints as in all the functions that apply constraints on the fly while add cell contributions into vectors and matrices.

Definition at line 2170 of file affine_constraints.h.

◆ locally_owned_dofs

template<typename number = double>
IndexSet AffineConstraints< number >::locally_owned_dofs
private

This IndexSet denotes the set of locally owned DoFs (or, more correctly, the locally owned vector elements when operating on parallel vectors) and is used in the distribute() function when determining for which degrees of freedom to import dependent DoFs.

Definition at line 2178 of file affine_constraints.h.

◆ local_lines

template<typename number = double>
IndexSet AffineConstraints< number >::local_lines
private

This IndexSet is used to limit the lines to save in the AffineConstraints to a subset. This is necessary, because the lines_cache vector would become too big in a distributed calculation.

Definition at line 2185 of file affine_constraints.h.

◆ needed_elements_for_distribute

template<typename number = double>
IndexSet AffineConstraints< number >::needed_elements_for_distribute
private

The set of elements that are necessary to call distribute(). Because we need to set locally owned constrained elements of a vector, the needed elements include all vector entries that these constrained elements depend on – which may be owned by other processes.

This variable is set in close().

Definition at line 2195 of file affine_constraints.h.

◆ sorted

template<typename number = double>
bool AffineConstraints< number >::sorted
private

Store whether the arrays are sorted. If so, no new entries can be added.

Definition at line 2200 of file affine_constraints.h.

◆ scratch_data

template<typename number = double>
Threads::ThreadLocalStorage< internal::AffineConstraints::ScratchData<number> > AffineConstraints< number >::scratch_data
mutableprivate

Definition at line 2204 of file affine_constraints.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 218 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 240 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 248 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when accessing the mutable members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers().

Definition at line 271 of file subscriptor.h.


The documentation for this class was generated from the following file: