dune-spgrid  2.7
capabilities.hh
Go to the documentation of this file.
1 #ifndef DUNE_SPGRID_CAPABILITIES_HH
2 #define DUNE_SPGRID_CAPABILITIES_HH
3 
4 #if HAVE_MPI
5 #include <mpi.h>
6 #endif
7 
8 #include <dune/geometry/type.hh>
9 
10 #include <dune/grid/common/capabilities.hh>
11 
13 
15 
21 namespace Dune
22 {
23 
24  // Capabilities
25  // ------------
26 
28  namespace Capabilities
29  {
30 
35  template< class ct, int dim, template< int > class Ref, class Comm >
36  struct hasSingleGeometryType< SPGrid< ct, dim, Ref, Comm > >
37  {
40  static const bool v = true;
42  static const unsigned int topologyId = Impl::CubeTopology< dim >::type::id;
43  };
44 
58  template< class ct, int dim, template< int > class Ref, class Comm >
59  struct isCartesian< SPGrid< ct, dim, Ref, Comm > >
60  {
62  static const bool v = true;
63  };
64 
70  template< class ct, int dim, template< int > class Ref, class Comm, int codim >
71  struct hasEntity< SPGrid< ct, dim, Ref, Comm >, codim >
72  {
75  static const bool v = ((codim >= 0) && (codim <= dim));
76  };
77 
78 #if HAVE_MPI
87  template< class ct, int dim, template< int > class Ref, int codim >
88  struct canCommunicate< SPGrid< ct, dim, Ref, MPI_Comm >, codim >
89  {
92  static const bool v = ((codim >= 0) && (codim <= dim));
93  };
94 #endif // #if HAVE_MPI
95 
100  template< class ct, int dim, template< int > class Ref, class Comm >
101  struct isLevelwiseConforming< SPGrid< ct, dim, Ref, Comm > >
102  {
104  static const bool v = true;
105  };
106 
111  template< class ct, int dim, template< int > class Ref, class Comm >
112  struct isLeafwiseConforming< SPGrid< ct, dim, Ref, Comm > >
113  {
115  static const bool v = true;
116  };
117 
122  template< class ct, int dim, template< int > class Ref, class Comm >
123  struct hasBackupRestoreFacilities< SPGrid< ct, dim, Ref, Comm > >
124  {
126  static const bool v = true;
127  };
128 
133  template< class ct, int dim, template< int > class Ref, class Comm >
134  struct threadSafe< SPGrid< ct, dim, Ref, Comm > >
135  {
137  static const bool v = false;
138  };
139 
144  template< class ct, int dim, template< int > class Ref, class Comm >
145  struct viewThreadSafe< SPGrid< ct, dim, Ref, Comm > >
146  {
148  static const bool v = false;
149  };
150 
151 
152 
153  // non-standard capabilities (see dune-fem)
154  // ----------------------------------------
155 
156  template< class Grid >
158 
159  template< class ct, int dim, template< int > class Ref, class Comm >
160  struct hasHierarchicIndexSet< SPGrid< ct, dim, Ref, Comm > >
161  {
162  static const bool v = true;
163  };
164 
165 
166  template< class Grid >
168 
175  template< class ct, int dim, template< int > class Ref, class Comm >
176  struct supportsCallbackAdaptation< SPGrid< ct, dim, Ref, Comm > >
177  {
179  static const bool v = true;
180  };
181 
182  } // namespace Capabilities
183 
184 
185 
186  // Extensions
187  // ----------
188 
189  namespace Extensions
190  {
191 
197  template< class ct, int dim, template< int > class Ref, class Comm, int codim >
198  struct SuperEntityIterator< SPGrid< ct, dim, Ref, Comm >, codim >
199  {
202  static const bool v = ((codim >= 0) && (codim <= dim));
203  };
204 
205  } // namespace Extensions
206 
207 } // namespace Dune
208 
209 #endif // #ifndef DUNE_SPGRID_CAPABILITIES_HH
Definition: iostream.hh:7
Does a grid support superentity iterators of a codimension?
Definition: extensions/superentityiterator.hh:82
static const bool v
by default, a grid does not support superentity iterators
Definition: extensions/superentityiterator.hh:84
Definition: capabilities.hh:157
Definition: capabilities.hh:167
structured, parallel DUNE grid
Definition: grid.hh:136
interface classes for superentity iterators