DUNE PDELab (2.7)

uggrid.hh
Go to the documentation of this file.
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4#ifndef DUNE_UGGRID_HH
5#define DUNE_UGGRID_HH
6
11#include <memory>
12
18
22
23#if HAVE_UG || DOXYGEN
24
25#ifdef ModelP
27#endif
28
29/* [Before reading the following: the macros UG_DIM_2 and UG_DIM_3 where named
30 * _2 and _3, respectively, up until ug-3.12.0.]
31 *
32 * The following lines including the necessary UG headers are somewhat
33 tricky. Here's what's happening:
34 UG can support two- and three-dimensional grids. You choose by setting
35 either UG_DIM_2 or UG_DIM_3 while compiling. This changes all sorts of stuff, in
36 particular data structures in the headers.
37 UG was never supposed to provide 2d and 3d grids at the same time.
38 However, when compiling it as c++, the dimension-dependent parts are
39 wrapped up cleanly in the namespaces UG::D2 and UG::D3, respectively. That
40 way it is possible to link together the UG lib for 2d and the one for 3d.
41 But we also need the headers twice! Once with UG_DIM_2 set and once with UG_DIM_3!
42 So here we go:*/
43
44/* The following define tells the UG headers that we want access to a few
45 special fields, for example the extra index fields in the element data structures.
46 This define remains only for backwards compatibility with older version of UG.
47 All dune-uggrid versions since 2016-08-05 do not need this #define or the #undef
48 further below. */
49#define FOR_DUNE
50
51// Set UG's space-dimension flag to 2d
52#ifdef UG_USE_NEW_DIMENSION_DEFINES
53#define UG_DIM_2
54#else
55#define _2
56#endif
57// And include all necessary UG headers
58#include "uggrid/ugincludes.hh"
59
60// Wrap a few large UG macros by functions before they get undef'ed away.
61// Here: The 2d-version of the macros
62#define UG_DIM 2
63#include "uggrid/ugwrapper.hh"
64#undef UG_DIM
65
66// UG defines a whole load of preprocessor macros. ug_undefs.hh undefines
67// them all, so we don't get name clashes.
68#include "uggrid/ug_undefs.hh"
69#ifdef UG_USE_NEW_DIMENSION_DEFINES
70#undef UG_DIM_2
71#else
72#undef _2
73#endif
74
75/* Now we're done with 2d, and we can do the whole thing over again for 3d */
76
77/* All macros set by UG have been unset. This includes the macros that ensure
78 single inclusion of headers. We can thus include them again. However, we
79 only want to include those headers again that contain dimension-dependent stuff.
80 Therefore, we set a few single-inclusion defines manually before including
81 ugincludes.hh again.
82 */
83#define UGTYPES_H
84#define __HEAPS__
85#define __UGENV__
86#define __DEVICESH__
87#ifdef ModelP
88#define __PPIF__
89#endif
90
91#ifdef UG_USE_NEW_DIMENSION_DEFINES
92#define UG_DIM_3
93#else
94#define _3
95#endif
96
97#include "uggrid/ugincludes.hh"
98
99// Wrap a few large UG macros by functions before they get undef'ed away.
100// This time it's the 3d-versions.
101#define UG_DIM 3
102#include "uggrid/ugwrapper.hh"
103#undef UG_DIM
104
105// undef all macros defined by UG
106#include "uggrid/ug_undefs.hh"
107
108#ifdef UG_USE_NEW_DIMENSION_DEFINES
109#undef UG_DIM_3
110#else
111#undef _3
112#endif
113#undef FOR_DUNE
114
115// The components of the UGGrid interface
116#include "uggrid/uggridgeometry.hh"
117#include "uggrid/uggridlocalgeometry.hh"
118#include "uggrid/uggridentity.hh"
119#include "uggrid/uggridentityseed.hh"
120#include "uggrid/uggridintersections.hh"
121#include "uggrid/uggridintersectioniterators.hh"
122#include "uggrid/uggridleveliterator.hh"
123#include "uggrid/uggridleafiterator.hh"
124#include "uggrid/uggridhieriterator.hh"
125#include "uggrid/uggridindexsets.hh"
126#include <dune/grid/uggrid/uggridviews.hh>
127#ifdef ModelP
128#include "uggrid/ugmessagebuffer.hh"
129#include "uggrid/uglbgatherscatter.hh"
130#endif
131
132// Not needed here, but included for user convenience
134
135#ifdef ModelP
136template <class DataHandle, int GridDim, int codim>
137const Dune::UGGrid<GridDim>* Dune::UGMessageBuffer<DataHandle, GridDim, codim>::grid_;
138
139template <class DataHandle, int GridDim, int codim>
140DataHandle *Dune::UGMessageBuffer<DataHandle,GridDim,codim>::duneDataHandle_ = nullptr;
141
142template <class DataHandle, int GridDim, int codim>
143int Dune::UGMessageBuffer<DataHandle,GridDim,codim>::level = -1;
144#endif // ModelP
145
146namespace Dune {
147
148#ifdef ModelP
149 using UGCollectiveCommunication = CollectiveCommunication<MPI_Comm>;
150#else
151 using UGCollectiveCommunication = CollectiveCommunication<No_Comm>;
152#endif
153
154 template<int dim>
155 struct UGGridFamily
156 {
157 typedef GridTraits<dim,dim,Dune::UGGrid<dim>,
158 UGGridGeometry,
159 UGGridEntity,
160 UGGridLevelIterator,
161 UGGridLeafIntersection,
162 UGGridLevelIntersection,
163 UGGridLeafIntersectionIterator,
164 UGGridLevelIntersectionIterator,
165 UGGridHierarchicIterator,
166 UGGridLeafIterator,
167 UGGridLevelIndexSet< const UGGrid<dim> >,
168 UGGridLeafIndexSet< const UGGrid<dim> >,
169 UGGridIdSet< const UGGrid<dim> >,
170 typename UG_NS<dim>::UG_ID_TYPE,
171 UGGridIdSet< const UGGrid<dim> >,
172 typename UG_NS<dim>::UG_ID_TYPE,
173 UGCollectiveCommunication,
174 UGGridLevelGridViewTraits,
175 UGGridLeafGridViewTraits,
176 UGGridEntitySeed,
177 UGGridLocalGeometry>
178 Traits;
179 };
180
181
182 //**********************************************************************
183 //
184 // --UGGrid
185 //
186 //**********************************************************************
187
225 template <int dim>
226 class UGGrid : public GridDefaultImplementation <dim, dim, double, UGGridFamily<dim> >
227 {
229
230 friend class UGGridGeometry<0,dim,const UGGrid<dim> >;
231 friend class UGGridGeometry<dim,dim,const UGGrid<dim> >;
232 friend class UGGridGeometry<1,2,const UGGrid<dim> >;
233 friend class UGGridGeometry<2,3,const UGGrid<dim> >;
234
235 friend class UGGridEntity <0,dim,const UGGrid<dim> >;
236 friend class UGGridEntity <1,dim,const UGGrid<dim> >;
237 friend class UGGridEntity <dim,dim,const UGGrid<dim> >;
238 friend class UGEdgeEntity <dim,const UGGrid<dim> >;
239 friend class UGGridHierarchicIterator<const UGGrid<dim> >;
240 friend class UGGridLeafIntersection<const UGGrid<dim> >;
241 friend class UGGridLevelIntersection<const UGGrid<dim> >;
242 friend class UGGridLeafIntersectionIterator<const UGGrid<dim> >;
243 friend class UGGridLevelIntersectionIterator<const UGGrid<dim> >;
244
245 friend class UGGridLevelIndexSet<const UGGrid<dim> >;
246 friend class UGGridLeafIndexSet<const UGGrid<dim> >;
247 friend class UGGridIdSet<const UGGrid<dim> >;
248 template <class GridImp_>
249 friend class UGGridLeafGridView;
250 template <class GridImp_>
251 friend class UGGridLevelGridView;
252
253 friend class GridFactory<UGGrid<dim> >;
254
255#ifdef ModelP
256 friend class UGLBGatherScatter;
257#endif
258
259 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
260 friend class UGGridLeafIterator;
261 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
262 friend class UGGridLevelIterator;
263
265 static_assert(dim==2 || dim==3, "Use UGGrid only for 2d and 3d!");
266
267 // The different instantiations are mutual friends so they can access
268 // each others numOfUGGrids field
269 friend class UGGrid<2>;
270 friend class UGGrid<3>;
271
272 //**********************************************************
273 // The Interface Methods
274 //**********************************************************
275 public:
277 typedef UGGridFamily<dim> GridFamily;
278
279 // the Traits
280 typedef typename UGGridFamily<dim>::Traits Traits;
281
283 typedef UG::DOUBLE ctype;
284
286 typedef unsigned int Rank;
287
291
293 ~UGGrid() noexcept(false);
294
297 int maxLevel() const;
298
300 template <typename Seed>
301 typename Traits::template Codim<Seed::codimension>::Entity
302 entity(const Seed& seed) const
303 {
304 const int codim = Seed::codimension;
305 return typename Traits::template Codim<codim>::Entity(UGGridEntity<codim,dim,const UGGrid<dim> >(seed.impl().target(),this));
306 }
307
310 int size (int level, int codim) const;
311
313 int size (int codim) const
314 {
315 return leafIndexSet().size(codim);
316 }
317
319 int size (int level, GeometryType type) const
320 {
321 return this->levelIndexSet(level).size(type);
322 }
323
325 int size (GeometryType type) const
326 {
327 return this->leafIndexSet().size(type);
328 }
329
331 size_t numBoundarySegments() const {
332 // The number is stored as a member of UGGrid upon grid creation.
333 // The corresponding data structure is not exported by UG. (It is in ug/dom/std/std_internal.h)
334 return numBoundarySegments_;
335 }
336
338 const typename Traits::GlobalIdSet& globalIdSet() const
339 {
340 return idSet_;
341 }
342
344 const typename Traits::LocalIdSet& localIdSet() const
345 {
346 return idSet_;
347 }
348
350 const typename Traits::LevelIndexSet& levelIndexSet(int level) const
351 {
352 if (level<0 || level>maxLevel())
353 DUNE_THROW(GridError, "levelIndexSet of nonexisting level " << level << " requested!");
354 return *levelIndexSets_[level];
355 }
356
358 const typename Traits::LeafIndexSet& leafIndexSet() const
359 {
360 return leafIndexSet_;
361 }
362
365
378 bool mark(int refCount, const typename Traits::template Codim<0>::Entity & e );
379
435 bool mark(const typename Traits::template Codim<0>::Entity & e,
436 typename UG_NS<dim>::RefinementRule rule,
437 int side=0);
438
440 int getMark(const typename Traits::template Codim<0>::Entity& e) const;
441
444 bool preAdapt();
445
447 bool adapt();
448
450 void postAdapt();
460 template<class DataHandle>
461 bool loadBalance (DataHandle& dataHandle)
462 {
463#ifdef ModelP
464 // gather element data
465 if (dataHandle.contains(dim, 0))
466 UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
467
468 // gather node data
469 if (dataHandle.contains(dim,dim))
470 UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
471#endif
472
473 // the load balancing step now also attaches
474 // the data to the entities and distributes it
475 loadBalance();
476
477#ifdef ModelP
478 // scatter element data
479 if (dataHandle.contains(dim, 0))
480 UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
481
482 // scatter node data
483 if (dataHandle.contains(dim,dim))
484 UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
485#endif
486
487 return true;
488 }
489
496 bool loadBalance(int minlevel=0);
497
528 bool loadBalance(const std::vector<Rank>& targetProcessors, unsigned int fromLevel);
529
539 template<class DataHandle>
540 bool loadBalance (const std::vector<Rank>& targetProcessors, unsigned int fromLevel, DataHandle& dataHandle)
541 {
542#ifdef ModelP
543 // gather element data
544 if (dataHandle.contains(dim, 0))
545 UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
546
547 // gather node data
548 if (dataHandle.contains(dim,dim))
549 UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
550#endif
551
552 // the load balancing step now also attaches
553 // the data to the entities and distributes it
554 loadBalance(targetProcessors,fromLevel);
555
556#ifdef ModelP
557 // scatter element data
558 if (dataHandle.contains(dim, 0))
559 UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
560
561 // scatter node data
562 if (dataHandle.contains(dim,dim))
563 UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
564#endif
565
566 return true;
567 }
568
571 {
572 return ccobj_;
573 }
574
575 protected:
576#ifdef ModelP
577 template <int codim, class GridView, class DataHandle>
578 void communicateUG_(const GridView& gv, int level,
579 DataHandle &dataHandle,
580 InterfaceType iftype,
581 CommunicationDirection dir) const
582 {
583 typename UG_NS<dim>::DDD_IF_DIR ugIfDir;
584 // Translate the communication direction from Dune-Speak to UG-Speak
585 if (dir==ForwardCommunication)
586 ugIfDir = UG_NS<dim>::IF_FORWARD();
587 else
588 ugIfDir = UG_NS<dim>::IF_BACKWARD();
589
590 typedef UGMessageBuffer<DataHandle,dim,codim> UGMsgBuf;
591 UGMsgBuf::duneDataHandle_ = &dataHandle;
592
593 UGMsgBuf::level = level;
594
595 std::vector<typename UG_NS<dim>::DDD_IF> ugIfs;
596 findDDDInterfaces_(ugIfs, iftype, codim);
597
598 unsigned bufSize = UGMsgBuf::ugBufferSize(gv);
599 if (!bufSize)
600 return; // we don't need to communicate if we don't have any data!
601 UGMsgBuf::grid_ = this;
602 for (unsigned i=0; i < ugIfs.size(); ++i)
603 UG_NS<dim>::DDD_IFOneway(
604#if DUNE_UGGRID_HAVE_DDDCONTEXT
605 multigrid_->dddContext(),
606#endif
607 ugIfs[i],
608 ugIfDir,
609 bufSize,
610 &UGMsgBuf::ugGather_,
611 &UGMsgBuf::ugScatter_);
612 }
613
614 void findDDDInterfaces_(std::vector<typename UG_NS<dim>::DDD_IF > &dddIfaces,
615 InterfaceType iftype,
616 int codim) const
617 {
618#if DUNE_UGGRID_HAVE_DDDCONTEXT
619# define DDD_CONTEXT multigrid_->dddContext()
620#else
621# define DDD_CONTEXT
622#endif
623 dddIfaces.clear();
624 if (codim == 0)
625 {
626 switch (iftype) {
628 // do not communicate anything: Elements can not be in
629 // the interior of two processes at the same time
630 return;
632 // The communicated elements are in the sender's
633 // interior and it does not matter what they are for
634 // the receiver
635 dddIfaces.push_back(UG_NS<dim>::ElementVHIF(DDD_CONTEXT));
636 return;
637 case All_All_Interface :
638 // It does neither matter what the communicated
639 // elements are for sender nor for the receiver. If
640 // they are seen by these two processes, data is send
641 // and received.
642 dddIfaces.push_back(UG_NS<dim>::ElementSymmVHIF(DDD_CONTEXT));
643 return;
644 default :
645 DUNE_THROW(GridError,
646 "Element communication not supported for "
647 "interfaces of type "
648 << iftype);
649 }
650 }
651 else if (codim == dim)
652 {
653 switch (iftype)
654 {
656 dddIfaces.push_back(UG_NS<dim>::BorderNodeSymmIF(DDD_CONTEXT));
657 return;
659 dddIfaces.push_back(UG_NS<dim>::NodeInteriorBorderAllIF(DDD_CONTEXT));
660 return;
661 case All_All_Interface :
662 dddIfaces.push_back(UG_NS<dim>::NodeAllIF(DDD_CONTEXT));
663 return;
664 default :
665 DUNE_THROW(GridError,
666 "Node communication not supported for "
667 "interfaces of type "
668 << iftype);
669 }
670 }
671 else if (codim == dim-1)
672 {
673 switch (iftype)
674 {
676 dddIfaces.push_back(UG_NS<dim>::BorderEdgeSymmIF(DDD_CONTEXT));
677 return;
679 dddIfaces.push_back(UG_NS<dim>::EdgeVHIF(DDD_CONTEXT));
680 return;
681 case All_All_Interface :
682 dddIfaces.push_back(UG_NS<dim>::EdgeSymmVHIF(DDD_CONTEXT));
683 return;
684 default :
685 DUNE_THROW(GridError,
686 "Edge communication not supported for "
687 "interfaces of type "
688 << iftype);
689 }
690 }
691 else if (codim == 1)
692 {
693 switch (iftype)
694 {
697 dddIfaces.push_back(UG_NS<dim>::FacetInteriorBorderAllIF(DDD_CONTEXT));
698 return;
699 default :
700 DUNE_THROW(GridError,
701 "Face communication not supported for "
702 "interfaces of type "
703 << iftype);
704 }
705 }
706 else
707 {
708 DUNE_THROW(GridError,
709 "Communication for codim "
710 << codim
711 << " entities is not yet supported "
712 << " by the DUNE UGGrid interface!");
713 }
714#undef DDD_CONTEXT
715 };
716#endif // ModelP
717
718 public:
719 // **********************************************************
720 // End of Interface Methods
721 // **********************************************************
722
730 void getChildrenOfSubface(const typename Traits::template Codim<0>::Entity & e,
731 int elementSide,
732 int maxl,
733 std::vector<typename Traits::template Codim<0>::Entity>& childElements,
734 std::vector<unsigned char>& childElementSides) const;
735
741 COPY
742 };
743
749 NONE
750 };
751
754 refinementType_ = type;
755 }
756
759 closureType_ = type;
760 }
761
768 static void setDefaultHeapSize(unsigned size)
769 DUNE_DEPRECATED_MSG("Do not set the UGGrid default heap size---it is ignored anyway!")
770 {}
771
775 void setPosition(const typename Traits::template Codim<dim>::Entity& e,
776 const FieldVector<double, dim>& pos);
777
782 void globalRefine(int n);
783
788 void saveState(const std::string& filename) const;
789
794 void loadState(const std::string& filename);
795
796 private:
798 typename UG_NS<dim>::MultiGrid* multigrid_;
799
802
808 void setIndices(bool setLevelZero,
809 std::vector<unsigned int>* nodePermutation);
810
811 // Each UGGrid object has a unique name to identify it in the
812 // UG environment structure
813 std::string name_;
814
815 // Our set of level indices
816 std::vector<std::shared_ptr<UGGridLevelIndexSet<const UGGrid<dim> > > > levelIndexSets_;
817
818 UGGridLeafIndexSet<const UGGrid<dim> > leafIndexSet_;
819
820 // One id set implementation
821 // Used for both the local and the global UGGrid id sets
822 UGGridIdSet<const UGGrid<dim> > idSet_;
823
825 RefinementType refinementType_;
826
828 ClosureType closureType_;
829
837 static int numOfUGGrids;
838
844 bool someElementHasBeenMarkedForRefinement_;
845
851 bool someElementHasBeenMarkedForCoarsening_;
852
854 std::vector<std::shared_ptr<BoundarySegment<dim> > > boundarySegments_;
855
861 unsigned int numBoundarySegments_;
862
863 }; // end Class UGGrid
864
865 namespace Capabilities
866 {
882 template<int dim>
883 struct hasEntity< UGGrid<dim>, 0>
884 {
885 static const bool v = true;
886 };
887
891 template<int dim>
892 struct hasEntity< UGGrid<dim>, dim>
893 {
894 static const bool v = true;
895 };
896
901 template<int dim>
902 struct hasEntityIterator<UGGrid<dim>, 0>
903 {
904 static const bool v = true;
905 };
906
911 template<int dim>
912 struct hasEntityIterator<UGGrid<dim>, dim>
913 {
914 static const bool v = true;
915 };
916
920 template<int dim>
922 {
923 static const bool v = true;
924 };
925
929 template<int dim>
931 {
932 static const bool v = false;
933 };
934
935 }
936
937} // namespace Dune
938
939#endif // HAVE_UG || DOXYGEN
940#endif // DUNE_UGGRID_HH
Base class for grid boundary segments of arbitrary geometry.
Wrapper class for entities.
Definition: entity.hh:64
vector space out of a tensor product of fields.
Definition: fvector.hh:96
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:280
Definition: grid.hh:855
Traits::LeafGridView leafGridView() const
View for the leaf grid for All_Partition.
Definition: grid.hh:875
bool loadBalance()
default implementation of load balance does nothing and returns false
Definition: grid.hh:941
Base class for exceptions in Dune grid modules.
Definition: exceptions.hh:18
Provide a generic factory class for unstructured grids.
Definition: gridfactory.hh:269
Grid view abstract base class.
Definition: gridview.hh:60
Id Set Interface.
Definition: indexidset.hh:441
Index Set Interface base class.
Definition: indexidset.hh:76
IndexType size(GeometryType type) const
Return total number of entities of given geometry type in entity set .
Definition: indexidset.hh:220
Front-end for the grid manager of the finite element toolbox UG.
Definition: uggrid.hh:227
void postAdapt()
Clean up refinement markers.
size_t numBoundarySegments() const
Return the number of boundary segments.
Definition: uggrid.hh:331
void setRefinementType(RefinementType type)
Sets the type of grid refinement.
Definition: uggrid.hh:753
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.
bool loadBalance(int minlevel=0)
Distributes this grid over the available nodes in a distributed machine.
~UGGrid() noexcept(false)
Destructor.
const UGCollectiveCommunication & comm() const
Definition: uggrid.hh:570
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:540
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.
UGGridFamily< dim > GridFamily
type of the used GridFamily for this grid
Definition: uggrid.hh:277
bool loadBalance(const std::vector< Rank > &targetProcessors, unsigned int fromLevel)
Distribute this grid over a distributed machine.
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.
UGGrid(UGCollectiveCommunication comm={})
Default constructor.
int maxLevel() const
const Traits::LevelIndexSet & levelIndexSet(int level) const
Access to the LevelIndexSets.
Definition: uggrid.hh:350
void setClosureType(ClosureType type)
Sets the type of grid refinement closure.
Definition: uggrid.hh:758
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition: uggrid.hh:344
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition: uggrid.hh:319
RefinementType
The different forms of grid refinement that UG supports.
Definition: uggrid.hh:737
@ COPY
New level consists of the refined elements and the unrefined ones, too.
Definition: uggrid.hh:741
@ LOCAL
New level consists only of the refined elements and the closure.
Definition: uggrid.hh:739
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: uggrid.hh:325
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:283
bool adapt()
Triggers the grid refinement process.
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition: uggrid.hh:338
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:461
Traits::template Codim< Seed::codimension >::Entity entity(const Seed &seed) const
Create an Entity from an EntitySeed.
Definition: uggrid.hh:302
unsigned int Rank
The type used for process ranks.
Definition: uggrid.hh:286
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.
ClosureType
Decide whether to add a green closure to locally refined grid sections or not.
Definition: uggrid.hh:745
@ GREEN
Standard red/green refinement.
Definition: uggrid.hh:747
@ NONE
No closure, results in nonconforming meshes.
Definition: uggrid.hh:749
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition: uggrid.hh:358
int size(int codim) const
number of leaf entities per codim in this process
Definition: uggrid.hh:313
static void setDefaultHeapSize(unsigned size)
Sets the default heap size.
Definition: uggrid.hh:768
A free function to provide the demangled class name of a given object or type as a string.
A set of traits classes to store static information about grid implementation.
Different resources needed by all grid implementations.
Implements an utility class that provides collective communication methods for sequential programs.
Definition of the DUNE_DEPRECATED macro for the case that config.h is not available.
A few common exception classes.
#define DUNE_DEPRECATED_MSG(text)
Mark some entity as deprecated.
Definition: deprecated.hh:169
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
CommunicationDirection
Define a type for communication direction parameter.
Definition: gridenums.hh:168
InterfaceType
Parameter to be used for the communication functions.
Definition: gridenums.hh:84
@ ForwardCommunication
communicate as given in InterfaceType
Definition: gridenums.hh:169
@ InteriorBorder_All_Interface
send interior and border, receive all entities
Definition: gridenums.hh:86
@ All_All_Interface
send all and receive all entities
Definition: gridenums.hh:89
@ InteriorBorder_InteriorBorder_Interface
send/receive interior and border entities
Definition: gridenums.hh:85
Implements an utility class that provides MPI's collective communication methods.
Helpers for dealing with MPI.
Dune namespace.
Definition: alignedallocator.hh:14
specialize with 'true' for all codims that a grid provides an iterator for (default=false)
Definition: capabilities.hh:72
Specialize with 'true' for all codims that a grid implements entities for. (default=false)
Definition: capabilities.hh:56
Specialize with 'true' if implementation guarantees a conforming leaf grid. (default=false)
Definition: capabilities.hh:113
Specialize with 'true' if implementation guarantees conforming level grids. (default=false)
Definition: capabilities.hh:104
Static tag representing a codimension.
Definition: dimension.hh:22
A traits struct that collects all associated types of one grid model.
Definition: grid.hh:1065
The specialization of the generic GridFactory for UGGrid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)