KASKADE 7 development version
Public Types | Public Member Functions | Static Public Attributes | List of all members
Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer > Class Template Reference

A class representing small KKT systems on a single subdomain. More...

#include <bddc.hh>

Detailed Description

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
mthe number of vectorial components in \( x \) (usually 1 for scalar systems, and 2 or 3 for vectorial ones)
Scalarthe field type \( K \), usually double
CoarseScalarthe field type \( K_c \) to be used for representing the coarse solver factorization, usually double, but could also be float
Transfera 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.

Public Types

using Scalar = Scalar_
 The scalar field type used in coefficients. More...
 
using XEntry = Dune::FieldVector< Scalar, components >
 
using XVector = Dune::BlockVector< XEntry >
 
using Vector = Dune::BlockVector< Dune::FieldVector< Scalar, 1 > >
 
using TransmissionType = typename Transfer::TransmissionType
 A container type used for exchanging coefficients between subdomains. More...
 
using AMatrix = NumaBCRSMatrix< Dune::FieldMatrix< Scalar, components, components > >
 
using SMatrix = DynamicMatrix< Dune::FieldMatrix< Scalar, 1, 1 > >
 A matrix type with scalar entries, used for the Schur complement matrix. More...
 

Public Member Functions

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_)
 
Construction.
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...
 
Restriction support.

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...
 
Coarse solve support.

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...
 
Prolongation support.

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...
 

Static Public Attributes

static const int components = m
 Number of (vectorial) components of \( u \). More...
 

Member Typedef Documentation

◆ AMatrix

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
using Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::AMatrix = NumaBCRSMatrix<Dune::FieldMatrix<Scalar,components,components> >

Definition at line 202 of file bddc.hh.

◆ Scalar

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
using Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::Scalar = Scalar_

The scalar field type used in coefficients.

Definition at line 190 of file bddc.hh.

◆ SMatrix

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
using Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::SMatrix = DynamicMatrix<Dune::FieldMatrix<Scalar,1,1> >

A matrix type with scalar entries, used for the Schur complement matrix.

Definition at line 205 of file bddc.hh.

◆ TransmissionType

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
using Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::TransmissionType = typename Transfer::TransmissionType

A container type used for exchanging coefficients between subdomains.

Definition at line 199 of file bddc.hh.

◆ Vector

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
using Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::Vector = Dune::BlockVector<Dune::FieldVector<Scalar,1> >

Definition at line 194 of file bddc.hh.

◆ XEntry

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
using Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::XEntry = Dune::FieldVector<Scalar,components>

Definition at line 192 of file bddc.hh.

◆ XVector

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
using Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::XVector = Dune::BlockVector<XEntry>

Definition at line 193 of file bddc.hh.

Constructor & Destructor Documentation

◆ Subdomain() [1/2]

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
template<int cComponents>
Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::Subdomain ( int  id,
AMatrix const &  A,
Interfaces const &  interfaces,
NumaBCRSMatrix< Dune::FieldMatrix< Scalar, cComponents, components > > const &  C 
)

Constructor.

This initializes the solution approximation to zero, which implies an admissible state over all subdomains.

Parameters
idthe subdomain number
Athe subdomain-local stiffness matrix
fthe subdomain-local right hand side
interfacesa list of pairs (s,ids) of neighboring subdomain number and the degrees of freedom shared with this subdomain (in the given order)
Cthe constraint matrix

◆ Subdomain() [2/2]

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
template<class Interface >
Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::Subdomain ( int  id,
AMatrix const &  A,
Interface const &  ifs 
)

Member Function Documentation

◆ getActiveSubdomains()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
std::vector< bool > Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::getActiveSubdomains ( )
inline

Definition at line 499 of file bddc.hh.

◆ getCorrection()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::getCorrection ( XVector dx) const

Provides the current correction vector \( dx \).

Parameters
[out]dxon exit, contains the current correction

◆ getInterfaceAveragingWeights()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::getInterfaceAveragingWeights ( int  s,
XVector rho 
) const
inline

Reports own interface averaging weights.

Parameters
snumber of adjacent subdomain
rhoweights on interface nodes (resized if needed)

Definition at line 248 of file bddc.hh.

◆ getInterfaceResidualUpdate()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::getInterfaceResidualUpdate ( int  s,
TransmissionType dres 
) const
inline

Gets interface residual that shall be shared with neighboring subdomain.

Parameters
sthe number of the neighboring subdomain
[out]resthe residual vector

Definition at line 321 of file bddc.hh.

◆ getInterfaceValues()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::getInterfaceValues ( int  s,
TransmissionType values 
) const
inline

Gets boundary values that shall be shared with neighboring subdomain.

Parameters
sthe number of the neighboring subdomain

Definition at line 425 of file bddc.hh.

◆ getRawCorrection()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::getRawCorrection ( XVector dx) const

◆ getResidual()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::getResidual ( XVector rr) const

◆ getRestrictedResidual()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
Vector const & Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::getRestrictedResidual ( ) const
inline

Definition at line 496 of file bddc.hh.

◆ getS()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::getS ( SMatrix S) const

Compute \( S \).

This computes

\[ S = \begin{bmatrix} 0 & I \end{bmatrix} \begin{bmatrix} A & C^T \\ C \end{bmatrix}^{-1} \begin{bmatrix} 0 \\ I \end{bmatrix} \]

◆ getSchurResidual()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::getSchurResidual ( Vector lambda) const

Provides the Schur residual.

This solves the system

\[ \begin{bmatrix} A & C^T \\ C \end{bmatrix} \begin{bmatrix} x+\delta x \\ \lambda +\delta\lambda \end{bmatrix} = \begin{bmatrix} f\\ 0 \end{bmatrix} \]

and returns the solution component \( \delta\lambda \) in the provided vector lambda. The vectors \( x \) and \( \lambda \) are maintained internally as solution approximation.

Parameters
[out]lambda

◆ getSolution() [1/2]

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
XVector Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::getSolution ( ) const

Provides the current solution vector \( x \).

Returns
the current solution approximation

Note that the result is returned by value, which may be less efficient than providing a result vector to be filled.

◆ getSolution() [2/2]

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::getSolution ( XVector x) const

Provides the current solution vector \( x \).

Parameters
[out]xon exit, contains the current solution approximation

◆ imposeDC_inActiveInterface()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::imposeDC_inActiveInterface ( Vector Aq,
Vector q 
)

◆ imposeDCq()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::imposeDCq ( Vector q)

◆ neighbors()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
std::vector< int > Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::neighbors ( ) const

Reports a list of its neighbor subdomains.

◆ preProlongate()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::preProlongate ( )
inline

Prepares the exchange of interface values.

Has to be called before getInterfaceValues() or setInterfaceValues().

Definition at line 416 of file bddc.hh.

◆ preRestrict()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::preRestrict ( )

First step of the restriction: Dirichlet solve.

This computes the transpose of the prolongation, i.e., with a suitable sorting of interior and bounday degrees of freedom, the

◆ prolongate()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
Dune::FieldVector< double, 2 > Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::prolongate ( )

Perform averaging on the shared interfaces, projecting the solution to the subspace of functions continuous across the interfaces.

Returns
a pair \( (\delta x ^T r, \delta x^T A \delta x) \) for the correction \( \delta x \) and the residual \( r \).

◆ restrict()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::restrict ( )

Performs the restriction of the residual vector.

◆ set_solver_type()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::set_solver_type ( bool  cg_solver)

Set solver type.

◆ setActiveId()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::setActiveId ( bool  active)
inline

Definition at line 440 of file bddc.hh.

◆ setActiveSubdomains_sub()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::setActiveSubdomains_sub ( std::vector< bool > &  activeSubdomains_)
inline

Definition at line 504 of file bddc.hh.

◆ setCorrection()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::setCorrection ( Vector const &  c)

Computes the correction vector \( [\delta x, \delta\lambda] \).

◆ setInterfaceAveragingWeights()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::setInterfaceAveragingWeights ( int  s,
XVector const &  rho 
)
inline

Accumulates interface weights of adjacent subdomains.

Parameters
snumber of adjacent subdomain
rhointerface node weights of neighboring subdomain

Definition at line 258 of file bddc.hh.

◆ setInterfaceResidualUpdate()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::setInterfaceResidualUpdate ( int  s,
TransmissionType const &  dres 
)
inline

Sets interface residual of the correction from neighboring subdomains.

Parameters
sthe number of the neighboring subdomain

Definition at line 331 of file bddc.hh.

◆ setInterfaceValues()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::setInterfaceValues ( int  s,
TransmissionType const &  values 
)
inline

Sets boundary values of the correction from neighboring subdomains.

Parameters
sthe number of the neighboring subdomain

Definition at line 435 of file bddc.hh.

◆ update_rhs()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::update_rhs ( XVector const &  f)

update rhs for the current subdomain

Parameters
frhs

◆ updateIteration()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::updateIteration ( int  cg_itr)

◆ updateSearchDirection()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::updateSearchDirection ( double  beta)

Update search direction q for cg solver.

◆ updateSearchDirectionNorm()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
Dune::FieldVector< double, 1 > Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::updateSearchDirectionNorm ( )

◆ updateSolution()

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
void Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::updateSolution ( Scalar  alpha)

Updates the solution approximation.

Member Data Documentation

◆ components

template<int m, class Scalar_ = double, class CoarseScalar = Scalar_, class Transfer = SpaceTransfer<m,Scalar_>>
const int Kaskade::BDDC::Subdomain< m, Scalar_, CoarseScalar, Transfer >::components = m
static

Number of (vectorial) components of \( u \).

Definition at line 185 of file bddc.hh.


The documentation for this class was generated from the following file: