3 #ifndef DUNE_PDELAB_GRIDFUNCTIONSPACE_GENERICDATAHANDLE_HH 4 #define DUNE_PDELAB_GRIDFUNCTIONSPACE_GENERICDATAHANDLE_HH 10 #include<dune/common/exceptions.hh> 11 #include <dune/common/parallel/mpihelper.hh> 13 #include <dune/grid/common/datahandleif.hh> 14 #include <dune/grid/common/gridenums.hh> 38 template<
typename GFS>
41 return gfs.dataHandleContains(codim);
44 template<
typename GFS>
47 return gfs.dataHandleFixedSize(codim);
50 template<
typename GFS,
typename Entity>
51 std::size_t
size(
const GFS& gfs,
const Entity&
e)
const 54 return gfs.dataHandleSize(e) *
sizeof(E) + (gfs.sendLeafSizes() ? TypeTree::TreeInfo<typename GFS::Ordering>::leafCount *
sizeof(
size_type) : 0);
70 template<
typename GFS>
73 return gfs.dataHandleContains(codim);
76 template<
typename GFS>
82 template<
typename GFS,
typename Entity>
83 std::size_t
size(
const GFS& gfs,
const Entity&
e)
const 85 return gfs.dataHandleContains(Entity::codimension) && gfs.entitySet().contains(e) ? _count : 0;
94 const std::size_t _count;
105 template<
typename GFS,
typename V,
typename GatherScatter,
typename CommunicationDescriptor = DOFDataCommunicationDescriptor<
typename V::ElementType> >
107 :
public Dune::CommDataHandleIF<GFSDataHandle<GFS,V,GatherScatter,CommunicationDescriptor>,typename CommunicationDescriptor::DataType>
112 typedef typename CommunicationDescriptor::DataType
DataType;
115 static const size_type leaf_count = TypeTree::TreeInfo<typename GFS::Ordering>::leafCount;
117 GFSDataHandle(
const GFS& gfs, V& v, GatherScatter gather_scatter = GatherScatter(), CommunicationDescriptor communication_descriptor = CommunicationDescriptor())
121 , _gather_scatter(gather_scatter)
122 , _communication_descriptor(communication_descriptor)
128 return _communication_descriptor.contains(_gfs,dim,codim);
134 return _communication_descriptor.fixedSize(_gfs,dim,codim);
141 template<
typename Entity>
142 size_type
size(
const Entity&
e)
const 144 return _communication_descriptor.size(_gfs,e);
148 template<
typename MessageBuffer,
typename Entity>
149 typename std::enable_if<
152 gather(MessageBuffer& buff,
const Entity&
e)
const 155 _index_cache.update(e);
156 _local_view.bind(_index_cache);
157 if (_gfs.sendLeafSizes())
160 for (
auto it = _index_cache.offsets().begin() + 1,
161 end_it = _index_cache.offsets().end();
165 buf_wrapper.
write(static_cast<typename CommunicationDescriptor::size_type>(*it));
169 if (_gather_scatter.gather(buf_wrapper,e,_local_view))
170 _local_view.commit();
171 _local_view.unbind();
175 template<
typename MessageBuffer,
typename Entity>
176 typename std::enable_if<
179 gather(MessageBuffer& buff,
const Entity&
e)
const 181 _index_cache.update(e);
182 _local_view.bind(_index_cache);
183 if (_gather_scatter.gather(buff,e,_local_view))
184 _local_view.commit();
185 _local_view.unbind();
194 template<
typename MessageBuffer,
typename Entity>
195 typename std::enable_if<
198 scatter(MessageBuffer& buff,
const Entity&
e, size_type n)
201 _index_cache.update(e);
202 _local_view.bind(_index_cache);
203 bool needs_commit =
false;
204 if (_gfs.sendLeafSizes())
208 for (
auto it = remote_offsets.begin() + 1,
209 end_it = remote_offsets.end();
213 typename CommunicationDescriptor::size_type data = 0;
214 buf_wrapper.
read(data);
218 needs_commit = _gather_scatter.scatter(buf_wrapper,remote_offsets,_index_cache.offsets(),
e,_local_view);
223 needs_commit = _gather_scatter.scatter(buf_wrapper,n /
sizeof(
typename CommunicationDescriptor::OriginalDataType),e,_local_view);
227 _local_view.commit();
229 _local_view.unbind();
238 template<
typename MessageBuffer,
typename Entity>
239 typename std::enable_if<
242 scatter(MessageBuffer& buff,
const Entity&
e, size_type n)
244 _index_cache.update(e);
245 _local_view.bind(_index_cache);
247 if (_gather_scatter.scatter(buff,n,e,_local_view))
248 _local_view.commit();
250 _local_view.unbind();
257 typedef typename V::template LocalView<IndexCache> LocalView;
260 mutable IndexCache _index_cache;
261 mutable LocalView _local_view;
262 mutable GatherScatter _gather_scatter;
263 CommunicationDescriptor _communication_descriptor;
268 template<
typename GatherScatter>
276 template<
typename MessageBuffer,
typename Entity,
typename LocalView>
277 bool gather(MessageBuffer& buff,
const Entity&
e,
const LocalView& local_view)
const 279 for (std::size_t i = 0; i < local_view.size(); ++i)
280 _gather_scatter.gather(buff,local_view[i]);
285 template<
typename MessageBuffer,
typename Entity,
typename LocalView>
286 bool scatter(MessageBuffer& buff, size_type n,
const Entity&
e, LocalView& local_view)
const 288 if (local_view.cache().gridFunctionSpace().entitySet().partitions().contains(e.partitionType()))
290 if (local_view.size() != n)
291 DUNE_THROW(
Exception,
"size mismatch in GridFunctionSpace data handle, have " << local_view.size() <<
"DOFs, but received " << n);
293 for (std::size_t i = 0; i < local_view.size(); ++i)
294 _gather_scatter.scatter(buff,local_view[i]);
299 if (local_view.size() != 0)
300 DUNE_THROW(
Exception,
"expected no DOFs in partition '" << e.partitionType() <<
"', but have " << local_view.size());
302 for (std::size_t i = 0; i < local_view.size(); ++i)
304 typename LocalView::ElementType dummy;
312 template<
typename MessageBuffer,
typename Offsets,
typename Entity,
typename LocalView>
313 bool scatter(MessageBuffer& buff,
const Offsets& remote_offsets,
const Offsets& local_offsets,
const Entity&
e, LocalView& local_view)
const 315 if (local_view.cache().gridFunctionSpace().entitySet().partitions().contains(e.partitionType()))
324 size_type remote_i = 0;
325 size_type local_i = 0;
326 bool needs_commit =
false;
327 for (size_type block = 1; block < local_offsets.size(); ++block)
330 if (remote_offsets[block] == remote_i)
332 local_i = local_offsets[block];
337 if (local_offsets[block] == local_i)
339 for (; remote_i < remote_offsets[block]; ++remote_i)
341 typename LocalView::ElementType dummy;
347 if (remote_offsets[block] - remote_i != local_offsets[block] - local_i)
348 DUNE_THROW(
Exception,
"size mismatch in GridFunctionSpace data handle block " << block <<
", have " << local_offsets[block] - local_i <<
"DOFs, but received " << remote_offsets[block] - remote_i);
350 for (; local_i < local_offsets[block]; ++local_i)
351 _gather_scatter.scatter(buff,local_view[local_i]);
353 remote_i = remote_offsets[block];
360 if (local_view.size() != 0)
361 DUNE_THROW(
Exception,
"expected no DOFs in partition '" << e.partitionType() <<
"', but have " << local_view.size());
363 for (std::size_t i = 0; i < remote_offsets.back(); ++i)
365 typename LocalView::ElementType dummy;
373 : _gather_scatter(gather_scatter)
378 GatherScatter _gather_scatter;
383 template<
typename GatherScatter>
391 template<
typename MessageBuffer,
typename Entity,
typename LocalView>
392 bool gather(MessageBuffer& buff,
const Entity&
e,
const LocalView& local_view)
const 394 for (std::size_t i = 0; i < local_view.size(); ++i)
395 _gather_scatter.gather(buff,e,local_view[i]);
400 template<
typename MessageBuffer,
typename Entity,
typename LocalView>
401 bool scatter(MessageBuffer& buff, size_type n,
const Entity&
e, LocalView& local_view)
const 403 if (local_view.cache().gridFunctionSpace().entitySet().partitions().contains(e.partitionType()))
405 if (local_view.size() != n)
406 DUNE_THROW(
Exception,
"size mismatch in GridFunctionSpace data handle, have " << local_view.size() <<
"DOFs, but received " << n);
408 for (std::size_t i = 0; i < local_view.size(); ++i)
409 _gather_scatter.scatter(buff,e,local_view[i]);
414 if (local_view.size() != 0)
415 DUNE_THROW(
Exception,
"expected no DOFs in partition '" << e.partitionType() <<
"', but have " << local_view.size());
417 for (std::size_t i = 0; i < local_view.size(); ++i)
419 typename LocalView::ElementType dummy;
427 template<
typename MessageBuffer,
typename Offsets,
typename Entity,
typename LocalView>
428 bool scatter(MessageBuffer& buff,
const Offsets& remote_offsets,
const Offsets& local_offsets,
const Entity&
e, LocalView& local_view)
const 430 if (local_view.cache().gridFunctionSpace().entitySet().partitions().contains(e.partitionType()))
432 size_type remote_i = 0;
433 size_type local_i = 0;
434 bool needs_commit =
false;
435 for (size_type block = 1; block < local_offsets.size(); ++block)
439 if (remote_offsets[block] == remote_i)
441 local_i = local_offsets[block];
446 if (local_offsets[block] == local_i)
448 for (; remote_i < remote_offsets[block]; ++remote_i)
450 typename LocalView::ElementType dummy;
456 if (remote_offsets[block] - remote_i != local_offsets[block] - local_i)
457 DUNE_THROW(
Exception,
"size mismatch in GridFunctionSpace data handle block " << block <<
", have " << local_offsets[block] - local_i <<
"DOFs, but received " << remote_offsets[block] - remote_i);
459 for (; local_i < local_offsets[block]; ++local_i)
460 _gather_scatter.scatter(buff,e,local_view[local_i]);
462 remote_i = remote_offsets[block];
469 if (local_view.size() != 0)
470 DUNE_THROW(
Exception,
"expected no DOFs in partition '" << e.partitionType() <<
"', but have " << local_view.size());
472 for (std::size_t i = 0; i < remote_offsets.back(); ++i)
474 typename LocalView::ElementType dummy;
482 : _gather_scatter(gather_scatter)
487 GatherScatter _gather_scatter;
492 template<
typename GatherScatter>
500 template<
typename MessageBuffer,
typename Entity,
typename LocalView>
501 bool gather(MessageBuffer& buff,
const Entity&
e,
const LocalView& local_view)
const 503 for (std::size_t i = 0; i < local_view.size(); ++i)
504 _gather_scatter.gather(buff,local_view.cache().containerIndex(i),local_view[i]);
509 template<
typename MessageBuffer,
typename Entity,
typename LocalView>
510 bool scatter(MessageBuffer& buff, size_type n,
const Entity&
e, LocalView& local_view)
const 512 if (local_view.cache().gridFunctionSpace().entitySet().partitions().contains(e.partitionType()))
514 if (local_view.size() != n)
515 DUNE_THROW(
Exception,
"size mismatch in GridFunctionSpace data handle, have " << local_view.size() <<
"DOFs, but received " << n);
517 for (std::size_t i = 0; i < local_view.size(); ++i)
518 _gather_scatter.scatter(buff,local_view.cache().containerIndex(i),local_view[i]);
524 if (local_view.size() != 0)
525 DUNE_THROW(
Exception,
"expected no DOFs in partition '" << e.partitionType() <<
"', but have " << local_view.size());
527 for (std::size_t i = 0; i < local_view.size(); ++i)
529 typename LocalView::ElementType dummy;
537 template<
typename MessageBuffer,
typename Offsets,
typename Entity,
typename LocalView>
538 bool scatter(MessageBuffer& buff,
const Offsets& remote_offsets,
const Offsets& local_offsets,
const Entity&
e, LocalView& local_view)
const 540 if (local_view.cache().gridFunctionSpace().entitySet().partitions().contains(e.partitionType()))
542 size_type remote_i = 0;
543 size_type local_i = 0;
544 bool needs_commit =
false;
545 for (size_type block = 1; block < local_offsets.size(); ++block)
549 if (remote_offsets[block] == remote_i)
551 local_i = local_offsets[block];
556 if (local_offsets[block] == local_i)
558 for (; remote_i < remote_offsets[block]; ++remote_i)
560 typename LocalView::ElementType dummy;
566 if (remote_offsets[block] - remote_i != local_offsets[block] - local_i)
567 DUNE_THROW(
Exception,
"size mismatch in GridFunctionSpace data handle block " << block <<
", have " << local_offsets[block] - local_i <<
"DOFs, but received " << remote_offsets[block] - remote_i);
569 for (; local_i < local_offsets[block]; ++local_i)
570 _gather_scatter.scatter(buff,local_view.cache().containerIndex(local_i),local_view[local_i]);
572 remote_i = remote_offsets[block];
579 if (local_view.size() != 0)
580 DUNE_THROW(
Exception,
"expected no DOFs in partition '" << e.partitionType() <<
"', but have " << local_view.size());
582 for (std::size_t i = 0; i < remote_offsets.back(); ++i)
584 typename LocalView::ElementType dummy;
593 : _gather_scatter(gather_scatter)
598 GatherScatter _gather_scatter;
606 template<
class MessageBuffer,
class DataType>
607 void gather (MessageBuffer& buff, DataType& data)
const 612 template<
class MessageBuffer,
class DataType>
613 void scatter (MessageBuffer& buff, DataType& data)
const 621 template<
class GFS,
class V>
623 :
public GFSDataHandle<GFS,V,DataGatherScatter<AddGatherScatter> >
637 template<
class MessageBuffer,
class DataType>
638 void gather (MessageBuffer& buff, DataType& data)
const 644 template<
class MessageBuffer,
class DataType>
645 void scatter (MessageBuffer& buff, DataType& data)
const 653 template<
class GFS,
class V>
655 :
public GFSDataHandle<GFS,V,DataGatherScatter<AddClearGatherScatter> >
669 template<
class MessageBuffer,
class DataType>
670 void gather (MessageBuffer& buff, DataType& data)
const 675 template<
class MessageBuffer,
class DataType>
676 void scatter (MessageBuffer& buff, DataType& data)
const 684 template<
class GFS,
class V>
686 :
public GFSDataHandle<GFS,V,DataGatherScatter<CopyGatherScatter> >
700 template<
class MessageBuffer,
class DataType>
701 void gather (MessageBuffer& buff, DataType& data)
const 706 template<
class MessageBuffer,
class DataType>
707 void scatter (MessageBuffer& buff, DataType& data)
const 711 data = std::min(data,x);
715 template<
class GFS,
class V>
717 :
public GFSDataHandle<GFS,V,DataGatherScatter<MinGatherScatter> >
731 template<
class MessageBuffer,
class DataType>
732 void gather (MessageBuffer& buff, DataType& data)
const 737 template<
class MessageBuffer,
class DataType>
738 void scatter (MessageBuffer& buff, DataType& data)
const 742 data = std::max(data,x);
746 template<
class GFS,
class V>
748 :
public GFSDataHandle<GFS,V,DataGatherScatter<MaxGatherScatter> >
772 template<
typename MessageBuffer,
typename Entity,
typename LocalView>
773 bool gather(MessageBuffer& buff,
const Entity&
e, LocalView& local_view)
const 776 const bool ghost = e.partitionType()!=Dune::InteriorEntity && e.partitionType()!=Dune::BorderEntity;
784 template<
typename MessageBuffer,
typename Entity,
typename LocalView>
785 bool scatter(MessageBuffer& buff, std::size_t n,
const Entity&
e, LocalView& local_view)
const 789 const bool ghost = e.partitionType()!=Dune::InteriorEntity && e.partitionType()!=Dune::BorderEntity;
795 for (std::size_t i = 0; i < local_view.size(); ++i)
796 local_view[i] = ghost;
811 template<
class GFS,
class V>
816 EntityDataCommunicationDescriptor<bool> >
826 "GhostDataHandle expects a vector of bool values");
858 template<
typename RankIndex>
864 template<
typename MessageBuffer,
typename Entity,
typename LocalView>
865 bool gather(MessageBuffer& buff,
const Entity&
e, LocalView& local_view)
const 874 template<
typename MessageBuffer,
typename Entity,
typename LocalView>
875 bool scatter(MessageBuffer& buff, std::size_t n,
const Entity&
e, LocalView& local_view)
const 878 const RankIndex unknown_rank = std::numeric_limits<RankIndex>::max();
881 const bool is_interior_or_border = (e.partitionType()==Dune::InteriorEntity || e.partitionType()==Dune::BorderEntity);
884 RankIndex received_rank;
885 buff.read(received_rank);
887 for (std::size_t i = 0; i < local_view.size(); ++i)
890 RankIndex current_rank = local_view[i];
897 if (!is_interior_or_border && current_rank == _rank)
898 current_rank = unknown_rank;
901 local_view[i] = std::min(current_rank,received_rank);
916 const RankIndex _rank;
929 template<
class GFS,
class V>
933 DisjointPartitioningGatherScatter<
934 typename V::ElementType
936 EntityDataCommunicationDescriptor<
937 typename V::ElementType
945 typename V::ElementType
948 typename V::ElementType
968 v_ = gfs_.gridView().comm().rank();
983 template<
typename MessageBuffer,
typename Entity,
typename LocalView>
984 bool gather(MessageBuffer& buff,
const Entity&
e, LocalView& local_view)
const 986 buff.write(local_view.size() > 0);
990 template<
typename MessageBuffer,
typename Entity,
typename LocalView>
991 bool scatter(MessageBuffer& buff, std::size_t n,
const Entity&
e, LocalView& local_view)
const 993 bool remote_entity_has_dofs;
994 buff.read(remote_entity_has_dofs);
996 for (std::size_t i = 0; i < local_view.size(); ++i)
998 local_view[i] |= remote_entity_has_dofs;
1013 template<
class GFS,
class V>
1017 SharedDOFGatherScatter,
1018 EntityDataCommunicationDescriptor<bool> >
1028 "SharedDOFDataHandle expects a vector of bool values");
1059 template<
typename GFS,
typename RankIndex>
1061 :
public Dune::CommDataHandleIF<GFSNeighborDataHandle<GFS,RankIndex>,RankIndex>
1075 , _neighbors(neighbors)
1081 return _gfs.dataHandleContains(codim);
1090 template<
typename Entity>
1096 template<
typename MessageBuffer,
typename Entity>
1097 void gather(MessageBuffer& buff,
const Entity&
e)
const 1102 template<
typename MessageBuffer,
typename Entity>
1103 void scatter(MessageBuffer& buff,
const Entity&
e, size_type n)
1107 _neighbors.insert(rank);
1113 const RankIndex _rank;
1114 std::set<RankIndex>& _neighbors;
1122 #endif // DUNE_PDELAB_GRIDFUNCTIONSPACE_GENERICDATAHANDLE_HH Definition: genericdatahandle.hh:697
std::size_t size_type
Definition: genericdatahandle.hh:389
AddDataHandle(const GFS &gfs_, V &v_)
Definition: genericdatahandle.hh:629
Definition: genericdatahandle.hh:654
void scatter(MessageBuffer &buff, DataType &data) const
Definition: genericdatahandle.hh:676
static const int dim
Definition: adaptivity.hh:83
Definition: genericdatahandle.hh:603
std::size_t size(const GFS &gfs, const Entity &e) const
Definition: genericdatahandle.hh:51
const Entity & e
Definition: localfunctionspace.hh:111
bool scatter(MessageBuffer &buff, const Offsets &remote_offsets, const Offsets &local_offsets, const Entity &e, LocalView &local_view) const
Definition: genericdatahandle.hh:538
void scatter(MessageBuffer &buff, DataType &data) const
Definition: genericdatahandle.hh:738
bool contains(int dim, int codim) const
Definition: genericdatahandle.hh:1078
GhostDataHandle(const GFS &gfs_, V &v_, bool init_vector=true)
Creates a new GhostDataHandle.
Definition: genericdatahandle.hh:840
Definition: genericdatahandle.hh:685
void gather(MessageBuffer &buff, DataType &data) const
Definition: genericdatahandle.hh:607
std::enable_if< !CommunicationDescriptor::wrap_buffer &&AlwaysTrue< Entity >::value >::type gather(MessageBuffer &buff, const Entity &e) const
pack data from user to message buffer - version without support for sending leaf ordering sizes ...
Definition: genericdatahandle.hh:179
std::enable_if< CommunicationDescriptor::wrap_buffer &&AlwaysTrue< Entity >::value >::type gather(MessageBuffer &buff, const Entity &e) const
pack data from user to message buffer - version with support for sending leaf ordering sizes ...
Definition: genericdatahandle.hh:152
bool scatter(MessageBuffer &buff, std::size_t n, const Entity &e, LocalView &local_view) const
Definition: genericdatahandle.hh:785
Definition: genericdatahandle.hh:716
bool fixedSize(const GFS &gfs, int dim, int codim) const
Definition: genericdatahandle.hh:45
MaxDataHandle(const GFS &gfs_, V &v_)
Definition: genericdatahandle.hh:754
GFSDataHandle(const GFS &gfs, V &v, GatherScatter gather_scatter=GatherScatter(), CommunicationDescriptor communication_descriptor=CommunicationDescriptor())
Definition: genericdatahandle.hh:117
E OriginalDataType
Definition: genericdatahandle.hh:36
bool scatter(MessageBuffer &buff, size_type n, const Entity &e, LocalView &local_view) const
Definition: genericdatahandle.hh:401
E DataType
Definition: genericdatahandle.hh:64
GatherScatter functor for creating a disjoint DOF partitioning.
Definition: genericdatahandle.hh:859
std::size_t size_type
Definition: genericdatahandle.hh:274
RankIndex DataType
Definition: genericdatahandle.hh:1069
bool gather(MessageBuffer &buff, const Entity &e, LocalView &local_view) const
Definition: genericdatahandle.hh:865
bool contains(int dim, int codim) const
returns true if data for this codim should be communicated
Definition: genericdatahandle.hh:126
std::enable_if< CommunicationDescriptor::wrap_buffer &&AlwaysTrue< Entity >::value >::type scatter(MessageBuffer &buff, const Entity &e, size_type n)
unpack data from message buffer to user
Definition: genericdatahandle.hh:198
void gather(MessageBuffer &buff, const Entity &e) const
Definition: genericdatahandle.hh:1097
void scatter(MessageBuffer &buff, DataType &data) const
Definition: genericdatahandle.hh:645
Base class for all PDELab exceptions.
Definition: exceptions.hh:17
bool gather(MessageBuffer &buff, const Entity &e, LocalView &local_view) const
Definition: genericdatahandle.hh:773
void gather(MessageBuffer &buff, DataType &data) const
Definition: genericdatahandle.hh:701
bool scatter(MessageBuffer &buff, const Offsets &remote_offsets, const Offsets &local_offsets, const Entity &e, LocalView &local_view) const
Definition: genericdatahandle.hh:313
For backward compatibility – Do not use this!
Definition: adaptivity.hh:27
void gather(MessageBuffer &buff, DataType &data) const
Definition: genericdatahandle.hh:732
DataContainerIndexGatherScatter(GatherScatter gather_scatter=GatherScatter())
Definition: genericdatahandle.hh:592
Definition: genericdatahandle.hh:747
Definition: genericdatahandle.hh:269
bool gather(MessageBuffer &buff, const Entity &e, const LocalView &local_view) const
Definition: genericdatahandle.hh:277
bool fixedSize(const GFS &gfs, int dim, int codim) const
Definition: genericdatahandle.hh:77
GFS::Traits::SizeType size_type
Definition: genericdatahandle.hh:113
DataGatherScatter(GatherScatter gather_scatter=GatherScatter())
Definition: genericdatahandle.hh:372
void scatter(MessageBuffer &buff, DataType &data) const
Definition: genericdatahandle.hh:613
GatherScatter functor for marking shared DOFs.
Definition: genericdatahandle.hh:980
Definition: genericdatahandle.hh:634
std::size_t size_type
size type to use if communicating leaf ordering sizes
Definition: genericdatahandle.hh:30
bool fixedsize(int dim, int codim) const
returns true if size per entity of given dim and codim is a constant
Definition: genericdatahandle.hh:132
Data handle for marking ghost DOFs.
Definition: genericdatahandle.hh:812
GFS::Traits::SizeType size_type
Definition: genericdatahandle.hh:1070
Definition: genericdatahandle.hh:728
bool gather(MessageBuffer &buff, const Entity &e, const LocalView &local_view) const
Definition: genericdatahandle.hh:392
SharedDOFDataHandle(const GFS &gfs_, V &v_, bool init_vector=true)
Creates a new SharedDOFDataHandle.
Definition: genericdatahandle.hh:1042
std::array< size_type, leaf_count+1 > Offsets
Definition: entityindexcache.hh:38
bool gather(MessageBuffer &buff, const Entity &e, const LocalView &local_view) const
Definition: genericdatahandle.hh:501
Wrapper for message buffers of grid DataHandles that allows for sending different types of data...
Definition: polymorphicbufferwrapper.hh:26
void read(T &data)
Definition: polymorphicbufferwrapper.hh:40
DisjointPartitioningDataHandle(const GFS &gfs_, V &v_, bool init_vector=true)
Creates a new DisjointPartitioningDataHandle.
Definition: genericdatahandle.hh:964
Communication descriptor for sending one item of type E per DOF.
Definition: genericdatahandle.hh:24
void gather(MessageBuffer &buff, DataType &data) const
Definition: genericdatahandle.hh:638
CommunicationDescriptor::DataType DataType
Definition: genericdatahandle.hh:112
Data handle for collecting set of neighboring MPI ranks.
Definition: genericdatahandle.hh:1060
bool scatter(MessageBuffer &buff, const Offsets &remote_offsets, const Offsets &local_offsets, const Entity &e, LocalView &local_view) const
Definition: genericdatahandle.hh:428
size_type size(Entity &e) const
Definition: genericdatahandle.hh:1091
char DataType
Definition: genericdatahandle.hh:27
static const bool wrap_buffer
Definition: genericdatahandle.hh:33
GatherScatter data handle for creating a disjoint DOF partitioning.
Definition: genericdatahandle.hh:930
void scatter(MessageBuffer &buff, const Entity &e, size_type n)
Definition: genericdatahandle.hh:1103
GatherScatter functor for marking ghost DOFs.
Definition: genericdatahandle.hh:768
bool scatter(MessageBuffer &buff, size_type n, const Entity &e, LocalView &local_view) const
Definition: genericdatahandle.hh:510
DataEntityGatherScatter(GatherScatter gather_scatter=GatherScatter())
Definition: genericdatahandle.hh:481
bool scatter(MessageBuffer &buff, size_type n, const Entity &e, LocalView &local_view) const
Definition: genericdatahandle.hh:286
Definition: genericdatahandle.hh:493
bool scatter(MessageBuffer &buff, std::size_t n, const Entity &e, LocalView &local_view) const
Definition: genericdatahandle.hh:991
Implement a data handle with a grid function space.
Definition: genericdatahandle.hh:106
bool gather(MessageBuffer &buff, const Entity &e, LocalView &local_view) const
Definition: genericdatahandle.hh:984
bool scatter(MessageBuffer &buff, std::size_t n, const Entity &e, LocalView &local_view) const
Definition: genericdatahandle.hh:875
GFSNeighborDataHandle(const GFS &gfs, RankIndex rank, std::set< RankIndex > &neighbors)
Definition: genericdatahandle.hh:1072
Data handle for marking shared DOFs.
Definition: genericdatahandle.hh:1014
bool contains(const GFS &gfs, int dim, int codim) const
Definition: genericdatahandle.hh:39
DisjointPartitioningGatherScatter(RankIndex rank)
Create a DisjointPartitioningGatherScatter object.
Definition: genericdatahandle.hh:910
Definition: genericdatahandle.hh:666
EntityDataCommunicationDescriptor(std::size_t count=1)
Definition: genericdatahandle.hh:88
std::size_t size_type
Definition: genericdatahandle.hh:498
void write(const T &data)
Definition: polymorphicbufferwrapper.hh:32
CopyDataHandle(const GFS &gfs_, V &v_)
Definition: genericdatahandle.hh:692
Definition: genericdatahandle.hh:622
static const unsigned int value
Definition: gridfunctionspace/tags.hh:139
std::size_t size(const GFS &gfs, const Entity &e) const
Definition: genericdatahandle.hh:83
void gather(MessageBuffer &buff, DataType &data) const
Definition: genericdatahandle.hh:670
Definition: genericdatahandle.hh:384
std::enable_if< !CommunicationDescriptor::wrap_buffer &&AlwaysTrue< Entity >::value >::type scatter(MessageBuffer &buff, const Entity &e, size_type n)
unpack data from message buffer to user
Definition: genericdatahandle.hh:242
size_type size(const Entity &e) const
how many objects of type DataType have to be sent for a given entity
Definition: genericdatahandle.hh:142
AddClearDataHandle(const GFS &gfs_, V &v_)
Definition: genericdatahandle.hh:661
MinDataHandle(const GFS &gfs_, V &v_)
Definition: genericdatahandle.hh:723
Communication descriptor for sending count items of type E per entity with attached DOFs...
Definition: genericdatahandle.hh:61
bool contains(const GFS &gfs, int dim, int codim) const
Definition: genericdatahandle.hh:71
bool fixedsize(int dim, int codim) const
Definition: genericdatahandle.hh:1084
void scatter(MessageBuffer &buff, DataType &data) const
Definition: genericdatahandle.hh:707