DUNE-FEM (unstable)

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