Dune Core Modules (2.8.0)

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
17
21
22#if HAVE_UG || DOXYGEN
23
24#ifdef ModelP
26#endif
27
28/* [Before reading the following: the macros UG_DIM_2 and UG_DIM_3 where named
29 * _2 and _3, respectively, up until ug-3.12.0.]
30 *
31 * The following lines including the necessary UG headers are somewhat
32 tricky. Here's what's happening:
33 UG can support two- and three-dimensional grids. You choose by setting
34 either UG_DIM_2 or UG_DIM_3 while compiling. This changes all sorts of stuff, in
35 particular data structures in the headers.
36 UG was never supposed to provide 2d and 3d grids at the same time.
37 However, when compiling it as c++, the dimension-dependent parts are
38 wrapped up cleanly in the namespaces UG::D2 and UG::D3, respectively. That
39 way it is possible to link together the UG lib for 2d and the one for 3d.
40 But we also need the headers twice! Once with UG_DIM_2 set and once with UG_DIM_3!
41 So here we go:*/
42
43/* The following define tells the UG headers that we want access to a few
44 special fields, for example the extra index fields in the element data structures.
45 This define remains only for backwards compatibility with older version of UG.
46 All dune-uggrid versions since 2016-08-05 do not need this #define or the #undef
47 further below. */
48#define FOR_DUNE
49
50// Set UG's space-dimension flag to 2d
51#define UG_DIM_2
52// And include all necessary UG headers
53#include "uggrid/ugincludes.hh"
54#undef DUNE_UGINCLUDES_HH
55
56// Wrap a few large UG macros by functions before they get undef'ed away.
57// Here: The 2d-version of the macros
58#define UG_DIM 2
59#include "uggrid/ugwrapper.hh"
60#undef UG_DIM
61
62// UG defines a whole load of preprocessor macros. ug_undefs.hh undefines
63// them all, so we don't get name clashes.
64#include "uggrid/ug_undefs.hh"
65#undef UG_DIM_2
66
67/* Now we're done with 2d, and we can do the whole thing over again for 3d */
68
69/* All macros set by UG have been unset. This includes the macros that ensure
70 single inclusion of headers. We can thus include them again. However, we
71 only want to include those headers again that contain dimension-dependent stuff.
72 Therefore, we set a few single-inclusion defines manually before including
73 ugincludes.hh again.
74 */
75#define UGTYPES_H
76#define __HEAPS__
77#define __UGENV__
78#define __DEVICESH__
79#ifdef ModelP
80#define __PPIF__
81#endif
82
83#define UG_DIM_3
84#include "uggrid/ugincludes.hh"
85#undef DUNE_UGINCLUDES_HH
86
87// Wrap a few large UG macros by functions before they get undef'ed away.
88// This time it's the 3d-versions.
89#define UG_DIM 3
90#include "uggrid/ugwrapper.hh"
91#undef UG_DIM
92
93// undef all macros defined by UG
94#include "uggrid/ug_undefs.hh"
95
96#undef UG_DIM_3
97#undef FOR_DUNE
98
99// The components of the UGGrid interface
100#include "uggrid/uggridgeometry.hh"
101#include "uggrid/uggridlocalgeometry.hh"
102#include "uggrid/uggridentity.hh"
103#include "uggrid/uggridentityseed.hh"
104#include "uggrid/uggridintersections.hh"
105#include "uggrid/uggridintersectioniterators.hh"
106#include "uggrid/uggridleveliterator.hh"
107#include "uggrid/uggridleafiterator.hh"
108#include "uggrid/uggridhieriterator.hh"
109#include "uggrid/uggridindexsets.hh"
110#include <dune/grid/uggrid/uggridviews.hh>
111#ifdef ModelP
112#include "uggrid/ugmessagebuffer.hh"
113#include "uggrid/uglbgatherscatter.hh"
114#endif
115
116// Not needed here, but included for user convenience
118
119#ifdef ModelP
120template <class DataHandle, int GridDim, int codim>
121const Dune::UGGrid<GridDim>* Dune::UGMessageBuffer<DataHandle, GridDim, codim>::grid_;
122
123template <class DataHandle, int GridDim, int codim>
124DataHandle *Dune::UGMessageBuffer<DataHandle,GridDim,codim>::duneDataHandle_ = nullptr;
125
126template <class DataHandle, int GridDim, int codim>
127int Dune::UGMessageBuffer<DataHandle,GridDim,codim>::level = -1;
128#endif // ModelP
129
130namespace Dune {
131
132#ifdef ModelP
133 using UGCollectiveCommunication = CollectiveCommunication<MPI_Comm>;
134#else
135 using UGCollectiveCommunication = CollectiveCommunication<No_Comm>;
136#endif
137
138 template<int dim>
139 struct UGGridFamily
140 {
141 typedef GridTraits<dim,dim,Dune::UGGrid<dim>,
142 UGGridGeometry,
143 UGGridEntity,
144 UGGridLevelIterator,
145 UGGridLeafIntersection,
146 UGGridLevelIntersection,
147 UGGridLeafIntersectionIterator,
148 UGGridLevelIntersectionIterator,
149 UGGridHierarchicIterator,
150 UGGridLeafIterator,
151 UGGridLevelIndexSet< const UGGrid<dim> >,
152 UGGridLeafIndexSet< const UGGrid<dim> >,
153 UGGridIdSet< const UGGrid<dim> >,
154 typename UG_NS<dim>::UG_ID_TYPE,
155 UGGridIdSet< const UGGrid<dim> >,
156 typename UG_NS<dim>::UG_ID_TYPE,
157 UGCollectiveCommunication,
158 UGGridLevelGridViewTraits,
159 UGGridLeafGridViewTraits,
160 UGGridEntitySeed,
161 UGGridLocalGeometry>
162 Traits;
163 };
164
165
166 //**********************************************************************
167 //
168 // --UGGrid
169 //
170 //**********************************************************************
171
203 template <int dim>
204 class UGGrid : public GridDefaultImplementation <dim, dim, double, UGGridFamily<dim> >
205 {
207
208 friend class UGGridGeometry<0,dim,const UGGrid<dim> >;
209 friend class UGGridGeometry<dim,dim,const UGGrid<dim> >;
210 friend class UGGridGeometry<1,2,const UGGrid<dim> >;
211 friend class UGGridGeometry<2,3,const UGGrid<dim> >;
212
213 friend class UGGridEntity <0,dim,const UGGrid<dim> >;
214 friend class UGGridEntity <1,dim,const UGGrid<dim> >;
215 friend class UGGridEntity <2,dim,const UGGrid<dim> >;
216 friend class UGGridEntity <dim,dim,const UGGrid<dim> >;
217 friend class UGGridHierarchicIterator<const UGGrid<dim> >;
218 friend class UGGridLeafIntersection<const UGGrid<dim> >;
219 friend class UGGridLevelIntersection<const UGGrid<dim> >;
220 friend class UGGridLeafIntersectionIterator<const UGGrid<dim> >;
221 friend class UGGridLevelIntersectionIterator<const UGGrid<dim> >;
222
223 friend class UGGridLevelIndexSet<const UGGrid<dim> >;
224 friend class UGGridLeafIndexSet<const UGGrid<dim> >;
225 friend class UGGridIdSet<const UGGrid<dim> >;
226 template <class GridImp_>
227 friend class UGGridLeafGridView;
228 template <class GridImp_>
229 friend class UGGridLevelGridView;
230
231 friend class GridFactory<UGGrid<dim> >;
232
233#ifdef ModelP
234 friend class UGLBGatherScatter;
235#endif
236
237 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
238 friend class UGGridLeafIterator;
239 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
240 friend class UGGridLevelIterator;
241
243 static_assert(dim==2 || dim==3, "Use UGGrid only for 2d and 3d!");
244
245 // The different instantiations are mutual friends so they can access
246 // each others numOfUGGrids field
247 friend class UGGrid<2>;
248 friend class UGGrid<3>;
249
250 //**********************************************************
251 // The Interface Methods
252 //**********************************************************
253 public:
255 typedef UGGridFamily<dim> GridFamily;
256
257 // the Traits
258 typedef typename UGGridFamily<dim>::Traits Traits;
259
261 typedef UG::DOUBLE ctype;
262
264 typedef unsigned int Rank;
265
269
271 ~UGGrid() noexcept(false);
272
275 int maxLevel() const;
276
278 template <typename Seed>
279 typename Traits::template Codim<Seed::codimension>::Entity
280 entity(const Seed& seed) const
281 {
282 const int codim = Seed::codimension;
283 return typename Traits::template Codim<codim>::Entity(UGGridEntity<codim,dim,const UGGrid<dim> >(seed.impl().target(),this));
284 }
285
288 int size (int level, int codim) const;
289
291 int size (int codim) const
292 {
293 return leafIndexSet().size(codim);
294 }
295
297 int size (int level, GeometryType type) const
298 {
299 return this->levelIndexSet(level).size(type);
300 }
301
303 int size (GeometryType type) const
304 {
305 return this->leafIndexSet().size(type);
306 }
307
309 size_t numBoundarySegments() const {
310 // The number is stored as a member of UGGrid upon grid creation.
311 // The corresponding data structure is not exported by UG. (It is in ug/dom/std/std_internal.h)
312 return numBoundarySegments_;
313 }
314
316 const typename Traits::GlobalIdSet& globalIdSet() const
317 {
318 return idSet_;
319 }
320
322 const typename Traits::LocalIdSet& localIdSet() const
323 {
324 return idSet_;
325 }
326
328 const typename Traits::LevelIndexSet& levelIndexSet(int level) const
329 {
330 if (level<0 || level>maxLevel())
331 DUNE_THROW(GridError, "levelIndexSet of nonexisting level " << level << " requested!");
332 return *levelIndexSets_[level];
333 }
334
336 const typename Traits::LeafIndexSet& leafIndexSet() const
337 {
338 return leafIndexSet_;
339 }
340
343
356 bool mark(int refCount, const typename Traits::template Codim<0>::Entity & e );
357
413 bool mark(const typename Traits::template Codim<0>::Entity & e,
414 typename UG_NS<dim>::RefinementRule rule,
415 int side=0);
416
418 int getMark(const typename Traits::template Codim<0>::Entity& e) const;
419
422 bool preAdapt();
423
425 bool adapt();
426
428 void postAdapt();
438 template<class DataHandle>
439 bool loadBalance (DataHandle& dataHandle)
440 {
441#ifdef ModelP
442 // gather element data
443 if (dataHandle.contains(dim, 0))
444 UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
445
446 // gather node data
447 if (dataHandle.contains(dim,dim))
448 UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
449#endif
450
451 // the load balancing step now also attaches
452 // the data to the entities and distributes it
453 loadBalance();
454
455#ifdef ModelP
456 // scatter element data
457 if (dataHandle.contains(dim, 0))
458 UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
459
460 // scatter node data
461 if (dataHandle.contains(dim,dim))
462 UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
463#endif
464
465 return true;
466 }
467
474 bool loadBalance(int minlevel=0);
475
506 bool loadBalance(const std::vector<Rank>& targetProcessors, unsigned int fromLevel);
507
517 template<class DataHandle>
518 bool loadBalance (const std::vector<Rank>& targetProcessors, unsigned int fromLevel, DataHandle& dataHandle)
519 {
520#ifdef ModelP
521 // gather element data
522 if (dataHandle.contains(dim, 0))
523 UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
524
525 // gather node data
526 if (dataHandle.contains(dim,dim))
527 UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
528#endif
529
530 // the load balancing step now also attaches
531 // the data to the entities and distributes it
532 loadBalance(targetProcessors,fromLevel);
533
534#ifdef ModelP
535 // scatter element data
536 if (dataHandle.contains(dim, 0))
537 UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
538
539 // scatter node data
540 if (dataHandle.contains(dim,dim))
541 UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
542#endif
543
544 return true;
545 }
546
549 {
550 return ccobj_;
551 }
552
553 protected:
554#ifdef ModelP
555 template <int codim, class GridView, class DataHandle>
556 void communicateUG_(const GridView& gv, int level,
557 DataHandle &dataHandle,
558 InterfaceType iftype,
559 CommunicationDirection dir) const
560 {
561 typename UG_NS<dim>::DDD_IF_DIR ugIfDir;
562 // Translate the communication direction from Dune-Speak to UG-Speak
563 if (dir==ForwardCommunication)
564 ugIfDir = UG_NS<dim>::IF_FORWARD();
565 else
566 ugIfDir = UG_NS<dim>::IF_BACKWARD();
567
568 typedef UGMessageBuffer<DataHandle,dim,codim> UGMsgBuf;
569 UGMsgBuf::duneDataHandle_ = &dataHandle;
570
571 UGMsgBuf::level = level;
572
573 std::vector<typename UG_NS<dim>::DDD_IF> ugIfs = findDDDInterfaces(iftype, codim);
574
575 unsigned bufSize = UGMsgBuf::ugBufferSize(gv);
576 if (!bufSize)
577 return; // we don't need to communicate if we don't have any data!
578 UGMsgBuf::grid_ = this;
579 for (unsigned i=0; i < ugIfs.size(); ++i)
580 UG_NS<dim>::DDD_IFOneway(multigrid_->dddContext(),
581 ugIfs[i],
582 ugIfDir,
583 bufSize,
584 &UGMsgBuf::ugGather_,
585 &UGMsgBuf::ugScatter_);
586 }
587
589 std::vector<typename UG_NS<dim>::DDD_IF> findDDDInterfaces(InterfaceType iftype,
590 int codim) const;
591#endif
592 public:
593 // **********************************************************
594 // End of Interface Methods
595 // **********************************************************
596
604 void getChildrenOfSubface(const typename Traits::template Codim<0>::Entity & e,
605 int elementSide,
606 int maxl,
607 std::vector<typename Traits::template Codim<0>::Entity>& childElements,
608 std::vector<unsigned char>& childElementSides) const;
609
615 COPY
616 };
617
623 NONE
624 };
625
628 refinementType_ = type;
629 }
630
633 closureType_ = type;
634 }
635
639 void setPosition(const typename Traits::template Codim<dim>::Entity& e,
640 const FieldVector<double, dim>& pos);
641
646 void globalRefine(int n);
647
652 void saveState(const std::string& filename) const;
653
658 void loadState(const std::string& filename);
659
660 private:
662 typename UG_NS<dim>::MultiGrid* multigrid_;
663
666
672 void setIndices(bool setLevelZero,
673 std::vector<unsigned int>* nodePermutation);
674
675 // Each UGGrid object has a unique name to identify it in the
676 // UG environment structure
677 std::string name_;
678
679 // Our set of level indices
680 std::vector<std::shared_ptr<UGGridLevelIndexSet<const UGGrid<dim> > > > levelIndexSets_;
681
682 UGGridLeafIndexSet<const UGGrid<dim> > leafIndexSet_;
683
684 // One id set implementation
685 // Used for both the local and the global UGGrid id sets
686 UGGridIdSet<const UGGrid<dim> > idSet_;
687
689 RefinementType refinementType_;
690
692 ClosureType closureType_;
693
701 static int numOfUGGrids;
702
708 bool someElementHasBeenMarkedForRefinement_;
709
715 bool someElementHasBeenMarkedForCoarsening_;
716
718 std::vector<std::shared_ptr<BoundarySegment<dim> > > boundarySegments_;
719
725 unsigned int numBoundarySegments_;
726
727 }; // end Class UGGrid
728
729 namespace Capabilities
730 {
746 template<int dim, int codim>
747 struct hasEntity< UGGrid<dim>, codim>
748 {
749 static const bool v = true;
750 };
751
757 template<int dim, int codim>
758 struct hasEntityIterator<UGGrid<dim>, codim>
759 {
760 static const bool v = false;
761 };
762
767 template<int dim>
768 struct hasEntityIterator<UGGrid<dim>, 0>
769 {
770 static const bool v = true;
771 };
772
777 template<int dim>
778 struct hasEntityIterator<UGGrid<dim>, dim>
779 {
780 static const bool v = true;
781 };
782
786 template<int dim, int codim>
787 struct canCommunicate<UGGrid<dim>, codim>
788 {
789 static const bool v = (codim>=0 && codim<=dim);
790 };
791
795 template<int dim>
797 {
798 static const bool v = true;
799 };
800
804 template<int dim>
806 {
807 static const bool v = false;
808 };
809
810 }
811
812} // namespace Dune
813
814#endif // HAVE_UG || DOXYGEN
815#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:95
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:123
Definition: grid.hh:851
Traits::LeafGridView leafGridView() const
View for the leaf grid for All_Partition.
Definition: grid.hh:871
bool loadBalance()
default implementation of load balance does nothing and returns false
Definition: grid.hh:937
Base class for exceptions in Dune grid modules.
Definition: exceptions.hh:18
Provide a generic factory class for unstructured grids.
Definition: gridfactory.hh:312
Grid view abstract base class.
Definition: gridview.hh:63
Id Set Interface.
Definition: indexidset.hh:450
Index Set Interface base class.
Definition: indexidset.hh:76
auto size(GeometryType type) const
Return total number of entities of given geometry type in entity set .
Definition: indexidset.hh:221
Front-end for the grid manager of the finite element toolbox UG3.
Definition: uggrid.hh:205
void postAdapt()
Clean up refinement markers.
size_t numBoundarySegments() const
Return the number of boundary segments.
Definition: uggrid.hh:309
void setRefinementType(RefinementType type)
Sets the type of grid refinement.
Definition: uggrid.hh:627
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:548
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:518
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:255
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:328
void setClosureType(ClosureType type)
Sets the type of grid refinement closure.
Definition: uggrid.hh:632
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition: uggrid.hh:322
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition: uggrid.hh:297
RefinementType
The different forms of grid refinement that UG supports.
Definition: uggrid.hh:611
@ COPY
New level consists of the refined elements and the unrefined ones, too.
Definition: uggrid.hh:615
@ LOCAL
New level consists only of the refined elements and the closure.
Definition: uggrid.hh:613
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: uggrid.hh:303
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:261
bool adapt()
Triggers the grid refinement process.
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition: uggrid.hh:316
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:439
Traits::template Codim< Seed::codimension >::Entity entity(const Seed &seed) const
Create an Entity from an EntitySeed.
Definition: uggrid.hh:280
unsigned int Rank
The type used for process ranks.
Definition: uggrid.hh:264
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:619
@ GREEN
Standard red/green refinement.
Definition: uggrid.hh:621
@ NONE
No closure, results in nonconforming meshes.
Definition: uggrid.hh:623
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition: uggrid.hh:336
int size(int codim) const
number of leaf entities per codim in this process
Definition: uggrid.hh:291
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.
A few common exception classes.
#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
Implements an utility class that provides MPI's collective communication methods.
Helpers for dealing with MPI.
Dune namespace.
Definition: alignedallocator.hh:11
specialize with 'true' for all codims that a grid can communicate data on (default=false)
Definition: capabilities.hh:95
specialize with 'true' for all codims that a grid provides an iterator for (default=hasEntity<codim>:...
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:984
The specialization of the generic GridFactory for UGGrid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)