7 #include <dune/istl/owneroverlapcopy.hh> 8 #include <dune/istl/solvercategory.hh> 9 #include <dune/istl/paamg/pinfo.hh> 11 #include <amdis/Environment.hpp> 12 #include <amdis/functions/GlobalIdSet.hpp> 13 #include <amdis/linearalgebra/Communication.hpp> 14 #include <amdis/linearalgebra/ParallelIndexSet.hpp> 19 template <
class G,
class L>
22 using Impl = Dune::Amg::SequentialInformation;
29 template <
class Basis>
32 assert(cat == Dune::SolverCategory::sequential);
35 typename Dune::SolverCategory::Category category()
const 37 return Dune::SolverCategory::sequential;
40 Impl const&
get()
const 50 template <
class Basis>
51 void update(Basis
const&) { }
59 template <
class GlobalId,
class LocalIndex>
63 using Impl = Dune::OwnerOverlapCopyCommunication<GlobalId, LocalIndex>;
67 using IndexSet =
typename Impl::ParallelIndexSet;
68 using RemoteIndices =
typename Impl::RemoteIndices;
71 template <
class Basis>
78 IndexSet
const& indexSet()
const 81 return impl_->indexSet();
84 RemoteIndices
const& remoteIndices()
const 87 return impl_->remoteIndices();
90 typename Dune::SolverCategory::Category category()
const 95 Impl const&
get()
const 106 template <
class Basis>
107 void update(Basis
const& basis)
111 auto& pis = impl_->indexSet();
114 auto& ris = impl_->remoteIndices();
115 ris.setIndexSets(pis, pis, impl_->communicator());
116 ris.template rebuild<true>();
120 typename Dune::SolverCategory::Category cat_;
121 std::unique_ptr<Impl> impl_ =
nullptr;
126 template <
class G,
class L>
136 template <
class G,
class L>
149 template <
class Basis>
150 static Communication create(Basis
const& basis, std::string
const& prefix =
"");
155 #include <amdis/linearalgebra/istl/Communication.inc.hpp>
Dummy implementation for ISTL-specific communication when no MPI is found.
Definition: Communication.hpp:20
Contains all classes needed for solving linear and non linear equation systems.
Definition: AdaptBase.hpp:6
auto cat(Dune::TypeTree::HybridTreePath< S... > const &tp0, Dune::TypeTree::HybridTreePath< T... > const &tp1)
Concatenate two treepaths.
Definition: TreePath.hpp:209
Implementation of a creator pattern for Communication types.
Definition: Communication.hpp:47
void buildParallelIndexSet(Basis const &basis, PIS ¶llelIndexSet)
Fills a parallelIndexSet with indices from a basis.
Definition: ParallelIndexSet.hpp:24
static Dune::MPIHelper::MPICommunicator comm()
Return the MPI_Comm object (or a fake communicator)
Definition: Environment.hpp:86