G+Smo  24.08.0
Geometry + Simulation Modules
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
gsIetiSystem< T > Class Template Reference

Detailed Description

template<typename T>
class gismo::gsIetiSystem< T >

This class represents a IETI system.

The IETI saddle point system is a system of the form:

\[ \begin{pmatrix} \tilde A_1 & & & & \tilde B_1^\top \\ & \tilde A_2 & & & \tilde B_2^\top \\ & & \ddots & & \vdots \\ & & & \tilde A_N & \tilde B_N^\top \\ \tilde B_1 & \tilde B_2 & \cdots & \tilde B_N & 0 \\ \end{pmatrix} \]

The corresponding Schur complement is

\[ \sum_{k=1}^N \tilde B_k \tilde A_k^{-1} \tilde B_k^\top \]

For a standard IETI-DP setup, \( \tilde A_k \) and \( \tilde B_k \) are obtained from the original matrices \( A_k \) and \( B_k \) by eliminating the primal dofs or by incorporating a constraint that sets them to zero.

This class does not have any special treatment for the primal problem of a IETI-DP solver. Thus, the primal problem is just another subdomain and in case of IETI-DP, we have N=K+1, where K is the number of patches.

The matrices \( \tilde A_k \) are stored in a vector that can be accessed (read and write) via localMatrixOp. To allow certain matrix-free variants, each of them is stored as gsLinearOperator.

The inverses \( \tilde A_k^{-1} \) are stored as vector accessable via localSolverOp. If the matrices \( \tilde A_k \) are stored as gsMatrixOp< gsSparseMatrix<T> >, LU solvers will be automatically created on the fly if needed. Otherwise or if the caller wants other local solvers (like inexact ones), the vector can be populated by the caller.

The matrices \( \tilde B_k \) are stored in a vector accessible via jumpMatrix.

The right-hand sides are stored in a vector accessible via localRhs.

Public Member Functions

void addSubdomain (JumpMatrixPtr jumpMatrix, OpPtr localMatrixOp, Matrix localRhs, OpPtr localSolverOp=OpPtr())
 
std::vector< MatrixconstructSolutionFromLagrangeMultipliers (const Matrix &multipliers) const
 Returns the local solutions for the individual subdomains. More...
 
std::vector< MatrixconstructSolutionFromSaddlePoint (const Matrix &x) const
 Returns the local solutions for the individual subdomains. More...
 
JumpMatrixPtrjumpMatrix (index_t k)
 Access the jump matrix.
 
OpPtrlocalMatrixOp (index_t k)
 Access the local stiffness matrix (as gsLinearOperator)
 
MatrixlocalRhs (index_t k)
 Access the local right-hand side.
 
OpPtrlocalSolverOp (index_t k)
 Access the local solver operator.
 
index_t nLagrangeMultipliers () const
 Returns the number of Lagrange multipliers. More...
 
void reserve (index_t n)
 Reserves the memory required to store the number of subdomains. More...
 
Matrix rhsForSaddlePoint () const
 Returns the right-hand-side that is required for the saddle point formulation of the IETI problem.
 
Matrix rhsForSchurComplement () const
 Returns the right-hand-side that is required for the Schur complement formulation of the IETI problem. More...
 
OpPtr saddlePointProblem () const
 Returns gsLinearOperator that represents the IETI problem as saddle point problem.
 
OpPtr schurComplement () const
 Returns gsLinearOperator that represents the Schur complement for the IETI problem. More...
 

Private Types

typedef gsSparseMatrix< T,
RowMajor > 
JumpMatrix
 Sparse matrix type for jumps.
 
typedef memory::shared_ptr
< JumpMatrix
JumpMatrixPtr
 Shared pointer to sparse matrix type for jumps.
 
typedef gsMatrix< T > Matrix
 Matrix type.
 
typedef gsLinearOperator< T > Op
 Linear operator.
 
typedef memory::shared_ptr< OpOpPtr
 Shared pointer to linear operator.
 
typedef gsSparseMatrix< T > SparseMatrix
 Sparse matrix type.
 
typedef gsMatrixOp< SparseMatrixSparseMatrixOp
 Sparse matrix type as operatpr.
 

Private Member Functions

void setupSparseLUSolvers () const
 Setup solvers if not provided by user.
 

Private Attributes

std::vector< JumpMatrixPtrm_jumpMatrices
 Stores the jump matrices.
 
std::vector< OpPtrm_localMatrixOps
 Stores the local matrix ops \( \tilde A_k \).
 
std::vector< Matrixm_localRhs
 Stores the local right-hand sides.
 
std::vector< OpPtrm_localSolverOps
 Stores the local solvers.
 

Member Function Documentation

void addSubdomain ( JumpMatrixPtr  jumpMatrix,
OpPtr  localMatrixOp,
Matrix  localRhs,
OpPtr  localSolverOp = OpPtr() 
)

Adds a new subdomain

Subdomain might be, e.g., a patch-local problem or the primal problem

Parameters
jumpMatrixThe associated jump matrix
localMatrixOpThe operator that represents the local stiffness matrix
localRhsThe contribution to the right-hand side
localSolverOpThe operator that represents a solver for the local problem. This parameter is optional; the solver is created automatically if needed.
std::vector< gsMatrix< T > > constructSolutionFromLagrangeMultipliers ( const Matrix multipliers) const

Returns the local solutions for the individual subdomains.

Parameters
multipliersThe Lagrange multipliers previously computed (based on the Schur complement form)

If the local solvers have not been provided, this function will generate them from the localMatrixOp if they are gsMatrixOp<gsSparseMatrix<T>>

std::vector< gsMatrix< T > > constructSolutionFromSaddlePoint ( const Matrix x) const

Returns the local solutions for the individual subdomains.

Parameters
xThe solution computed using the saddle point system
index_t nLagrangeMultipliers ( ) const
inline

Returns the number of Lagrange multipliers.

This requires that at least one jump matrix has been set.

void reserve ( index_t  n)

Reserves the memory required to store the number of subdomains.

Parameters
nNumber of subdomains
gsMatrix< T > rhsForSchurComplement ( ) const

Returns the right-hand-side that is required for the Schur complement formulation of the IETI problem.

If the local solvers have not been provided, this function will generate them from the localMatrixOp if they are gsMatrixOp<gsSparseMatrix<T>>

gsIetiSystem< T >::OpPtr schurComplement ( ) const

Returns gsLinearOperator that represents the Schur complement for the IETI problem.

If the local solvers have not been provided, this function will generate them from the localMatrixOp if they are gsMatrixOp<gsSparseMatrix<T>>