1 #ifndef DUNE_ALU3DGRIDGRID_HH
2 #define DUNE_ALU3DGRIDGRID_HH
3 
4 //- System includes
5 #include <memory>
6 #include <vector>
7 
8 //- Dune includes
9 #include <dune/grid/common/capabilities.hh>
10 #include <dune/alugrid/common/interfaces.hh>
11 #include <dune/common/bigunsignedint.hh>
12 #include <dune/common/version.hh>
13 
14 #include <dune/geometry/referenceelements.hh>
15 
16 #include <dune/grid/common/grid.hh>
17 #include <dune/alugrid/common/defaultindexsets.hh>
18 #include <dune/grid/common/sizecache.hh>
19 #include <dune/alugrid/common/intersectioniteratorwrapper.hh>
20 #include <dune/grid/common/datahandleif.hh>
21 
22 #include <dune/alugrid/common/typetraits.hh>
23 
24 // bnd projection stuff
25 #include <dune/grid/common/boundaryprojection.hh>
26 #include <dune/alugrid/common/bndprojection.hh>
27 #include <dune/alugrid/common/backuprestore.hh>
28 #include <dune/alugrid/common/macrogridview.hh>
29 #include <dune/alugrid/common/twists.hh>
30 
31 //- Local includes
32 #include "alu3dinclude.hh"
33 #include "topology.hh"
34 #include "indexsets.hh"
35 #include "datahandle.hh"
36 
37 #include <dune/alugrid/3d/communication.hh>
38 #include <dune/alugrid/3d/gridview.hh>
39 
40 #include <dune/common/parallel/mpihelper.hh>
41 
42 #if ALU3DGRID_PARALLEL
43 #include <dune/common/parallel/mpicommunication.hh>
44 #else
45 #include <dune/common/parallel/communication.hh>
46 #endif
47 
48 namespace Dune
49 {
50   // Forward declarations
51   template<int cd, int dim, class GridImp>
52   class ALU3dGridEntity;
53   template<int cd, PartitionIteratorType pitype, class GridImp >
54   class ALU3dGridLevelIterator;
55   template<int cd, class GridImp >
56   class ALU3dGridEntityPointerBase;
57   template<int cd, class GridImp >
58   class ALU3dGridEntitySeed;
59   template<int cd, class GridImp >
60   class ALU3dGridEntityPointer;
61   template<int mydim, int coorddim, class GridImp>
62   class ALU3dGridGeometry;
63   template<class GridImp>
64   class ALU3dGridHierarchicIterator;
65   template<class GridImp>
66   class ALU3dGridIntersectionIterator;
67   template<class GridImp>
68   class ALU3dGridLevelIntersectionIterator;
69   template<int codim, PartitionIteratorType pitype, class GridImp>
70   class ALU3dGridLeafIterator;
71   template <int mydim, int coorddim, class GridImp>
72   class ALU3dGridMakeableEntity;
73   template <class GridImp>
74   class ALU3dGridFaceGeometryInfo;
75   template< int, int, ALU3dGridElementType, class >
76   class ALU3dGridGlobalIdSet;
77   template< int, int, ALU3dGridElementType, class >
78   class ALU3dGridLocalIdSet;
79   template< int, int, ALU3dGridElementType, class >
80   class ALU3dGridHierarchicIndexSet;
81   template< class >
82   class ALU3dGridFactory;
83   template <class GridImp, class GeometryImp, int nChild>
84   class ALULocalGeometryStorage;
85 
86 
87 
88   // Internal Forward Declarations
89   // -----------------------------
90 
91 #if ALU3DGRID_PARALLEL
92   template<int dim, int dimworld,  ALU3dGridElementType elType, class Comm = ALUGridMPIComm >
93   class ALU3dGrid;
94 #else // #if ALU3DGRID_PARALLEL
95   template< int dim, int dimworld, ALU3dGridElementType elType, class Comm = ALUGridNoComm >
96   class ALU3dGrid;
97 #endif // #else // #if ALU3DGRID_PARALLEL
98 
99 
100   // Internal Forward Declarations
101   // -----------------------------
102 
103   template < int dim, int dimw, class Comm >
104   struct ALUGridBaseGrid< dim, dimw, cube, Comm >
105   {
106     typedef ALU3dGrid< dim, dimw, hexa, Comm >  BaseGrid ;
107   };
108 
109   template < int dim, int dimw, class Comm >
110   struct ALUGridBaseGrid< dim, dimw, simplex, Comm >
111   {
112     typedef ALU3dGrid< dim, dimw, tetra, Comm >  BaseGrid ;
113   };
114 
115 
116   template< int dim, int dimworld, ALU3dGridElementType elType, class Comm >
117   struct ALU3dGridCommunications;
118 
119   template< int dim, int dimworld, ALU3dGridElementType elType >
120   struct ALU3dGridCommunications< dim, dimworld, elType, ALUGridNoComm >
121   {
122     typedef ALU3dGridLocalIdSet< dim, dimworld, elType, ALUGridNoComm > GlobalIdSet;
123     typedef int GlobalId;
124 
125     typedef ALU3DSPACE GitterDuneImpl GitterImplType;
126 
127     typedef Dune::CollectiveCommunication< No_Comm > CollectiveCommunication;
128 
ALU3dGridCommunicationsDune::ALU3dGridCommunications129     explicit ALU3dGridCommunications ( ALUGridNoComm comm ) {}
130 
nlinksDune::ALU3dGridCommunications131     int nlinks () const { return 0; }
132 
createALUGridDune::ALU3dGridCommunications133     GitterImplType *createALUGrid ( const std::string &macroName, const ALU3DSPACE ProjectVertexPtrPair& projections,
134                                     const bool conformingRefinement )
135     {
136       GitterImplType* grid = ( macroName.empty() ) ?
137         new GitterImplType( dim, conformingRefinement ) : new GitterImplType ( dim, conformingRefinement, macroName.c_str(), projections );
138       return grid ;
139     }
140 
createALUGridDune::ALU3dGridCommunications141     GitterImplType *createALUGrid ( std::istream& stream, const ALU3DSPACE ProjectVertexPtrPair& projection,
142                                     const bool conformingRefinement )
143     {
144       return new GitterImplType ( dim, conformingRefinement, stream, projection );
145     }
146 
147     // ALUGridNoComm casts into No_Comm and MPI_Comm and here the default is MPI_COMM_SELF
defaultCommDune::ALU3dGridCommunications148     static ALUGridNoComm defaultComm () { return ALUGridNoComm(); }
149 
getRankDune::ALU3dGridCommunications150     static int getRank ( ALUGridNoComm comm ) { return 0; }
151 
getBuilderDune::ALU3dGridCommunications152     static typename ALU3DSPACE Gitter::Geometric::BuilderIF &getBuilder ( GitterImplType &grid )
153     {
154       ALU3DSPACE Gitter::Geometric::BuilderIF* builder =
155         dynamic_cast< ALU3DSPACE Gitter::Geometric::BuilderIF* >( &grid.container() );
156       if( ! builder )
157         DUNE_THROW(InvalidStateException,"dynamic_cast of ALUGrid builder failed");
158       return *builder;
159     }
160 
completeGridDune::ALU3dGridCommunications161     static void completeGrid ( GitterImplType &grid ) {}
162 
printDune::ALU3dGridCommunications163     void print( std::ostream& out ) const
164     {}
165 
166     CollectiveCommunication ccobj_;
167   };
168 
169 #if ALU3DGRID_PARALLEL
170   template< int dim, int dimworld, ALU3dGridElementType elType >
171   struct ALU3dGridCommunications< dim, dimworld, elType, ALUGridMPIComm >
172   {
173     typedef ALU3dGridGlobalIdSet< dim, dimworld, elType, ALUGridMPIComm > GlobalIdSet;
174     typedef ALUGridId< ALUMacroKey > GlobalId;
175 
176     typedef ALU3DSPACE GitterDunePll GitterImplType;
177 
178     typedef Dune::CollectiveCommunication< MPI_Comm > CollectiveCommunication;
179 
ALU3dGridCommunicationsDune::ALU3dGridCommunications180     explicit ALU3dGridCommunications ( MPI_Comm comm )
181     : ccobj_( comm ), mpAccess_( comm )
182     {}
183 
nlinksDune::ALU3dGridCommunications184     int nlinks () const { return mpAccess_.sendLinks(); }
185 
createALUGridDune::ALU3dGridCommunications186     GitterImplType *createALUGrid ( const std::string &macroName, const ALU3DSPACE ProjectVertexPtrPair& projections,
187                                     const bool conformingRefinement )
188     {
189       return new GitterImplType( dim, conformingRefinement, macroName.c_str(), mpAccess_, projections );
190     }
191 
createALUGridDune::ALU3dGridCommunications192     GitterImplType *createALUGrid ( std::istream& stream, const ALU3DSPACE ProjectVertexPtrPair& projections,
193                                     const bool conformingRefinement )
194     {
195       return new GitterImplType ( dim, conformingRefinement, stream, mpAccess_, projections );
196     }
197 
198     // ALUGridMPIComm casts into MPI_Comm and the default is MPI_COMM_WORLD
defaultCommDune::ALU3dGridCommunications199     static ALUGridMPIComm defaultComm () { return ALUGridMPIComm(); }
200 
getRankDune::ALU3dGridCommunications201     static int getRank ( MPI_Comm comm )
202     {
203       int rank = 0;
204       MPI_Comm_rank( comm, &rank );
205       return rank;
206     }
207 
printDune::ALU3dGridCommunications208     void print( std::ostream& out ) const
209     {
210       mpAccess_.printLinkage( out );
211     }
212 
getBuilderDune::ALU3dGridCommunications213     static typename ALU3DSPACE Gitter::Geometric::BuilderIF &getBuilder ( GitterImplType &grid )
214     {
215       ALU3DSPACE Gitter::Geometric::BuilderIF* builder =
216         dynamic_cast< ALU3DSPACE Gitter::Geometric::BuilderIF* >( &grid.containerPll() );
217       if( ! builder )
218         DUNE_THROW(InvalidStateException,"dynamic_cast of ALUGrid builder failed");
219       return *builder;
220     }
221 
completeGridDune::ALU3dGridCommunications222     static void completeGrid ( GitterImplType &grid )
223     {
224       // setup communication patterns
225       grid.notifyMacroGridChanges();
226       // rebuild ghost cells
227       grid.rebuildGhostCells();
228     }
229 
230     CollectiveCommunication ccobj_;
231     ALU3DSPACE MpAccessMPI mpAccess_;
232   };
233 #endif // #if ALU3DGRID_PARALLEL
234 
235 
236 
237   // ALU3dGridTwist
238   // --------------
239 
240   template< int dim, ALU3dGridElementType elType, int codim >
241   struct ALU3dGridTwists;
242 
243   template<int dim>
244   struct ALU3dGridTwists< dim, tetra, 0 >
245   {
246     static const unsigned int topoId = GeometryTypes::simplex(dim).id();
247     typedef TrivialTwists< topoId, dim > Type;
248   };
249 
250   template<int dim>
251   struct ALU3dGridTwists< dim, hexa, 0 >
252   {
253     static const unsigned int topoId = GeometryTypes::cube(dim).id();
254     typedef TrivialTwists< topoId, dim > Type;
255   };
256 
257   template< int dim, ALU3dGridElementType elType >
258   struct ALU3dGridTwists< dim, elType, 1 >
259   {
260     typedef ALUTwists< dim == 2 ? 2 : ElementTopologyMapping< elType >::numVerticesPerFace, dim-1 > Type;
261   };
262 
263   template< ALU3dGridElementType elType >
264   struct ALU3dGridTwists< 3, elType, 2 >
265   {
266     typedef ALUTwists< 2, 1 > Type;
267   };
268 
269   template< ALU3dGridElementType elType >
270   struct ALU3dGridTwists< 2, elType, 2 >
271   {
272     typedef TrivialTwists< 0u, 0 > Type;
273   };
274 
275   template< int dim, ALU3dGridElementType elType >
276   struct ALU3dGridTwists< dim, elType, 3 >
277   {
278     typedef TrivialTwists< 0u, 0 > Type;
279   };
280 
281 
282 
283   // ALU3dGridFamily
284   // ---------------
285 
286   template< int dimG, int dimW, ALU3dGridElementType elType, class Comm >
287   struct ALU3dGridFamily
288   {
289     static const int dim = dimG;
290     static const int dimworld = dimW;
291 
292     typedef ALU3dGrid< dim, dimworld, elType, Comm > GridImp;
293     typedef ALU3dGridFamily< dim, dimworld, elType, Comm > GridFamily;
294 
295     //! Type of the local id set
296     typedef ALU3dGridLocalIdSet< dim, dimworld, elType, Comm > LocalIdSetImp;
297 
298     //! Type of the global id set
299     typedef typename ALU3dGridCommunications< dim, dimworld, elType, Comm >::GlobalIdSet GlobalIdSetImp;
300 
301     //! type of ALU3dGrids global id
302     typedef typename ALU3dGridCommunications< dim, dimworld, elType, Comm >::GlobalId GlobalIdType;
303 
304     //! type of ALU3dGrids local id
305     typedef int LocalIdType;
306 
307     struct Traits
308     {
309       //! type of ALU3dGrids local id
310       typedef typename GridFamily::LocalIdType LocalIdType;
311 
312       //! type of ALU3dGrids global id
313       typedef typename GridFamily::GlobalIdType GlobalIdType;
314 
315       typedef typename GridFamily::GridImp Grid;
316 
317       typedef Dune::Intersection< const Grid, LeafIntersectionWrapper< const Grid > > LeafIntersection;
318       typedef Dune::Intersection< const Grid, LevelIntersectionWrapper< const Grid > > LevelIntersection;
319 
320       typedef Dune::IntersectionIterator< const Grid, LeafIntersectionIteratorWrapper< const Grid >, LeafIntersectionWrapper< const Grid > > IntersectionIterator;
321 
322       typedef Dune::IntersectionIterator< const Grid, LeafIntersectionIteratorWrapper< const Grid >, LeafIntersectionWrapper< const Grid > > LeafIntersectionIterator;
323       typedef Dune::IntersectionIterator< const Grid, LevelIntersectionIteratorWrapper< const Grid >, LevelIntersectionWrapper< const Grid > > LevelIntersectionIterator;
324 
325       typedef Dune::EntityIterator< 0, const Grid, ALU3dGridHierarchicIterator< const Grid > > HierarchicIterator;
326 
327       typedef DuneBoundaryProjection< dimworld > DuneBoundaryProjectionType;
328       typedef std::vector< const DuneBoundaryProjectionType * > DuneBoundaryProjectionVector;
329 
330       template< int cd >
331       struct Codim
332       {
333         typedef typename ALU3dGridTwists< dim, elType, cd >::Type Twists;
334         typedef typename Twists::Twist Twist;
335 
336         // IMPORTANT: Codim<codim>::Geometry == Geometry<dim-codim,dimw>
337         typedef ALU3dGridGeometry< dim-cd, dimworld, const Grid > GeometryImpl;
338         typedef ALU3dGridGeometry< dim-cd, dim, const Grid > LocalGeometryImpl;
339         typedef Dune::Geometry< dim-cd, dimworld, const Grid, ALU3dGridGeometry > Geometry;
340         typedef Dune::Geometry< dim-cd, dim, const Grid, ALU3dGridGeometry > LocalGeometry;
341 
342         typedef ALU3dGridEntity< cd, dim, const Grid > EntityImp;
343         typedef Dune::Entity< cd, dim, const Grid, ALU3dGridEntity > Entity;
344 
345         // minimal information to generate entities
346         typedef ALU3dGridEntitySeed< cd , const Grid> EntitySeed ;
347 
348         template< PartitionIteratorType pitype >
349         struct Partition
350         {
351           typedef Dune::EntityIterator< cd, const Grid, ALU3dGridLevelIterator< cd, pitype, const Grid > > LevelIterator;
352           typedef Dune::EntityIterator< cd, const Grid, ALU3dGridLeafIterator< cd, pitype, const Grid > > LeafIterator;
353         }; // struct Partition
354 
355         typedef typename Partition< All_Partition >::LevelIterator LevelIterator;
356         typedef typename Partition< All_Partition >::LeafIterator LeafIterator;
357       }; // struct Codim
358 
359       template< PartitionIteratorType pitype >
360       struct Partition
361       {
362         typedef Dune::GridView< ALU3dLevelGridViewTraits< const Grid, pitype > > LevelGridView;
363         typedef Dune::GridView< ALU3dLeafGridViewTraits< const Grid, pitype > >  LeafGridView;
364         typedef Dune::MacroGridView<const Grid, pitype> MacroGridView;
365       }; // struct Partition
366 
367       typedef typename Partition< All_Partition > :: MacroGridView   MacroGridView;
368       typedef typename Partition< All_Partition > :: LeafGridView    LeafGridView;
369       typedef typename Partition< All_Partition > :: LevelGridView   LevelGridView;
370 
371       //! Type of the level index set
372       typedef DefaultIndexSet< Grid, typename Codim< 0 > :: LevelIterator > LevelIndexSetImp;
373 
374       //! Type of the leaf index set
375       typedef DefaultIndexSet< Grid, typename Codim< 0 > :: LeafIterator > LeafIndexSetImp;
376 
377       typedef IndexSet< Grid, LevelIndexSetImp > LevelIndexSet;
378       typedef IndexSet< Grid, LeafIndexSetImp > LeafIndexSet;
379       typedef IdSet< Grid, LocalIdSetImp, LocalIdType > LocalIdSet;
380       typedef IdSet< Grid, GlobalIdSetImp, GlobalIdType > GlobalIdSet;
381 
382       //! Type of the communication class
383       typedef typename ALU3dGridCommunications< dim, dimworld, elType, Comm >::CollectiveCommunication CollectiveCommunication;
384     }; // struct Traits
385 
386     //! Type of the level index set implementation
387     typedef typename Traits :: LevelIndexSetImp  LevelIndexSetImp;
388 
389     //! Type of the leaf index set implementation
390     typedef typename Traits :: LeafIndexSetImp   LeafIndexSetImp;
391 
392   }; // struct ALU3dGridFamily
393 
394 
395 
396   //**********************************************************************
397   //
398   // --ALU3dGrid
399   // --Grid
400   //
401   //**********************************************************************
402 
403   /**
404      \brief [<em> provides \ref Dune::Grid </em>]
405      \brief 3D grid with support for hexahedrons and tetrahedrons.
406      The ALU3dGrid implements the Dune GridInterface for 3d tetrahedral and
407      hexahedral meshes. This grid can be locally adapted and used in parallel
408      computations using dynamic load balancing.
409 
410      @note
411      Adaptive parallel grid supporting dynamic load balancing, written
412      mainly by Bernard Schupp. This grid supports hexahedrons and tetrahedrons.
413 
414      (see ALUGrid homepage: http://www.mathematik.uni-freiburg.de/IAM/Research/alugrid/)
415 
416      Two tools are available for partitioning :
417      \li Metis ( version 4.0 and higher, see http://glaros.dtc.umn.edu/gkhome/views/metis/metis/ )
418      \li ParMETIS ( http://glaros.dtc.umn.edu/gkhome/metis/parmetis/overview )
419 
420      For installation instructions see http://www.dune-project.org/external_libraries/install_alugrid.html .
421      @author Robert Kloefkorn
422   */
423   template< int dim, int dimworld,  ALU3dGridElementType elType, class Comm >
424   class ALU3dGrid
425   : public GridDefaultImplementation< dim, dimworld, alu3d_ctype,
426                                       ALU3dGridFamily< dim, dimworld, elType, Comm > >,
427     public HasObjectStream,
428     public HasHierarchicIndexSet
429   {
430     typedef ALU3dGrid< dim, dimworld, elType, Comm > ThisType;
431     typedef GridDefaultImplementation< dim, dimworld, alu3d_ctype, ALU3dGridFamily< dim, dimworld, elType, Comm > > BaseType;
432 
433     // for compatibility: MyType := ThisType
434     typedef ThisType MyType;
435 
436     // friend declarations
437     friend class ALU3dGridEntity< 0, dim, const ThisType>;
438     friend class ALU3dGridEntity< 1, dim, const ThisType>;
439     friend class ALU3dGridEntity< 2, dim, const ThisType>;
440     friend class ALU3dGridEntity< dim, dim, const ThisType>;
441 
442     friend class ALU3dGridIntersectionIterator< ThisType >;
443 
444     friend class ALU3dGridEntityPointerBase< 0, const ThisType >;
445     friend class ALU3dGridEntityPointerBase< 1, const ThisType >;
446     friend class ALU3dGridEntityPointerBase< 2, const ThisType >;
447     friend class ALU3dGridEntityPointerBase< dim, const ThisType >;
448 
449     friend class ALU3dGridEntityPointer< 0, const ThisType >;
450     friend class ALU3dGridEntityPointer< 1, const ThisType >;
451     friend class ALU3dGridEntityPointer< 2, const ThisType >;
452     friend class ALU3dGridEntityPointer< dim, const ThisType >;
453 
454     friend class ALU3dGridIntersectionIterator< const ThisType >;
455     friend class ALU3dGridHierarchicIterator< const ThisType >;
456 
457     friend class ALU3dGridHierarchicIndexSet< dim, dimworld, elType, Comm >;
458     friend class ALU3dGridGlobalIdSet< dim, dimworld, elType, Comm >;
459     friend class ALU3dGridLocalIdSet< dim, dimworld, elType, Comm >;
460 
461     // new intersection iterator is a wrapper which get itersectioniteratoimp as pointers
462   public:
463     typedef ALU3dGridIntersectionIterator<const ThisType>
464       IntersectionIteratorImp;
465     typedef ALU3dGridIntersectionIterator<const ThisType>
466       LeafIntersectionIteratorImp;
467     typedef ALU3dGridLevelIntersectionIterator<const ThisType>
468       LevelIntersectionIteratorImp;
469 
470     friend class IntersectionIteratorWrapper < const ThisType, LeafIntersectionIteratorImp > ;
471     friend class IntersectionIteratorWrapper < const ThisType, LevelIntersectionIteratorImp > ;
472     friend class LeafIntersectionIteratorWrapper < const ThisType > ;
473     friend class LevelIntersectionIteratorWrapper< const ThisType > ;
474 
475     //**********************************************************
476     // The Interface Methods
477     //**********************************************************
478   public:
479     enum { refineStepsForHalf = 1 };
480 
481     static const ALU3dGridElementType elementType = elType;
482 
483     typedef typename ALU3DSPACE GatherScatterType::ObjectStreamType ObjectStreamType;
484     typedef ObjectStreamType  InStreamType ;
485     typedef ObjectStreamType  OutStreamType ;
486 
487     typedef ALU3dGridFamily< dim, dimworld, elType, Comm > GridFamily;
488     typedef typename GridFamily::Traits Traits;
489 
490     static const int dimension =      BaseType::dimension;
491     static const int dimensionworld = BaseType::dimensionworld;
492 
493     template< int codim >
494     struct Codim
495       : public BaseType::template Codim< codim >
496     {
497       typedef typename Traits::template Codim< codim >::EntityImp   EntityImp;
498       typedef typename Traits::template Codim< codim >::Twists      Twists;
499       typedef typename Twists::Twist                                Twist;
500     };
501 
502   protected:
503     typedef MakeableInterfaceObject< typename Traits::template Codim< 0 >::Geometry > GeometryObject;
504     friend class ALULocalGeometryStorage< const ThisType, GeometryObject, 8 >;
505 
506   public:
507     /** \brief Types for GridView */
508     template <PartitionIteratorType pitype>
509     struct Partition
510     {
511       typedef typename GridFamily::Traits::template Partition<pitype>::LevelGridView
512          LevelGridView;
513       typedef typename GridFamily::Traits::template Partition<pitype>::LeafGridView
514          LeafGridView;
515       typedef typename GridFamily::Traits::template Partition<pitype>::MacroGridView
516          MacroGridView;
517     };
518     /** \brief View types for All_Partition */
519     typedef typename Partition< All_Partition > :: LevelGridView LevelGridView;
520     typedef typename Partition< All_Partition > :: LeafGridView LeafGridView;
521     typedef typename Partition< All_Partition > :: MacroGridView MacroGridView;
522 
523     //! Type of the hierarchic index set
524     typedef ALU3dGridHierarchicIndexSet< dim, dimworld, elType, Comm > HierarchicIndexSet;
525 
526     //! Type of the level index set, needed by data handle
527     typedef typename GridFamily::LevelIndexSetImp LevelIndexSetImp;
528     //! Type of the leaf index set, needed by data handle
529     typedef typename GridFamily::LeafIndexSetImp LeafIndexSetImp;
530 
531     // type of container for reference elements
532     typedef ReferenceElements< alu3d_ctype, dim > ReferenceElementContainerType;
533     // type of container for reference faces
534     typedef ReferenceElements< alu3d_ctype, dim-1 > ReferenceFaceContainerType;
535 
536     // type of reference element
537     typedef std::decay_t< decltype( ReferenceElementContainerType::general( std::declval< const Dune::GeometryType & >() ) ) > ReferenceElementType;
538     // type of reference face
539     typedef std::decay_t< decltype( ReferenceFaceContainerType::general( std::declval< const Dune::GeometryType & >() ) ) > ReferenceFaceType;
540 
541     //! \brief boundary projection type
542     typedef typename Traits::DuneBoundaryProjectionType DuneBoundaryProjectionType;
543 
544     //! type of vertex projection
545     typedef ALU3DSPACE ProjectVertex  ALUGridVertexProjectionType;
546 
547     //! type of ALUGrid Vertex Projection Interface (shared_ptr)
548     typedef ALU3DSPACE ProjectVertexPtr       ALUGridVertexProjectionPointerType;
549     typedef ALU3DSPACE ProjectVertexPtrPair   ALUGridVertexProjectionPairType;
550 
551     //! type of collective communication object
552     typedef typename Traits::CollectiveCommunication CollectiveCommunication;
553 
554     typedef ALULeafCommunication< dim, dimworld, elType, Comm > LeafCommunication;
555     typedef ALULevelCommunication< dim, dimworld, elType, Comm > LevelCommunication;
556 
557   protected:
558     //! Type of the local id set
559     typedef typename GridFamily::LocalIdSetImp LocalIdSetImp;
560 
561     typedef typename GridFamily::GlobalIdSetImp GlobalIdSetImp;
562 
563   public:
564     //! Type of the global id set
565     typedef typename Traits::GlobalIdSet GlobalIdSet;
566 
567     //! Type of the local id set
568     typedef typename Traits::LocalIdSet LocalIdSet;
569 
570   protected:
571     typedef ALU3dGridLeafIterator< 0, All_Partition, const ThisType > LeafIteratorImp;
572     typedef typename Traits::template Codim< 0 >::LeafIterator LeafIteratorType;
573     typedef typename Traits::template Codim< 0 >::LeafIterator LeafIterator;
574 
575     typedef ALU3dGridHierarchicIterator< const ThisType > HierarchicIteratorImp;
576 
577     typedef typename ALU3dImplTraits< elType, Comm >::GitterImplType GitterImplType;
578 
579     //! element chunk for refinement
580     enum {
581       //! \brief normal default number of new elements for new adapt method
582       newElementsChunk_ = 128 };
583 
584     //! upper estimate on number of elements that could be created when a new element is created
585     enum {
586       /** \brief if one element is refined then it
587           causes apporximately not more than
588           this number of new elements  */
589       refineEstimate_ = 8 };
590 
591   public:
592     typedef Comm MPICommunicatorType;
593 
594     typedef ALU3dGridCommunications< dim, dimworld, elType, Comm > Communications;
595 
596   protected:
597     typedef ALU3dGridVertexList< Comm >     VertexListType;
598     typedef ALU3dGridLeafVertexList< Comm > LeafVertexListType;
599 
600     typedef DefaultBoundarySegmentIndexSet< ThisType > BoundarySegmentIndexSetType;
601 
602   public:
603     //! Constructor which reads an ALU3dGrid Macro Triang file
604     //! or given GridFile
605     ALU3dGrid ( const std::string &macroTriangFilename,
606                 const MPICommunicatorType mpiComm,
607                 const ALUGridVertexProjectionPairType& bndPrj,
608                 const ALUGridRefinementType refinementType );
609 
610     //! \brief Desctructor
~ALU3dGrid()611     virtual ~ALU3dGrid() {}
612 
613     //! \brief for grid identification
614     static inline std::string name ();
615 
616     /** \brief  Return maximum level defined in this grid. Levels are numbered
617         maxLevel with 0 the coarsest level.
618       */
619     int maxLevel() const;
620 
621     //! Iterator to first entity of given codim on level
622     template<int cd, PartitionIteratorType pitype>
623     typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator
624     lbegin (int level) const;
625 
626     //! one past the end on this level
627     template<int cd, PartitionIteratorType pitype>
628     typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator
629     lend (int level) const;
630 
631     //! Iterator to first entity of given codim on level
632     template<int cd>
633     typename Traits::template Codim<cd>::
634     template Partition<All_Partition>::LevelIterator
635     lbegin (int level) const;
636 
637     //! one past the end on this level
638     template<int cd>
639     typename Traits::template Codim<cd>::
640     template Partition<All_Partition>::LevelIterator
641     lend (int level) const;
642 
643     typedef LeafIntersectionIteratorWrapper  < const ThisType > LefInterItWrapperType;
644     typedef LevelIntersectionIteratorWrapper < const ThisType > LvlInterItWrapperType;
645 
646     typename Traits::LeafIntersectionIterator
ileafbegin(const typename Traits::template Codim<0>::Entity & entity) const647     ileafbegin( const typename Traits::template Codim< 0 >::Entity& entity ) const
648     {
649       return LefInterItWrapperType( *this,
650                                     entity.impl(),
651                                     entity.level(), false );
652     }
653 
654     typename Traits::LeafIntersectionIterator
ileafend(const typename Traits::template Codim<0>::Entity & entity) const655     ileafend( const typename Traits::template Codim< 0 >::Entity& entity ) const
656     {
657       return LefInterItWrapperType( *this,
658                                     entity.impl(),
659                                     entity.level(), true );
660     }
661 
662     typename Traits::LevelIntersectionIterator
ilevelbegin(const typename Traits::template Codim<0>::Entity & entity) const663     ilevelbegin( const typename Traits::template Codim< 0 >::Entity& entity ) const
664     {
665       return LvlInterItWrapperType( *this,
666                                     entity.impl(),
667                                     entity.level(), false );
668     }
669 
670     typename Traits::LevelIntersectionIterator
ilevelend(const typename Traits::template Codim<0>::Entity & entity) const671     ilevelend( const typename Traits::template Codim< 0 >::Entity& entity ) const
672     {
673       return LvlInterItWrapperType( *this,
674                                     entity.impl(),
675                                     entity.level(), true );
676     }
677 
678   public:
679     //! General definiton for a leaf iterator
680     template <int codim, PartitionIteratorType pitype>
681     typename Traits::template Codim<codim>::template Partition<pitype>::LeafIterator
682     leafbegin() const;
683 
684     //! General definition for an end iterator on leaf level
685     template <int codim, PartitionIteratorType pitype>
686     typename Traits::template Codim<codim>::template Partition<pitype>::LeafIterator
687     leafend() const;
688 
689     //! General definiton for a leaf iterator
690     template <int codim>
691     typename Traits::template Codim<codim>::LeafIterator
692     leafbegin() const;
693 
694     //! General definition for an end iterator on leaf level
695     template <int codim>
696     typename Traits::template Codim<codim>::LeafIterator
697     leafend() const;
698 
699   public:
700     //! number of grid entities per level and codim
701     int size (int level, int cd) const;
702 
703     //! number of leaf entities per codim in this process
704     int size (int codim) const;
705 
706     //! number of entities per level and geometry type in this process
707     int size (int level, GeometryType type) const;
708 
709     //! number of boundary segments
710     size_t numBoundarySegments() const;
711 
712     //! number of leaf entities per geometry type in this process
713     int size (GeometryType type) const;
714 
715     //! number of grid entities on all levels for given codim
716     int global_size (int cd) const ;
717 
718     // (no interface method) number of grid entities in the entire grid for given codim
719     int hierSetSize (int cd) const;
720 
721     //! get global id set of grid
globalIdSet() const722     const GlobalIdSet &globalIdSet () const
723     {
724       if( !globalIdSet_ )
725       {
726         globalIdSet_.reset( new GlobalIdSetImp( *this ) );
727       }
728       return *globalIdSet_;
729     }
730 
731     //! View for te macro grid with some alu specific methods
732     template<PartitionIteratorType pitype>
macroGridView() const733     typename Partition<pitype>::MacroGridView macroGridView() const
734     {
735       typedef typename Traits::template Partition<pitype>::MacroGridView View;
736       return View(*this);
737     }
738 
739     //! View for te macro grid with some alu specific methods (All_Partition)
macroGridView() const740     MacroGridView macroGridView() const
741     {
742       typedef MacroGridView View;
743       return View(*this);
744     }
745 
746     //! get global id set of grid
localIdSet() const747     const LocalIdSet & localIdSet () const { return localIdSet_; }
748 
749     //! get leaf index set of the grid
750     const typename Traits :: LeafIndexSet & leafIndexSet () const;
751 
752     //! get level index set of the grid
753     const typename Traits :: LevelIndexSet & levelIndexSet (int level) const;
754 
755     /** \brief return instance of level index set
756         \note if index set for this level has not been created then this
757         instance will be deleted once the shared_ptr goes out of scope.
758     */
759     std::shared_ptr< LevelIndexSetImp > accessLevelIndexSet ( int level ) const;
760 
761   protected:
762     std::shared_ptr< LevelIndexSetImp > createLevelIndexSet ( int level ) const;
763 
764   public:
765     template< int cd >
twists(GeometryType type) const766     typename Codim< cd >::Twists twists ( GeometryType type ) const
767     {
768       assert( type.dim() == dimension - cd );
769       assert( elType == tetra ? type.isSimplex() : type.isCube() );
770       return typename Traits::template Codim< cd >::Twists();
771     }
772 
773   protected:
774     typedef ALU3DSPACE GatherScatter GatherScatterType;
775 
776     /** \brief Calculates load of each process and repartition the grid if neccessary.
777         For parameters of the load balancing process see the README file
778         of the ALUGrid package.
779        \param data the data handler class that must implement three methods:
780           \code
781           // calls data inline on macro element. From there the data of
782           // all children can be written to the message buffer.
783           // MessageBufferImp implements the MessageBufferIF interface.
784           template<class MessageBufferImp>
785           void inlineData ( MessageBufferImp& buff, Dune::Entity<0> & e);
786 
787           // calls data xtract on macro element. From there the data of
788           // all children can be restored from the message buffer.
789           // numChildren is the number of all children underneath the
790           // macro element e.
791           // MessageBufferImp implements the MessageBufferIF interface.
792           template<class MessageBufferImp>
793           void xtractData ( MessageBufferImp& buff, Dune::Entity<0> & e, size_t numChildren );
794 
795           // This method is called at the end of the load balancing process
796           // before adaptation markers are removed. Here the user can apply
797           // a data compression or other features. This method can be
798           // empty if nothing should be done.
799           void compress ();
800           \endcode
801 
802          \return true if the grid has changed
803     */
804     bool loadBalance ( GatherScatterType* lbData );
805 
806   public:
807     /** \brief Set load balancing method and lower and upper bound for decision
808                on whether to load balance or not.
809 
810         \note  Possible choices of load balancing methods are:
811           \code
812             // no load balancing
813             NONE = 0
814 
815             // collect all to rank 0
816             COLLECT = 1,
817 
818             // assuming the elements to be ordered by a
819             // space filling curve approach
820             // here, the edges in the graph are neglected
821             // parallel version
822             ALUGRID_SpaceFillingCurveLinkage = 4,
823             // serial version that requires the whole graph to be avaiable
824             ALUGRID_SpaceFillingCurveSerialLinkage = 5,
825 
826             // METIS method for graph partitioning (with linkage storage)
827             //METIS_PartGraphKwayLinkage      = 6,
828             //METIS_PartGraphRecursiveLinkage = 7,
829 
830             // ALU sfc without linkage
831             ALUGRID_SpaceFillingCurve       = 9,
832             ALUGRID_SpaceFillingCurveSerial = 10,
833 
834             // METIS method for graph partitioning
835             METIS_PartGraphKway = 11,
836             METIS_PartGraphRecursive = 12,
837 
838             // ZOLTAN partitioning
839             ZOLTAN_LB_HSFC = 13 ,
840             ZOLTAN_LB_GraphPartitioning = 14 ,
841             ZOLTAN_LB_PARMETIS = 15
842          \endcode
843     */
setLoadBalanceMethod(const int mthd,const double ldbUnder=0.0,const double ldbOver=1.2)844     static void setLoadBalanceMethod( const int mthd,
845                                       const double ldbUnder = 0.0,
846                                       const double ldbOver  = 1.2 )
847     {
848       using DataBase = ALU3DSPACE LoadBalancer::DataBase ;
849       if( mthd < int( DataBase::NONE ) && mthd > DataBase::ZOLTAN_LB_PARMETIS )
850       {
851         DUNE_THROW(InvalidStateException,"ALUGrid::setLoadBalanceMethod: wrong method passed, check documentation for correect values");
852       }
853 
854       ALU3DSPACE ALUGridExternalParameters::setLoadBalanceParameters( mthd, ldbUnder, ldbOver );
855     }
856 
857 
858 
859     /** \brief Calculates load of each process and repartition by using ALUGrid's default partitioning method.
860                The specific load balancing algorithm is selected from a file alugrid.cfg.
861         \return true if grid has changed
862     */
loadBalance()863     bool loadBalance ()
864     {
865       return loadBalance( (GatherScatterType* ) 0 );
866     }
867 
868     /** \brief Calculates load of each process and repartition by using ALUGrid's default partitioning method.
869                The specific load balancing algorithm is selected from a file alugrid.cfg.
870         \param  optional dataHandleIF data handle that implements the Dune::CommDataHandleIF interface to include
871                 user data during load balancing
872         \return true if grid has changed
873     */
874     template< class DataHandleImpl, class Data >
loadBalance(CommDataHandleIF<DataHandleImpl,Data> & dataHandleIF)875     bool loadBalance ( CommDataHandleIF< DataHandleImpl, Data > &dataHandleIF )
876     {
877       typedef ALU3DSPACE GatherScatterLoadBalanceDataHandle
878             < ThisType, GatherScatterType, DataHandleImpl, Data, false > DataHandleType;
879       DataHandleType dataHandle( *this, dataHandleIF );
880 
881       // call the above loadBalance method with general GatherScatterType
882       return loadBalance( &dataHandle );
883     }
884 
885     /** \brief Calculates load of each process and repartition by using ALUGrid's default partitioning method,
886                the partitioning can be optimized by providing weights for each element on the macro grid.
887                The specific load balancing algorithm is selected from a file alugrid.cfg.
888         \param  weights class with double operator()(const Entity<0>&) returning a weight for each element
889                 which the includes in its internal loadbalancing process - for ALUGrid these are all macro elements.
890         \param  dataHandleIF data handle that implements the Dune::CommDataHandleIF interface to include
891                 user data during load balancing
892         \return true if grid has changed
893     */
894     template< class LBWeights, class DataHandleImpl, class Data >
loadBalance(LBWeights & weights,CommDataHandleIF<DataHandleImpl,Data> & dataHandleIF)895     bool loadBalance ( LBWeights &weights,
896                        CommDataHandleIF< DataHandleImpl, Data > &dataHandleIF )
897     {
898       typedef ALU3DSPACE GatherScatterLoadBalanceDataHandle
899             < ThisType, LBWeights, DataHandleImpl, Data, false > DataHandleType;
900       DataHandleType dataHandle( *this, dataHandleIF, weights );
901 
902       // call the above loadBalance method with general GatherScatterType
903       return loadBalance( &dataHandle );
904     }
905 
906     /** \brief Calculates load of each process and repartition by using ALUGrid's default partitioning method,
907                the partitioning can be optimized by providing weights for each element on the macro grid.
908                The specific load balancing algorithm is selected from a file alugrid.cfg.
909         \param  weights class with double operator()(const Entity<0>&) returning a weight for each element
910                 which the includes in its internal loadbalancing process - for ALUGrid these are all macro elements.
911         \return true if grid has changed
912     */
913     template< class LBWeights >
loadBalance(LBWeights & weights)914     typename std::enable_if< !IsDataHandle< LBWeights >::value, bool >::type loadBalance ( LBWeights &weights )
915     {
916       typedef ALU3DSPACE GatherScatterLoadBalance < ThisType, LBWeights, false > LoadBalanceHandleType;
917       LoadBalanceHandleType loadBalanceHandle( *this, weights );
918       return loadBalance( &loadBalanceHandle );
919     }
920 
921     /** \brief Distribute the grid based on a user defined partitioning.
922         \param  destinations class with int operator()(const Entity<0>&) returning the new owner process
923                 of this element. A destination has to be provided for all elements in the grid hierarchy
924                 but depending on the grid implementation it is possibly called on a subset only.
925                 The elements for which the method is called will be moved to the new processor together
926                 with all children. ALUGrid requires destinations for all macro elements.
927         \return true if grid has changed
928     */
929     template< class LBDestinations >
repartition(LBDestinations & destinations)930     bool repartition ( LBDestinations &destinations )
931     {
932       typedef ALU3DSPACE GatherScatterLoadBalance< ThisType, LBDestinations, true > LoadBalanceHandleType ;
933       LoadBalanceHandleType loadBalanceHandle( *this, destinations );
934       return loadBalance( &loadBalanceHandle );
935     }
936 
937     /** \brief Distribute the grid based on a user defined partitioning.
938         \param  destinations class with int operator()(const Entity<0>&) returning the new owner process
939                 of this element. A destination has to be provided for all elements in the grid hierarchy
940                 but depending on the grid implementation it is possibly called on a subset only.
941                 The elements for which the method is called will be moved to the new processor together
942                 with all children. ALUGrid requires destinations for all macro elements.
943         \param  dataHandleIF data handle that implements the Dune::CommDataHandleIF interface to include
944                 user data during load balancing
945         \return true if grid has changed
946     */
947     template< class LBDestinations, class DataHandleImpl, class Data >
repartition(LBDestinations & destinations,CommDataHandleIF<DataHandleImpl,Data> & dataHandleIF)948     bool repartition ( LBDestinations &destinations,
949                        CommDataHandleIF< DataHandleImpl, Data > &dataHandleIF )
950     {
951       typedef ALU3DSPACE GatherScatterLoadBalanceDataHandle< ThisType, LBDestinations, DataHandleImpl, Data, true > DataHandleType;
952       DataHandleType dataHandle( *this, dataHandleIF, destinations );
953 
954       // call the above loadBalance method with general GatherScatterType
955       return loadBalance( &dataHandle );
956     }
957 
958 
959     /** \brief ghostSize is one for codim 0 and zero otherwise for this grid  */
960     int ghostSize (int level, int codim) const;
961 
962     /** \brief overlapSize is zero for this grid  */
overlapSize(int level,int codim) const963     int overlapSize (int level, int codim) const { return 0; }
964 
965     /** \brief ghostSize is one for codim 0 and zero otherwise for this grid  */
966     int ghostSize (int codim) const;
967 
968     /** \brief overlapSize is zero for this grid  */
overlapSize(int codim) const969     int overlapSize (int codim) const { return 0; }
970 
971     /** \brief @copydoc Dune::Grid::communicate */
972     template< class DataHandle, class Data >
communicate(CommDataHandleIF<DataHandle,Data> & data,InterfaceType iftype,CommunicationDirection dir,int level) const973     LevelCommunication communicate ( CommDataHandleIF< DataHandle, Data > &data,
974                                      InterfaceType iftype,
975                                      CommunicationDirection dir,
976                                      int level ) const
977     {
978       return LevelCommunication( *this, data, iftype, dir, level );
979     }
980 
981     /** \brief Communicate information on distributed entities on the leaf grid.
982        Template parameter is a model of Dune::CommDataHandleIF.
983      */
984     template< class DataHandle, class Data >
communicate(CommDataHandleIF<DataHandle,Data> & data,InterfaceType iftype,CommunicationDirection dir) const985     LeafCommunication communicate ( CommDataHandleIF< DataHandle, Data > &data,
986                                     InterfaceType iftype,
987                                     CommunicationDirection dir ) const
988     {
989       return LeafCommunication( *this, data, iftype, dir );
990     }
991 
992   protected:
993     // load balance and compress memory if possible
994     void finalizeGridCreation();
995 
996     //! clear all entity new markers
997     void clearIsNewMarkers( );
998 
999   public:
1000     /** \brief @copydoc Dune::Grid::comm() */
comm() const1001     const CollectiveCommunication &comm () const { return communications().ccobj_; }
1002 
1003     //! returns if a least one entity was marked for coarsening
1004     bool preAdapt ( );
1005 
1006     //! clear all entity new markers if lockPostAdapt_ is set
1007     void postAdapt ( );
1008 
1009     /** \brief  @copydoc Dune::Grid::adapt() */
1010     bool adapt ();
1011 
1012     /** \brief  @copydoc Dune::Grid::adapt()
1013         \param handle handler for restriction and prolongation operations
1014         which is a Model of the AdaptDataHandleInterface class.
1015     */
1016     template< class GridImp, class DataHandle >
1017     bool adapt ( AdaptDataHandleInterface< GridImp, DataHandle > &handle );
1018 
1019     //! uses the interface, mark on entity and refineLocal
1020     void globalRefine ( int refCount );
1021 
1022     template< class GridImp, class DataHandle >
1023     void globalRefine ( int refCount, AdaptDataHandleInterface< GridImp, DataHandle > &handle );
1024 
1025     //**********************************************************
1026     // End of Interface Methods
1027     //**********************************************************
1028 
1029     /** \brief write macro grid in ALUGrid macro format to path/filename.rank */
1030     bool writeMacroGrid( const std::string path, const std::string filename,
1031                          const ALU3DSPACE MacroFileHeader::Format format = ALU3DSPACE MacroFileHeader::defaultFormat ) const ;
1032 
1033     /** \brief backup to ostream */
1034     void backup( std::ostream&, const ALU3DSPACE MacroFileHeader::Format format ) const ;
1035 
1036     /** \brief restore from istream */
1037     void restore( std::istream& ) ;
1038 
1039     // (no interface method) get hierarchic index set of the grid
hierarchicIndexSet() const1040     const HierarchicIndexSet & hierarchicIndexSet () const { return hIndexSet_; }
1041 
1042     // no interface method, but has to be public
1043     void updateStatus ();
1044 
1045     //! @copydoc Dune::Grid::mark
1046     bool mark( int refCount , const typename Traits::template Codim<0>::Entity & e);
1047 
1048     //! @copydoc Dune::Grid::getMark
1049     int getMark( const typename Traits::template Codim<0>::Entity & e) const;
1050 
1051   public:
defaultCommunicator()1052     static MPICommunicatorType defaultCommunicator ()
1053     {
1054       return Communications::defaultComm();
1055     }
1056 
1057     //! deliver all geometry types used in this grid
geomTypes(int codim) const1058     const std::vector<GeometryType>& geomTypes (int codim) const { return geomTypes_[codim]; }
1059 
1060     // return reference to org ALU3dGrid
1061     // private method, but otherwise we have to friend class all possible
1062     // types of LevelIterator ==> later
1063     GitterImplType &myGrid () const;
1064 
createALUGrid(const std::string & macroName)1065     virtual GitterImplType *createALUGrid ( const std::string &macroName )
1066     {
1067       alugrid_assert ( communications_ );
1068       return communications_->createALUGrid( macroName, vertexProjections(), conformingRefinement() );
1069     }
1070 
createALUGrid(std::istream & stream)1071     virtual GitterImplType *createALUGrid ( std::istream& stream )
1072     {
1073       alugrid_assert ( communications_ );
1074       return communications_->createALUGrid( stream, vertexProjections(), conformingRefinement() );
1075     }
1076 
vertexProjections() const1077     ALUGridVertexProjectionPairType vertexProjections() const
1078     {
1079       return  vertexProjections_ ;
1080     }
1081 
1082     // return appropriate ALUGrid builder
getBuilder() const1083     virtual typename ALU3DSPACE Gitter::Geometric::BuilderIF &getBuilder () const
1084     {
1085       return Communications::getBuilder( myGrid() );
1086     }
1087 
1088     // helper function for factory
completeGrid()1089     virtual void completeGrid ()
1090     {
1091       Communications::completeGrid( myGrid() );
1092       clearIsNewMarkers();
1093       // update macro boundary segment index
1094       macroBoundarySegmentIndexSet_.invalidate();
1095     }
1096 
1097     //! return reference to Dune reference element according to elType
referenceElement()1098     static const ReferenceElementType& referenceElement()
1099     {
1100       static const auto& refElem = ( elType == tetra ) ?
1101           Dune::ReferenceElements< alu3d_ctype, dimension >::simplex() :
1102           Dune::ReferenceElements< alu3d_ctype, dimension >::cube();
1103       return refElem ;
1104     }
1105 
1106     //! return reference to Dune face reference element according to elType
faceReferenceElement()1107     static const ReferenceFaceType& faceReferenceElement()
1108     {
1109       static const auto& refElem = ( elType == tetra ) ?
1110           Dune::ReferenceElements< alu3d_ctype, dimension-1 >::simplex() :
1111           Dune::ReferenceElements< alu3d_ctype, dimension-1 >::cube();
1112       return refElem ;
1113     }
1114 
1115     template < class EntitySeed >
1116     typename Traits :: template Codim< EntitySeed :: codimension > :: EntityPointer
entityPointer(const EntitySeed & seed) const1117     entityPointer( const EntitySeed& seed ) const
1118     {
1119       enum { codim = EntitySeed :: codimension };
1120       return typename Traits :: template Codim< codim > :: EntityPointerImpl( seed );
1121     }
1122 
1123     template < class EntitySeed >
1124     typename Traits :: template Codim< EntitySeed :: codimension > :: Entity
entity(const EntitySeed & seed) const1125     entity( const EntitySeed& seed ) const
1126     {
1127       typedef typename Traits :: template Codim< EntitySeed :: codimension > :: Entity  Entity;
1128       return Entity( typename Traits :: template Codim< EntitySeed :: codimension > :: EntityImp( seed ) );
1129     }
1130 
1131     // number of links to other processors, for internal use only
nlinks() const1132     int nlinks () const { return communications().nlinks(); }
1133 
getLeafVertexList() const1134     LeafVertexListType & getLeafVertexList() const
1135     {
1136       if( !leafVertexList_.up2Date() ) leafVertexList_.setupVxList(*this);
1137       return leafVertexList_;
1138     }
1139 
getLevelOfLeafVertex(const typename ALU3dImplTraits<elType,Comm>::VertexType & vertex) const1140     int getLevelOfLeafVertex ( const typename ALU3dImplTraits< elType, Comm >::VertexType &vertex ) const
1141     {
1142       alugrid_assert ( leafVertexList_.up2Date() );
1143       return leafVertexList_.getLevel(vertex);
1144     }
1145 
getVertexList(int level) const1146     VertexListType & getVertexList(int level) const
1147     {
1148       alugrid_assert ( level >= 0 );
1149       alugrid_assert ( level <= maxLevel() );
1150       VertexListType & vxList = vertexList_[level];
1151       if(!vxList.up2Date()) vxList.setupVxList(*this,level);
1152       return vxList;
1153     }
1154 
getGhostLeafList(int codim) const1155     ALU3dGridItemListType & getGhostLeafList(int codim) const
1156     {
1157       alugrid_assert ( codim >= 1 );
1158       alugrid_assert ( codim <= 3 );
1159       return ghostLeafList_[codim-1];
1160     }
1161 
getGhostLevelList(int codim,int level) const1162     ALU3dGridItemListType & getGhostLevelList(int codim, int level) const
1163     {
1164       alugrid_assert ( codim >= 1 );
1165       alugrid_assert ( codim <= 3 );
1166 
1167       alugrid_assert ( level >= 0 );
1168       alugrid_assert ( level <= maxLevel() );
1169       alugrid_assert ( level < int(ghostLevelList_[codim-1].size()) );
1170       return ghostLevelList_[codim-1][level];
1171     }
1172 
getEdgeList(int level) const1173     ALU3dGridItemListType & getEdgeList(int level) const
1174     {
1175       alugrid_assert ( level >= 0 );
1176       alugrid_assert ( level <= maxLevel() );
1177       return levelEdgeList_[level];
1178     }
1179 
1180   protected:
1181     //! Copy constructor should not be used
1182     ALU3dGrid( const ThisType & );
1183 
1184     //! assignment operator should not be used
1185     const ThisType &operator= ( const ThisType & );
1186 
1187     //! reset size and global size, update Level- and LeafIndexSet, if they exist
1188     void calcExtras();
1189 
1190     //! calculate maxlevel
1191     void calcMaxLevel();
1192 
1193     //! make grid walkthrough and calc global size
1194     void recalcGlobalSize();
1195 
1196     //! check whether macro grid format is of our type
1197     void checkMacroGridFile (const std::string filename);
1198 
1199     //! check whether macro grid has the right element type
1200     void checkMacroGrid ();
1201 
communications() const1202     const Communications &communications () const
1203     {
1204       alugrid_assert ( communications_ );
1205       return *communications_;
1206     }
1207 
1208     // initialize geometry types and return correct geometryInFather storage
1209     void makeGeometries();
1210 
1211   public:
1212     // return true if conforming refinement is enabled
conformingRefinement() const1213     bool conformingRefinement() const
1214     {
1215       return (refinementType_ == conforming) ;
1216     }
1217 
1218     // return true if ghost cells are available
ghostCellsEnabled() const1219     bool ghostCellsEnabled () const
1220     {
1221       return comm().size() > 1 && myGrid().ghostCellsEnabled();
1222     }
1223 
macroBoundarySegmentIndexSet() const1224     const BoundarySegmentIndexSetType& macroBoundarySegmentIndexSet() const
1225     {
1226       if( ! macroBoundarySegmentIndexSet_.valid() )
1227       {
1228         macroBoundarySegmentIndexSet_.update( macroGridView() );
1229       }
1230       alugrid_assert( macroBoundarySegmentIndexSet_.valid() );
1231       return macroBoundarySegmentIndexSet_;
1232     }
1233 
1234   protected:
1235     /////////////////////////////////////////////////////////////////
1236     //
1237     // Internal variables
1238     //
1239     /////////////////////////////////////////////////////////////////
1240 
1241     // the real ALU grid
1242     mutable std::unique_ptr< GitterImplType > mygrid_;
1243 
1244     // max level of grid
1245     int maxlevel_;
1246 
1247     // count how much elements where marked
1248     mutable int coarsenMarked_;
1249     mutable int refineMarked_;
1250 
1251     // at the moment the number of different geom types is 1
1252     enum { numberOfGeomTypes = 1 };
1253     std::vector< std::vector<GeometryType> > geomTypes_;
1254 
1255     // our hierarchic index set
1256     HierarchicIndexSet hIndexSet_;
1257 
1258     // out global id set
1259     mutable std::unique_ptr< GlobalIdSetImp > globalIdSet_;
1260 
1261     // out global id set
1262     LocalIdSetImp localIdSet_;
1263 
1264     // the level index set ( default type )
1265     mutable std::vector < std::shared_ptr< LevelIndexSetImp > > levelIndexVec_;
1266 
1267     // the leaf index set
1268     mutable std::unique_ptr< LeafIndexSetImp > leafIndexSet_;
1269 
1270     mutable std::vector< VertexListType > vertexList_;
1271 
1272     //the ghostleaf list is used in alu3diterators, where we use the internal aluIterators
1273     // the vertex codim there is 3, so the list has to fulfill that
1274     mutable ALU3dGridItemListType ghostLeafList_[ 3 ];
1275     mutable std::vector< ALU3dGridItemListType > ghostLevelList_[ 3 ];
1276 
1277     mutable std::vector< ALU3dGridItemListType > levelEdgeList_;
1278 
1279     mutable LeafVertexListType leafVertexList_;
1280 
1281     // the type of our size cache
1282     typedef SizeCache<MyType> SizeCacheType;
1283     std::unique_ptr< SizeCacheType > sizeCache_;
1284 
1285     // macro boundary segment index
1286     mutable BoundarySegmentIndexSetType macroBoundarySegmentIndexSet_;
1287 
1288     // variable to ensure that postAdapt ist called after adapt
1289     bool lockPostAdapt_;
1290 
1291     // boundary projection for vertices
1292     // pair: first is globalProjection_ for boundaries
1293     // second is surfaceProjection_ for manifolds
1294     ALUGridVertexProjectionPairType vertexProjections_ ;
1295 
1296     // pointer to communications object
1297     std::unique_ptr< Communications > communications_;
1298 
1299     // refinement type (nonconforming or conforming)
1300     const ALUGridRefinementType refinementType_ ;
1301   }; // end class ALU3dGrid
1302 
1303 
1304     bool checkMacroGrid ( ALU3dGridElementType elType ,
1305                           const std::string filename );
1306     const char* elType2Name( ALU3dGridElementType elType );
1307 
1308   namespace Capabilities
1309   {
1310 
1311     template< int dim, int dimworld, ALU3dGridElementType elType, class Comm, int cdim >
1312     struct hasEntity< Dune::ALU3dGrid< dim, dimworld, elType, Comm >, cdim >
1313     {
1314       static const bool v = true;
1315     };
1316 
1317     template< int dim, int dimworld,  ALU3dGridElementType elType, class Comm >
1318     struct isLevelwiseConforming< ALU3dGrid< dim, dimworld, elType, Comm > >
1319     {
1320       static const bool v = true;
1321     };
1322 
1323     template< int dim, int dimworld, ALU3dGridElementType elType, class Comm >
1324     struct hasBackupRestoreFacilities< ALU3dGrid< dim, dimworld, elType, Comm > >
1325     {
1326       static const bool v = true;
1327     };
1328 
1329   } // end namespace Capabilities
1330 
1331 } // end namespace Dune
1332 
1333 #include "grid_inline.hh"
1334 #if COMPILE_ALUGRID_INLINE
1335   #include "grid_imp.cc"
1336 #endif
1337 #endif
1338