DUNE-FEM (unstable)

datahandle.hh
1#ifndef DUNE_FEM_GRIDPART_GEOGRIDPART_DATAHANDLE_HH
2#define DUNE_FEM_GRIDPART_GEOGRIDPART_DATAHANDLE_HH
3
4#include <type_traits>
5
7
8#include <dune/fem/gridpart/geogridpart/entity.hh>
9
10namespace Dune
11{
12
13 namespace Fem
14 {
15
16 // GeoDataHandle
17 // -------------
18
19 template< class GridFamily, class WrappedHandle >
20 class GeoDataHandle
21 : public CommDataHandleIF< GeoDataHandle< GridFamily, WrappedHandle >, typename WrappedHandle::DataType >
22 {
23 typedef typename std::remove_const< GridFamily >::type::Traits Traits;
24
25 template< class HostEntity >
26 struct EntityProxy;
27
28 public:
29 typedef typename Traits::CoordFunctionType CoordFunctionType;
30
31 GeoDataHandle ( const CoordFunctionType &coordFunction, WrappedHandle &handle )
32 : coordFunction_( coordFunction ),
33 wrappedHandle_( handle )
34 {}
35
36 bool contains ( int dim, int codim ) const
37 {
38 return wrappedHandle_.contains( dim, codim );
39 }
40
41 bool fixedSize ( int dim, int codim ) const
42 {
43 return wrappedHandle_.fixedSize( dim, codim );
44 }
45
46 template< class HostEntity >
47 size_t size ( const HostEntity &hostEntity ) const
48 {
49 EntityProxy< HostEntity > proxy( coordFunction_, hostEntity );
50 return wrappedHandle_.size( *proxy );
51 }
52
53 template< class MessageBuffer, class HostEntity >
54 void gather ( MessageBuffer &buffer, const HostEntity &hostEntity ) const
55 {
56 EntityProxy< HostEntity > proxy( coordFunction_, hostEntity );
57 wrappedHandle_.gather( buffer, *proxy );
58 }
59
60 template< class MessageBuffer, class HostEntity >
61 void scatter ( MessageBuffer &buffer, const HostEntity &hostEntity, size_t size )
62 {
63 EntityProxy< HostEntity > proxy( coordFunction_, hostEntity );
64 wrappedHandle_.scatter( buffer, *proxy, size );
65 }
66
67 private:
68 const CoordFunctionType &coordFunction_;
69 WrappedHandle &wrappedHandle_;
70 };
71
72
73
74 // GeoDataHandle::EntityProxy
75 // --------------------------
76
77 template< class GridFamily, class WrappedHandle >
78 template< class HostEntity >
79 struct GeoDataHandle< GridFamily, WrappedHandle >::EntityProxy
80 {
81 static const int dimension = HostEntity::dimension;
82 static const int codimension = HostEntity::codimension;
83
84 typedef typename GridFamily::template Codim< codimension > :: Entity Entity;
85
86 protected:
87 typedef typename Entity::Implementation EntityImpl;
88
89 public:
90 EntityProxy ( const CoordFunctionType &coordFunction, const HostEntity &hostEntity )
91 : entity_( EntityImpl( coordFunction, hostEntity ) )
92 {}
93
94 const Entity &operator* () const
95 {
96 return entity_;
97 }
98
99 private:
100 Entity entity_;
101 };
102
103 } // namespace Fem
104
105} // namespace Dune
106
107#endif // #ifndef DUNE_FEM_GRIDPART_GEOGRIDPART_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
EntityImp< cd, dim, GridImp > Implementation
type of underlying implementation
Definition: entity.hh:73
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 21, 23:30, 2024)