dune-grid 2.9.0
uggrid.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5
6#ifndef DUNE_UGGRID_HH
7#define DUNE_UGGRID_HH
8
13#include <memory>
14
15#include <dune/common/classname.hh>
16#include <dune/common/parallel/communication.hh>
17#include <dune/common/exceptions.hh>
18#include <dune/common/parallel/mpihelper.hh>
19
23
24#if HAVE_DUNE_UGGRID || DOXYGEN
25
26#ifdef ModelP
27#include <dune/common/parallel/mpicommunication.hh>
28#endif
29
30/* [Before reading the following: the macros UG_DIM_2 and UG_DIM_3 where named
31 * _2 and _3, respectively, up until ug-3.12.0.]
32 *
33 * The following lines including the necessary UG headers are somewhat
34 tricky. Here's what's happening:
35 UG can support two- and three-dimensional grids. You choose by setting
36 either UG_DIM_2 or UG_DIM_3 while compiling. This changes all sorts of stuff, in
37 particular data structures in the headers.
38 UG was never supposed to provide 2d and 3d grids at the same time.
39 However, when compiling it as c++, the dimension-dependent parts are
40 wrapped up cleanly in the namespaces UG::D2 and UG::D3, respectively. That
41 way it is possible to link together the UG lib for 2d and the one for 3d.
42 But we also need the headers twice! Once with UG_DIM_2 set and once with UG_DIM_3!
43 So here we go:*/
44
45/* The following define tells the UG headers that we want access to a few
46 special fields, for example the extra index fields in the element data structures.
47 This define remains only for backwards compatibility with older version of UG.
48 All dune-uggrid versions since 2016-08-05 do not need this #define or the #undef
49 further below. */
50#define FOR_DUNE
51
52// Set UG's space-dimension flag to 2d
53#define UG_DIM_2
54// And include all necessary UG headers
55#include "uggrid/ugincludes.hh"
56#undef DUNE_UGINCLUDES_HH
57
58// Wrap a few large UG macros by functions before they get undef'ed away.
59// Here: The 2d-version of the macros
60#define UG_DIM 2
61#include "uggrid/ugwrapper.hh"
62#undef UG_DIM
63
64// UG defines a whole load of preprocessor macros. ug_undefs.hh undefines
65// them all, so we don't get name clashes.
66#include "uggrid/ug_undefs.hh"
67#undef UG_DIM_2
68
69/* Now we're done with 2d, and we can do the whole thing over again for 3d */
70
71/* All macros set by UG have been unset. This includes the macros that ensure
72 single inclusion of headers. We can thus include them again. However, we
73 only want to include those headers again that contain dimension-dependent stuff.
74 Therefore, we set a few single-inclusion defines manually before including
75 ugincludes.hh again.
76 */
77#define UGTYPES_H
78#define __HEAPS__
79#define __UGENV__
80#define __DEVICESH__
81#ifdef ModelP
82#define __PPIF__
83#endif
84
85#define UG_DIM_3
86#include "uggrid/ugincludes.hh"
87#undef DUNE_UGINCLUDES_HH
88
89// Wrap a few large UG macros by functions before they get undef'ed away.
90// This time it's the 3d-versions.
91#define UG_DIM 3
92#include "uggrid/ugwrapper.hh"
93#undef UG_DIM
94
95// undef all macros defined by UG
96#include "uggrid/ug_undefs.hh"
97
98#undef UG_DIM_3
99#undef FOR_DUNE
100
101// The components of the UGGrid interface
102#include "uggrid/uggridgeometry.hh"
103#include "uggrid/uggridlocalgeometry.hh"
104#include "uggrid/uggridentity.hh"
105#include "uggrid/uggridentityseed.hh"
106#include "uggrid/uggridintersections.hh"
107#include "uggrid/uggridintersectioniterators.hh"
108#include "uggrid/uggridleveliterator.hh"
109#include "uggrid/uggridleafiterator.hh"
110#include "uggrid/uggridhieriterator.hh"
111#include "uggrid/uggridindexsets.hh"
112#include <dune/grid/uggrid/uggridviews.hh>
113#ifdef ModelP
114#include "uggrid/ugmessagebuffer.hh"
115#include "uggrid/uglbgatherscatter.hh"
116#endif
117
118// Not needed here, but included for user convenience
120
121#ifdef ModelP
122template <class DataHandle, int GridDim, int codim>
123const Dune::UGGrid<GridDim>* Dune::UGMessageBuffer<DataHandle, GridDim, codim>::grid_;
124
125template <class DataHandle, int GridDim, int codim>
126DataHandle *Dune::UGMessageBuffer<DataHandle,GridDim,codim>::duneDataHandle_ = nullptr;
127
128template <class DataHandle, int GridDim, int codim>
129int Dune::UGMessageBuffer<DataHandle,GridDim,codim>::level = -1;
130#endif // ModelP
131
132namespace Dune {
133
134#ifdef ModelP
135 using UGCommunication = Communication<MPI_Comm>;
136#else
137 using UGCommunication = Communication<No_Comm>;
138#endif
139
140 template<int dim>
142 {
144 UGGridGeometry,
145 UGGridEntity,
146 UGGridLevelIterator,
147 UGGridLeafIntersection,
148 UGGridLevelIntersection,
149 UGGridLeafIntersectionIterator,
150 UGGridLevelIntersectionIterator,
151 UGGridHierarchicIterator,
152 UGGridLeafIterator,
153 UGGridLevelIndexSet< const UGGrid<dim> >,
154 UGGridLeafIndexSet< const UGGrid<dim> >,
155 UGGridIdSet< const UGGrid<dim> >,
156 typename UG_NS<dim>::UG_ID_TYPE,
157 UGGridIdSet< const UGGrid<dim> >,
158 typename UG_NS<dim>::UG_ID_TYPE,
160 UGGridLevelGridViewTraits,
161 UGGridLeafGridViewTraits,
162 UGGridEntitySeed,
163 UGGridLocalGeometry>
165 };
166
167
168 //**********************************************************************
169 //
170 // --UGGrid
171 //
172 //**********************************************************************
173
205 template <int dim>
206 class UGGrid : public GridDefaultImplementation <dim, dim, double, UGGridFamily<dim> >
207 {
209
210 friend class UGGridGeometry<0,dim,const UGGrid<dim> >;
211 friend class UGGridGeometry<dim,dim,const UGGrid<dim> >;
212 friend class UGGridGeometry<1,2,const UGGrid<dim> >;
213 friend class UGGridGeometry<2,3,const UGGrid<dim> >;
214
215 friend class UGGridEntity <0,dim,const UGGrid<dim> >;
216 friend class UGGridEntity <1,dim,const UGGrid<dim> >;
217 friend class UGGridEntity <2,dim,const UGGrid<dim> >;
218 friend class UGGridEntity <dim,dim,const UGGrid<dim> >;
219 friend class UGGridHierarchicIterator<const UGGrid<dim> >;
220 friend class UGGridLeafIntersection<const UGGrid<dim> >;
221 friend class UGGridLevelIntersection<const UGGrid<dim> >;
222 friend class UGGridLeafIntersectionIterator<const UGGrid<dim> >;
223 friend class UGGridLevelIntersectionIterator<const UGGrid<dim> >;
224
225 friend class UGGridLevelIndexSet<const UGGrid<dim> >;
226 friend class UGGridLeafIndexSet<const UGGrid<dim> >;
227 friend class UGGridIdSet<const UGGrid<dim> >;
228 template <class GridImp_>
229 friend class UGGridLeafGridView;
230 template <class GridImp_>
232
233 friend class GridFactory<UGGrid<dim> >;
234
235#ifdef ModelP
236 friend class UGLBGatherScatter;
237#endif
238
239 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
240 friend class UGGridLeafIterator;
241 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
243
245 static_assert(dim==2 || dim==3, "Use UGGrid only for 2d and 3d!");
246
247 // The different instantiations are mutual friends so they can access
248 // each others numOfUGGrids field
249 friend class UGGrid<2>;
250 friend class UGGrid<3>;
251
252 //**********************************************************
253 // The Interface Methods
254 //**********************************************************
255 public:
258
259 // the Traits
261
263 typedef UG::DOUBLE ctype;
264
266 typedef unsigned int Rank;
267
271
273 ~UGGrid() noexcept(false);
274
277 int maxLevel() const;
278
280 template <typename Seed>
281 typename Traits::template Codim<Seed::codimension>::Entity
282 entity(const Seed& seed) const
283 {
284 const int codim = Seed::codimension;
285 return typename Traits::template Codim<codim>::Entity(UGGridEntity<codim,dim,const UGGrid<dim> >(seed.impl().target(),this));
286 }
287
290 int size (int level, int codim) const;
291
293 int size (int codim) const
294 {
295 return leafIndexSet().size(codim);
296 }
297
299 int size (int level, GeometryType type) const
300 {
301 return this->levelIndexSet(level).size(type);
302 }
303
305 int size (GeometryType type) const
306 {
307 return this->leafIndexSet().size(type);
308 }
309
311 size_t numBoundarySegments() const {
312 // The number is stored as a member of UGGrid upon grid creation.
313 // The corresponding data structure is not exported by UG. (It is in ug/dom/std/std_internal.h)
314 return numBoundarySegments_;
315 }
316
318 const typename Traits::GlobalIdSet& globalIdSet() const
319 {
320 return idSet_;
321 }
322
324 const typename Traits::LocalIdSet& localIdSet() const
325 {
326 return idSet_;
327 }
328
330 const typename Traits::LevelIndexSet& levelIndexSet(int level) const
331 {
332 if (level<0 || level>maxLevel())
333 DUNE_THROW(GridError, "levelIndexSet of nonexisting level " << level << " requested!");
334 return *levelIndexSets_[level];
335 }
336
338 const typename Traits::LeafIndexSet& leafIndexSet() const
339 {
340 return leafIndexSet_;
341 }
342
345
358 bool mark(int refCount, const typename Traits::template Codim<0>::Entity & e );
359
415 bool mark(const typename Traits::template Codim<0>::Entity & e,
416 typename UG_NS<dim>::RefinementRule rule,
417 int side=0);
418
420 int getMark(const typename Traits::template Codim<0>::Entity& e) const;
421
424 bool preAdapt();
425
427 bool adapt();
428
430 void postAdapt();
440 template<class DataHandle>
441 bool loadBalance (DataHandle& dataHandle)
442 {
443#ifdef ModelP
444 // gather element data
445 if (dataHandle.contains(dim, 0))
446 UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
447
448 // gather node data
449 if (dataHandle.contains(dim,dim))
450 UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
451#endif
452
453 // the load balancing step now also attaches
454 // the data to the entities and distributes it
455 loadBalance();
456
457#ifdef ModelP
458 // scatter element data
459 if (dataHandle.contains(dim, 0))
460 UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
461
462 // scatter node data
463 if (dataHandle.contains(dim,dim))
464 UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
465#endif
466
467 return true;
468 }
469
476 bool loadBalance(int minlevel=0);
477
508 bool loadBalance(const std::vector<Rank>& targetProcessors, unsigned int fromLevel);
509
519 template<class DataHandle>
520 bool loadBalance (const std::vector<Rank>& targetProcessors, unsigned int fromLevel, DataHandle& dataHandle)
521 {
522#ifdef ModelP
523 // gather element data
524 if (dataHandle.contains(dim, 0))
525 UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
526
527 // gather node data
528 if (dataHandle.contains(dim,dim))
529 UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
530#endif
531
532 // the load balancing step now also attaches
533 // the data to the entities and distributes it
534 loadBalance(targetProcessors,fromLevel);
535
536#ifdef ModelP
537 // scatter element data
538 if (dataHandle.contains(dim, 0))
539 UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
540
541 // scatter node data
542 if (dataHandle.contains(dim,dim))
543 UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
544#endif
545
546 return true;
547 }
548
550 const UGCommunication& comm () const
551 {
552 return ccobj_;
553 }
554
555 protected:
556#ifdef ModelP
557 template <int codim, class GridView, class DataHandle>
558 void communicateUG_(const GridView& gv, int level,
559 DataHandle &dataHandle,
560 InterfaceType iftype,
561 CommunicationDirection dir) const
562 {
563 typename UG_NS<dim>::DDD_IF_DIR ugIfDir;
564 // Translate the communication direction from Dune-Speak to UG-Speak
565 if (dir==ForwardCommunication)
566 ugIfDir = UG_NS<dim>::IF_FORWARD();
567 else
568 ugIfDir = UG_NS<dim>::IF_BACKWARD();
569
570 typedef UGMessageBuffer<DataHandle,dim,codim> UGMsgBuf;
571 UGMsgBuf::duneDataHandle_ = &dataHandle;
572
573 UGMsgBuf::level = level;
574
575 std::vector<typename UG_NS<dim>::DDD_IF> ugIfs = findDDDInterfaces(iftype, codim);
576
577 unsigned bufSize = UGMsgBuf::ugBufferSize(gv);
578 if (!bufSize)
579 return; // we don't need to communicate if we don't have any data!
580 UGMsgBuf::grid_ = this;
581 for (unsigned i=0; i < ugIfs.size(); ++i)
582 UG_NS<dim>::DDD_IFOneway(multigrid_->dddContext(),
583 ugIfs[i],
584 ugIfDir,
585 bufSize,
586 &UGMsgBuf::ugGather_,
587 &UGMsgBuf::ugScatter_);
588 }
589
591 std::vector<typename UG_NS<dim>::DDD_IF> findDDDInterfaces(InterfaceType iftype,
592 int codim) const;
593#endif
594 public:
595 // **********************************************************
596 // End of Interface Methods
597 // **********************************************************
598
606 void getChildrenOfSubface(const typename Traits::template Codim<0>::Entity & e,
607 int elementSide,
608 int maxl,
609 std::vector<typename Traits::template Codim<0>::Entity>& childElements,
610 std::vector<unsigned char>& childElementSides) const;
611
617 COPY
618 };
619
625 NONE
626 };
627
630 refinementType_ = type;
631 }
632
635 closureType_ = type;
636 }
637
641 void setPosition(const typename Traits::template Codim<dim>::Entity& e,
642 const FieldVector<double, dim>& pos);
643
648 void globalRefine(int n);
649
654 void saveState(const std::string& filename) const;
655
660 void loadState(const std::string& filename);
661
662 private:
664 typename UG_NS<dim>::MultiGrid* multigrid_;
665
667 UGCommunication ccobj_;
668
674 void setIndices(bool setLevelZero,
675 std::vector<unsigned int>* nodePermutation);
676
677 // Each UGGrid object has a unique name to identify it in the
678 // UG environment structure
679 std::string name_;
680
681 // Our set of level indices
682 std::vector<std::shared_ptr<UGGridLevelIndexSet<const UGGrid<dim> > > > levelIndexSets_;
683
684 UGGridLeafIndexSet<const UGGrid<dim> > leafIndexSet_;
685
686 // One id set implementation
687 // Used for both the local and the global UGGrid id sets
688 UGGridIdSet<const UGGrid<dim> > idSet_;
689
691 RefinementType refinementType_;
692
694 ClosureType closureType_;
695
703 static int numOfUGGrids;
704
710 bool someElementHasBeenMarkedForRefinement_;
711
717 bool someElementHasBeenMarkedForCoarsening_;
718
720 std::vector<std::shared_ptr<BoundarySegment<dim> > > boundarySegments_;
721
727 unsigned int numBoundarySegments_;
728
729 }; // end Class UGGrid
730
731 namespace Capabilities
732 {
748 template<int dim, int codim>
749 struct hasEntity< UGGrid<dim>, codim>
750 {
751 static const bool v = true;
752 };
753
759 template<int dim, int codim>
760 struct hasEntityIterator<UGGrid<dim>, codim>
761 {
762 static const bool v = false;
763 };
764
769 template<int dim>
770 struct hasEntityIterator<UGGrid<dim>, 0>
771 {
772 static const bool v = true;
773 };
774
779 template<int dim>
780 struct hasEntityIterator<UGGrid<dim>, dim>
781 {
782 static const bool v = true;
783 };
784
788 template<int dim, int codim>
789 struct canCommunicate<UGGrid<dim>, codim>
790 {
791 static const bool v = (codim>=0 && codim<=dim);
792 };
793
797 template<int dim>
799 {
800 static const bool v = true;
801 };
802
806 template<int dim>
808 {
809 static const bool v = false;
810 };
811
812 }
813
814} // namespace Dune
815
816#endif // HAVE_DUNE_UGGRID || DOXYGEN
817#endif // DUNE_UGGRID_HH
Base class for grid boundary segments of arbitrary geometry.
The specialization of the generic GridFactory for UGGrid.
CommunicationDirection
Define a type for communication direction parameter.
Definition: gridenums.hh:170
InterfaceType
Parameter to be used for the communication functions.
Definition: gridenums.hh:86
@ ForwardCommunication
communicate as given in InterfaceType
Definition: gridenums.hh:171
Include standard header files.
Definition: agrid.hh:60
Communication< No_Comm > UGCommunication
Definition: uggrid.hh:137
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
Specialize with 'true' for all codims that a grid implements entities for. (default=false)
Definition: common/capabilities.hh:58
static const bool v
Definition: common/capabilities.hh:59
specialize with 'true' for all codims that a grid provides an iterator for (default=hasEntity<codim>:...
Definition: common/capabilities.hh:74
static const bool v
Definition: common/capabilities.hh:75
specialize with 'true' for all codims that a grid can communicate data on (default=false)
Definition: common/capabilities.hh:97
static const bool v
Definition: common/capabilities.hh:98
Specialize with 'true' if implementation guarantees conforming level grids. (default=false)
Definition: common/capabilities.hh:106
static const bool v
Definition: common/capabilities.hh:107
Specialize with 'true' if implementation guarantees a conforming leaf grid. (default=false)
Definition: common/capabilities.hh:115
static const bool v
Definition: common/capabilities.hh:116
Wrapper class for entities.
Definition: common/entity.hh:66
Base class for exceptions in Dune grid modules.
Definition: exceptions.hh:20
Definition: common/grid.hh:862
Traits::LeafGridView leafGridView() const
View for the leaf grid for All_Partition.
Definition: common/grid.hh:882
bool loadBalance()
default implementation of load balance does nothing and returns false
Definition: common/grid.hh:948
Index Set Interface base class.
Definition: indexidset.hh:78
auto size(GeometryType type) const
Return total number of entities of given geometry type in entity set .
Definition: indexidset.hh:223
Id Set Interface.
Definition: indexidset.hh:452
A traits struct that collects all associated types of one grid model.
Definition: common/grid.hh:995
Provide a generic factory class for unstructured grids.
Definition: common/gridfactory.hh:314
Grid view abstract base class.
Definition: common/gridview.hh:66
Definition: uggrid.hh:142
GridTraits< dim, dim, Dune::UGGrid< dim >, UGGridGeometry, UGGridEntity, UGGridLevelIterator, UGGridLeafIntersection, UGGridLevelIntersection, UGGridLeafIntersectionIterator, UGGridLevelIntersectionIterator, UGGridHierarchicIterator, UGGridLeafIterator, UGGridLevelIndexSet< const UGGrid< dim > >, UGGridLeafIndexSet< const UGGrid< dim > >, UGGridIdSet< const UGGrid< dim > >, typename UG_NS< dim >::UG_ID_TYPE, UGGridIdSet< const UGGrid< dim > >, typename UG_NS< dim >::UG_ID_TYPE, UGCommunication, UGGridLevelGridViewTraits, UGGridLeafGridViewTraits, UGGridEntitySeed, UGGridLocalGeometry > Traits
Definition: uggrid.hh:164
Front-end for the grid manager of the finite element toolbox UG3.
Definition: uggrid.hh:207
void postAdapt()
Clean up refinement markers.
size_t numBoundarySegments() const
Return the number of boundary segments.
Definition: uggrid.hh:311
void setRefinementType(RefinementType type)
Sets the type of grid refinement.
Definition: uggrid.hh:629
int getMark(const typename Traits::template Codim< 0 >::Entity &e) const
Query whether element is marked for refinement.
bool mark(int refCount, const typename Traits::template Codim< 0 >::Entity &e)
Mark element for refinement.
friend class UGGridLeafIterator
Definition: uggrid.hh:240
bool loadBalance(int minlevel=0)
Distributes this grid over the available nodes in a distributed machine.
~UGGrid() noexcept(false)
Destructor.
UGGridFamily< dim >::Traits Traits
Definition: uggrid.hh:260
friend class UGGridLevelIterator
Definition: uggrid.hh:242
bool loadBalance(const std::vector< Rank > &targetProcessors, unsigned int fromLevel, DataHandle &dataHandle)
Distributes the grid over the processes of a parallel machine, and sends data along with it.
Definition: uggrid.hh:520
void getChildrenOfSubface(const typename Traits::template Codim< 0 >::Entity &e, int elementSide, int maxl, std::vector< typename Traits::template Codim< 0 >::Entity > &childElements, std::vector< unsigned char > &childElementSides) const
Rudimentary substitute for a hierarchic iterator on faces.
friend class UGGridLevelGridView
Definition: uggrid.hh:231
UGGridFamily< dim > GridFamily
type of the used GridFamily for this grid
Definition: uggrid.hh:257
bool loadBalance(const std::vector< Rank > &targetProcessors, unsigned int fromLevel)
Distribute this grid over a distributed machine.
const UGCommunication & comm() const
Definition: uggrid.hh:550
bool mark(const typename Traits::template Codim< 0 >::Entity &e, typename UG_NS< dim >::RefinementRule rule, int side=0)
Mark method accepting a UG refinement rule.
int maxLevel() const
const Traits::LevelIndexSet & levelIndexSet(int level) const
Access to the LevelIndexSets.
Definition: uggrid.hh:330
void setClosureType(ClosureType type)
Sets the type of grid refinement closure.
Definition: uggrid.hh:634
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition: uggrid.hh:324
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition: uggrid.hh:299
RefinementType
The different forms of grid refinement that UG supports.
Definition: uggrid.hh:613
@ COPY
New level consists of the refined elements and the unrefined ones, too.
Definition: uggrid.hh:617
@ LOCAL
New level consists only of the refined elements and the closure.
Definition: uggrid.hh:615
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: uggrid.hh:305
void setPosition(const typename Traits::template Codim< dim >::Entity &e, const FieldVector< double, dim > &pos)
Sets a vertex to a new position.
int size(int level, int codim) const
Number of grid entities per level and codim.
UG::DOUBLE ctype
The type used to store coordinates.
Definition: uggrid.hh:263
bool adapt()
Triggers the grid refinement process.
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition: uggrid.hh:318
friend class UGGridLeafGridView
Definition: uggrid.hh:229
void loadState(const std::string &filename)
Read entire grid hierarchy from disk.
bool loadBalance(DataHandle &dataHandle)
Distributes the grid and some data over the available nodes in a distributed machine.
Definition: uggrid.hh:441
Traits::template Codim< Seed::codimension >::Entity entity(const Seed &seed) const
Create an Entity from an EntitySeed.
Definition: uggrid.hh:282
unsigned int Rank
The type used for process ranks.
Definition: uggrid.hh:266
void saveState(const std::string &filename) const
Save entire grid hierarchy to disk.
bool preAdapt()
returns true, if some elements might be coarsend during grid adaption, here always returns true
void globalRefine(int n)
Does uniform refinement.
UGGrid(UGCommunication comm={})
Default constructor.
ClosureType
Decide whether to add a green closure to locally refined grid sections or not.
Definition: uggrid.hh:621
@ GREEN
Standard red/green refinement.
Definition: uggrid.hh:623
@ NONE
No closure, results in nonconforming meshes.
Definition: uggrid.hh:625
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition: uggrid.hh:338
int size(int codim) const
number of leaf entities per codim in this process
Definition: uggrid.hh:293
A set of traits classes to store static information about grid implementation.
Different resources needed by all grid implementations.