DUNE PDELab (git)

Dune::Amg::LevelTransferPolicy< FO, CO > Class Template Referenceabstract

Abstract base class for transfer between levels and creation of the coarse level system. More...

#include <dune/istl/paamg/twolevelmethod.hh>

Public Types

typedef FO FineOperatorType
 The linear operator of the finel level system. Has to be derived from AssembledLinearOperator.
 
typedef FineOperatorType::range_type FineRangeType
 The type of the range of the fine level operator.
 
typedef FineOperatorType::domain_type FineDomainType
 The type of the domain of the fine level operator.
 
typedef CO CoarseOperatorType
 The linear operator of the finel level system. Has to be derived from AssembledLinearOperator.
 
typedef CoarseOperatorType::range_type CoarseRangeType
 The type of the range of the coarse level operator.
 
typedef CoarseOperatorType::domain_type CoarseDomainType
 The type of the domain of the coarse level operator.
 

Public Member Functions

std::shared_ptr< CoarseOperatorType > & getCoarseLevelOperator ()
 Get the coarse level operator. More...
 
CoarseRangeTypegetCoarseLevelRhs ()
 Get the coarse level right hand side. More...
 
CoarseDomainTypegetCoarseLevelLhs ()
 Get the coarse level left hand side. More...
 
virtual void moveToCoarseLevel (const FineRangeType &fineRhs)=0
 Transfers the data to the coarse level. More...
 
virtual void moveToFineLevel (FineDomainType &fineLhs)=0
 Updates the fine level linear system after the correction of the coarse levels system. More...
 
virtual void createCoarseLevelSystem (const FineOperatorType &fineOperator)=0
 Algebraically creates the coarse level system. More...
 
virtual LevelTransferPolicyclone () const =0
 Clone the current object.
 
virtual ~LevelTransferPolicy ()
 Destructor.
 

Protected Attributes

CoarseRangeType rhs_
 The coarse level rhs.
 
CoarseDomainType lhs_
 The coarse level lhs.
 
std::shared_ptr< CoarseOperatorTypeoperator_
 the coarse level linear operator.
 

Detailed Description

template<class FO, class CO>
class Dune::Amg::LevelTransferPolicy< FO, CO >

Abstract base class for transfer between levels and creation of the coarse level system.

Template Parameters
FOThe type of the linear operator of the finel level system. Has to be derived from AssembledLinearOperator.
COThe type of the linear operator of the coarse level system. Has to be derived from AssembledLinearOperator.

Member Function Documentation

◆ createCoarseLevelSystem()

template<class FO , class CO >
virtual void Dune::Amg::LevelTransferPolicy< FO, CO >::createCoarseLevelSystem ( const FineOperatorType fineOperator)
pure virtual

Algebraically creates the coarse level system.

After returning from this function the coarse level operator can be accessed using getCoarseLevelOperator().

Parameters
fineOperatorThe operator of the fine level system.

Implemented in Dune::Amg::AggregationLevelTransferPolicy< O, C >.

Referenced by Dune::Amg::TwoLevelMethod< FO, CSP, S >::TwoLevelMethod().

◆ getCoarseLevelLhs()

template<class FO , class CO >
CoarseDomainType & Dune::Amg::LevelTransferPolicy< FO, CO >::getCoarseLevelLhs ( )
inline

Get the coarse level left hand side.

Returns
The coarse level leftt hand side.

References Dune::Amg::LevelTransferPolicy< FO, CO >::lhs_.

◆ getCoarseLevelOperator()

template<class FO , class CO >
std::shared_ptr< CoarseOperatorType > & Dune::Amg::LevelTransferPolicy< FO, CO >::getCoarseLevelOperator ( )
inline

Get the coarse level operator.

Returns
A shared pointer to the coarse level system.

References Dune::Amg::LevelTransferPolicy< FO, CO >::operator_.

◆ getCoarseLevelRhs()

template<class FO , class CO >
CoarseRangeType & Dune::Amg::LevelTransferPolicy< FO, CO >::getCoarseLevelRhs ( )
inline

Get the coarse level right hand side.

Returns
The coarse level right hand side.

References Dune::Amg::LevelTransferPolicy< FO, CO >::rhs_.

◆ moveToCoarseLevel()

template<class FO , class CO >
virtual void Dune::Amg::LevelTransferPolicy< FO, CO >::moveToCoarseLevel ( const FineRangeType fineRhs)
pure virtual

Transfers the data to the coarse level.

Restricts the residual to the right hand side of the coarse level system and initializes the left hand side of the coarse level system. These can afterwards be accessed using getCoarseLevelRhs() and getCoarseLevelLhs().

Parameters
fineRhsThe current right hand side of the fine level system.

◆ moveToFineLevel()

template<class FO , class CO >
virtual void Dune::Amg::LevelTransferPolicy< FO, CO >::moveToFineLevel ( FineDomainType fineLhs)
pure virtual

Updates the fine level linear system after the correction of the coarse levels system.

After returning from this function the coarse level correction will have been added to fine level system.

Parameters
[in,out]fineLhsThe left hand side of the fine level to update with the coarse level correction.

Implemented in Dune::Amg::AggregationLevelTransferPolicy< O, C >.


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)