Dune Core Modules (unstable)
A Parallel Algebraic Multigrid based on Agglomeration. More...
Files | |
file | aggregates.hh |
Provides classes for the Coloring process of AMG. | |
file | amg.hh |
The AMG preconditioner. | |
file | construction.hh |
Helper classes for the construction of classes without empty constructor. | |
file | dependency.hh |
Provides classes for initializing the link attributes of a matrix graph. | |
file | galerkin.hh |
Provides a class for building the galerkin product based on a aggregation scheme. | |
file | globalaggregates.hh |
Provdes class for identifying aggregates globally. | |
file | graph.hh |
Provides classes for building the matrix graph. | |
file | hierarchy.hh |
Provides a classes representing the hierarchies in AMG. | |
file | indicescoarsener.hh |
Provides a class for building the index set and remote indices on the coarse level. | |
file | kamg.hh |
Provides an algebraic multigrid using a Krylov cycle. | |
file | matrixhierarchy.hh |
Provides a classes representing the hierarchies in AMG. | |
file | parameters.hh |
Parameter classes for customizing AMG. | |
file | properties.hh |
Provides classes for handling internal properties in a graph. | |
file | smoother.hh |
Classes for the generic construction and application of the smoothers. | |
file | transfer.hh |
Prolongation and restriction for amg. | |
file | twolevelmethod.hh |
Algebraic twolevel methods. | |
Namespaces | |
namespace | Dune |
Dune namespace. | |
Classes | |
class | Dune::Amg::AggregationCriterion< T > |
Base class of all aggregation criterions. More... | |
class | Dune::Amg::SymmetricMatrixDependency< M, N > |
Dependency policy for symmetric matrices. More... | |
class | Dune::Amg::Dependency< M, N > |
Dependency policy for symmetric matrices. More... | |
class | Dune::Amg::SymmetricDependency< M, N > |
Dependency policy for symmetric matrices. More... | |
class | Dune::Amg::Diagonal< N > |
Norm that uses only the [N][N] entry of the block to determine couplings. More... | |
class | Dune::Amg::FirstDiagonal |
Norm that uses only the [0][0] entry of the block to determine couplings. More... | |
struct | Dune::Amg::RowSum |
Functor using the row sum (infinity) norm to determine strong couplings. More... | |
class | Dune::Amg::SymmetricCriterion< M, Norm > |
Criterion taking advantage of symmetric matrices. More... | |
class | Dune::Amg::UnSymmetricCriterion< M, Norm > |
Criterion suitable for unsymmetric matrices. More... | |
class | Dune::Amg::Aggregator< G > |
Class for building the aggregates. More... | |
class | Dune::Amg::AggregatesMap< V > |
Class providing information about the mapping of the vertices onto aggregates. More... | |
class | Dune::Amg::AggregatesMap< V >::DummyEdgeVisitor |
A Dummy visitor that does nothing for each visited edge. More... | |
class | Dune::Amg::Aggregate< G, S > |
A class for temporarily storing the vertices of an aggregate in. More... | |
class | Dune::Amg::KAMG< M, X, S, PI, K, A > |
an algebraic multigrid method using a Krylov-cycle. More... | |
class | Dune::Amg::KAmgTwoGrid< AMG > |
Two grid operator for AMG with Krylov cycle. More... | |
class | Dune::Amg::AMG< M, X, S, PI, A > |
Parallel algebraic multigrid based on agglomeration. More... | |
struct | Dune::Amg::ConstructionTraits< T > |
Traits class for generically constructing non default constructable types. More... | |
class | Dune::Amg::EdgeProperties |
Class representing the properties of an edge in the matrix graph. More... | |
class | Dune::Amg::VertexProperties |
Class representing a node in the matrix graph. More... | |
class | Dune::Amg::SparsityBuilder< M > |
Functor for building the sparsity pattern of the matrix using examineConnectivity. More... | |
class | Dune::Amg::BaseConnectivityConstructor::ConnectedBuilder< G, S, V > |
Visitor for identifying connected aggregates during a breadthFirstSearch. More... | |
class | Dune::Amg::MatrixGraph< M > |
The (undirected) graph of a matrix. More... | |
class | Dune::Amg::SubGraph< G, T > |
A subgraph of a graph. More... | |
class | Dune::Amg::VertexPropertiesGraph< G, VP, VM > |
Attaches properties to the vertices of a graph. More... | |
class | Dune::Amg::PropertiesGraph< G, VP, EP, VM, EM > |
Attaches properties to the edges and vertices of a graph. More... | |
class | Dune::Amg::GraphVertexPropertiesSelector< G > |
Wrapper to access the internal edge properties of a graph via operator[]() More... | |
class | Dune::Amg::GraphEdgePropertiesSelector< G > |
Wrapper to access the internal vertex properties of a graph via operator[]() More... | |
class | Dune::Amg::Hierarchy< T, A > |
A hierarchy of containers (e.g. matrices or vectors) More... | |
class | Dune::Amg::IndicesCoarsener< OwnerOverlapCopyCommunication< G, L >, E > |
Coarsen Indices in the parallel case. More... | |
class | Dune::Amg::IndicesCoarsener< SequentialInformation, E > |
Coarsen Indices in the sequential case. More... | |
class | Dune::Amg::MatrixHierarchy< M, PI, A > |
The hierarchies build by the coarsening process. More... | |
class | Dune::Amg::CoarsenCriterion< T > |
The criterion describing the stop criteria for the coarsening process. More... | |
class | Dune::Amg::DependencyParameters |
Parameters needed to check whether a node depends on another. More... | |
class | Dune::Amg::AggregationParameters |
Parameters needed for the aggregation process. More... | |
class | Dune::Amg::CoarseningParameters |
Parameters for the complete coarsening process. More... | |
class | Dune::Amg::Parameters |
All parameters for AMG. More... | |
struct | Dune::Amg::VertexVisitedTag |
Tag idnetifying the visited property of a vertex. More... | |
class | Dune::Amg::RandomAccessBundledPropertyMap< C, K, i, T, R > |
A property map that extracts one property out of a bundle using operator[]() More... | |
struct | Dune::Amg::DefaultSmootherArgs< T > |
The default class for the smoother arguments. More... | |
struct | Dune::Amg::SmootherTraits< T > |
Traits class for getting the attribute class of a smoother. More... | |
class | Dune::Amg::DefaultConstructionArgs< T > |
Construction Arguments for the default smoothers. More... | |
struct | Dune::Amg::ConstructionTraits< SeqSSOR< M, X, Y, l > > |
Policy for the construction of the SeqSSOR smoother. More... | |
struct | Dune::Amg::ConstructionTraits< SeqSOR< M, X, Y, l > > |
Policy for the construction of the SeqSOR smoother. More... | |
struct | Dune::Amg::ConstructionTraits< SeqJac< M, X, Y, l > > |
Policy for the construction of the SeqJac smoother. More... | |
struct | Dune::Amg::ConstructionTraits< Richardson< X, Y > > |
Policy for the construction of the Richardson smoother. More... | |
struct | Dune::Amg::ConstructionTraits< SeqILU< M, X, Y > > |
Policy for the construction of the SeqILU smoother. More... | |
struct | Dune::Amg::ConstructionTraits< ParSSOR< M, X, Y, C > > |
Policy for the construction of the ParSSOR smoother. More... | |
struct | Dune::Amg::SmootherApplier< T > |
Helper class for applying the smoothers. More... | |
Typedefs | |
typedef T | Dune::Amg::AggregationCriterion< T >::DependencyPolicy |
The policy for calculating the dependency graph. | |
typedef M | Dune::Amg::SymmetricMatrixDependency< M, N >::Matrix |
The matrix type we build the dependency of. | |
typedef N | Dune::Amg::SymmetricMatrixDependency< M, N >::Norm |
The norm to use for examining the matrix entries. | |
typedef Matrix::row_type | Dune::Amg::SymmetricMatrixDependency< M, N >::Row |
Constant Row iterator of the matrix. | |
typedef Matrix::ConstColIterator | Dune::Amg::SymmetricMatrixDependency< M, N >::ColIter |
Constant column iterator of the matrix. | |
typedef Matrix::field_type | Dune::Amg::SymmetricMatrixDependency< M, N >::field_type |
The current max value. | |
typedef M | Dune::Amg::Dependency< M, N >::Matrix |
The matrix type we build the dependency of. | |
typedef N | Dune::Amg::Dependency< M, N >::Norm |
The norm to use for examining the matrix entries. | |
typedef Matrix::row_type | Dune::Amg::Dependency< M, N >::Row |
Constant Row iterator of the matrix. | |
typedef Matrix::ConstColIterator | Dune::Amg::Dependency< M, N >::ColIter |
Constant column iterator of the matrix. | |
typedef Matrix::field_type | Dune::Amg::Dependency< M, N >::field_type |
The current max value. | |
typedef M | Dune::Amg::SymmetricDependency< M, N >::Matrix |
The matrix type we build the dependency of. | |
typedef N | Dune::Amg::SymmetricDependency< M, N >::Norm |
The norm to use for examining the matrix entries. | |
typedef Matrix::row_type | Dune::Amg::SymmetricDependency< M, N >::Row |
Constant Row iterator of the matrix. | |
typedef Matrix::ConstColIterator | Dune::Amg::SymmetricDependency< M, N >::ColIter |
Constant column iterator of the matrix. | |
typedef Matrix::field_type | Dune::Amg::SymmetricDependency< M, N >::field_type |
The current max value. | |
typedef V | Dune::Amg::AggregatesMap< V >::VertexDescriptor |
The vertex descriptor type. | |
typedef V | Dune::Amg::AggregatesMap< V >::AggregateDescriptor |
The aggregate descriptor type. | |
typedef PoolAllocator< VertexDescriptor, 100 > | Dune::Amg::AggregatesMap< V >::Allocator |
The allocator we use for our lists and the set. | |
typedef SLList< VertexDescriptor, Allocator > | Dune::Amg::AggregatesMap< V >::VertexList |
The type of a single linked list of vertex descriptors. | |
typedef MatrixGraph::VertexDescriptor | Dune::Amg::Aggregate< G, S >::Vertex |
The vertex descriptor type. | |
typedef PoolAllocator< Vertex, 100 > | Dune::Amg::Aggregate< G, S >::Allocator |
The allocator we use for our lists and the set. | |
typedef S | Dune::Amg::Aggregate< G, S >::VertexSet |
The type of a single linked list of vertex descriptors. | |
typedef VertexSet::const_iterator | Dune::Amg::Aggregate< G, S >::const_iterator |
Const iterator over a vertex list. | |
typedef std::size_t * | Dune::Amg::Aggregate< G, S >::SphereMap |
Type of the mapping of aggregate members onto distance spheres. | |
typedef G | Dune::Amg::Aggregator< G >::MatrixGraph |
The matrix graph type used. | |
typedef MatrixGraph::VertexDescriptor | Dune::Amg::Aggregator< G >::Vertex |
The vertex identifier. | |
typedef MatrixGraph::VertexDescriptor | Dune::Amg::Aggregator< G >::AggregateDescriptor |
The type of the aggregate descriptor. | |
typedef V | Dune::Amg::Aggregator< G >::AggregateVisitor< V >::Visitor |
The type of the adapted visitor. | |
typedef M | Dune::Amg::AMG< M, X, S, PI, A >::Operator |
The matrix operator type. | |
typedef PI | Dune::Amg::AMG< M, X, S, PI, A >::ParallelInformation |
The type of the parallel information. Either OwnerOverlapCommunication or another type describing the parallel data distribution and providing communication methods. | |
typedef MatrixHierarchy< M, ParallelInformation, A > | Dune::Amg::AMG< M, X, S, PI, A >::OperatorHierarchy |
The operator hierarchy type. | |
typedef OperatorHierarchy::ParallelInformationHierarchy | Dune::Amg::AMG< M, X, S, PI, A >::ParallelInformationHierarchy |
The parallal data distribution hierarchy type. | |
typedef X | Dune::Amg::AMG< M, X, S, PI, A >::Domain |
The domain type. | |
typedef X | Dune::Amg::AMG< M, X, S, PI, A >::Range |
The range type. | |
typedef InverseOperator< X, X > | Dune::Amg::AMG< M, X, S, PI, A >::CoarseSolver |
the type of the coarse solver. | |
typedef S | Dune::Amg::AMG< M, X, S, PI, A >::Smoother |
The type of the smoother. More... | |
typedef SmootherTraits< Smoother >::Arguments | Dune::Amg::AMG< M, X, S, PI, A >::SmootherArgs |
The argument type for the construction of the smoother. | |
typedef const void * | Dune::Amg::ConstructionTraits< T >::Arguments |
A type holding all the arguments needed to call the constructor. | |
typedef G | Dune::Amg::PropertyGraphVertexPropertyMap< G, i >::Graph |
The type of the graph with internal properties. | |
typedef std::bitset< VertexProperties::SIZE > | Dune::Amg::PropertyGraphVertexPropertyMap< G, i >::BitSet |
The type of the bitset. | |
typedef BitSet::reference | Dune::Amg::PropertyGraphVertexPropertyMap< G, i >::Reference |
The reference type. | |
typedef bool | Dune::Amg::PropertyGraphVertexPropertyMap< G, i >::ValueType |
The value type. | |
typedef G::VertexDescriptor | Dune::Amg::PropertyGraphVertexPropertyMap< G, i >::Vertex |
The vertex descriptor. | |
typedef T | Dune::Amg::OverlapVertex< T >::Aggregate |
The aggregate descriptor. | |
typedef T | Dune::Amg::OverlapVertex< T >::Vertex |
The vertex descriptor. | |
typedef G | Dune::Amg::BaseConnectivityConstructor::ConnectedBuilder< G, S, V >::Graph |
The type of the graph. | |
typedef Graph::ConstEdgeIterator | Dune::Amg::BaseConnectivityConstructor::ConnectedBuilder< G, S, V >::ConstEdgeIterator |
The constant edge iterator. | |
typedef S | Dune::Amg::BaseConnectivityConstructor::ConnectedBuilder< G, S, V >::Set |
The type of the connected set. | |
typedef V | Dune::Amg::BaseConnectivityConstructor::ConnectedBuilder< G, S, V >::VisitedMap |
The type of the map for marking vertices as visited. | |
typedef Graph::VertexDescriptor | Dune::Amg::BaseConnectivityConstructor::ConnectedBuilder< G, S, V >::Vertex |
The vertex descriptor of the graph. | |
typedef E | Dune::Amg::ParallelIndicesCoarsener< T, E >::ExcludedAttributes |
The set of excluded attributes. | |
typedef T | Dune::Amg::ParallelIndicesCoarsener< T, E >::ParallelInformation |
The type of the parallel information. | |
typedef ParallelIndexSet::GlobalIndex | Dune::Amg::ParallelIndicesCoarsener< T, E >::GlobalIndex |
The type of the global index. | |
typedef ParallelIndexSet::LocalIndex | Dune::Amg::ParallelIndicesCoarsener< T, E >::LocalIndex |
The type of the local index. | |
typedef LocalIndex::Attribute | Dune::Amg::ParallelIndicesCoarsener< T, E >::Attribute |
The type of the attribute. | |
typedef Dune::RemoteIndices< ParallelIndexSet > | Dune::Amg::ParallelIndicesCoarsener< T, E >::RemoteIndices |
The type of the remote indices. | |
typedef C | Dune::Amg::RandomAccessBundledPropertyMap< C, K, i, T, R >::Container |
The container that holds the properties. | |
typedef R | Dune::Amg::RandomAccessBundledPropertyMap< C, K, i, T, R >::Reference |
The reference type of the container. | |
typedef K | Dune::Amg::RandomAccessBundledPropertyMap< C, K, i, T, R >::Key |
The key of the property map. | |
typedef LvaluePropertyMapTag | Dune::Amg::RandomAccessBundledPropertyMap< C, K, i, T, R >::Category |
The category of the property map. | |
typedef FieldTraits< T >::real_type | Dune::Amg::DefaultSmootherArgs< T >::RelaxationFactor |
The type of the relaxation factor. | |
Enumerations | |
enum | |
Flags of the link. | |
enum | { Dune::Amg::PropertyGraphVertexPropertyMap< G, i >::index = i } |
enum | { Dune::Amg::MAX_PROCESSES = 72000 } |
enum | Dune::Amg::AccumulationMode { Dune::Amg::noAccu = 0 , Dune::Amg::atOnceAccu =1 , Dune::Amg::successiveAccu =2 } |
Identifiers for the different accumulation modes. More... | |
enum | { Dune::Amg::RandomAccessBundledPropertyMap< C, K, i, T, R >::index = i } |
Functions | |
Dune::Amg::AggregationCriterion< T >::AggregationCriterion () | |
Constructor. More... | |
void | Dune::Amg::AggregationCriterion< T >::setDefaultValuesIsotropic (std::size_t dim, std::size_t diameter=2) |
Sets reasonable default values for an isotropic problem. More... | |
void | Dune::Amg::AggregationCriterion< T >::setDefaultValuesAnisotropic (std::size_t dim, std::size_t diameter=2) |
Sets reasonable default values for an aisotropic problem. More... | |
template<class M > | |
FieldTraits< typenameM::field_type >::real_type | Dune::Amg::Diagonal< N >::operator() (const M &m, typename std::enable_if_t<!Dune::IsNumber< M >::value > *sfinae=nullptr) const |
compute the norm of a matrix. More... | |
template<class M > | |
auto | Dune::Amg::Diagonal< N >::operator() (const M &m, typename std::enable_if_t< Dune::IsNumber< M >::value > *sfinae=nullptr) const |
Compute the norm of a scalar. More... | |
template<class M > | |
FieldTraits< typenameM::field_type >::real_type | Dune::Amg::RowSum::operator() (const M &m) const |
compute the norm of a matrix. More... | |
template<class M > | |
FieldTraits< typenameM::field_type >::real_type | Dune::Amg::FrobeniusNorm::operator() (const M &m) const |
compute the norm of a matrix. More... | |
template<class M > | |
FieldTraits< typenameM::field_type >::real_type | Dune::Amg::AlwaysOneNorm::operator() (const M &m) const |
compute the norm of a matrix. More... | |
Dune::Amg::AggregatesMap< V >::AggregatesMap () | |
Constructs without allocating memory. | |
Dune::Amg::AggregatesMap< V >::AggregatesMap (std::size_t noVertices) | |
Constructs with allocating memory. More... | |
Dune::Amg::AggregatesMap< V >::~AggregatesMap () | |
Destructor. | |
template<class M , class G , class C > | |
std::tuple< int, int, int, int > | Dune::Amg::AggregatesMap< V >::buildAggregates (const M &matrix, G &graph, const C &criterion, bool finestLevel) |
Build the aggregates. More... | |
template<bool reset, class G , class F , class VM > | |
std::size_t | Dune::Amg::AggregatesMap< V >::breadthFirstSearch (const VertexDescriptor &start, const AggregateDescriptor &aggregate, const G &graph, F &aggregateVisitor, VM &visitedMap) const |
Breadth first search within an aggregate. More... | |
template<bool remove, bool reset, class G , class L , class F1 , class F2 , class VM > | |
std::size_t | Dune::Amg::AggregatesMap< V >::breadthFirstSearch (const VertexDescriptor &start, const AggregateDescriptor &aggregate, const G &graph, L &visited, F1 &aggregateVisitor, F2 &nonAggregateVisitor, VM &visitedMap) const |
Breadth first search within an aggregate. More... | |
void | Dune::Amg::AggregatesMap< V >::allocate (std::size_t noVertices) |
Allocate memory for holding the information. More... | |
std::size_t | Dune::Amg::AggregatesMap< V >::noVertices () const |
Get the number of vertices. | |
void | Dune::Amg::AggregatesMap< V >::free () |
Free the allocated memory. | |
AggregateDescriptor & | Dune::Amg::AggregatesMap< V >::operator[] (const VertexDescriptor &v) |
Get the aggregate a vertex belongs to. More... | |
const AggregateDescriptor & | Dune::Amg::AggregatesMap< V >::operator[] (const VertexDescriptor &v) const |
Get the aggregate a vertex belongs to. More... | |
template<class G , class C > | |
void | Dune::Amg::buildDependency (G &graph, const typename C::Matrix &matrix, C criterion, bool finestLevel) |
Build the dependency of the matrix graph. | |
Dune::Amg::Aggregate< G, S >::Aggregate (MatrixGraph &graph, AggregatesMap< Vertex > &aggregates, VertexSet &connectivity, std::vector< Vertex > &front_) | |
Constructor. More... | |
void | Dune::Amg::Aggregate< G, S >::reconstruct (const Vertex &vertex) |
Reconstruct the aggregat from an seed node. More... | |
void | Dune::Amg::Aggregate< G, S >::seed (const Vertex &vertex) |
Initialize the aggregate with one vertex. | |
void | Dune::Amg::Aggregate< G, S >::add (const Vertex &vertex) |
Add a vertex to the aggregate. | |
void | Dune::Amg::Aggregate< G, S >::clear () |
Clear the aggregate. | |
VertexSet::size_type | Dune::Amg::Aggregate< G, S >::size () |
Get the size of the aggregate. | |
VertexSet::size_type | Dune::Amg::Aggregate< G, S >::connectSize () |
Get the number of connections to other aggregates. | |
int | Dune::Amg::Aggregate< G, S >::id () |
Get the id identifying the aggregate. | |
const_iterator | Dune::Amg::Aggregate< G, S >::begin () const |
get an iterator over the vertices of the aggregate. | |
const_iterator | Dune::Amg::Aggregate< G, S >::end () const |
get an iterator over the vertices of the aggregate. | |
Dune::Amg::Aggregator< G >::Aggregator () | |
Constructor. | |
Dune::Amg::Aggregator< G >::~Aggregator () | |
Destructor. | |
template<class M , class C > | |
std::tuple< int, int, int, int > | Dune::Amg::Aggregator< G >::build (const M &m, G &graph, AggregatesMap< Vertex > &aggregates, const C &c, bool finestLevel) |
Build the aggregates. More... | |
Dune::Amg::Aggregator< G >::AggregateVisitor< V >::AggregateVisitor (const AggregatesMap< Vertex > &aggregates, const AggregateDescriptor &aggregate, Visitor &visitor) | |
Constructor. More... | |
void | Dune::Amg::Aggregator< G >::AggregateVisitor< V >::operator() (const typename MatrixGraph::ConstEdgeIterator &edge) |
Examine an edge. More... | |
Dune::Amg::Aggregator< G >::Counter::Counter () | |
Constructor. | |
int | Dune::Amg::Aggregator< G >::Counter::value () |
Access the current count. | |
void | Dune::Amg::Aggregator< G >::Counter::increment () |
Increment counter. | |
void | Dune::Amg::Aggregator< G >::Counter::decrement () |
Decrement counter. | |
Dune::Amg::Aggregator< G >::FrontNeighbourCounter::FrontNeighbourCounter (const MatrixGraph &front) | |
Constructor. More... | |
Dune::Amg::Aggregator< G >::ConnectivityCounter::ConnectivityCounter (const VertexSet &connected, const AggregatesMap< Vertex > &aggregates) | |
Constructor. More... | |
Dune::Amg::Aggregator< G >::DependencyCounter::DependencyCounter () | |
Constructor. | |
Dune::Amg::Aggregator< G >::FrontMarker::FrontMarker (std::vector< Vertex > &front, MatrixGraph &graph) | |
Constructor. More... | |
Dune::Amg::AMG< M, X, S, PI, A >::AMG (OperatorHierarchy &matrices, CoarseSolver &coarseSolver, const SmootherArgs &smootherArgs, const Parameters &parms) | |
Construct a new amg with a specific coarse solver. More... | |
template<class C > | |
Dune::Amg::AMG< M, X, S, PI, A >::AMG (const Operator &fineOperator, const C &criterion, const SmootherArgs &smootherArgs=SmootherArgs(), const ParallelInformation &pinfo=ParallelInformation()) | |
Construct an AMG with an inexact coarse solver based on the smoother. More... | |
Dune::Amg::AMG< M, X, S, PI, A >::AMG (std::shared_ptr< const Operator > fineOperator, const ParameterTree &configuration, const ParallelInformation &pinfo=ParallelInformation()) | |
Constructor an AMG via ParameterTree. More... | |
Dune::Amg::AMG< M, X, S, PI, A >::AMG (const AMG &amg) | |
Copy constructor. | |
void | Dune::Amg::AMG< M, X, S, PI, A >::pre (Domain &x, Range &b) |
Prepare the preconditioner. More... | |
void | Dune::Amg::AMG< M, X, S, PI, A >::apply (Domain &v, const Range &d) |
Apply one step of the preconditioner to the system A(v)=d. More... | |
virtual SolverCategory::Category | Dune::Amg::AMG< M, X, S, PI, A >::category () const |
Category of the preconditioner (see SolverCategory::Category) | |
void | Dune::Amg::AMG< M, X, S, PI, A >::post (Domain &x) |
Clean up. More... | |
template<class A1 > | |
void | Dune::Amg::AMG< M, X, S, PI, A >::getCoarsestAggregateNumbers (std::vector< std::size_t, A1 > &cont) |
Get the aggregate number of each unknown on the coarsest level. More... | |
void | Dune::Amg::AMG< M, X, S, PI, A >::recalculateHierarchy () |
Recalculate the matrix hierarchy. More... | |
bool | Dune::Amg::AMG< M, X, S, PI, A >::usesDirectCoarseLevelSolver () const |
Check whether the coarse solver used is a direct solver. More... | |
static std::shared_ptr< T > | Dune::Amg::ConstructionTraits< T >::construct (Arguments &args) |
Construct an object with the specified arguments. More... | |
Dune::Amg::EdgeProperties::EdgeProperties () | |
Constructor. | |
std::bitset< SIZE >::reference | Dune::Amg::EdgeProperties::operator[] (std::size_t v) |
Access the bits directly. | |
bool | Dune::Amg::EdgeProperties::operator[] (std::size_t v) const |
Access the bits directly. | |
bool | Dune::Amg::EdgeProperties::depends () const |
Checks whether the vertex the edge points to depends on the vertex the edge starts. More... | |
void | Dune::Amg::EdgeProperties::setDepends () |
Marks the edge as one of which the end point depends on the starting point. | |
void | Dune::Amg::EdgeProperties::resetDepends () |
Resets the depends flag. | |
bool | Dune::Amg::EdgeProperties::influences () const |
Checks whether the start vertex is influenced by the end vertex. More... | |
void | Dune::Amg::EdgeProperties::setInfluences () |
Marks the edge as one of which the start vertex by the end vertex. | |
void | Dune::Amg::EdgeProperties::resetInfluences () |
Resets the influence flag. | |
bool | Dune::Amg::EdgeProperties::isOneWay () const |
Checks whether the edge is one way. I.e. either the influence or the depends flag but is set. | |
bool | Dune::Amg::EdgeProperties::isTwoWay () const |
Checks whether the edge is two way. I.e. both the influence flag and the depends flag are that. | |
bool | Dune::Amg::EdgeProperties::isStrong () const |
Checks whether the edge is strong. I.e. the influence or depends flag is set. | |
void | Dune::Amg::EdgeProperties::reset () |
Reset all flags. | |
void | Dune::Amg::EdgeProperties::printFlags () const |
Prints the attributes of the edge for debugging. | |
Dune::Amg::VertexProperties::VertexProperties () | |
Constructor. | |
std::bitset< SIZE >::reference | Dune::Amg::VertexProperties::operator[] (std::size_t v) |
Access the bits directly. | |
bool | Dune::Amg::VertexProperties::operator[] (std::size_t v) const |
Access the bits directly. | |
void | Dune::Amg::VertexProperties::setIsolated () |
Marks that node as being isolated. More... | |
bool | Dune::Amg::VertexProperties::isolated () const |
Checks whether the node is isolated. | |
void | Dune::Amg::VertexProperties::resetIsolated () |
Resets the isolated flag. | |
void | Dune::Amg::VertexProperties::setVisited () |
Mark the node as already visited. | |
bool | Dune::Amg::VertexProperties::visited () const |
Checks whether the node is marked as visited. | |
void | Dune::Amg::VertexProperties::resetVisited () |
Resets the visited flag. | |
void | Dune::Amg::VertexProperties::setFront () |
Marks the node as belonging to the current clusters front. | |
bool | Dune::Amg::VertexProperties::front () const |
Checks whether the node is marked as a front node. | |
void | Dune::Amg::VertexProperties::resetFront () |
Resets the front node flag. | |
void | Dune::Amg::VertexProperties::setExcludedBorder () |
Marks the vertex as excluded from the aggregation. | |
bool | Dune::Amg::VertexProperties::excludedBorder () const |
Tests whether the vertex is excluded from the aggregation. More... | |
void | Dune::Amg::VertexProperties::resetExcludedBorder () |
Marks the vertex as included in the aggregation. | |
void | Dune::Amg::VertexProperties::reset () |
Reset all flags. | |
Dune::Amg::PropertyGraphVertexPropertyMap< G, i >::PropertyGraphVertexPropertyMap (G &g) | |
Constructor. More... | |
Dune::Amg::PropertyGraphVertexPropertyMap< G, i >::PropertyGraphVertexPropertyMap () | |
Default constructor. | |
Reference | Dune::Amg::PropertyGraphVertexPropertyMap< G, i >::operator[] (const Vertex &vertex) const |
Get the properties associated to a vertex. More... | |
Dune::Amg::SparsityBuilder< M >::SparsityBuilder (M &matrix) | |
Constructor. More... | |
template<class M , class V , class I , class O > | |
void | Dune::Amg::BaseGalerkinProduct::calculate (const M &fine, const AggregatesMap< V > &aggregates, M &coarse, const I &pinfo, const O ©) |
Calculate the galerkin product. More... | |
template<class G , class V , class Set > | |
G::MutableMatrix * | Dune::Amg::GalerkinProduct< T >::build (G &fineGraph, V &visitedMap, const ParallelInformation &pinfo, AggregatesMap< typename G::VertexDescriptor > &aggregates, const typename G::Matrix::size_type &size, const Set ©) |
Calculates the coarse matrix via a Galerkin product. More... | |
template<class G , class V , class Set > | |
G::MutableMatrix * | Dune::Amg::GalerkinProduct< SequentialInformation >::build (G &fineGraph, V &visitedMap, const SequentialInformation &pinfo, const AggregatesMap< typename G::VertexDescriptor > &aggregates, const typename G::Matrix::size_type &size, const Set ©) |
Calculates the coarse matrix via a Galerkin product. More... | |
template<class R , class G , class V > | |
static void | Dune::Amg::BaseConnectivityConstructor::constructNonOverlapConnectivity (R &row, G &graph, V &visitedMap, const AggregatesMap< typename G::VertexDescriptor > &aggregates, const typename G::VertexDescriptor &seed) |
Construct the connectivity of an aggregate in the overlap. | |
Dune::Amg::BaseConnectivityConstructor::ConnectedBuilder< G, S, V >::ConnectedBuilder (const AggregatesMap< Vertex > &aggregates, Graph &graph, VisitedMap &visitedMap, Set &connected) | |
Constructor. More... | |
void | Dune::Amg::BaseConnectivityConstructor::ConnectedBuilder< G, S, V >::operator() (const ConstEdgeIterator &edge) |
Process an edge pointing to another aggregate. More... | |
template<class G , class V > | |
int | Dune::Amg::visitNeighbours (const G &graph, const typename G::VertexDescriptor &vertex, V &visitor) |
Visit all neighbour vertices of a vertex in a graph. More... | |
template<typename Graph , typename VM > | |
static Graph::VertexDescriptor | Dune::Amg::ParallelIndicesCoarsener< T, E >::coarsen (ParallelInformation &fineInfo, Graph &fineGraph, VM &visitedMap, AggregatesMap< typename Graph::VertexDescriptor > &aggregates, ParallelInformation &coarseInfo, typename Graph::VertexDescriptor noAggregates, bool useFixedOrder=false) |
Build the coarse index set after the aggregatio. More... | |
Dune::Amg::DependencyParameters::DependencyParameters () | |
Constructor. | |
void | Dune::Amg::DependencyParameters::setBeta (double b) |
Set threshold for marking nodes as isolated. The default value is 1.0E-5. | |
double | Dune::Amg::DependencyParameters::beta () const |
Get the threshold for marking nodes as isolated. The default value is 1.0E-5. More... | |
void | Dune::Amg::DependencyParameters::setAlpha (double a) |
Set the scaling value for marking connections as strong. Default value is 1/3. | |
double | Dune::Amg::DependencyParameters::alpha () const |
Get the scaling value for marking connections as strong. Default value is 1/3. | |
Dune::Amg::AggregationParameters::AggregationParameters () | |
Constructor. More... | |
void | Dune::Amg::AggregationParameters::setDefaultValuesIsotropic (std::size_t dim, std::size_t diameter=2) |
Sets reasonable default values for an isotropic problem. More... | |
void | Dune::Amg::AggregationParameters::setDefaultValuesAnisotropic (std::size_t dim, std::size_t diameter=2) |
Sets reasonable default values for an anisotropic problem. More... | |
std::size_t | Dune::Amg::AggregationParameters::maxDistance () const |
Get the maximal distance allowed between two nodes in a aggregate. More... | |
void | Dune::Amg::AggregationParameters::setMaxDistance (std::size_t distance) |
Set the maximal distance allowed between two nodes in a aggregate. More... | |
bool | Dune::Amg::AggregationParameters::skipIsolated () const |
Whether isolated aggregates will not be represented on the coarse level. More... | |
void | Dune::Amg::AggregationParameters::setSkipIsolated (bool skip) |
Set whether isolated aggregates will not be represented on the coarse level. More... | |
std::size_t | Dune::Amg::AggregationParameters::minAggregateSize () const |
Get the minimum number of nodes a aggregate has to consist of. More... | |
void | Dune::Amg::AggregationParameters::setMinAggregateSize (std::size_t size) |
Set the minimum number of nodes a aggregate has to consist of. More... | |
std::size_t | Dune::Amg::AggregationParameters::maxAggregateSize () const |
Get the maximum number of nodes a aggregate is allowed to have. More... | |
void | Dune::Amg::AggregationParameters::setMaxAggregateSize (std::size_t size) |
Set the maximum number of nodes a aggregate is allowed to have. More... | |
std::size_t | Dune::Amg::AggregationParameters::maxConnectivity () const |
Get the maximum number of connections a aggregate is allowed to have. More... | |
void | Dune::Amg::AggregationParameters::setMaxConnectivity (std::size_t connectivity) |
Set the maximum number of connections a aggregate is allowed to have. More... | |
void | Dune::Amg::CoarseningParameters::setMaxLevel (int l) |
Set the maximum number of levels allowed in the hierarchy. | |
int | Dune::Amg::CoarseningParameters::maxLevel () const |
Get the maximum number of levels allowed in the hierarchy. | |
void | Dune::Amg::CoarseningParameters::setCoarsenTarget (int nodes) |
Set the maximum number of unknowns allowed on the coarsest level. | |
int | Dune::Amg::CoarseningParameters::coarsenTarget () const |
Get the maximum number of unknowns allowed on the coarsest level. | |
void | Dune::Amg::CoarseningParameters::setMinCoarsenRate (double rate) |
Set the minimum coarsening rate to be achieved in each coarsening. More... | |
double | Dune::Amg::CoarseningParameters::minCoarsenRate () const |
Get the minimum coarsening rate to be achieved. | |
AccumulationMode | Dune::Amg::CoarseningParameters::accumulate () const |
Whether the data should be accumulated on fewer processes on coarser levels. | |
void | Dune::Amg::CoarseningParameters::setAccumulate (AccumulationMode accu) |
Set whether the data should be accumulated on fewer processes on coarser levels. | |
bool | Dune::Amg::CoarseningParameters::useFixedOrder () const |
Check if the indices for the coarser levels should be created in a fixed order. | |
void | Dune::Amg::CoarseningParameters::setProlongationDampingFactor (double d) |
Set the damping factor for the prolongation. More... | |
double | Dune::Amg::CoarseningParameters::getProlongationDampingFactor () const |
Get the damping factor for the prolongation. More... | |
Dune::Amg::CoarseningParameters::CoarseningParameters (int maxLevel=100, int coarsenTarget=1000, double minCoarsenRate=1.2, double prolongDamp=1.6, AccumulationMode accumulate=successiveAccu, bool useFixedOrder=false) | |
Constructor. More... | |
void | Dune::Amg::Parameters::setDebugLevel (int level) |
Set the debugging level. More... | |
int | Dune::Amg::Parameters::debugLevel () const |
Get the debugging Level. More... | |
void | Dune::Amg::Parameters::setNoPreSmoothSteps (std::size_t steps) |
Set the number of presmoothing steps to apply. More... | |
std::size_t | Dune::Amg::Parameters::getNoPreSmoothSteps () const |
Get the number of presmoothing steps to apply. More... | |
void | Dune::Amg::Parameters::setNoPostSmoothSteps (std::size_t steps) |
Set the number of postsmoothing steps to apply. More... | |
std::size_t | Dune::Amg::Parameters::getNoPostSmoothSteps () const |
Get the number of postsmoothing steps to apply. More... | |
void | Dune::Amg::Parameters::setGamma (std::size_t gamma) |
Set the value of gamma; 1 for V-cycle, 2 for W-cycle. | |
std::size_t | Dune::Amg::Parameters::getGamma () const |
Get the value of gamma; 1 for V-cycle, 2 for W-cycle. | |
void | Dune::Amg::Parameters::setAdditive (bool additive) |
Set whether to use additive multigrid. More... | |
bool | Dune::Amg::Parameters::getAdditive () const |
Get whether to use additive multigrid. More... | |
Dune::Amg::Parameters::Parameters (int maxLevel=100, int coarsenTarget=1000, double minCoarsenRate=1.2, double prolongDamp=1.6, AccumulationMode accumulate=successiveAccu, bool useFixedOrder=false) | |
Constructor. More... | |
Reference | Dune::Amg::RandomAccessBundledPropertyMap< C, K, i, T, R >::operator[] (const Key &key) const |
Get the property for a key. More... | |
Dune::Amg::RandomAccessBundledPropertyMap< C, K, i, T, R >::RandomAccessBundledPropertyMap (Container &container) | |
Constructor. More... | |
Dune::Amg::RandomAccessBundledPropertyMap< C, K, i, T, R >::RandomAccessBundledPropertyMap () | |
The default constructor. | |
Dune::Amg::DefaultSmootherArgs< T >::DefaultSmootherArgs () | |
Default constructor. | |
static void | Dune::Amg::SmootherApplier< T >::preSmooth (Smoother &smoother, Domain &v, const Range &d) |
apply pre smoothing in forward direction More... | |
static void | Dune::Amg::SmootherApplier< T >::postSmooth (Smoother &smoother, Domain &v, const Range &d) |
apply post smoothing in forward direction More... | |
template<typename LevelContext > | |
void | Dune::Amg::presmooth (LevelContext &levelContext, size_t steps) |
Apply pre smoothing on the current level. More... | |
template<typename LevelContext > | |
void | Dune::Amg::postsmooth (LevelContext &levelContext, size_t steps) |
Apply post smoothing on the current level. More... | |
Dune::Amg::Hierarchy< T, A >::Hierarchy (const std::shared_ptr< MemberType > &first) | |
Construct a new hierarchy. More... | |
Dune::Amg::Hierarchy< T, A >::Hierarchy (const Hierarchy &other) | |
Copy constructor (deep copy!). More... | |
std::size_t | Dune::Amg::Hierarchy< T, A >::levels () const |
Get the number of levels in the hierarchy. More... | |
void | Dune::Amg::Hierarchy< T, A >::addCoarser (Arguments &args) |
Add an element on a coarser level. More... | |
void | Dune::Amg::Hierarchy< T, A >::addFiner (Arguments &args) |
Add an element on a finer level. More... | |
Iterator | Dune::Amg::Hierarchy< T, A >::finest () |
Get an iterator positioned at the finest level. More... | |
Iterator | Dune::Amg::Hierarchy< T, A >::coarsest () |
Get an iterator positioned at the coarsest level. More... | |
ConstIterator | Dune::Amg::Hierarchy< T, A >::finest () const |
Get an iterator positioned at the finest level. More... | |
ConstIterator | Dune::Amg::Hierarchy< T, A >::coarsest () const |
Get an iterator positioned at the coarsest level. More... | |
Dune::Amg::KAMG< M, X, S, PI, K, A >::KAMG (OperatorHierarchy &matrices, CoarseSolver &coarseSolver, const SmootherArgs &smootherArgs, const Parameters &parms, std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1) | |
Construct a new amg with a specific coarse solver. More... | |
template<class C > | |
Dune::Amg::KAMG< M, X, S, PI, K, A >::KAMG (const Operator &fineOperator, const C &criterion, const SmootherArgs &smootherArgs=SmootherArgs(), std::size_t maxLevelKrylovSteps=3, double minDefectReduction=1e-1, const ParallelInformation &pinfo=ParallelInformation()) | |
Construct an AMG with an inexact coarse solver based on the smoother. More... | |
void | Dune::Amg::KAMG< M, X, S, PI, K, A >::pre (Domain &x, Range &b) |
Prepare the preconditioner. More... | |
void | Dune::Amg::KAMG< M, X, S, PI, K, A >::post (Domain &x) |
Clean up. More... | |
void | Dune::Amg::KAMG< M, X, S, PI, K, A >::apply (Domain &v, const Range &d) |
Apply one step of the preconditioner to the system A(v)=d. More... | |
Dune::Amg::MatrixHierarchy< M, PI, A >::MatrixHierarchy (std::shared_ptr< MatrixOperator > fineMatrix, std::shared_ptr< ParallelInformation > pinfo=std::make_shared< ParallelInformation >()) | |
Constructor. More... | |
template<typename O , typename T > | |
void | Dune::Amg::MatrixHierarchy< M, PI, A >::build (const T &criterion) |
Build the matrix hierarchy using aggregation. More... | |
const ParallelMatrixHierarchy & | Dune::Amg::MatrixHierarchy< M, PI, A >::matrices () const |
Get the matrix hierarchy. More... | |
const ParallelInformationHierarchy & | Dune::Amg::MatrixHierarchy< M, PI, A >::parallelInformation () const |
Get the hierarchy of the parallel data distribution information. More... | |
void | Dune::Amg::MatrixHierarchy< M, PI, A >::getCoarsestAggregatesOnFinest (std::vector< std::size_t > &data) const |
Get the mapping of fine level unknowns to coarse level aggregates. More... | |
const AggregatesMapList & | Dune::Amg::MatrixHierarchy< M, PI, A >::aggregatesMaps () const |
Get the hierarchy of the mappings of the nodes onto aggregates. More... | |
const RedistributeInfoList & | Dune::Amg::MatrixHierarchy< M, PI, A >::redistributeInformation () const |
Get the hierarchy of the information about redistributions,. More... | |
template<class V , class BA , class TA > | |
void | Dune::Amg::MatrixHierarchy< M, PI, A >::coarsenVector (Hierarchy< BlockVector< V, BA >, TA > &hierarchy) const |
Coarsen the vector hierarchy according to the matrix hierarchy. More... | |
template<class S , class TA > | |
void | Dune::Amg::MatrixHierarchy< M, PI, A >::coarsenSmoother (Hierarchy< S, TA > &smoothers, const typename SmootherTraits< S >::Arguments &args) const |
Coarsen the smoother hierarchy according to the matrix hierarchy. More... | |
template<class F > | |
void | Dune::Amg::MatrixHierarchy< M, PI, A >::recalculateGalerkin (const F ©Flags) |
Recalculate the galerkin products. More... | |
std::size_t | Dune::Amg::MatrixHierarchy< M, PI, A >::levels () const |
Get the number of levels in the hierarchy. More... | |
std::size_t | Dune::Amg::MatrixHierarchy< M, PI, A >::maxlevels () const |
Get the max number of levels in the hierarchy of processors. More... | |
bool | Dune::Amg::MatrixHierarchy< M, PI, A >::isBuilt () const |
Whether the hierarchy was built. More... | |
Variables | |
const Matrix * | Dune::Amg::SymmetricMatrixDependency< M, N >::matrix_ |
The matrix we work on. | |
Norm | Dune::Amg::SymmetricMatrixDependency< M, N >::norm_ |
The functor for calculating the norm. | |
int | Dune::Amg::SymmetricMatrixDependency< M, N >::row_ |
index of the currently evaluated row. | |
real_type | Dune::Amg::SymmetricMatrixDependency< M, N >::diagonal_ |
The norm of the current diagonal. | |
const Matrix * | Dune::Amg::Dependency< M, N >::matrix_ |
The matrix we work on. | |
Norm | Dune::Amg::Dependency< M, N >::norm_ |
The functor for calculating the norm. | |
int | Dune::Amg::Dependency< M, N >::row_ |
index of the currently evaluated row. | |
real_type | Dune::Amg::Dependency< M, N >::diagonal_ |
The norm of the current diagonal. | |
const Matrix * | Dune::Amg::SymmetricDependency< M, N >::matrix_ |
The matrix we work on. | |
Norm | Dune::Amg::SymmetricDependency< M, N >::norm_ |
The functor for calculating the norm. | |
int | Dune::Amg::SymmetricDependency< M, N >::row_ |
index of the currently evaluated row. | |
real_type | Dune::Amg::SymmetricDependency< M, N >::diagonal_ |
The norm of the current diagonal. | |
static const V | Dune::Amg::AggregatesMap< V >::UNAGGREGATED |
Identifier of not yet aggregated vertices. | |
static const V | Dune::Amg::AggregatesMap< V >::ISOLATED |
Identifier of isolated vertices. | |
Hierarchy< Smoother, A >::Iterator | Dune::Amg::AMG< M, X, S, PI, A >::LevelContext::smoother |
The iterator over the smoothers. | |
OperatorHierarchy::ParallelMatrixHierarchy::ConstIterator | Dune::Amg::AMG< M, X, S, PI, A >::LevelContext::matrix |
The iterator over the matrices. | |
ParallelInformationHierarchy::Iterator | Dune::Amg::AMG< M, X, S, PI, A >::LevelContext::pinfo |
The iterator over the parallel information. | |
OperatorHierarchy::RedistributeInfoList::const_iterator | Dune::Amg::AMG< M, X, S, PI, A >::LevelContext::redist |
The iterator over the redistribution information. | |
OperatorHierarchy::AggregatesMapList::const_iterator | Dune::Amg::AMG< M, X, S, PI, A >::LevelContext::aggregates |
The iterator over the aggregates maps. | |
Hierarchy< Domain, A >::Iterator | Dune::Amg::AMG< M, X, S, PI, A >::LevelContext::lhs |
The iterator over the left hand side. | |
Hierarchy< Domain, A >::Iterator | Dune::Amg::AMG< M, X, S, PI, A >::LevelContext::update |
The iterator over the updates. | |
Hierarchy< Range, A >::Iterator | Dune::Amg::AMG< M, X, S, PI, A >::LevelContext::rhs |
The iterator over the right hand sided. | |
std::size_t | Dune::Amg::AMG< M, X, S, PI, A >::LevelContext::level |
The level index. | |
Aggregate * | Dune::Amg::OverlapVertex< T >::aggregate |
The aggregate the vertex belongs to. | |
Vertex | Dune::Amg::OverlapVertex< T >::vertex |
The vertex descriptor. | |
int | Dune::Amg::DefaultSmootherArgs< T >::iterations |
The number of iterations to perform. | |
RelaxationFactor | Dune::Amg::DefaultSmootherArgs< T >::relaxationFactor |
The relaxation factor to use. | |
Detailed Description
A Parallel Algebraic Multigrid based on Agglomeration.
Typedef Documentation
◆ Smoother
typedef S Dune::Amg::AMG< M, X, S, PI, A >::Smoother |
The type of the smoother.
One of the preconditioners implementing the Preconditioner interface. Note that the smoother has to fit the ParallelInformation.
Enumeration Type Documentation
◆ anonymous enum
anonymous enum |
◆ anonymous enum
anonymous enum |
◆ anonymous enum
anonymous enum |
◆ AccumulationMode
Identifiers for the different accumulation modes.
Function Documentation
◆ addCoarser()
void Dune::Amg::Hierarchy< T, A >::addCoarser | ( | Arguments & | args | ) |
Add an element on a coarser level.
- Parameters
-
args The arguments needed for the construction.
References Dune::Amg::ConstructionTraits< T >::construct().
Referenced by Dune::Amg::MatrixHierarchy< M, PI, A >::coarsenSmoother().
◆ addFiner()
void Dune::Amg::Hierarchy< T, A >::addFiner | ( | Arguments & | args | ) |
Add an element on a finer level.
- Parameters
-
args The arguments needed for the construction.
References Dune::Amg::ConstructionTraits< T >::construct().
◆ Aggregate()
Dune::Amg::Aggregate< G, S >::Aggregate | ( | MatrixGraph & | graph, |
AggregatesMap< Vertex > & | aggregates, | ||
VertexSet & | connectivity, | ||
std::vector< Vertex > & | front_ | ||
) |
Constructor.
- Parameters
-
graph The matrix graph we work on. aggregates The mapping of vertices onto aggregates. connectivity The set of vertices connected to the aggregate. distance spheres. front_ The vertices of the current aggregate front.
◆ AggregatesMap()
Dune::Amg::AggregatesMap< V >::AggregatesMap | ( | std::size_t | noVertices | ) |
Constructs with allocating memory.
- Parameters
-
noVertices The number of vertices we will hold information for.
◆ aggregatesMaps()
const MatrixHierarchy< M, IS, A >::AggregatesMapList & Dune::Amg::MatrixHierarchy< M, IS, A >::aggregatesMaps |
Get the hierarchy of the mappings of the nodes onto aggregates.
- Returns
- The hierarchy of the mappings of the nodes onto aggregates.
◆ AggregateVisitor()
Dune::Amg::Aggregator< G >::AggregateVisitor< V >::AggregateVisitor | ( | const AggregatesMap< Vertex > & | aggregates, |
const AggregateDescriptor & | aggregate, | ||
Visitor & | visitor | ||
) |
Constructor.
- Parameters
-
aggregates The aggregate numbers of the vertices. aggregate The id of the aggregate to visit. visitor The visitor.
◆ AggregationCriterion()
|
inline |
Constructor.
The parameters will be initialized with default values suitable for 2D isotropic problems.
If that does not fit your needs either use setDefaultValuesIsotropic setDefaultValuesAnisotropic or setup the values by hand
◆ AggregationParameters()
|
inline |
Constructor.
The parameters will be initialized with default values suitable for 2D isotropic problems.
If that does not fit your needs either use setDefaultValuesIsotropic setDefaultValuesAnisotropic or setup the values by hand
◆ allocate()
void Dune::Amg::AggregatesMap< V >::allocate | ( | std::size_t | noVertices | ) |
Allocate memory for holding the information.
- Parameters
-
noVertices The total number of vertices to be mapped.
◆ AMG() [1/3]
Dune::Amg::AMG< M, X, S, PI, A >::AMG | ( | const Operator & | fineOperator, |
const C & | criterion, | ||
const SmootherArgs & | smootherArgs = SmootherArgs() , |
||
const ParallelInformation & | pinfo = ParallelInformation() |
||
) |
Construct an AMG with an inexact coarse solver based on the smoother.
As coarse solver a preconditioned CG method with the smoother as preconditioner will be used. The matrix hierarchy is built automatically.
- Parameters
-
fineOperator The operator on the fine level. criterion The criterion describing the coarsening strategy. E. g. SymmetricCriterion or UnsymmetricCriterion, and providing the parameters. smootherArgs The arguments for constructing the smoothers. pinfo The information about the parallel distribution of the data.
◆ AMG() [2/3]
Dune::Amg::AMG< M, X, S, PI, A >::AMG | ( | OperatorHierarchy & | matrices, |
CoarseSolver & | coarseSolver, | ||
const SmootherArgs & | smootherArgs, | ||
const Parameters & | parms | ||
) |
Construct a new amg with a specific coarse solver.
- Parameters
-
matrices The already set-up matrix hierarchy. coarseSolver The set up solver to use on the coarse grid, must match the coarse matrix in the matrix hierarchy. smootherArgs The arguments needed for thesmoother to use for pre and post smoothing. parms The parameters for the AMG.
◆ AMG() [3/3]
Dune::Amg::AMG< M, X, S, PI, A >::AMG | ( | std::shared_ptr< const Operator > | fineOperator, |
const ParameterTree & | configuration, | ||
const ParallelInformation & | pinfo = ParallelInformation() |
||
) |
Constructor an AMG via ParameterTree.
- Parameters
-
fineOperator The operator on the fine level. configuration ParameterTree containing AMG parameters. pinfo Optionally, specify ParallelInformation
ParameterTree Key | Meaning |
---|---|
smootherIterations | The number of iterations to perform. |
smootherRelaxation | The relaxation factor |
maxLevel | Maximum number of levels allowed in the hierarchy. |
coarsenTarget | Maximum number of unknowns on the coarsest level. |
minCoarseningRate | Coarsening will stop if the rate is below this threshold. |
accumulationMode | If and how data is agglomerated on coarser level to |
| fewer processors. ("atOnce": do agglomeration once and | to one process; "successive": Multiple agglomerations to | fewer processes until all data is on one process; | "none": Do no agglomeration at all and solve coarse level | iteratively). prolongationDampingFactor | Damping factor for the prolongation. alpha | Scaling value for marking connections as strong. beta | Threshold for marking nodes as isolated. additive | Whether to use additive multigrid. gamma | 1 for V-cycle, 2 for W-cycle. preSteps | Number of presmoothing steps. postSteps | Number of postsmoothing steps. verbosity | Output verbosity. default=2. criterionSymmetric | If true use SymmetricCriterion (default), else UnSymmetricCriterion strengthMeasure | What conversion to use to convert a matrix block | to a scalar when determining strength of connection: | diagonal (use a diagonal of row diagonalRowIndex, class Diagonal, default); | rowSum (rowSum norm), frobenius (Frobenius norm); | one (use always one and neglect the actual entries). diagonalRowIndex | The index to use for the diagonal strength (default 0) | if this is i and strengthMeasure is "diagonal", then | block[i][i] will be used when determining strength of | connection. defaultAggregationSizeMode | Whether to set default values depending on isotropy of | problem uses parameters "defaultAggregationDimension" and | "maxAggregateDistance" (isotropic: For and isotropic problem; | anisotropic: for an anisotropic problem). defaultAggregationDimension | Dimension of the problem (used for setting default aggregate size). maxAggregateDistance | Maximum distance in an aggregte (in term of minimum edges needed to travel. | one vertex to another within the aggregate). minAggregateSize | Minimum number of vertices an aggregate should consist of. maxAggregateSize | Maximum number of vertices an aggregate should consist of.
See ISTL_Factory for the ParameterTree layout and examples.
◆ apply() [1/2]
void Dune::Amg::AMG< M, X, S, PI, A >::apply | ( | Domain & | v, |
const Range & | d | ||
) |
Apply one step of the preconditioner to the system A(v)=d.
On entry v=0 and d=b-A(x) (although this might not be computed in that way. On exit v contains the update, i.e one step computes \( v = M^{-1} d \) where \( M \) is the approximate inverse of the operator \( A \) characterizing the preconditioner.
- Parameters
-
[out] v The update to be computed d The current defect.
◆ apply() [2/2]
void Dune::Amg::KAMG< M, X, S, P, K, A >::apply | ( | Domain & | v, |
const Range & | d | ||
) |
Apply one step of the preconditioner to the system A(v)=d.
On entry v=0 and d=b-A(x) (although this might not be computed in that way. On exit v contains the update, i.e one step computes \( v = M^{-1} d \) where \( M \) is the approximate inverse of the operator \( A \) characterizing the preconditioner.
- Parameters
-
[out] v The update to be computed d The current defect.
◆ beta()
|
inline |
Get the threshold for marking nodes as isolated. The default value is 1.0E-5.
- Returns
- beta
◆ breadthFirstSearch() [1/2]
std::size_t Dune::Amg::AggregatesMap< V >::breadthFirstSearch | ( | const VertexDescriptor & | start, |
const AggregateDescriptor & | aggregate, | ||
const G & | graph, | ||
F & | aggregateVisitor, | ||
VM & | visitedMap | ||
) | const |
Breadth first search within an aggregate.
- Template Parameters
-
reset If true the visited flags of the vertices will be reset after the search G The type of the graph we perform the search on F The type of the visitor to operate on the vertices
- Parameters
-
start The vertex where the search should start from. This does not need to belong to the aggregate. aggregate The aggregate id. graph The matrix graph to perform the search on. visitedMap A map to mark the already visited vertices aggregateVisitor A functor that is called with each G::ConstEdgeIterator with an edge pointing to the aggregate. Use DummyVisitor if these are of no interest.
◆ breadthFirstSearch() [2/2]
std::size_t Dune::Amg::AggregatesMap< V >::breadthFirstSearch | ( | const VertexDescriptor & | start, |
const AggregateDescriptor & | aggregate, | ||
const G & | graph, | ||
L & | visited, | ||
F1 & | aggregateVisitor, | ||
F2 & | nonAggregateVisitor, | ||
VM & | visitedMap | ||
) | const |
Breadth first search within an aggregate.
- Template Parameters
-
L A container type providing push_back(Vertex), and pop_front() in case remove is true remove If true the entries in the visited list will be removed. reset If true the visited flag will be reset after the search
- Parameters
-
start The vertex where the search should start from. This does not need to belong to the aggregate. aggregate The aggregate id. graph The matrix graph to perform the search on. visited A list to store the visited vertices in. aggregateVisitor A functor that is called with each G::ConstEdgeIterator with an edge pointing to the aggregate. Use DummyVisitor these are of no interest. nonAggregateVisitor A functor that is called with each G::ConstEdgeIterator with an edge pointing to another aggregate. Use DummyVisitor these are of no interest. visitedMap A map to mark the already visited vertices
◆ build() [1/4]
std::tuple< int, int, int, int > Dune::Amg::Aggregator< G >::build | ( | const M & | m, |
G & | graph, | ||
AggregatesMap< Vertex > & | aggregates, | ||
const C & | c, | ||
bool | finestLevel | ||
) |
Build the aggregates.
- Template Parameters
-
C The type of the coarsening Criterion to use
- Parameters
-
m The matrix to build the aggregates accordingly. graph A (sub) graph of the matrix. aggregates Aggregate map we will build. All entries should be initialized to UNAGGREGATED! c The coarsening criterion to use. finestLevel Whether this the finest level. In that case rows representing Dirichlet boundaries will be detected and ignored during aggregation.
- Returns
- A tuple of the total number of aggregates, the number of isolated aggregates, the number of isolated aggregates, the number of aggregates consisting only of one vertex, and the number of skipped aggregates built.
◆ build() [2/4]
void Dune::Amg::MatrixHierarchy< M, IS, A >::build | ( | const T & | criterion | ) |
Build the matrix hierarchy using aggregation.
criterion The criterion describing the aggregation process.
References Dune::Amg::atOnceAccu, Dune::countNonZeros(), Dune::Amg::MAX_PROCESSES, Dune::size(), and Dune::Amg::successiveAccu.
◆ build() [3/4]
G::MutableMatrix * Dune::Amg::GalerkinProduct< T >::build | ( | G & | fineGraph, |
V & | visitedMap, | ||
const ParallelInformation & | pinfo, | ||
AggregatesMap< typename G::VertexDescriptor > & | aggregates, | ||
const typename G::Matrix::size_type & | size, | ||
const Set & | copy | ||
) |
Calculates the coarse matrix via a Galerkin product.
- Parameters
-
fineGraph The graph of the fine matrix. visitedMap Map for marking vertices as visited. pinfo Parallel information about the fine level. aggregates The mapping of the fine level unknowns onto aggregates. size The number of columns and rows of the coarse matrix. copy The attribute set identifying the copy nodes of the graph.
◆ build() [4/4]
G::MutableMatrix * Dune::Amg::GalerkinProduct< SequentialInformation >::build | ( | G & | fineGraph, |
V & | visitedMap, | ||
const SequentialInformation & | pinfo, | ||
const AggregatesMap< typename G::VertexDescriptor > & | aggregates, | ||
const typename G::Matrix::size_type & | size, | ||
const Set & | copy | ||
) |
Calculates the coarse matrix via a Galerkin product.
- Parameters
-
fineGraph The graph of the fine matrix. visitedMap Map for marking vertices as visited. pinfo Parallel information about the fine level. aggregates The mapping of the fine level unknowns onto aggregates. size The number of columns and rows of the coarse matrix. copy The attribute set identifying the copy nodes of the graph.
References Dune::dinfo, Dune::size(), and Dune::GeometryTypes::vertex.
◆ buildAggregates()
std::tuple< int, int, int, int > Dune::Amg::AggregatesMap< V >::buildAggregates | ( | const M & | matrix, |
G & | graph, | ||
const C & | criterion, | ||
bool | finestLevel | ||
) |
Build the aggregates.
- Parameters
-
matrix The matrix describing the dependency. graph The graph corresponding to the matrix. criterion The aggregation criterion. finestLevel Whether this the finest level. In that case rows representing Dirichlet boundaries will be detected and ignored during aggregation.
- Returns
- A tuple of the total number of aggregates, the number of isolated aggregates, the number of isolated aggregates, the number of aggregates consisting only of one vertex, and the number of skipped aggregates built.
◆ calculate()
void Dune::Amg::BaseGalerkinProduct::calculate | ( | const M & | fine, |
const AggregatesMap< V > & | aggregates, | ||
M & | coarse, | ||
const I & | pinfo, | ||
const O & | copy | ||
) |
Calculate the galerkin product.
- Parameters
-
fine The fine matrix. aggregates The aggregate mapping. coarse The coarse Matrix. pinfo Parallel information about the fine level. copy The attribute set identifying the copy nodes of the graph.
◆ coarsen()
|
inlinestatic |
Build the coarse index set after the aggregatio.
- Parameters
-
fineInfo The parallel information at the fine level. fineGraph The graph of the fine lecel, visitedMap Map for marking vertices as visited. aggregates The mapping of unknowns onto aggregates. coarseInfo The information about the parallel data decomposition on the coarse level. noAggregates useFixedOrder Flag indicating if creating indices for the coarser level should be done in a fixed order, i.e., the order in which the rows were sent. If set to true, this makes the runs reproducible but it might slow down performance.
- Returns
- The number of unknowns on the coarse level.
◆ CoarseningParameters()
|
inline |
Constructor.
- Parameters
-
maxLevel The maximum number of levels allowed in the matrix hierarchy (default: 100). coarsenTarget If the number of nodes in the matrix is below this threshold the coarsening will stop (default: 1000). minCoarsenRate If the coarsening rate falls below this threshold the coarsening will stop (default: 1.2) prolongDamp The damping factor to apply to the prolongated update (default: 1.6) accumulate Whether to accumulate the data onto fewer processors on coarser levels. useFixedOrder Flag indicating if creating indices for the coarser level should be done in a fixed order, i.e., the order in which the rows were sent. If set to true, this makes the runs reproducible but it might slow down performance.
◆ coarsenSmoother()
void Dune::Amg::MatrixHierarchy< M, IS, A >::coarsenSmoother | ( | Hierarchy< S, TA > & | smoothers, |
const typename SmootherTraits< S >::Arguments & | args | ||
) | const |
Coarsen the smoother hierarchy according to the matrix hierarchy.
- Parameters
-
smoothers The smoother hierarchy to coarsen. args The arguments for the construction of the coarse level smoothers.
References Dune::Amg::Hierarchy< T, A >::addCoarser(), and Dune::Amg::Hierarchy< T, A >::levels().
◆ coarsenVector()
void Dune::Amg::MatrixHierarchy< M, IS, A >::coarsenVector | ( | Hierarchy< BlockVector< V, BA >, TA > & | hierarchy | ) | const |
Coarsen the vector hierarchy according to the matrix hierarchy.
- Parameters
-
hierarchy The vector hierarchy to coarsen.
References Dune::dvverb.
◆ coarsest() [1/2]
Hierarchy< T, A >::Iterator Dune::Amg::Hierarchy< T, A >::coarsest |
Get an iterator positioned at the coarsest level.
- Returns
- An iterator positioned at the coarsest level.
Referenced by Dune::Amg::AMG< M, X, S, PI, A >::post().
◆ coarsest() [2/2]
Hierarchy< T, A >::ConstIterator Dune::Amg::Hierarchy< T, A >::coarsest |
Get an iterator positioned at the coarsest level.
- Returns
- An iterator positioned at the coarsest level.
◆ ConnectedBuilder()
Dune::Amg::BaseConnectivityConstructor::ConnectedBuilder< G, S, V >::ConnectedBuilder | ( | const AggregatesMap< Vertex > & | aggregates, |
Graph & | graph, | ||
VisitedMap & | visitedMap, | ||
Set & | connected | ||
) |
Constructor.
- Parameters
-
aggregates The mapping of the vertices onto the aggregates. graph The graph to work on. visitedMap The map for marking vertices as visited connected The set to added the connected aggregates to.
◆ ConnectivityCounter()
Dune::Amg::Aggregator< G >::ConnectivityCounter::ConnectivityCounter | ( | const VertexSet & | connected, |
const AggregatesMap< Vertex > & | aggregates | ||
) |
Constructor.
- Parameters
-
connected The set of connected aggregates. aggregates Mapping of the vertices onto the aggregates. aggregates The mapping of aggregates to vertices.
◆ construct()
|
inlinestatic |
Construct an object with the specified arguments.
In the default implementation the copy constructor is called.
- Parameters
-
args The arguments for the construction.
Referenced by Dune::Amg::Hierarchy< T, A >::addCoarser(), and Dune::Amg::Hierarchy< T, A >::addFiner().
◆ debugLevel()
|
inline |
Get the debugging Level.
- Returns
- 0 if no debugging output will be generated.
◆ depends()
|
inline |
Checks whether the vertex the edge points to depends on the vertex the edge starts.
- Returns
- True if it depends on the starting point.
◆ excludedBorder()
|
inline |
Tests whether the vertex is excluded from the aggregation.
- Returns
- True if the vertex is excluded from the aggregation process.
◆ finest() [1/2]
Hierarchy< T, A >::Iterator Dune::Amg::Hierarchy< T, A >::finest |
Get an iterator positioned at the finest level.
- Returns
- An iterator positioned at the finest level.
◆ finest() [2/2]
Hierarchy< T, A >::ConstIterator Dune::Amg::Hierarchy< T, A >::finest |
Get an iterator positioned at the finest level.
- Returns
- An iterator positioned at the finest level.
◆ FrontMarker()
Dune::Amg::Aggregator< G >::FrontMarker::FrontMarker | ( | std::vector< Vertex > & | front, |
MatrixGraph & | graph | ||
) |
Constructor.
- Parameters
-
front The list to store the front vertices in. graph The matrix graph we work on.
◆ FrontNeighbourCounter()
Dune::Amg::Aggregator< G >::FrontNeighbourCounter::FrontNeighbourCounter | ( | const MatrixGraph & | front | ) |
Constructor.
- Parameters
-
front The vertices of the front.
◆ getAdditive()
|
inline |
Get whether to use additive multigrid.
- Returns
- True if multigrid should be additive.
◆ getCoarsestAggregateNumbers()
void Dune::Amg::AMG< M, X, S, PI, A >::getCoarsestAggregateNumbers | ( | std::vector< std::size_t, A1 > & | cont | ) |
Get the aggregate number of each unknown on the coarsest level.
- Parameters
-
cont The random access container to store the numbers in.
◆ getCoarsestAggregatesOnFinest()
void Dune::Amg::MatrixHierarchy< M, IS, A >::getCoarsestAggregatesOnFinest | ( | std::vector< std::size_t > & | data | ) | const |
Get the mapping of fine level unknowns to coarse level aggregates.
For each fine level unknown i the corresponding data[i] is the aggregate it belongs to on the coarsest level.
- Parameters
-
[out] data The mapping of fine level unknowns to coarse level aggregates.
References Dune::Amg::AggregatesMap< V >::ISOLATED, Dune::Hybrid::max, and Dune::size().
◆ getNoPostSmoothSteps()
|
inline |
Get the number of postsmoothing steps to apply.
- Returns
- The number of steps:
◆ getNoPreSmoothSteps()
|
inline |
Get the number of presmoothing steps to apply.
- Returns
- The number of steps:
◆ getProlongationDampingFactor()
|
inline |
Get the damping factor for the prolongation.
- Returns
- d The damping factor.
◆ Hierarchy() [1/2]
Dune::Amg::Hierarchy< T, A >::Hierarchy | ( | const Hierarchy< T, A > & | other | ) |
Copy constructor (deep copy!).
deep copy of a given hierarchy
◆ Hierarchy() [2/2]
Dune::Amg::Hierarchy< T, A >::Hierarchy | ( | const std::shared_ptr< MemberType > & | first | ) |
Construct a new hierarchy.
- Parameters
-
first std::shared_ptr to the first element in the hierarchy.
◆ influences()
|
inline |
Checks whether the start vertex is influenced by the end vertex.
- Returns
- True if it is influenced.
◆ isBuilt()
bool Dune::Amg::MatrixHierarchy< M, IS, A >::isBuilt |
Whether the hierarchy was built.
- Returns
- true if the MatrixHierarchy::build method was called.
◆ KAMG() [1/2]
Dune::Amg::KAMG< M, X, S, P, K, A >::KAMG | ( | const Operator & | fineOperator, |
const C & | criterion, | ||
const SmootherArgs & | smootherArgs = SmootherArgs() , |
||
std::size_t | maxLevelKrylovSteps = 3 , |
||
double | minDefectReduction = 1e-1 , |
||
const ParallelInformation & | pinfo = ParallelInformation() |
||
) |
Construct an AMG with an inexact coarse solver based on the smoother.
As coarse solver a preconditioned CG method with the smoother as preconditioner will be used. The matrix hierarchy is built automatically.
- Parameters
-
fineOperator The operator on the fine level. criterion The criterion describing the coarsening strategy. E. g. SymmetricCriterion or UnsymmetricCriterion, and providing the parameters. smootherArgs The arguments for constructing the smoothers. maxLevelKrylovSteps maximum of krylov iterations on a particular level (default=3) minDefectReduction minimal defect reduction during the krylov iterations on a particular level (default=1e-1) pinfo The information about the parallel distribution of the data.
◆ KAMG() [2/2]
Dune::Amg::KAMG< M, X, S, P, K, A >::KAMG | ( | OperatorHierarchy & | matrices, |
CoarseSolver & | coarseSolver, | ||
const SmootherArgs & | smootherArgs, | ||
const Parameters & | parms, | ||
std::size_t | maxLevelKrylovSteps = 3 , |
||
double | minDefectReduction = 1e-1 |
||
) |
Construct a new amg with a specific coarse solver.
- Parameters
-
matrices The already set-up matrix hierarchy. coarseSolver The set up solver to use on the coarse grid, must match the coarse matrix in the matrix hierarchy. smootherArgs The arguments needed for thesmoother to use for pre and post smoothing. parms The parameters for the AMG. maxLevelKrylovSteps maximum of krylov iterations on a particular level (default=3) minDefectReduction minimal defect reduction during the krylov iterations on a particular level (default=1e-1)
◆ levels() [1/2]
std::size_t Dune::Amg::Hierarchy< T, A >::levels |
Get the number of levels in the hierarchy.
- Returns
- The number of levels.
Referenced by Dune::Amg::MatrixHierarchy< M, PI, A >::coarsenSmoother().
◆ levels() [2/2]
std::size_t Dune::Amg::MatrixHierarchy< M, IS, A >::levels |
Get the number of levels in the hierarchy.
- Returns
- The number of levels.
◆ matrices()
const MatrixHierarchy< M, IS, A >::ParallelMatrixHierarchy & Dune::Amg::MatrixHierarchy< M, IS, A >::matrices |
Get the matrix hierarchy.
- Returns
- The matrix hierarchy.
◆ MatrixHierarchy()
Dune::Amg::MatrixHierarchy< M, IS, A >::MatrixHierarchy | ( | std::shared_ptr< MatrixOperator > | fineMatrix, |
std::shared_ptr< ParallelInformation > | pinfo = std::make_shared<ParallelInformation>() |
||
) |
Constructor.
- Parameters
-
fineMatrix The matrix to coarsen. pinfo The information about the parallel data decomposition at the first level.
References Dune::SolverCategory::category(), and DUNE_THROW.
◆ maxAggregateSize()
|
inline |
Get the maximum number of nodes a aggregate is allowed to have.
- Returns
- The maximum number of nodes.
◆ maxConnectivity()
|
inline |
Get the maximum number of connections a aggregate is allowed to have.
This limit exists to achieve sparsity of the coarse matrix. the default value is 15.
- Returns
- The maximum number of connections a aggregate is allowed to have.
◆ maxDistance()
|
inline |
Get the maximal distance allowed between two nodes in a aggregate.
The distance between two nodes in a aggregate is the minimal number of edges it takes to travel from one node to the other without leaving the aggregate.
- Returns
- The maximum distance allowed.
◆ maxlevels()
std::size_t Dune::Amg::MatrixHierarchy< M, IS, A >::maxlevels |
Get the max number of levels in the hierarchy of processors.
- Returns
- The maximum number of levels.
◆ minAggregateSize()
|
inline |
Get the minimum number of nodes a aggregate has to consist of.
- Returns
- The minimum number of nodes.
◆ operator()() [1/7]
void Dune::Amg::BaseConnectivityConstructor::ConnectedBuilder< G, S, V >::operator() | ( | const ConstEdgeIterator & | edge | ) |
Process an edge pointing to another aggregate.
- Parameters
-
edge The iterator positioned at the edge.
References Dune::GeometryTypes::vertex.
◆ operator()() [2/7]
|
inline |
compute the norm of a matrix.
- Parameters
-
m The matrix row to compute the norm of.
◆ operator()() [3/7]
|
inline |
compute the norm of a matrix.
- Parameters
-
m The matrix row to compute the norm of.
◆ operator()() [4/7]
|
inline |
compute the norm of a matrix.
- Parameters
-
m The matrix row to compute the norm of.
◆ operator()() [5/7]
|
inline |
Compute the norm of a scalar.
- Parameters
-
m The scalar to compute the norm of
◆ operator()() [6/7]
|
inline |
compute the norm of a matrix.
- Parameters
-
m The matrix to compute the norm of
◆ operator()() [7/7]
void Dune::Amg::Aggregator< G >::AggregateVisitor< V >::operator() | ( | const typename MatrixGraph::ConstEdgeIterator & | edge | ) |
Examine an edge.
The edge will be examined by the adapted visitor if it belongs to the right aggregate.
◆ operator[]() [1/4]
|
inline |
Get the property for a key.
- Parameters
-
key The key.
- Returns
- The corresponding property.
References Dune::Amg::RandomAccessBundledPropertyMap< C, K, i, T, R >::index.
◆ operator[]() [2/4]
|
inline |
Get the properties associated to a vertex.
- Parameters
-
vertex The vertex whose Properties we want.
References Dune::GeometryTypes::vertex.
◆ operator[]() [3/4]
AggregateDescriptor & Dune::Amg::AggregatesMap< V >::operator[] | ( | const VertexDescriptor & | v | ) |
Get the aggregate a vertex belongs to.
- Parameters
-
v The vertex we want to know the aggregate of.
- Returns
- The aggregate the vertex is mapped to.
◆ operator[]() [4/4]
const AggregateDescriptor & Dune::Amg::AggregatesMap< V >::operator[] | ( | const VertexDescriptor & | v | ) | const |
Get the aggregate a vertex belongs to.
- Parameters
-
v The vertex we want to know the aggregate of.
- Returns
- The aggregate the vertex is mapped to.
◆ parallelInformation()
const MatrixHierarchy< M, IS, A >::ParallelInformationHierarchy & Dune::Amg::MatrixHierarchy< M, IS, A >::parallelInformation |
Get the hierarchy of the parallel data distribution information.
- Returns
- The hierarchy of the parallel data distribution information.
◆ Parameters()
|
inline |
Constructor.
- Parameters
-
maxLevel The maximum number of levels allowed in the matrix hierarchy (default: 100). coarsenTarget If the number of nodes in the matrix is below this threshold the coarsening will stop (default: 1000). minCoarsenRate If the coarsening rate falls below this threshold the coarsening will stop (default: 1.2) prolongDamp The damping factor to apply to the prolongated update (default: 1.6) accumulate Whether to accumulate the data onto fewer processors on coarser levels.
◆ post() [1/2]
|
virtual |
Clean up.
This method is called after the last apply call for the linear system to be solved. Memory may be deallocated safely here. x is the solution of the linear equation.
- Parameters
-
x The right hand side of the equation.
Implements Dune::Preconditioner< X, X >.
References Dune::Amg::Hierarchy< T, A >::coarsest().
◆ post() [2/2]
|
virtual |
Clean up.
This method is called after the last apply call for the linear system to be solved. Memory may be deallocated safely here. x is the solution of the linear equation.
- Parameters
-
x The right hand side of the equation.
Implements Dune::Preconditioner< X, X >.
◆ postsmooth()
void Dune::Amg::postsmooth | ( | LevelContext & | levelContext, |
size_t | steps | ||
) |
Apply post smoothing on the current level.
- Parameters
-
levelContext the iterators of the current level. steps The number of smoothing steps to apply.
References Dune::Amg::SmootherApplier< T >::postSmooth().
Referenced by Dune::Amg::KAmgTwoGrid< AMG >::apply().
◆ postSmooth()
|
inlinestatic |
apply post smoothing in forward direction
- Parameters
-
smoother The smoother to use. d The current defect. v handle to store the update in.
Referenced by Dune::Amg::postsmooth().
◆ pre() [1/2]
void Dune::Amg::AMG< M, X, S, PI, A >::pre | ( | Domain & | x, |
Range & | b | ||
) |
Prepare the preconditioner.
A solver solves a linear operator equation A(x)=b by applying one or several steps of the preconditioner. The method pre() is called before the first apply operation. b and x are right hand side and solution vector of the linear system respectively. It may. e.g., scale the system, allocate memory or compute a (I)LU decomposition. Note: The ILU decomposition could also be computed in the constructor or with a separate method of the derived method if several linear systems with the same matrix are to be solved.
- Note
- if a preconditioner is copied (e.g. for a second thread) again the pre() method has to be called to ensure proper memory management.
- Parameters
-
x The left hand side of the equation. b The right hand side of the equation.
◆ pre() [2/2]
void Dune::Amg::KAMG< M, X, S, P, K, A >::pre | ( | Domain & | x, |
Range & | b | ||
) |
Prepare the preconditioner.
A solver solves a linear operator equation A(x)=b by applying one or several steps of the preconditioner. The method pre() is called before the first apply operation. b and x are right hand side and solution vector of the linear system respectively. It may. e.g., scale the system, allocate memory or compute a (I)LU decomposition. Note: The ILU decomposition could also be computed in the constructor or with a separate method of the derived method if several linear systems with the same matrix are to be solved.
- Note
- if a preconditioner is copied (e.g. for a second thread) again the pre() method has to be called to ensure proper memory management.
- Parameters
-
x The left hand side of the equation. b The right hand side of the equation.
◆ presmooth()
void Dune::Amg::presmooth | ( | LevelContext & | levelContext, |
size_t | steps | ||
) |
Apply pre smoothing on the current level.
- Parameters
-
levelContext the iterators of the current level. steps The number of smoothing steps to apply.
References Dune::Amg::SmootherApplier< T >::preSmooth().
Referenced by Dune::Amg::KAmgTwoGrid< AMG >::apply().
◆ preSmooth()
|
inlinestatic |
apply pre smoothing in forward direction
- Parameters
-
smoother The smoother to use. d The current defect. v handle to store the update in.
Referenced by Dune::Amg::presmooth().
◆ PropertyGraphVertexPropertyMap()
|
inline |
Constructor.
- Parameters
-
g The graph whose properties we access.
◆ RandomAccessBundledPropertyMap()
|
inline |
Constructor.
- Parameters
-
container The container with the property bundle.
◆ recalculateGalerkin()
void Dune::Amg::MatrixHierarchy< M, IS, A >::recalculateGalerkin | ( | const F & | copyFlags | ) |
Recalculate the galerkin products.
If the data of the fine matrix changes but not its sparsity pattern this will recalculate all coarser levels without starting the expensive aggregation process all over again.
◆ recalculateHierarchy()
|
inline |
Recalculate the matrix hierarchy.
It is assumed that the coarsening for the changed fine level matrix would yield the same aggregates. In this case it suffices to recalculate all the Galerkin products for the matrices of the coarser levels.
◆ reconstruct()
void Dune::Amg::Aggregate< G, S >::reconstruct | ( | const Vertex & | vertex | ) |
Reconstruct the aggregat from an seed node.
Will determine all vertices of the same aggregate and reference those.
◆ redistributeInformation()
const MatrixHierarchy< M, IS, A >::RedistributeInfoList & Dune::Amg::MatrixHierarchy< M, IS, A >::redistributeInformation |
Get the hierarchy of the information about redistributions,.
- Returns
- The hierarchy of the information about redistributions of the data to fewer processes.
◆ setAdditive()
|
inline |
Set whether to use additive multigrid.
- Parameters
-
additive True if multigrid should be additive.
◆ setDebugLevel()
|
inline |
Set the debugging level.
- Parameters
-
level If 0 no debugging output will be generated.
- Warning
- In parallel the level has to be consistent over all procceses.
◆ setDefaultValuesAnisotropic() [1/2]
|
inline |
Sets reasonable default values for an aisotropic problem.
Reasonable means that we should end up with cube aggregates with sides of diameter 2 and sides in one dimension that are longer (e.g. for 3D: 2x2x3).
- Parameters
-
dim The dimension of the problem. diameter The preferred diameter for the aggregation.
References Dune::Amg::AggregationCriterion< T >::setDefaultValuesIsotropic().
◆ setDefaultValuesAnisotropic() [2/2]
|
inline |
Sets reasonable default values for an anisotropic problem.
Reasonable means that we should end up with cube aggregates with sides of diameter 2 and sides in one dimension that are longer (e.g. for 3D: 2x2x3).
- Parameters
-
dim The dimension of the problem. diameter The preferred diameter for the aggregation.
References Dune::Amg::AggregationParameters::setDefaultValuesIsotropic().
◆ setDefaultValuesIsotropic() [1/2]
|
inline |
Sets reasonable default values for an isotropic problem.
Reasonable means that we should end up with cube aggregates of diameter 2.
- Parameters
-
dim The dimension of the problem. diameter The preferred diameter for the aggregation.
Referenced by Dune::Amg::AggregationCriterion< T >::setDefaultValuesAnisotropic().
◆ setDefaultValuesIsotropic() [2/2]
|
inline |
Sets reasonable default values for an isotropic problem.
Reasonable means that we should end up with cube aggregates of diameter 2.
- Parameters
-
dim The dimension of the problem. diameter The preferred diameter for the aggregation.
Referenced by Dune::Amg::AggregationParameters::setDefaultValuesAnisotropic().
◆ setIsolated()
|
inline |
Marks that node as being isolated.
A node is isolated if it ha not got any strong connections to other nodes in the matrix graph.
◆ setMaxAggregateSize()
|
inline |
Set the maximum number of nodes a aggregate is allowed to have.
The default values is 6.
- Parameters
-
size The maximum number of nodes.
References Dune::size().
◆ setMaxConnectivity()
|
inline |
Set the maximum number of connections a aggregate is allowed to have.
This limit exists to achieve sparsity of the coarse matrix. the default value is 15.
- Parameters
-
connectivity The maximum number of connections a aggregate is allowed to have.
◆ setMaxDistance()
|
inline |
Set the maximal distance allowed between two nodes in a aggregate.
The distance between two nodes in a aggregate is the minimal number of edges it takes to travel from one node to the other without leaving the aggregate. The default value is 2.
- Parameters
-
distance The maximum distance allowed.
◆ setMinAggregateSize()
|
inline |
Set the minimum number of nodes a aggregate has to consist of.
the default value is 4.
References Dune::size().
◆ setMinCoarsenRate()
|
inline |
Set the minimum coarsening rate to be achieved in each coarsening.
The default value is 1.2
◆ setNoPostSmoothSteps()
|
inline |
Set the number of postsmoothing steps to apply.
- Parameters
-
steps The number of steps:
◆ setNoPreSmoothSteps()
|
inline |
Set the number of presmoothing steps to apply.
- Parameters
-
steps The number of steps:
◆ setProlongationDampingFactor()
|
inline |
Set the damping factor for the prolongation.
- Parameters
-
d The new damping factor.
◆ setSkipIsolated()
|
inline |
Set whether isolated aggregates will not be represented on the coarse level.
- Parameters
-
skip True if these aggregates will be skipped.
◆ skipIsolated()
|
inline |
Whether isolated aggregates will not be represented on the coarse level.
- Returns
- True if these aggregates will be skipped.
◆ SparsityBuilder()
Dune::Amg::SparsityBuilder< M >::SparsityBuilder | ( | M & | matrix | ) |
Constructor.
- Parameters
-
matrix The matrix whose sparsity pattern we should set up.
◆ usesDirectCoarseLevelSolver()
bool Dune::Amg::AMG< M, X, S, PI, A >::usesDirectCoarseLevelSolver |
Check whether the coarse solver used is a direct solver.
- Returns
- True if the coarse level solver is a direct solver.
◆ visitNeighbours()
int Dune::Amg::visitNeighbours | ( | const G & | graph, |
const typename G::VertexDescriptor & | vertex, | ||
V & | visitor | ||
) |
Visit all neighbour vertices of a vertex in a graph.
- Parameters
-
graph The graph whose vertices to visit. vertex The vertex whose neighbours we want to visit. visitor The visitor evaluated for each EdgeIterator (by its method operator()(const ConstEdgeIterator& edge)
- Returns
- The number of neighbours of the vertex.