DUNE-FEM (unstable)

datahandle.hh
1#ifndef DUNE_IDGRID_DATAHANDLE_HH
2#define DUNE_IDGRID_DATAHANDLE_HH
3
4//- C++ include
5#include <type_traits>
6
7//- dune-grid includes
9
10#if HAVE_DUNE_ALUGRID
11#include <dune/alugrid/common/ldbhandleif.hh>
12#endif
13
14namespace Dune
15{
16
17 // IdDataHandle
18 // ------------
19
20 template< class WrappedHandle, class GridFamily >
21 class IdDataHandle
22 : public CommDataHandleIF< IdDataHandle< WrappedHandle, GridFamily >, typename WrappedHandle::DataType >
23 {
24 protected:
25 typedef IdDataHandle< WrappedHandle, GridFamily > ThisType;
26
27 // type of traits
28 typedef typename std::remove_const< GridFamily >::type::Traits Traits;
29
30 typedef typename Traits :: ExtraData ExtraData ;
31
32 template< int codim >
33 struct Codim
34 {
35 // type of entity
36 typedef typename Traits::template Codim< codim >::Entity EntityType;
37 };
38
39 public:
40 // type of data to be communicated
41 typedef typename WrappedHandle::DataType DataType;
42
43 typedef CommDataHandleIF< ThisType, DataType > DataHandleIFType;
44
45 IdDataHandle ( const ThisType & ) = delete;
46
47 IdDataHandle ( ExtraData data, WrappedHandle &wrappedHandle )
48 : wrappedHandle_( wrappedHandle ),
49 data_( data )
50 {}
51
52 bool contains ( int dim, int codim ) const
53 {
54 return wrappedHandle_.contains( dim, codim );
55 }
56
57 bool fixedSize ( int dim, int codim ) const
58 {
59 return wrappedHandle_.fixedSize( dim, codim );
60 }
61
62 template< class HostEntity >
63 size_t size ( const HostEntity &hostEntity ) const
64 {
65 typedef typename Codim< HostEntity::codimension >::EntityType EntityType;
66 const EntityType entity( typename EntityType::Implementation( data(), hostEntity ) );
67 return wrappedHandle_.size( entity );
68 }
69
70 template< class MessageBuffer, class HostEntity >
71 void gather ( MessageBuffer &buffer, const HostEntity &hostEntity ) const
72 {
73 typedef typename Codim< HostEntity::codimension >::EntityType EntityType;
74 const EntityType entity( typename EntityType::Implementation( data(), hostEntity ) );
75 wrappedHandle_.gather( buffer, entity );
76 }
77
78 template< class MessageBuffer, class HostEntity >
79 void scatter ( MessageBuffer &buffer, const HostEntity &hostEntity, size_t size )
80 {
81 typedef typename Codim< HostEntity::codimension >::EntityType EntityType;
82 const EntityType entity( typename EntityType::Implementation( data(), hostEntity ) );
83 wrappedHandle_.scatter( buffer, entity, size );
84 }
85
86 ExtraData data() const { return data_; }
87
88 protected:
89 WrappedHandle &wrappedHandle_;
90 ExtraData data_;
91 };
92
93} // namespace Dune
94
95#endif // #ifndef DUNE_IDGRID_DATAHANDLE_HH
void scatter(MessageBufferImp &buff, const EntityType &e, size_t n)
unpack data from message buffer to user.
Definition: datahandleif.hh:143
void gather(MessageBufferImp &buff, const EntityType &e) const
pack data from user to message buffer
Definition: datahandleif.hh:129
Describes the parallel communication interface class for MessageBuffers and DataHandles.
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
constexpr std::bool_constant<((II==value)||...)> contains(std::integer_sequence< T, II... >, std::integral_constant< T, value >)
Checks whether or not a given sequence contains a value.
Definition: integersequence.hh:137
std::size_t fixedSize
The number of data items per index if it is fixed, 0 otherwise.
Definition: variablesizecommunicator.hh:264
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 13, 23:29, 2024)