dune-fem 2.8.0
Loading...
Searching...
No Matches
gridpart/filteredgridpart/datahandle.hh
Go to the documentation of this file.
1#ifndef DUNE_FEM_GRIDPART_FILTEREDGRIDPART_DATAHANDLE_HH
2#define DUNE_FEM_GRIDPART_FILTEREDGRIDPART_DATAHANDLE_HH
3
4// C++
5#include <type_traits>
6
7//- dune-grid includes
8#include <dune/grid/common/datahandleif.hh>
9
10
11namespace Dune
12{
13
14 namespace Fem
15 {
16
17 // Forward declaration
18 // -------------------
19
20 template< class HostGridPartImp, class FilterImp, bool useFilteredIndexSet >
21 class FilteredGridPart;
22
23
24
25 // FilteredGridPartDataHandle
26 // --------------------------
27
28 template< class WrappedHandle, class GridPart >
30 : public CommDataHandleIF< FilteredGridPartDataHandle< WrappedHandle, GridPart >, typename WrappedHandle::DataType >
31 {
32 typedef CommDataHandleIF< FilteredGridPartDataHandle< WrappedHandle, GridPart >,
33 typename WrappedHandle::DataType > BaseType;
34 typedef GridPart GridPartType;
35 typedef typename std::remove_const< GridPartType >::type::Traits Traits;
36
37 public:
38 FilteredGridPartDataHandle ( WrappedHandle &dataHandle, const GridPart &gridPart )
39 : gridPart_( gridPart ),
40 wrappedHandle_( dataHandle )
41 { }
42
43 bool contains ( int dim, int codim ) const
44 {
45 return wrappedHandle_.contains( dim, codim );
46 }
47
48 bool fixedSize ( int dim, int codim ) const
49 {
50 return false;
51 }
52
53 template< class HostEntity >
54 size_t size ( const HostEntity &hostEntity ) const
55 {
56 if( gridPart().contains( hostEntity ) )
57 return wrappedHandle_.size( hostEntity );
58 else
59 return 0;
60 }
61
62 template< class MessageBuffer, class HostEntity >
63 void gather ( MessageBuffer &buffer, const HostEntity &hostEntity ) const
64 {
65 if( gridPart().contains( hostEntity ) )
66 wrappedHandle_.gather( buffer, hostEntity );
67 }
68
69 template< class MessageBuffer, class HostEntity >
70 void scatter ( MessageBuffer &buffer, const HostEntity &hostEntity, size_t size )
71 {
72 if( gridPart().contains( hostEntity ) )
73 wrappedHandle_.scatter( buffer, hostEntity, size );
74 else
75 {
76 typename BaseType::DataType tmp;
77 for (size_t i=0;i<size;++i)
78 buffer.read(tmp);
79 }
80 }
81
82 protected:
83 const GridPart &gridPart () const
84 {
85 return gridPart_;
86 }
87
88 private:
89 const GridPart &gridPart_;
90 WrappedHandle &wrappedHandle_;
91 };
92
93 } // namespace Fem
94
95} // namespace Dune
96
97#endif // #ifndef DUNE_FEM_GRIDPART_FILTEREDGRIDPART_DATAHANDLE_HH
Definition: bindguard.hh:11
Definition: gridpart/filteredgridpart/datahandle.hh:31
size_t size(const HostEntity &hostEntity) const
Definition: gridpart/filteredgridpart/datahandle.hh:54
void gather(MessageBuffer &buffer, const HostEntity &hostEntity) const
Definition: gridpart/filteredgridpart/datahandle.hh:63
bool contains(int dim, int codim) const
Definition: gridpart/filteredgridpart/datahandle.hh:43
bool fixedSize(int dim, int codim) const
Definition: gridpart/filteredgridpart/datahandle.hh:48
void scatter(MessageBuffer &buffer, const HostEntity &hostEntity, size_t size)
Definition: gridpart/filteredgridpart/datahandle.hh:70
FilteredGridPartDataHandle(WrappedHandle &dataHandle, const GridPart &gridPart)
Definition: gridpart/filteredgridpart/datahandle.hh:38
const GridPart & gridPart() const
Definition: gridpart/filteredgridpart/datahandle.hh:83