DUNE PDELab (git)

Dune::PDELab::SubdomainProjectedCoarseSpace< GFS, M, X, PIH > Class Template Reference

This constructs a coarse space from a per-subdomain local basis. More...

#include <dune/pdelab/backend/istl/geneo/subdomainprojectedcoarsespace.hh>

Public Member Functions

 SubdomainProjectedCoarseSpace (const GFS &gfs, const M &AF_exterior_, std::shared_ptr< SubdomainBasis< X > > subdomainbasis, const PIH &parallelhelper, int verbosity=1)
 Constructor. More...
 
void restrict (const X &fine, COARSE_V &restricted) const override
 Restricts a vector defined on a subdomain to the coarse space. More...
 
void prolongate (const COARSE_V &coarse, X &prolongated) const override
 Prolongates a vector defined on the coarse space to the subdomain. More...
 
std::shared_ptr< COARSE_Mget_coarse_system () override
 Returns the matrix representing the coarse basis. More...
 
rank_type basis_size () override
 Returns the size of the coarse basis. More...
 

Detailed Description

template<class GFS, class M, class X, class PIH>
class Dune::PDELab::SubdomainProjectedCoarseSpace< GFS, M, X, PIH >

This constructs a coarse space from a per-subdomain local basis.

The per-subdomain coarse basis is communicated to each subdomain's neighbors, a global coarse system based on those is constructed and distributed to across all processes. In the process, the per-subdomain basis functions are extended by zeros, resulting in a sparse system.

Constructor & Destructor Documentation

◆ SubdomainProjectedCoarseSpace()

template<class GFS , class M , class X , class PIH >
Dune::PDELab::SubdomainProjectedCoarseSpace< GFS, M, X, PIH >::SubdomainProjectedCoarseSpace ( const GFS &  gfs,
const M &  AF_exterior_,
std::shared_ptr< SubdomainBasis< X > >  subdomainbasis,
const PIH &  parallelhelper,
int  verbosity = 1 
)
inline

Constructor.

Parameters
gfsGrid function space.
AF_exterior_Stiffness matrix of the problem to be solved.
subdomainbasisPer-subdomain coarse basis.
verbosityVerbosity.

Member Function Documentation

◆ basis_size()

template<class GFS , class M , class X , class PIH >
rank_type Dune::PDELab::SubdomainProjectedCoarseSpace< GFS, M, X, PIH >::basis_size ( )
inlineoverridevirtual

Returns the size of the coarse basis.

Returns
Size of the basis

Implements CoarseSpace< X >.

◆ get_coarse_system()

template<class GFS , class M , class X , class PIH >
std::shared_ptr< COARSE_M > Dune::PDELab::SubdomainProjectedCoarseSpace< GFS, M, X, PIH >::get_coarse_system ( )
inlineoverridevirtual

Returns the matrix representing the coarse basis.

Returns
The coarse matrix

Implements CoarseSpace< X >.

◆ prolongate()

template<class GFS , class M , class X , class PIH >
void Dune::PDELab::SubdomainProjectedCoarseSpace< GFS, M, X, PIH >::prolongate ( const COARSE_V coarse,
X &  prolongated 
) const
inlineoverridevirtual

Prolongates a vector defined on the coarse space to the subdomain.

Parameters
[in]vThe coarse space vector to be prolongated
[out]prolongatedThe prolongation in subdomain space.

Implements CoarseSpace< X >.

◆ restrict()

template<class GFS , class M , class X , class PIH >
void Dune::PDELab::SubdomainProjectedCoarseSpace< GFS, M, X, PIH >::restrict ( const X &  fine,
COARSE_V restricted 
) const
inlineoverridevirtual

Restricts a vector defined on a subdomain to the coarse space.

Parameters
[in]dThe subdomain space vector to be restricted
[out]restrictedResulting restriction in coarse space. Must be of size given by basis_size().

Implements CoarseSpace< X >.


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 24, 23:30, 2024)