dune-grid 2.11
Loading...
Searching...
No Matches
uggrid.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © 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
12
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
20#include <dune-grid-config.hh> // HAVE_DUNE_UGGRID
24
25#if HAVE_DUNE_UGGRID || DOXYGEN
26
27#ifdef ModelP
28#include <dune/common/parallel/mpicommunication.hh>
29#endif
30
31/* [Before reading the following: the macros UG_DIM_2 and UG_DIM_3 where named
32 * _2 and _3, respectively, up until ug-3.12.0.]
33 *
34 * The following lines including the necessary UG headers are somewhat
35 tricky. Here's what's happening:
36 UG can support two- and three-dimensional grids. You choose by setting
37 either UG_DIM_2 or UG_DIM_3 while compiling. This changes all sorts of stuff, in
38 particular data structures in the headers.
39 UG was never supposed to provide 2d and 3d grids at the same time.
40 However, when compiling it as c++, the dimension-dependent parts are
41 wrapped up cleanly in the namespaces UG::D2 and UG::D3, respectively. That
42 way it is possible to link together the UG lib for 2d and the one for 3d.
43 But we also need the headers twice! Once with UG_DIM_2 set and once with UG_DIM_3!
44 So here we go:*/
45
46// Set UG's space-dimension flag to 2d
47#define UG_DIM_2
48// And include all necessary UG headers
49#include "uggrid/ugincludes.hh"
50#undef DUNE_UGINCLUDES_HH
51
52// Wrap a few large UG macros by functions before they get undef'ed away.
53// Here: The 2d-version of the macros
54#define UG_DIM 2
55#include "uggrid/ugwrapper.hh"
56#undef UG_DIM
57
58// UG defines a whole load of preprocessor macros. ug_undefs.hh undefines
59// them all, so we don't get name clashes.
60#include "uggrid/ug_undefs.hh"
61#undef UG_DIM_2
62
63/* Now we're done with 2d, and we can do the whole thing over again for 3d */
64
65/* All macros set by UG have been unset. This includes the macros that ensure
66 single inclusion of headers. We can thus include them again. However, we
67 only want to include those headers again that contain dimension-dependent stuff.
68 Therefore, we set a few single-inclusion defines manually before including
69 ugincludes.hh again.
70 */
71#define UGTYPES_H
72#define __HEAPS__
73#define __UGENV__
74#define __DEVICESH__
75#ifdef ModelP
76#define __PPIF__
77#endif
78
79#define UG_DIM_3
80#include "uggrid/ugincludes.hh"
81#undef DUNE_UGINCLUDES_HH
82
83// Wrap a few large UG macros by functions before they get undef'ed away.
84// This time it's the 3d-versions.
85#define UG_DIM 3
86#include "uggrid/ugwrapper.hh"
87#undef UG_DIM
88
89// undef all macros defined by UG
90#include "uggrid/ug_undefs.hh"
91
92#undef UG_DIM_3
93
94// The components of the UGGrid interface
95#include "uggrid/uggridgeometry.hh"
96#include "uggrid/uggridlocalgeometry.hh"
97#include "uggrid/uggridentity.hh"
98#include "uggrid/uggridentityseed.hh"
99#include "uggrid/uggridintersections.hh"
100#include "uggrid/uggridintersectioniterators.hh"
101#include "uggrid/uggridleveliterator.hh"
102#include "uggrid/uggridleafiterator.hh"
103#include "uggrid/uggridhieriterator.hh"
104#include "uggrid/uggridindexsets.hh"
105#include <dune/grid/uggrid/uggridviews.hh>
106#ifdef ModelP
107#include "uggrid/ugmessagebuffer.hh"
108#include "uggrid/uglbgatherscatter.hh"
109#endif
110
111// Not needed here, but included for user convenience
113
114#ifdef ModelP
115template <class DataHandle, int GridDim, int codim>
116const Dune::UGGrid<GridDim>* Dune::UGMessageBuffer<DataHandle, GridDim, codim>::grid_;
117
118template <class DataHandle, int GridDim, int codim>
119DataHandle *Dune::UGMessageBuffer<DataHandle,GridDim,codim>::duneDataHandle_ = nullptr;
120
121template <class DataHandle, int GridDim, int codim>
122int Dune::UGMessageBuffer<DataHandle,GridDim,codim>::level = -1;
123#endif // ModelP
124
125namespace Dune {
126
127#ifdef ModelP
128 using UGCommunication = Communication<MPI_Comm>;
129#else
130 using UGCommunication = Communication<No_Comm>;
131#endif
132
133 template<int dim>
135 {
137 UGGridGeometry,
138 UGGridEntity,
139 UGGridLevelIterator,
140 UGGridLeafIntersection,
141 UGGridLevelIntersection,
142 UGGridLeafIntersectionIterator,
143 UGGridLevelIntersectionIterator,
144 UGGridHierarchicIterator,
145 UGGridLeafIterator,
146 UGGridLevelIndexSet< const UGGrid<dim> >,
147 UGGridLeafIndexSet< const UGGrid<dim> >,
148 UGGridIdSet< const UGGrid<dim> >,
149 typename UG_NS<dim>::UG_ID_TYPE,
150 UGGridIdSet< const UGGrid<dim> >,
151 typename UG_NS<dim>::UG_ID_TYPE,
153 UGGridLevelGridViewTraits,
154 UGGridLeafGridViewTraits,
155 UGGridEntitySeed,
156 UGGridLocalGeometry>
158 };
159
160
161 //**********************************************************************
162 //
163 // --UGGrid
164 //
165 //**********************************************************************
166
198 template <int dim>
199 class UGGrid : public GridDefaultImplementation <dim, dim, double, UGGridFamily<dim> >
200 {
201 typedef GridDefaultImplementation <dim, dim, double, UGGridFamily<dim> > Base;
202
203 friend class UGGridGeometry<0,dim,const UGGrid<dim> >;
204 friend class UGGridGeometry<dim,dim,const UGGrid<dim> >;
205 friend class UGGridGeometry<1,2,const UGGrid<dim> >;
206 friend class UGGridGeometry<2,3,const UGGrid<dim> >;
207
208 friend class UGGridEntity <0,dim,const UGGrid<dim> >;
209 friend class UGGridEntity <1,dim,const UGGrid<dim> >;
210 friend class UGGridEntity <2,dim,const UGGrid<dim> >;
211 friend class UGGridEntity <dim,dim,const UGGrid<dim> >;
212 friend class UGGridHierarchicIterator<const UGGrid<dim> >;
213 friend class UGGridLeafIntersection<const UGGrid<dim> >;
214 friend class UGGridLevelIntersection<const UGGrid<dim> >;
215 friend class UGGridLeafIntersectionIterator<const UGGrid<dim> >;
216 friend class UGGridLevelIntersectionIterator<const UGGrid<dim> >;
217
218 friend class UGGridLevelIndexSet<const UGGrid<dim> >;
219 friend class UGGridLeafIndexSet<const UGGrid<dim> >;
220 friend class UGGridIdSet<const UGGrid<dim> >;
221 template <class GridImp_>
222 friend class UGGridLeafGridView;
223 template <class GridImp_>
225
226 friend class GridFactory<UGGrid<dim> >;
227
228#ifdef ModelP
229 friend class UGLBGatherScatter;
230#endif
231
232 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
233 friend class UGGridLeafIterator;
234 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
236
238 static_assert(dim==2 || dim==3, "Use UGGrid only for 2d and 3d!");
239
240 // The different instantiations are mutual friends so they can access
241 // each others numOfUGGrids field
242 friend class UGGrid<2>;
243 friend class UGGrid<3>;
244
245 //**********************************************************
246 // The Interface Methods
247 //**********************************************************
248 public:
251
252 // the Traits
254
256 typedef UG::DOUBLE ctype;
257
259 typedef unsigned int Rank;
260
264
266 ~UGGrid() noexcept(false);
267
270 int maxLevel() const;
271
273 template <typename Seed>
274 typename Traits::template Codim<Seed::codimension>::Entity
275 entity(const Seed& seed) const
276 {
277 const int codim = Seed::codimension;
278 return typename Traits::template Codim<codim>::Entity(UGGridEntity<codim,dim,const UGGrid<dim> >(seed.impl().target(),this));
279 }
280
283 int size (int level, int codim) const;
284
286 int size (int codim) const
287 {
288 return leafIndexSet().size(codim);
289 }
290
292 int size (int level, GeometryType type) const
293 {
294 return this->levelIndexSet(level).size(type);
295 }
296
298 int size (GeometryType type) const
299 {
300 return this->leafIndexSet().size(type);
301 }
302
304 size_t numBoundarySegments() const {
305 // The number is stored as a member of UGGrid upon grid creation.
306 // The corresponding data structure is not exported by UG. (It is in ug/dom/std/std_internal.h)
307 return numBoundarySegments_;
308 }
309
311 const typename Traits::GlobalIdSet& globalIdSet() const
312 {
313 return idSet_;
314 }
315
317 const typename Traits::LocalIdSet& localIdSet() const
318 {
319 return idSet_;
320 }
321
323 const typename Traits::LevelIndexSet& levelIndexSet(int level) const
324 {
325 if (level<0 || level>maxLevel())
326 DUNE_THROW(GridError, "levelIndexSet of nonexisting level " << level << " requested!");
327 return *levelIndexSets_[level];
328 }
329
331 const typename Traits::LeafIndexSet& leafIndexSet() const
332 {
333 return leafIndexSet_;
334 }
335
338
351 bool mark(int refCount, const typename Traits::template Codim<0>::Entity & e );
352
408 bool mark(const typename Traits::template Codim<0>::Entity & e,
409 typename UG_NS<dim>::RefinementRule rule,
410 int side=0);
411
413 int getMark(const typename Traits::template Codim<0>::Entity& e) const;
414
417 bool preAdapt();
418
420 bool adapt();
421
423 void postAdapt();
425
433 template<class DataHandle>
434 bool loadBalance (DataHandle& dataHandle)
435 {
436#ifdef ModelP
437 // gather element data
438 if (dataHandle.contains(dim, 0))
439 UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
440
441 // gather node data
442 if (dataHandle.contains(dim,dim))
443 UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
444#endif
445
446 // the load balancing step now also attaches
447 // the data to the entities and distributes it
448 loadBalance();
449
450#ifdef ModelP
451 // scatter element data
452 if (dataHandle.contains(dim, 0))
453 UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
454
455 // scatter node data
456 if (dataHandle.contains(dim,dim))
457 UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
458#endif
459
460 return true;
461 }
462
469 bool loadBalance(int minlevel=0);
470
501 bool loadBalance(const std::vector<Rank>& targetProcessors, unsigned int fromLevel);
502
512 template<class DataHandle>
513 bool loadBalance (const std::vector<Rank>& targetProcessors, unsigned int fromLevel, DataHandle& dataHandle)
514 {
515#ifdef ModelP
516 // gather element data
517 if (dataHandle.contains(dim, 0))
518 UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
519
520 // gather node data
521 if (dataHandle.contains(dim,dim))
522 UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
523#endif
524
525 // the load balancing step now also attaches
526 // the data to the entities and distributes it
527 loadBalance(targetProcessors,fromLevel);
528
529#ifdef ModelP
530 // scatter element data
531 if (dataHandle.contains(dim, 0))
532 UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
533
534 // scatter node data
535 if (dataHandle.contains(dim,dim))
536 UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
537#endif
538
539 return true;
540 }
541
543 const UGCommunication& comm () const
544 {
545 return ccobj_;
546 }
547
548 protected:
549#ifdef ModelP
550 template <int codim, class GridView, class DataHandle>
551 void communicateUG_(const GridView& gv, int level,
552 DataHandle &dataHandle,
553 InterfaceType iftype,
554 CommunicationDirection dir) const
555 {
556 // Translate the communication direction from Dune-Speak to UG-Speak
557 const auto ugIfDir = (dir==ForwardCommunication) ? UG_NS<dim>::IF_FORWARD() : UG_NS<dim>::IF_BACKWARD();
558
559 // Set up the message buffer
560 // No actual object is constructed, because the DDD_IFOneway method called below
561 // needs *static* methods. Therefore, everything is currently routed
562 // via static data members of the UGMessageBuffer class.
563 using UGMsgBuf = UGMessageBuffer<DataHandle,dim,codim>;
564 UGMsgBuf::duneDataHandle_ = &dataHandle;
565 UGMsgBuf::level = level;
566 UGMsgBuf::grid_ = this;
567
568 const std::vector<typename UG_NS<dim>::DDD_IF> ugIfs = findDDDInterfaces(iftype, codim);
569
570 const auto bufSize = UGMsgBuf::ugBufferSize(gv);
571 if (!bufSize)
572 return; // we don't need to communicate if we don't have any data!
573 for (auto&& dddInterface : ugIfs)
574 UG_NS<dim>::DDD_IFOneway(multigrid_->dddContext(),
575 dddInterface,
576 ugIfDir,
577 bufSize,
578 &UGMsgBuf::ugGather_,
579 &UGMsgBuf::ugScatter_);
580 }
581
583 std::vector<typename UG_NS<dim>::DDD_IF> findDDDInterfaces(InterfaceType iftype,
584 int codim) const;
585#endif
586 public:
587 // **********************************************************
588 // End of Interface Methods
589 // **********************************************************
590
598 void getChildrenOfSubface(const typename Traits::template Codim<0>::Entity & e,
599 int elementSide,
600 int maxl,
601 std::vector<typename Traits::template Codim<0>::Entity>& childElements,
602 std::vector<unsigned char>& childElementSides) const;
603
611
619
622 refinementType_ = type;
623 }
624
627 closureType_ = type;
628 }
629
633 void setPosition(const typename Traits::template Codim<dim>::Entity& e,
634 const FieldVector<double, dim>& pos);
635
640 void globalRefine(int n);
641
646 void saveState(const std::string& filename) const;
647
652 void loadState(const std::string& filename);
653
654 private:
656 typename UG_NS<dim>::MultiGrid* multigrid_;
657
663 typename UG_NS<dim>::STD_BVP* bvp_;
664
666 UGCommunication ccobj_;
667
673 void setIndices(bool setLevelZero,
674 std::vector<unsigned int>* nodePermutation);
675
676 // Each UGGrid object has a unique name to identify it in the
677 // UG environment structure
678 std::string name_;
679
680 // Our set of level indices
681 std::vector<std::shared_ptr<UGGridLevelIndexSet<const UGGrid<dim> > > > levelIndexSets_;
682
683 UGGridLeafIndexSet<const UGGrid<dim> > leafIndexSet_;
684
685 // One id set implementation
686 // Used for both the local and the global UGGrid id sets
687 UGGridIdSet<const UGGrid<dim> > idSet_;
688
690 RefinementType refinementType_;
691
693 ClosureType closureType_;
694
702 static int numOfUGGrids;
703
709 bool someElementHasBeenMarkedForRefinement_;
710
716 bool someElementHasBeenMarkedForCoarsening_;
717
719 std::vector<std::shared_ptr<BoundarySegment<dim> > > boundarySegments_;
720
726 unsigned int numBoundarySegments_;
727
728 }; // end Class UGGrid
729
730 namespace Capabilities
731 {
735
739
743
747 template<int dim, int codim>
748 struct hasEntity< UGGrid<dim>, codim>
749 {
750 static const bool v = true;
751 };
752
758 template<int dim, int codim>
759 struct hasEntityIterator<UGGrid<dim>, codim>
760 {
761 static const bool v = false;
762 };
763
768 template<int dim>
769 struct hasEntityIterator<UGGrid<dim>, 0>
770 {
771 static const bool v = true;
772 };
773
778 template<int dim>
779 struct hasEntityIterator<UGGrid<dim>, dim>
780 {
781 static const bool v = true;
782 };
783
787 template<int dim, int codim>
788 struct canCommunicate<UGGrid<dim>, codim>
789 {
790 static const bool v = (codim>=0 && codim<=dim);
791 };
792
796 template<int dim>
798 {
799 static const bool v = true;
800 };
801
805 template<int dim>
807 {
808 static const bool v = false;
809 };
810
814 template<int dim>
815 struct viewThreadSafe< UGGrid<dim> > {
816 static const bool v = true;
817 };
818
819 }
820
821} // namespace Dune
822
823#endif // HAVE_DUNE_UGGRID || DOXYGEN
824#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:130
Contains all capabilities classes.
Definition albertagrid/capabilities.hh:29
Specialize with 'true' for all codims that a grid implements entities for. (default=false).
Definition common/capabilities.hh:58
specialize with 'true' for all codims that a grid provides an iterator for (default=hasEntity<codim>:...
Definition common/capabilities.hh:74
specialize with 'true' for all codims that a grid can communicate data on (default=false)
Definition common/capabilities.hh:97
Specialize with 'true' if implementation guarantees conforming level grids. (default=false).
Definition common/capabilities.hh:106
Specialize with 'true' if implementation guarantees a conforming leaf grid. (default=false).
Definition common/capabilities.hh:115
Specialize with 'true' if the grid implementation is thread safe, while it is not modified....
Definition common/capabilities.hh:169
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:848
Traits::LeafGridView leafGridView() const
Definition common/grid.hh:868
auto size(GeometryType type) const
Return total number of entities of given geometry type in entity set .
Definition common/indexidset.hh:223
A Traits struct that collects all associated types of one implementation.
Definition common/grid.hh:411
A traits struct that collects all associated types of one grid model.
Definition common/grid.hh:1013
IndexSet< const GridImp, LeafIndexSetImp, LeafIndexType, LeafGeometryTypes > LeafIndexSet
The type of the leaf index set.
Definition common/grid.hh:1082
IndexSet< const GridImp, LevelIndexSetImp, LevelIndexType, LevelGeometryTypes > LevelIndexSet
The type of the level index set.
Definition common/grid.hh:1080
IdSet< const GridImp, LocalIdSetImp, LIDType > LocalIdSet
The type of the local id set.
Definition common/grid.hh:1086
IdSet< const GridImp, GlobalIdSetImp, GIDType > GlobalIdSet
The type of the global id set.
Definition common/grid.hh:1084
Provide a generic factory class for unstructured grids.
Definition common/gridfactory.hh:275
Grid view abstract base class.
Definition common/gridview.hh:66
Definition uggrid.hh:135
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:157
Front-end for the grid manager of the finite element toolbox UG3.
Definition uggrid.hh:200
void postAdapt()
Clean up refinement markers.
size_t numBoundarySegments() const
Return the number of boundary segments.
Definition uggrid.hh:304
void setRefinementType(RefinementType type)
Sets the type of grid refinement.
Definition uggrid.hh:621
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:233
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:253
friend class UGGridLevelIterator
Definition uggrid.hh:235
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:513
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:224
UGGridFamily< dim > GridFamily
type of the used GridFamily for this grid
Definition uggrid.hh:250
bool loadBalance(const std::vector< Rank > &targetProcessors, unsigned int fromLevel)
Distribute this grid over a distributed machine.
const UGCommunication & comm() const
Definition uggrid.hh:543
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:323
void setClosureType(ClosureType type)
Sets the type of grid refinement closure.
Definition uggrid.hh:626
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition uggrid.hh:317
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition uggrid.hh:292
RefinementType
The different forms of grid refinement that UG supports.
Definition uggrid.hh:605
@ COPY
Definition uggrid.hh:609
@ LOCAL
Definition uggrid.hh:607
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition uggrid.hh:298
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:256
bool adapt()
Triggers the grid refinement process.
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition uggrid.hh:311
friend class UGGridLeafGridView
Definition uggrid.hh:222
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:434
Traits::template Codim< Seed::codimension >::Entity entity(const Seed &seed) const
Definition uggrid.hh:275
unsigned int Rank
The type used for process ranks.
Definition uggrid.hh:259
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:613
@ GREEN
Definition uggrid.hh:615
@ NONE
Definition uggrid.hh:617
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition uggrid.hh:331
int size(int codim) const
number of leaf entities per codim in this process
Definition uggrid.hh:286
static const bool v
Definition uggrid.hh:750
static const bool v
Definition uggrid.hh:761
static const bool v
Definition uggrid.hh:771
static const bool v
Definition uggrid.hh:781
static const bool v
Definition uggrid.hh:790
static const bool v
Definition uggrid.hh:799
static const bool v
Definition uggrid.hh:808
static const bool v
Definition uggrid.hh:816
A set of traits classes to store static information about grid implementation.
Different resources needed by all grid implementations.