Dune Core Modules (2.5.2)
A LeveTransferPolicy that used aggregation to construct the coarse level system. More...
#include <dune/istl/paamg/twolevelmethod.hh>
Public Member Functions | |
void | createCoarseLevelSystem (const O &fineOperator) |
Algebraically creates the coarse level system. More... | |
void | moveToFineLevel (typename FatherType::FineDomainType &fineLhs) |
Updates the fine level linear system after the correction of the coarse levels system. More... | |
AggregationLevelTransferPolicy * | clone () const |
Clone the current object. | |
std::shared_ptr< CoarseOperatorType > & | getCoarseLevelOperator () |
Get the coarse level operator. More... | |
CoarseRangeType & | getCoarseLevelRhs () |
Get the coarse level right hand side. More... | |
CoarseDomainType & | getCoarseLevelLhs () |
Get the coarse level left hand side. More... | |
virtual void | moveToCoarseLevel (const FineRangeType &fineRhs)=0 |
Transfers the data to the coarse level. More... | |
Protected Attributes | |
CoarseRangeType | rhs_ |
The coarse level rhs. | |
CoarseDomainType | lhs_ |
The coarse level lhs. | |
std::shared_ptr< CoarseOperatorType > | operator_ |
the coarse level linear operator. | |
Detailed Description
class Dune::Amg::AggregationLevelTransferPolicy< O, C >
A LeveTransferPolicy that used aggregation to construct the coarse level system.
- Template Parameters
-
O The type of the fine and coarse level operator. C The criterion that describes the aggregation procedure.
Member Function Documentation
◆ createCoarseLevelSystem()
|
inlinevirtual |
Algebraically creates the coarse level system.
After returning from this function the coarse level operator can be accessed using getCoarseLevelOperator().
- Parameters
-
fineOperator The operator of the fine level system.
Implements Dune::Amg::LevelTransferPolicy< O, O >.
References Dune::Amg::LevelTransferPolicy< O, O >::lhs_, Dune::Amg::PropertiesGraph< G, VP, EP, VM, EM >::maxVertex(), Dune::Amg::LevelTransferPolicy< O, O >::operator_, and Dune::Amg::LevelTransferPolicy< O, O >::rhs_.
◆ getCoarseLevelLhs()
|
inlineinherited |
Get the coarse level left hand side.
- Returns
- The coarse level leftt hand side.
◆ getCoarseLevelOperator()
|
inlineinherited |
Get the coarse level operator.
- Returns
- A shared pointer to the coarse level system.
◆ getCoarseLevelRhs()
|
inlineinherited |
Get the coarse level right hand side.
- Returns
- The coarse level right hand side.
◆ moveToCoarseLevel()
|
pure virtualinherited |
Transfers the data to the coarse level.
Restricts the residual to the right hand side of the coarse level system and initialies the left hand side of the coarse level system. These can afterwards be accessed usinf getCoarseLevelRhs() and getCoarseLevelLhs().
- Parameters
-
fineDefect The current residual of the fine level system.
◆ moveToFineLevel()
|
inlinevirtual |
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] fineLhs The left hand side of the fine level to update with the coarse level correction.
Implements Dune::Amg::LevelTransferPolicy< O, O >.
References Dune::Amg::LevelTransferPolicy< O, O >::lhs_.
The documentation for this class was generated from the following file:
- dune/istl/paamg/twolevelmethod.hh