Reference documentation for deal.II version 8.5.1
|
#include <deal.II/lac/schur_matrix.h>
Public Member Functions | |
SchurMatrix (const MA_inverse &Ainv, const MB &B, const MDt &Dt, const MC &C, VectorMemory< BlockVector< double > > &mem, const std::vector< types::global_dof_index > &signature=std::vector< types::global_dof_index >(0)) | |
void | prepare_rhs (BlockVector< double > &dst, const BlockVector< double > &src) const |
void | vmult (BlockVector< double > &dst, const BlockVector< double > &src) const |
double | residual (BlockVector< double > &dst, const BlockVector< double > &src, const BlockVector< double > &rhs) const |
void | postprocess (BlockVector< double > &dst, const BlockVector< double > &src, const BlockVector< double > &rhs) const |
void | debug_level (unsigned int l) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) |
void | subscribe (const char *identifier=0) const |
void | unsubscribe (const char *identifier=0) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Member Functions | |
SchurMatrix (const SchurMatrix< MA_inverse, MB, MDt, MC > &) | |
SchurMatrix & | operator= (const SchurMatrix< MA_inverse, MB, MDt, MC > &) |
Private Attributes | |
const SmartPointer< const MA_inverse, SchurMatrix< MA_inverse, MB, MDt, MC > > | Ainv |
const SmartPointer< const MB, SchurMatrix< MA_inverse, MB, MDt, MC > > | B |
const SmartPointer< const MDt, SchurMatrix< MA_inverse, MB, MDt, MC > > | Dt |
const SmartPointer< const MC, SchurMatrix< MA_inverse, MB, MDt, MC > > | C |
VectorMemory< BlockVector< double > > & | mem |
std::vector< types::global_dof_index > | signature |
unsigned int | debug |
Additional Inherited Members | |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, char *arg2, std::string &arg3) |
static ::ExceptionBase & | ExcNoSubscriber (char *arg1, char *arg2) |
Schur complement of a block matrix.
Given a non-singular matrix A
(often positive definite) and a positive semi-definite matrix C
as well as matrices B
and Dt
of full rank, this class implements a new matrix, the Schur complement a the system of equations of the structure
* / \ / \ / \ * | A Dt | | u | - | f | * | -B C | | p | - | g | * \ / \ / \ / *
Multiplication with the Schur matrix S
is the operation
* S p = C p + B A-inverse Dt-transpose p, *
which is an operation within the space for p
.
The data handed to the Schur matrix are as follows:
A:
the inverse of A
is stored, instead of A
. This allows the application to use the most efficient form of inversion, iterative or direct.
B
, C:
these matrices are stored "as is".
Dt:
the computation of the Schur complement involves the function Tvmult
of the matrix Dt
, not vmult!
This way, it is sufficient to build only one matrix B
for the symmetric Schur complement and use it twice.
All matrices involved are of arbitrary type and vectors are BlockVectors. This way, SchurMatrix
can be coupled with any matrix classes providing vmult
and Tvmult
and can be even nested. Since SmartPointers of matrices are stored, the matrix blocks should be derived from Subscriptor.
Since the Schur complement of a matrix corresponds to a Gaussian block elimination, the right hand side of the condensed system must be preprocessed. Furthermore, the eliminated variable must be reconstructed after solving.
* g = g + B A-inverse f * u = A-inverse (f - D-transpose p) *
Applying these transformations, the solution of the system above by a SchurMatrix
schur
is coded as follows:
Definition at line 97 of file schur_matrix.h.
SchurMatrix< MA_inverse, MB, MDt, MC >::SchurMatrix | ( | const MA_inverse & | Ainv, |
const MB & | B, | ||
const MDt & | Dt, | ||
const MC & | C, | ||
VectorMemory< BlockVector< double > > & | mem, | ||
const std::vector< types::global_dof_index > & | signature = std::vector<types::global_dof_index>(0) |
||
) |
Constructor. This constructor receives all the matrices needed. Furthermore, it gets a reference to a memory structure for obtaining block vectors.
Optionally, the length of the u-vector
can be provided.
For the meaning of the matrices see the class documentation.
Definition at line 206 of file schur_matrix.h.
|
private |
No copy constructor.
void SchurMatrix< MA_inverse, MB, MDt, MC >::prepare_rhs | ( | BlockVector< double > & | dst, |
const BlockVector< double > & | src | ||
) | const |
Do block elimination of the right hand side. Given right hand sides for both components of the block system, this function provides the right hand side for the Schur complement.
The result is stored in the first argument, which is also part of the input data. If it is necessary to conserve the data, dst
must be copied before calling this function. This is reasonable, since in many cases, only the pre-processed right hand side is needed.
Definition at line 282 of file schur_matrix.h.
void SchurMatrix< MA_inverse, MB, MDt, MC >::vmult | ( | BlockVector< double > & | dst, |
const BlockVector< double > & | src | ||
) | const |
Multiplication with the Schur complement.
Definition at line 231 of file schur_matrix.h.
double SchurMatrix< MA_inverse, MB, MDt, MC >::residual | ( | BlockVector< double > & | dst, |
const BlockVector< double > & | src, | ||
const BlockVector< double > & | rhs | ||
) | const |
Computation of the residual of the Schur complement.
Definition at line 269 of file schur_matrix.h.
void SchurMatrix< MA_inverse, MB, MDt, MC >::postprocess | ( | BlockVector< double > & | dst, |
const BlockVector< double > & | src, | ||
const BlockVector< double > & | rhs | ||
) | const |
Compute the eliminated variable from the solution of the Schur complement problem.
Definition at line 311 of file schur_matrix.h.
void SchurMatrix< MA_inverse, MB, MDt, MC >::debug_level | ( | unsigned int | l | ) |
Select debugging information for log-file. Debug level 1 is defined and writes the norm of every vector before and after each operation. Debug level 0 turns off debugging information.
Definition at line 223 of file schur_matrix.h.
|
private |
No assignment.
|
private |
Pointer to inverse of upper left block.
Definition at line 172 of file schur_matrix.h.
|
private |
Pointer to lower left block.
Definition at line 176 of file schur_matrix.h.
|
private |
Pointer to transpose of upper right block.
Definition at line 180 of file schur_matrix.h.
|
private |
Pointer to lower right block.
Definition at line 184 of file schur_matrix.h.
|
private |
Auxiliary memory for vectors.
Definition at line 188 of file schur_matrix.h.
|
private |
Optional signature of the u-vector
.
Definition at line 193 of file schur_matrix.h.
|
private |
Switch for debugging information.
Definition at line 198 of file schur_matrix.h.