Reference documentation for deal.II version 9.1.1
|
#include <deal.II/lac/filtered_matrix.h>
Classes | |
class | Accessor |
class | const_iterator |
struct | PairComparison |
Public Types | |
using | size_type = types::global_dof_index |
using | IndexValuePair = std::pair< size_type, double > |
Public Member Functions | |
std::size_t | memory_consumption () const |
Constructors and initialization | |
FilteredMatrix () | |
FilteredMatrix (const FilteredMatrix &fm) | |
template<typename MatrixType > | |
FilteredMatrix (const MatrixType &matrix, const bool expect_constrained_source=false) | |
FilteredMatrix & | operator= (const FilteredMatrix &fm) |
template<typename MatrixType > | |
void | initialize (const MatrixType &m, const bool expect_constrained_source=false) |
void | clear () |
Managing constraints | |
void | add_constraint (const size_type i, const double v) |
template<class ConstraintList > | |
void | add_constraints (const ConstraintList &new_constraints) |
void | clear_constraints () |
void | apply_constraints (VectorType &v, const bool matrix_is_symmetric) const |
void | apply_constraints (VectorType &v) const |
void | vmult (VectorType &dst, const VectorType &src) const |
void | Tvmult (VectorType &dst, const VectorType &src) const |
void | vmult_add (VectorType &dst, const VectorType &src) const |
void | Tvmult_add (VectorType &dst, const VectorType &src) const |
Iterators | |
const_iterator | begin () const |
const_iterator | end () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
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 |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Types | |
using | const_index_value_iterator = typename std::vector< IndexValuePair >::const_iterator |
Private Member Functions | |
void | pre_filter (VectorType &v) const |
void | post_filter (const VectorType &in, VectorType &out) const |
Private Attributes | |
bool | expect_constrained_source |
std::shared_ptr< PointerMatrixBase< VectorType > > | matrix |
std::vector< IndexValuePair > | constraints |
Friends | |
class | FilteredMatrixBlock< VectorType > |
Additional Inherited Members | |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
This class is a wrapper for linear systems of equations with simple equality constraints fixing individual degrees of freedom to a certain value such as when using Dirichlet boundary values.
In order to accomplish this, the vmult(), Tvmult(), vmult_add() and Tvmult_add functions modify the same function of the original matrix such as if all constrained entries of the source vector were zero. Additionally, all constrained entries of the destination vector are set to zero.
Usage is simple: create an object of this type, point it to a matrix that shall be used for \(A\) above (either through the constructor, the copy constructor, or the set_referenced_matrix() function), specify the list of boundary values or other constraints (through the add_constraints() function), and then for each required solution modify the right hand side vector (through apply_constraints()) and use this object as matrix object in a linear solver. As linear solvers should only use vmult() and residual() functions of a matrix class, this class should be as good a matrix as any other for that purpose.
Furthermore, also the precondition_Jacobi() function is provided (since the computation of diagonal elements of the filtered matrix \(A_X\) is simple), so you can use this as a preconditioner. Some other functions useful for matrices are also available.
A typical code snippet showing the above steps is as follows:
The function MatrixTools::apply_boundary_values() does exactly the same that this class does, except for the fact that that function actually modifies the matrix. Consequently, it is only possible to solve with a matrix to which MatrixTools::apply_boundary_values() was applied for one right hand side and one set of boundary values since the modification of the right hand side depends on the original matrix.
While this is a feasible method in cases where only one solution of the linear system is required, for example in solving linear stationary systems, one would often like to have the ability to solve multiple times with the same matrix in nonlinear problems (where one often does not want to update the Hessian between Newton steps, despite having different right hand sides in subsequent steps) or time dependent problems, without having to re-assemble the matrix or copy it to temporary matrices with which one then can work. For these cases, this class is meant.
Mathematically speaking, it is used to represent a system of linear equations \(Ax=b\) with the constraint that \(B_D x = g_D\), where \(B_D\) is a rectangular matrix with exactly one \(1\) in each row, and these \(1\)s in those columns representing constrained degrees of freedom (e.g. for Dirichlet boundary nodes, thus the index \(D\)) and zeroes for all other diagonal entries, and \(g_D\) having the requested nodal values for these constrained nodes. Thus, the underdetermined equation \(B_D x = g_D\) fixes only the constrained nodes and does not impose any condition on the others. We note that \(B_D B_D^T = 1_D\), where \(1_D\) is the identity matrix with dimension as large as the number of constrained degrees of freedom. Likewise, \(B_D^T B_D\) is the diagonal matrix with diagonal entries \(0\) or \(1\) that, when applied to a vector, leaves all constrained nodes untouched and deletes all unconstrained ones.
For solving such a system of equations, we first write down the Lagrangian \(L=1/2 x^T A x - x^T b + l^T B_D x\), where \(l\) is a Lagrange multiplier for the constraints. The stationarity condition then reads
The first equation then reads \(B_D^T l = b-Ax\). On the other hand, if we left-multiply the first equation by \(B_D^T B_D\), we obtain \(B_D^T B_D A x + B_D^T l = B_D^T B_D b\) after equating \(B_D B_D^T\) to the identity matrix. Inserting the previous equality, this yields \((A - B_D^T B_D A) x = (1 - B_D^T B_D)b\). Since \(x=(1 - B_D^T B_D) x + B_D^T B_D x = (1 - B_D^T B_D) x + B_D^T g_D\), we can restate the linear system: \(A_D x = (1 - B_D^T B_D)b - (1 - B_D^T B_D) A B^T g_D\), where \(A_D = (1 - B_D^T B_D) A (1 - B_D^T B_D)\) is the matrix where all rows and columns corresponding to constrained nodes have been deleted.
The last system of equation only defines the value of the unconstrained nodes, while the constrained ones are determined by the equation \(B_D x = g_D\). We can combine these two linear systems by using the zeroed out rows of \(A_D\): if we set the diagonal to \(1\) and the corresponding zeroed out element of the right hand side to that of \(g_D\), then this fixes the constrained elements as well. We can write this as follows: \(A_X x = (1 - B_D^T B_D)b - (1 - B_D^T B_D) A B^T g_D + B_D^T g_D\), where \(A_X = A_D + B_D^T B_D\). Note that the two parts of the latter matrix operate on disjoint subspaces (the first on the unconstrained nodes, the latter on the constrained ones).
In iterative solvers, it is not actually necessary to compute \(A_X\) explicitly, since only matrix-vector operations need to be performed. This can be done in a three-step procedure that first clears all elements in the incoming vector that belong to constrained nodes, then performs the product with the matrix \(A\), then clears again. This class is a wrapper to this procedure, it takes a pointer to a matrix with which to perform matrix- vector products, and does the cleaning of constrained elements itself. This class therefore implements an overloaded vmult
function that does the matrix-vector product, as well as Tvmult
for transpose matrix-vector multiplication and residual
for residual computation, and can thus be used as a matrix replacement in linear solvers.
It also has the ability to generate the modification of the right hand side, through the apply_constraints() function.
This class takes as template arguments a matrix and a vector class. The former must provide vmult
, vmult_add
, Tvmult
, and residual
member function that operate on the vector type (the second template argument). The latter template parameter must provide access to individual elements through operator()
, assignment through operator=
.
The functions that operate as a matrix and do not change the internal state of this object are synchronized and thus threadsafe. Consequently, you do not need to serialize calls to vmult
or residual
.
Definition at line 199 of file filtered_matrix.h.
using FilteredMatrix< VectorType >::size_type = types::global_dof_index |
Declare the type of container size.
Definition at line 207 of file filtered_matrix.h.
using FilteredMatrix< VectorType >::IndexValuePair = std::pair<size_type, double> |
Typedef defining a type that represents a pair of degree of freedom index and the value it shall have.
Definition at line 332 of file filtered_matrix.h.
|
private |
Declare an abbreviation for an iterator into the array constraint pairs, since that data type is so often used and is rather awkward to write out each time.
Definition at line 524 of file filtered_matrix.h.
|
inline |
Default constructor. You will have to set the matrix to be used later using initialize().
Definition at line 723 of file filtered_matrix.h.
|
inline |
Copy constructor. Use the matrix and the constraints set in the given object for the present one as well.
Definition at line 730 of file filtered_matrix.h.
|
inline |
Constructor. Use the given matrix for future operations.
m:
The matrix being used in multiplications.expect_constrained_source:
See documentation of expect_constrained_source. Definition at line 741 of file filtered_matrix.h.
|
inline |
Copy operator. Take over matrix and constraints from the other object.
Definition at line 752 of file filtered_matrix.h.
|
inline |
Set the matrix to be used further on. You will probably also want to call the clear_constraints() function if constraints were previously added.
m:
The matrix being used in multiplications.expect_constrained_source:
See documentation of expect_constrained_source. Definition at line 713 of file filtered_matrix.h.
|
inline |
Delete all constraints and the matrix pointer.
Definition at line 808 of file filtered_matrix.h.
|
inline |
Add the constraint that the value with index i
should have the value v
.
Definition at line 764 of file filtered_matrix.h.
|
inline |
Add a list of constraints to the ones already managed by this object. The actual data type of this list must be so that dereferenced iterators are pairs of indices and the corresponding values to be enforced on the respective solution vector's entry. Thus, the data type might be, for example, a std::list
or std::vector
of IndexValuePair objects, but also a std::map<unsigned, double>
.
The second component of these pairs will only be used in apply_constraints(). The first is used to set values to zero in matrix vector multiplications.
It is an error if the argument contains an entry for a degree of freedom that has already been constrained previously.
Definition at line 776 of file filtered_matrix.h.
|
inline |
Delete the list of constraints presently in use.
Definition at line 797 of file filtered_matrix.h.
|
inline |
Vector operations Apply the constraints to a right hand side vector. This needs to be done before starting to solve with the filtered matrix. If the matrix is symmetric (i.e. the matrix itself, not only its sparsity pattern), set the second parameter to true
to use a faster algorithm. Note: This method is deprecated as matrix_is_symmetric parameter is no longer used.
Definition at line 818 of file filtered_matrix.h.
|
inline |
Apply the constraints to a right hand side vector. This needs to be done before starting to solve with the filtered matrix.
Definition at line 828 of file filtered_matrix.h.
|
inline |
Matrix-vector multiplication: this operation performs pre_filter(), multiplication with the stored matrix and post_filter() in that order.
Definition at line 889 of file filtered_matrix.h.
|
inline |
Matrix-vector multiplication: this operation performs pre_filter(), transposed multiplication with the stored matrix and post_filter() in that order.
Definition at line 916 of file filtered_matrix.h.
|
inline |
Adding matrix-vector multiplication.
dst
. We expect that in most cases this is the required behavior. Definition at line 943 of file filtered_matrix.h.
|
inline |
Adding transpose matrix-vector multiplication:
dst
. We expect that in most cases this is the required behavior. Definition at line 971 of file filtered_matrix.h.
|
inline |
Iterator to the first constraint.
Definition at line 686 of file filtered_matrix.h.
|
inline |
Final iterator.
Definition at line 694 of file filtered_matrix.h.
|
inline |
Determine an estimate for the memory consumption (in bytes) of this object. Since we are not the owner of the matrix referenced, its memory consumption is not included.
Definition at line 999 of file filtered_matrix.h.
|
inlineprivate |
Do the pre-filtering step, i.e. zero out those components that belong to constrained degrees of freedom.
Definition at line 857 of file filtered_matrix.h.
|
inlineprivate |
Do the postfiltering step, i.e. set constrained degrees of freedom to the value of the input vector, as the matrix contains only ones on the diagonal for these degrees of freedom.
Definition at line 871 of file filtered_matrix.h.
|
friend |
FilteredMatrixBlock accesses pre_filter() and post_filter().
Definition at line 569 of file filtered_matrix.h.
|
private |
Determine, whether multiplications can expect that the source vector has all constrained entries set to zero.
If so, the auxiliary vector can be avoided and memory as well as time can be saved.
We expect this for instance in Newton's method, where the residual already should be zero on constrained nodes. This is, because there is no test function in these nodes.
Definition at line 516 of file filtered_matrix.h.
|
private |
Pointer to the sparsity pattern used for this matrix.
Definition at line 542 of file filtered_matrix.h.
|
private |
Sorted list of pairs denoting the index of the variable and the value to which it shall be fixed.
Definition at line 548 of file filtered_matrix.h.