template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
class Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >
A class representing small KKT systems on a single subdomain.
This is intended to work together with BDDCSolver and KKTSolver, and represents subdomains, providing subdomain-local storage and processing capabilities. This is the basis for distributed domain decomposition solvers.
A Subdomain object represents the system
\[
B = \begin{bmatrix} A & C^T \\ C \end{bmatrix},
\]
and provides in particular access to the "Schur complement"
\[
S = \begin{bmatrix} 0 & I \end{bmatrix}
\begin{bmatrix} A & C^T \\ C \end{bmatrix}^{-1}
\begin{bmatrix} 0 \\ I \end{bmatrix}.
\]
Since \( A \) is only positive semidefinite in general, the true Schur complement \( - C A^{-1} C^T \) cannot be formed, but is equal in case of invertible \( A \).
The subdomain holds \( A_i \), \( f_i \), and the constraints block \( C \). It also holds a local solution vector \( u_i \). For the overall picture, we refer to KKTSolver.
- Template Parameters
-
| m | the number of vectorial components in \( x \) (usually 1 for scalar systems, and 2 or 3 for vectorial ones) |
| Scalar | the field type \( K \), usually double |
| CoarseScalar | the field type \( K_c \) to be used for representing the coarse solver factorization, usually double, but could also be float |
| Transfer | a space transfer policy class |
- Todo:
- UMFPACK does not provide
float factorization. Extract L and U factors and convert to float?
This is a reference implementation.
Definition at line 179 of file bddc.hh.
|
| std::vector< int > | neighbors () const |
| | Reports a list of its neighbor subdomains. More...
|
| |
| void | updateSolution (Scalar alpha) |
| | Updates the solution approximation. More...
|
| |
| void | getSolution (XVector &x) const |
| | Provides the current solution vector \( x \). More...
|
| |
| XVector | getSolution () const |
| | Provides the current solution vector \( x \). More...
|
| |
| void | set_solver_type (bool cg_solver) |
| | Set solver type. More...
|
| |
| void | updateIteration (int cg_itr) |
| |
| void | getCorrection (XVector &dx) const |
| | Provides the current correction vector \( dx \). More...
|
| |
| void | getRawCorrection (XVector &dx) const |
| |
| void | getResidual (XVector &rr) const |
| |
| Vector const & | getRestrictedResidual () const |
| |
| std::vector< bool > | getActiveSubdomains () |
| |
| void | setActiveSubdomains_sub (std::vector< bool > &activeSubdomains_) |
| |
|
| template<int cComponents> |
| | Subdomain (int id, AMatrix const &A, Interfaces const &interfaces, NumaBCRSMatrix< Dune::FieldMatrix< Scalar, cComponents, components > > const &C) |
| | Constructor. More...
|
| |
| template<class Interface > |
| | Subdomain (int id, AMatrix const &A, Interface const &ifs) |
| |
| void | imposeDC_inActiveInterface (Vector &Aq, Vector &q) |
| |
| void | imposeDCq (Vector &q) |
| |
| void | update_rhs (XVector const &f) |
| | update rhs for the current subdomain More...
|
| |
| void | getInterfaceAveragingWeights (int s, XVector &rho) const |
| | Reports own interface averaging weights. More...
|
| |
| void | setInterfaceAveragingWeights (int s, XVector const &rho) |
| | Accumulates interface weights of adjacent subdomains. More...
|
| |
|
These methods implement the "restriction" of the residual to the "coarse" space.
The restriction is defined as the transpose of the prolongation, i.e. \( \Pi^T = \pi^T E^T \) with averaging \( \pi \) and extension \( E \).
The extension operator \( E^T \) acts locally, realizing the dual harmonic extension
\[
\begin{bmatrix} r_i \\ r_D \\ \bar r_D \end{bmatrix}_{\rm new}
= \begin{bmatrix} I & \\
H^T & \\
-H^T & I \end{bmatrix}
\begin{bmatrix} r_i \\ r_D \end{bmatrix}
\]
working on interior residual \( r_i \), interface residual \( r_D \), and averaged interface residual \( \bar r_D \). Here, \( H = A_{ii}^{-1} A_{iD} \) denotes the harmonic extension of interface values into the interior. This local operation is performed by the preRestrict() method.
The averaging (for a single degree of freedom shared by \( k \) subdomains) is
\[ \pi^T = \frac{1}{k} \begin{bmatrix} 1 & \dots & 1 \\ \vdots & \ddots & \vdots \\
1 & \dots & 1 \end{bmatrix}, \]
therefore symmetric, and acts on the interface residual values \( r_D \) only. It involves communication between subdomains, exchanging their interface residual \( r_D \), which is realized by the getInterfaceResidual() and setInterfaceResidual() methods. Finally, restrict() must be called as a sentinel to denote the end of communication, and to perform the actual averaging.
|
| void | preRestrict () |
| | First step of the restriction: Dirichlet solve. More...
|
| |
| void | getInterfaceResidualUpdate (int s, TransmissionType &dres) const |
| | Gets interface residual that shall be shared with neighboring subdomain. More...
|
| |
| void | setInterfaceResidualUpdate (int s, TransmissionType const &dres) |
| | Sets interface residual of the correction from neighboring subdomains. More...
|
| |
| void | restrict () |
| | Performs the restriction of the residual vector. More...
|
| |
|
These methods implement the subdomain-local parts of the coarse space solver.
|
| void | getS (SMatrix &S) const |
| | Compute \( S \). More...
|
| |
| void | getSchurResidual (Vector &lambda) const |
| | Provides the Schur residual. More...
|
| |
| void | setCorrection (Vector const &c) |
| | Computes the correction vector \( [\delta x, \delta\lambda] \). More...
|
| |
|
These methods implement the "prolongation" of "coarse" space corrections to the "fine" space.
The prolongation is defined as \( \Pi = E\pi \) with averaging \( \pi \) and extension \( E \). The averaging (for a single degree of freedom shared by \( k \) subdomains) is
\[ \pi = \frac{1}{k} \begin{bmatrix} 1 & \dots & 1 \\ \vdots & \ddots & \vdots \\
1 & \dots & 1 \end{bmatrix} \]
and acts on the interface degrees of freedom \( x_D \) only. It involves communication between subdomains, exchanging their interface values \( x_D \), which is realized by the getInterfaceValues() and setInterfaceValues() methods.
The extension operator \( E \) acts locally, realizing the harmonic extension
\[
\begin{bmatrix} x_i \\ x_D \end{bmatrix}_{\rm new}
= \begin{bmatrix} I & H & - H \\
& & I \end{bmatrix}
\begin{bmatrix} x_i \\ x_D \\ \bar x_D \end{bmatrix}
\]
working on interior dofs \( x_i \), interface dofs \( x_D \), and averaged interface values \( \bar x_D \). Here, \( H = A_{ii}^{-1} A_{iD} \) denotes the harmonic extension of interface values into the interior. This local operation is performed by the prolongate() method.
|
| void | preProlongate () |
| | Prepares the exchange of interface values. More...
|
| |
| void | getInterfaceValues (int s, TransmissionType &values) const |
| | Gets boundary values that shall be shared with neighboring subdomain. More...
|
| |
| void | setInterfaceValues (int s, TransmissionType const &values) |
| | Sets boundary values of the correction from neighboring subdomains. More...
|
| |
| void | setActiveId (bool active) |
| |
| Dune::FieldVector< double, 2 > | prolongate () |
| | Perform averaging on the shared interfaces, projecting the solution to the subspace of functions continuous across the interfaces. More...
|
| |
| Dune::FieldVector< double, 1 > | updateSearchDirectionNorm () |
| |
| void | updateSearchDirection (double beta) |
| | Update search direction q for cg solver. More...
|
| |