3 #ifndef DUNE_PDELAB_BACKEND_ISTL_PARALLELHELPER_HH
4 #define DUNE_PDELAB_BACKEND_ISTL_PARALLELHELPER_HH
8 #include <dune/common/deprecated.hh>
9 #include <dune/common/parallel/mpihelper.hh>
10 #include <dune/common/stdstreams.hh>
12 #include <dune/istl/owneroverlapcopy.hh>
13 #include <dune/istl/solvercategory.hh>
14 #include <dune/istl/operators.hh>
15 #include <dune/istl/solvers.hh>
16 #include <dune/istl/preconditioners.hh>
17 #include <dune/istl/scalarproducts.hh>
18 #include <dune/istl/paamg/amg.hh>
19 #include <dune/istl/paamg/pinfo.hh>
20 #include <dune/istl/io.hh>
21 #include <dune/istl/superlu.hh>
44 template<
typename GFS>
49 typedef int RankIndex;
57 typedef typename GFS::Ordering::Traits::ContainerIndex ContainerIndex;
63 , _rank(gfs.gridView().comm().
rank())
71 if (gfs.ordering().containedPartitions() == NonOverlappingPartitionSelector::partition_mask)
75 _interiorBorder_all_interface = InteriorBorder_InteriorBorder_Interface;
76 _all_all_interface = InteriorBorder_InteriorBorder_Interface;
81 _interiorBorder_all_interface = InteriorBorder_All_Interface;
82 _all_all_interface = All_All_Interface;
85 if (_gfs.gridView().comm().size()>1)
91 gdh(_gfs,_ghosts,
false);
92 _gfs.gridView().communicate(gdh,_interiorBorder_all_interface,Dune::ForwardCommunication);
98 _gfs.gridView().communicate(pdh,_interiorBorder_all_interface,Dune::ForwardCommunication);
116 template<
typename X,
typename Mask>
119 typename Mask::const_iterator mask_it = mask.begin();
120 for (
typename X::iterator it = x.begin(),
128 template<
typename X,
typename Mask>
131 typename Mask::const_iterator mask_it = mask.begin();
132 for (
typename X::iterator it = x.begin(),
136 *it = (*mask_it == _rank ? *it :
typename X::field_type(0));
142 bool owned(
const ContainerIndex& i)
const
144 return _ranks[i] == _rank;
154 template<
typename X,
typename Y>
155 typename PromotionTraits<
156 typename X::field_type,
157 typename Y::field_type
173 template<
typename X,
typename Y,
typename Mask>
174 typename PromotionTraits<
175 typename X::field_type,
176 typename Y::field_type
180 typedef typename PromotionTraits<
181 typename X::field_type,
182 typename Y::field_type
183 >::PromotedType result_type;
187 typename Y::const_iterator y_it = y.begin();
188 typename Mask::const_iterator mask_it = mask.begin();
189 for (
typename X::const_iterator x_it = x.begin(),
192 ++x_it, ++y_it, ++mask_it)
200 template<
typename X,
typename Y,
typename Mask>
201 typename PromotionTraits<
202 typename X::field_type,
203 typename Y::field_type
207 typedef typename PromotionTraits<
208 typename X::field_type,
209 typename Y::field_type
210 >::PromotedType result_type;
214 typename Y::const_iterator y_it = y.begin();
215 typename Mask::const_iterator mask_it = mask.begin();
216 for (
typename X::const_iterator x_it = x.begin(),
219 ++x_it, ++y_it, ++mask_it)
220 r += (*mask_it == _rank ? Dune::dot(*x_it,*y_it) : result_type(0));
252 template<
typename MatrixType,
typename Comm>
260 bool owned_for_amg(std::size_t i)
const
268 bool is_ghost_for_amg(std::size_t i)
const
278 const RankIndex _rank;
284 InterfaceType _interiorBorder_all_interface;
287 InterfaceType _all_all_interface;
292 template<
typename GFS>
293 template<
typename M,
typename C>
299 const bool is_bcrs_matrix =
307 const bool block_type_is_field_matrix =
316 static_assert(is_bcrs_matrix && block_type_is_field_matrix,
"matrix structure not compatible with AMG");
325 typedef typename GFS::Traits::GridViewType GV;
326 typedef typename RankVector::size_type size_type;
327 const GV& gv = _gfs.gridView();
330 const bool need_communication = _gfs.gridView().comm().size() > 1;
334 BoolVector sharedDOF(_gfs,
false);
336 if (need_communication)
339 _gfs.gridView().communicate(data_handle,_all_all_interface,Dune::ForwardCommunication);
343 typedef typename C::ParallelIndexSet::GlobalIndex GlobalIndex;
344 GlobalIndex count = 0;
346 for (size_type i = 0; i < sharedDOF.N(); ++i)
347 if (owned_for_amg(i) &&
native(sharedDOF)[i][0])
350 dverb << gv.comm().rank() <<
": shared block count is " << count.touint() << std::endl;
353 std::vector<GlobalIndex> counts(_gfs.gridView().comm().size());
354 _gfs.gridView().comm().allgather(&count, 1, &(counts[0]));
357 GlobalIndex start = std::accumulate(counts.begin(),counts.begin() + _rank,GlobalIndex(0));
360 GIVector scalarIndices(_gfs, std::numeric_limits<GlobalIndex>::max());
362 for (size_type i = 0; i < sharedDOF.N(); ++i)
363 if (owned_for_amg(i) &&
native(sharedDOF)[i][0])
365 native(scalarIndices)[i][0] = start;
370 if (need_communication)
373 _gfs.gridView().communicate(data_handle,_interiorBorder_all_interface,Dune::ForwardCommunication);
377 c.indexSet().beginResize();
378 for (size_type i=0; i<scalarIndices.N(); ++i)
380 Dune::OwnerOverlapCopyAttributeSet::AttributeSet attr;
381 if(
native(scalarIndices)[i][0] != std::numeric_limits<GlobalIndex>::max())
384 if (owned_for_amg(i))
387 attr = Dune::OwnerOverlapCopyAttributeSet::owner;
389 else if (is_ghost_for_amg(i) && c.getSolverCategory() ==
static_cast<int>(SolverCategory::nonoverlapping))
392 attr = Dune::OwnerOverlapCopyAttributeSet::overlap;
396 attr = Dune::OwnerOverlapCopyAttributeSet::copy;
398 c.indexSet().add(
native(scalarIndices)[i][0],
typename C::ParallelIndexSet::LocalIndex(i,attr));
401 c.indexSet().endResize();
404 std::set<int> neighbors;
406 if (need_communication)
409 _gfs.gridView().communicate(data_handle,_all_all_interface,Dune::ForwardCommunication);
412 c.remoteIndices().setNeighbours(neighbors);
413 c.remoteIndices().template rebuild<false>();
418 template<
int s,
bool isFakeMPIHelper>
421 typedef Dune::Amg::SequentialInformation
type;
431 typedef OwnerOverlapCopyCommunication<bigunsignedint<s>,
int>
type;
442 #endif // DUNE_PDELAB_BACKEND_ISTL_PARALLELHELPER_HH
GatherScatter data handle for creating a disjoint DOF partitioning.
Definition: genericdatahandle.hh:930
bool owned(const ContainerIndex &i) const
Tests whether the given index is owned by this process.
Definition: parallelhelper.hh:142
tags::container< T >::type container_tag(const T &)
Gets instance of container tag associated with T.
Definition: backend/istl/tags.hh:246
Data handle for marking shared DOFs.
Definition: genericdatahandle.hh:1014
Tag describing a BlockVector.
Definition: backend/istl/tags.hh:23
bool isGhost(const ContainerIndex &i) const
Tests whether the given index belongs to a ghost DOF.
Definition: parallelhelper.hh:148
void maskForeignDOFs(X &x) const
Mask out all DOFs not owned by the current process with 0.
Definition: parallelhelper.hh:106
typename impl::BackendVectorSelector< GridFunctionSpace, FieldType >::Type Vector
alias of the return type of BackendVectorSelector
Definition: backend/interface.hh:112
Data handle for collecting set of neighboring MPI ranks.
Definition: genericdatahandle.hh:1060
Definition: parallelhelper.hh:45
static const unsigned int value
Definition: gridfunctionspace/tags.hh:175
Tag describing a BCRSMatrix.
Definition: backend/istl/tags.hh:60
ParallelHelper(const GFS &gfs, int verbose=1)
Definition: parallelhelper.hh:61
std::enable_if< std::is_base_of< impl::WrapperBase, T >::value, Native< T > & >::type native(T &t)
Definition: backend/interface.hh:198
Definition: genericdatahandle.hh:716
OwnerOverlapCopyCommunication< bigunsignedint< s >, int > type
Definition: parallelhelper.hh:431
Tag describing an arbitrary FieldVector.
Definition: backend/istl/tags.hh:43
Definition: adaptivity.hh:27
Dune::Amg::SequentialInformation type
Definition: parallelhelper.hh:421
typename native_type< T >::type Native
Alias of the native container type associated with T or T itself if it is not a backend wrapper...
Definition: backend/interface.hh:182
void createIndexSetAndProjectForAMG(MatrixType &m, Comm &c)
Makes the matrix consistent and creates the parallel information for AMG.
Tag describing an arbitrary FieldMatrix.
Definition: backend/istl/tags.hh:80
PromotionTraits< typename X::field_type, typename Y::field_type >::PromotedType disjointDot(const X &x, const Y &y) const
Calculates the (rank-local) dot product of x and y on the disjoint partition defined by the helper...
Definition: parallelhelper.hh:159
const std::string s
Definition: function.hh:1102
RankIndex rank() const
Returns the MPI rank of this process.
Definition: parallelhelper.hh:228
Definition: parallelhelper.hh:419
Data handle for marking ghost DOFs.
Definition: genericdatahandle.hh:812
Extracts the container tag from T.
Definition: backend/istl/tags.hh:142