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

Detailed Description

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

This class represents the scaled Dirichlet preconditioner for a IETI problem.

Its formal representation is

\[ \sum_{k=1}^K \hat B_k D_k^{-1} S_k D_k^{-1} \hat B_k^\top \]

It is a preconditioner for the Schur complement of the IETI system (as represented by gsIetiSystem)

\[ \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} \]

For a standard IETI-dp setup, we additionally have a primal problem, thus N=K+1. In this case, the matrices \( \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 is done by gsPrimalSystem.

The matrices \( S_k \) are stored in a vector accesible via localSchurOp. As usual, they are stored in form of a vector of gsLinearOperator s. These operators represent the Schur-complements of the matrices \( A_k \) with respect to the degrees of freedom on the skeleton.

The jump matrices \( \hat B_k \) are accessible via jumpMatrix. These matrices usually differ from the matrices \( \tilde B_k \) from the IETI-system since – for the preconditioner – the jump matrices have to be restricted to the skeleton.

If the matrices \( A_k \) and \( B_k \) are given, the function restrictToSkeleton allows to compute the corresponding matrices \( S_k \) and \( \hat B_k \). The degrees of freedom belonging to the skeleton can be specified by the caller. The caller can use the function getSkeletonDofs to extract this information from the jump matrices, i.e., skeleton dofs are those that are affected by a Lagrange multiplier. (Alternatively, the caller might use the corresponding function from the class gsIetiMapper, which uses the gsDofMapper s and might yield different results.)

The scaling matrices \( D_k \) are stored in a vector accessible via scalingMatrix. They can be provided by the caller or generated by calling setupMultiplicityScaling.

Classes

struct  Blocks
 

Public Member Functions

void addSubdomain (JumpMatrixPtr jumpMatrix, OpPtr localSchurOp)
 
JumpMatrixPtrjumpMatrix (index_t k)
 Access the jump matrix.
 
MatrixlocalScaling (index_t k)
 Access the local scaling matrix (as row vector)
 
OpPtrlocalSchurOps (index_t k)
 Access the local Schur complements operator.
 
index_t nLagrangeMultipliers () const
 Returns the number of Lagrange multipliers. More...
 
OpPtr preconditioner () const
 This returns the preconditioner as gsLinearOperator. More...
 
void reserve (index_t n)
 Reserves the memory required to store the given number of subdomain. More...
 
void setupMultiplicityScaling ()
 This sets up the member vector localScaling based on multiplicity scaling. More...
 

Static Public Member Functions

static Blocks matrixBlocks (const SparseMatrix &localMatrix, const std::vector< index_t > &dofs)
 Computes the matrix blocks with respect to the given dofs. More...
 
static JumpMatrix restrictJumpMatrix (const JumpMatrix &jumpMatrix, const std::vector< index_t > &dofs)
 Restricts the jump matrix to the given dofs. More...
 
static std::pair< JumpMatrix,
OpPtr
restrictToSkeleton (const JumpMatrix &jumpMatrix, const SparseMatrix &localMatrix, const std::vector< index_t > &dofs)
 
static OpPtr schurComplement (Blocks matrixBlocks, OpPtr solver)
 Computes the Schur complement with respect to the block A11 of matrixBlocks. More...
 
static OpPtr schurComplement (const SparseMatrix &localMatrix, const std::vector< index_t > &dofs)
 Computes the Schur complement of the matrix with respect to the given dofs using a sparse Cholesky solver. More...
 
static gsSortedVector< index_tskeletonDofs (const JumpMatrix &jumpMatrix)
 Extracts the skeleton dofs from the jump matrix. More...
 

Public Attributes

std::vector< JumpMatrixPtrm_jumpMatrices
 The jump matrices \( \hat B_k \).
 
std::vector< Matrixm_localScaling
 The diagonal entries of \( D_k \) as vectors.
 
std::vector< OpPtrm_localSchurOps
 The local Schur complements \( S_k \).
 

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.
 

Member Function Documentation

void addSubdomain ( JumpMatrixPtr  jumpMatrix,
OpPtr  localSchurOp 
)
inline

Adds a new subdomain

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

Parameters
jumpMatrixThe associated jump matrix
localSchurOpThe operator that represents the local Schur complement
Note
These two parameters can also be provided as std::pair as provided by restrictToSkeleton
gsScaledDirichletPrec< T >::Blocks matrixBlocks ( const SparseMatrix localMatrix,
const std::vector< index_t > &  dofs 
)
static

Computes the matrix blocks with respect to the given dofs.

Parameters
localMatrixThe local stiffness matrix
dofsThe corresponding degrees of freedom (usually skeleton dofs)

If 0 corresponds to the list of dofs and 1 remains to the others, this function returns the blocks A00, A10, A01, A11 of A

index_t nLagrangeMultipliers ( ) const
inline

Returns the number of Lagrange multipliers.

This requires that at least one subdomain was defined.

gsScaledDirichletPrec< T >::OpPtr preconditioner ( ) const

This returns the preconditioner as gsLinearOperator.

This requires that the subdomains have been defined first.

void reserve ( index_t  n)
inline

Reserves the memory required to store the given number of subdomain.

Parameters
nNumber of subdomains
gsScaledDirichletPrec< T >::JumpMatrix restrictJumpMatrix ( const JumpMatrix jumpMatrix,
const std::vector< index_t > &  dofs 
)
static

Restricts the jump matrix to the given dofs.

Parameters
jumpMatrixThe jump matrix
dofsThe corresponding degrees of freedom (usually skeleton dofs)
static std::pair<JumpMatrix,OpPtr> restrictToSkeleton ( const JumpMatrix jumpMatrix,
const SparseMatrix localMatrix,
const std::vector< index_t > &  dofs 
)
inlinestatic

Restricts the jump matrix and the local stiffness matrix to the skeleton

Parameters
jumpMatrixThe jump matrix
localMatrixThe local stiffness matrix
dofsThe degrees of freedom on the skeleton

The implementation is basically:

return std::pair<gsSparseMatrix<T,RowMajor>,typename gsLinearOperator<T>::Ptr>(
);
gsScaledDirichletPrec< T >::OpPtr schurComplement ( Blocks  matrixBlocks,
OpPtr  solver 
)
static

Computes the Schur complement with respect to the block A11 of matrixBlocks.

Parameters
matrixBlocksThe blocks of the stiffess matrix as returned by matrixBlocks
solverA gsLinearOperator that realizes the inverse of matrixBlocks.A11
static OpPtr schurComplement ( const SparseMatrix localMatrix,
const std::vector< index_t > &  dofs 
)
inlinestatic

Computes the Schur complement of the matrix with respect to the given dofs using a sparse Cholesky solver.

Parameters
localMatrixThe local stiffness matrix
dofsThe degrees of freedom for which the Schur complement is taken

The implementation is basically:

auto blocks = gsScaledDirichletPrec<T>::matrixBlocks(mat, dofs);
return gsScaledDirichletPrec<T>::schurComplement( blocks, makeSparseCholeskySolver(blocks.A11) );
void setupMultiplicityScaling ( )

This sets up the member vector localScaling based on multiplicity scaling.

This requires that the subdomains have been defined first.

gsSortedVector< index_t > skeletonDofs ( const JumpMatrix jumpMatrix)
static

Extracts the skeleton dofs from the jump matrix.

Parameters
jumpMatrixThe jump matrix

This means that a dof is considered to be on the skeleton iff at least one Lagrange multiplier acts on it. This might lead to other results than the function that is provided by gsIetiMapper.