DUNE PDELab (git)

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
13#include <memory>
14
19
20#include <dune-grid-config.hh> // HAVE_DUNE_UGGRID
24
25#if HAVE_DUNE_UGGRID || DOXYGEN
26
27#ifdef ModelP
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/* The following define tells the UG headers that we want access to a few
47 special fields, for example the extra index fields in the element data structures.
48 This define remains only for backwards compatibility with older version of UG.
49 All dune-uggrid versions since 2016-08-05 do not need this #define or the #undef
50 further below. */
51#define FOR_DUNE
52
53// Set UG's space-dimension flag to 2d
54#define UG_DIM_2
55// And include all necessary UG headers
56#include "uggrid/ugincludes.hh"
57#undef DUNE_UGINCLUDES_HH
58
59// Wrap a few large UG macros by functions before they get undef'ed away.
60// Here: The 2d-version of the macros
61#define UG_DIM 2
62#include "uggrid/ugwrapper.hh"
63#undef UG_DIM
64
65// UG defines a whole load of preprocessor macros. ug_undefs.hh undefines
66// them all, so we don't get name clashes.
67#include "uggrid/ug_undefs.hh"
68#undef UG_DIM_2
69
70/* Now we're done with 2d, and we can do the whole thing over again for 3d */
71
72/* All macros set by UG have been unset. This includes the macros that ensure
73 single inclusion of headers. We can thus include them again. However, we
74 only want to include those headers again that contain dimension-dependent stuff.
75 Therefore, we set a few single-inclusion defines manually before including
76 ugincludes.hh again.
77 */
78#define UGTYPES_H
79#define __HEAPS__
80#define __UGENV__
81#define __DEVICESH__
82#ifdef ModelP
83#define __PPIF__
84#endif
85
86#define UG_DIM_3
87#include "uggrid/ugincludes.hh"
88#undef DUNE_UGINCLUDES_HH
89
90// Wrap a few large UG macros by functions before they get undef'ed away.
91// This time it's the 3d-versions.
92#define UG_DIM 3
93#include "uggrid/ugwrapper.hh"
94#undef UG_DIM
95
96// undef all macros defined by UG
97#include "uggrid/ug_undefs.hh"
98
99#undef UG_DIM_3
100#undef FOR_DUNE
101
102// The components of the UGGrid interface
103#include "uggrid/uggridgeometry.hh"
104#include "uggrid/uggridlocalgeometry.hh"
105#include "uggrid/uggridentity.hh"
106#include "uggrid/uggridentityseed.hh"
107#include "uggrid/uggridintersections.hh"
108#include "uggrid/uggridintersectioniterators.hh"
109#include "uggrid/uggridleveliterator.hh"
110#include "uggrid/uggridleafiterator.hh"
111#include "uggrid/uggridhieriterator.hh"
112#include "uggrid/uggridindexsets.hh"
113#include <dune/grid/uggrid/uggridviews.hh>
114#ifdef ModelP
115#include "uggrid/ugmessagebuffer.hh"
116#include "uggrid/uglbgatherscatter.hh"
117#endif
118
119// Not needed here, but included for user convenience
121
122#ifdef ModelP
123template <class DataHandle, int GridDim, int codim>
124const Dune::UGGrid<GridDim>* Dune::UGMessageBuffer<DataHandle, GridDim, codim>::grid_;
125
126template <class DataHandle, int GridDim, int codim>
127DataHandle *Dune::UGMessageBuffer<DataHandle,GridDim,codim>::duneDataHandle_ = nullptr;
128
129template <class DataHandle, int GridDim, int codim>
130int Dune::UGMessageBuffer<DataHandle,GridDim,codim>::level = -1;
131#endif // ModelP
132
133namespace Dune {
134
135#ifdef ModelP
136 using UGCommunication = Communication<MPI_Comm>;
137#else
138 using UGCommunication = Communication<No_Comm>;
139#endif
140
141 template<int dim>
142 struct UGGridFamily
143 {
144 typedef GridTraits<dim,dim,Dune::UGGrid<dim>,
145 UGGridGeometry,
146 UGGridEntity,
147 UGGridLevelIterator,
148 UGGridLeafIntersection,
149 UGGridLevelIntersection,
150 UGGridLeafIntersectionIterator,
151 UGGridLevelIntersectionIterator,
152 UGGridHierarchicIterator,
153 UGGridLeafIterator,
154 UGGridLevelIndexSet< const UGGrid<dim> >,
155 UGGridLeafIndexSet< const UGGrid<dim> >,
156 UGGridIdSet< const UGGrid<dim> >,
157 typename UG_NS<dim>::UG_ID_TYPE,
158 UGGridIdSet< const UGGrid<dim> >,
159 typename UG_NS<dim>::UG_ID_TYPE,
160 UGCommunication,
161 UGGridLevelGridViewTraits,
162 UGGridLeafGridViewTraits,
163 UGGridEntitySeed,
164 UGGridLocalGeometry>
165 Traits;
166 };
167
168
169 //**********************************************************************
170 //
171 // --UGGrid
172 //
173 //**********************************************************************
174
206 template <int dim>
207 class UGGrid : public GridDefaultImplementation <dim, dim, double, UGGridFamily<dim> >
208 {
210
211 friend class UGGridGeometry<0,dim,const UGGrid<dim> >;
212 friend class UGGridGeometry<dim,dim,const UGGrid<dim> >;
213 friend class UGGridGeometry<1,2,const UGGrid<dim> >;
214 friend class UGGridGeometry<2,3,const UGGrid<dim> >;
215
216 friend class UGGridEntity <0,dim,const UGGrid<dim> >;
217 friend class UGGridEntity <1,dim,const UGGrid<dim> >;
218 friend class UGGridEntity <2,dim,const UGGrid<dim> >;
219 friend class UGGridEntity <dim,dim,const UGGrid<dim> >;
220 friend class UGGridHierarchicIterator<const UGGrid<dim> >;
221 friend class UGGridLeafIntersection<const UGGrid<dim> >;
222 friend class UGGridLevelIntersection<const UGGrid<dim> >;
223 friend class UGGridLeafIntersectionIterator<const UGGrid<dim> >;
224 friend class UGGridLevelIntersectionIterator<const UGGrid<dim> >;
225
226 friend class UGGridLevelIndexSet<const UGGrid<dim> >;
227 friend class UGGridLeafIndexSet<const UGGrid<dim> >;
228 friend class UGGridIdSet<const UGGrid<dim> >;
229 template <class GridImp_>
230 friend class UGGridLeafGridView;
231 template <class GridImp_>
232 friend class UGGridLevelGridView;
233
234 friend class GridFactory<UGGrid<dim> >;
235
236#ifdef ModelP
237 friend class UGLBGatherScatter;
238#endif
239
240 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
241 friend class UGGridLeafIterator;
242 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
243 friend class UGGridLevelIterator;
244
246 static_assert(dim==2 || dim==3, "Use UGGrid only for 2d and 3d!");
247
248 // The different instantiations are mutual friends so they can access
249 // each others numOfUGGrids field
250 friend class UGGrid<2>;
251 friend class UGGrid<3>;
252
253 //**********************************************************
254 // The Interface Methods
255 //**********************************************************
256 public:
258 typedef UGGridFamily<dim> GridFamily;
259
260 // the Traits
261 typedef typename UGGridFamily<dim>::Traits Traits;
262
264 typedef UG::DOUBLE ctype;
265
267 typedef unsigned int Rank;
268
272
274 ~UGGrid() noexcept(false);
275
278 int maxLevel() const;
279
281 template <typename Seed>
282 typename Traits::template Codim<Seed::codimension>::Entity
283 entity(const Seed& seed) const
284 {
285 const int codim = Seed::codimension;
286 return typename Traits::template Codim<codim>::Entity(UGGridEntity<codim,dim,const UGGrid<dim> >(seed.impl().target(),this));
287 }
288
291 int size (int level, int codim) const;
292
294 int size (int codim) const
295 {
296 return leafIndexSet().size(codim);
297 }
298
300 int size (int level, GeometryType type) const
301 {
302 return this->levelIndexSet(level).size(type);
303 }
304
306 int size (GeometryType type) const
307 {
308 return this->leafIndexSet().size(type);
309 }
310
312 size_t numBoundarySegments() const {
313 // The number is stored as a member of UGGrid upon grid creation.
314 // The corresponding data structure is not exported by UG. (It is in ug/dom/std/std_internal.h)
315 return numBoundarySegments_;
316 }
317
319 const typename Traits::GlobalIdSet& globalIdSet() const
320 {
321 return idSet_;
322 }
323
325 const typename Traits::LocalIdSet& localIdSet() const
326 {
327 return idSet_;
328 }
329
331 const typename Traits::LevelIndexSet& levelIndexSet(int level) const
332 {
333 if (level<0 || level>maxLevel())
334 DUNE_THROW(GridError, "levelIndexSet of nonexisting level " << level << " requested!");
335 return *levelIndexSets_[level];
336 }
337
339 const typename Traits::LeafIndexSet& leafIndexSet() const
340 {
341 return leafIndexSet_;
342 }
343
346
359 bool mark(int refCount, const typename Traits::template Codim<0>::Entity & e );
360
416 bool mark(const typename Traits::template Codim<0>::Entity & e,
417 typename UG_NS<dim>::RefinementRule rule,
418 int side=0);
419
421 int getMark(const typename Traits::template Codim<0>::Entity& e) const;
422
425 bool preAdapt();
426
428 bool adapt();
429
431 void postAdapt();
441 template<class DataHandle>
442 bool loadBalance (DataHandle& dataHandle)
443 {
444#ifdef ModelP
445 // gather element data
446 if (dataHandle.contains(dim, 0))
447 UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
448
449 // gather node data
450 if (dataHandle.contains(dim,dim))
451 UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
452#endif
453
454 // the load balancing step now also attaches
455 // the data to the entities and distributes it
456 loadBalance();
457
458#ifdef ModelP
459 // scatter element data
460 if (dataHandle.contains(dim, 0))
461 UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
462
463 // scatter node data
464 if (dataHandle.contains(dim,dim))
465 UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
466#endif
467
468 return true;
469 }
470
477 bool loadBalance(int minlevel=0);
478
509 bool loadBalance(const std::vector<Rank>& targetProcessors, unsigned int fromLevel);
510
520 template<class DataHandle>
521 bool loadBalance (const std::vector<Rank>& targetProcessors, unsigned int fromLevel, DataHandle& dataHandle)
522 {
523#ifdef ModelP
524 // gather element data
525 if (dataHandle.contains(dim, 0))
526 UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
527
528 // gather node data
529 if (dataHandle.contains(dim,dim))
530 UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
531#endif
532
533 // the load balancing step now also attaches
534 // the data to the entities and distributes it
535 loadBalance(targetProcessors,fromLevel);
536
537#ifdef ModelP
538 // scatter element data
539 if (dataHandle.contains(dim, 0))
540 UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
541
542 // scatter node data
543 if (dataHandle.contains(dim,dim))
544 UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
545#endif
546
547 return true;
548 }
549
551 const UGCommunication& comm () const
552 {
553 return ccobj_;
554 }
555
556 protected:
557#ifdef ModelP
558 template <int codim, class GridView, class DataHandle>
559 void communicateUG_(const GridView& gv, int level,
560 DataHandle &dataHandle,
561 InterfaceType iftype,
562 CommunicationDirection dir) const
563 {
564 typename UG_NS<dim>::DDD_IF_DIR ugIfDir;
565 // Translate the communication direction from Dune-Speak to UG-Speak
566 if (dir==ForwardCommunication)
567 ugIfDir = UG_NS<dim>::IF_FORWARD();
568 else
569 ugIfDir = UG_NS<dim>::IF_BACKWARD();
570
571 typedef UGMessageBuffer<DataHandle,dim,codim> UGMsgBuf;
572 UGMsgBuf::duneDataHandle_ = &dataHandle;
573
574 UGMsgBuf::level = level;
575
576 std::vector<typename UG_NS<dim>::DDD_IF> ugIfs = findDDDInterfaces(iftype, codim);
577
578 unsigned bufSize = UGMsgBuf::ugBufferSize(gv);
579 if (!bufSize)
580 return; // we don't need to communicate if we don't have any data!
581 UGMsgBuf::grid_ = this;
582 for (unsigned i=0; i < ugIfs.size(); ++i)
583 UG_NS<dim>::DDD_IFOneway(multigrid_->dddContext(),
584 ugIfs[i],
585 ugIfDir,
586 bufSize,
587 &UGMsgBuf::ugGather_,
588 &UGMsgBuf::ugScatter_);
589 }
590
592 std::vector<typename UG_NS<dim>::DDD_IF> findDDDInterfaces(InterfaceType iftype,
593 int codim) const;
594#endif
595 public:
596 // **********************************************************
597 // End of Interface Methods
598 // **********************************************************
599
607 void getChildrenOfSubface(const typename Traits::template Codim<0>::Entity & e,
608 int elementSide,
609 int maxl,
610 std::vector<typename Traits::template Codim<0>::Entity>& childElements,
611 std::vector<unsigned char>& childElementSides) const;
612
618 COPY
619 };
620
626 NONE
627 };
628
631 refinementType_ = type;
632 }
633
636 closureType_ = type;
637 }
638
642 void setPosition(const typename Traits::template Codim<dim>::Entity& e,
643 const FieldVector<double, dim>& pos);
644
649 void globalRefine(int n);
650
655 void saveState(const std::string& filename) const;
656
661 void loadState(const std::string& filename);
662
663 private:
665 typename UG_NS<dim>::MultiGrid* multigrid_;
666
668 UGCommunication ccobj_;
669
675 void setIndices(bool setLevelZero,
676 std::vector<unsigned int>* nodePermutation);
677
678 // Each UGGrid object has a unique name to identify it in the
679 // UG environment structure
680 std::string name_;
681
682 // Our set of level indices
683 std::vector<std::shared_ptr<UGGridLevelIndexSet<const UGGrid<dim> > > > levelIndexSets_;
684
685 UGGridLeafIndexSet<const UGGrid<dim> > leafIndexSet_;
686
687 // One id set implementation
688 // Used for both the local and the global UGGrid id sets
689 UGGridIdSet<const UGGrid<dim> > idSet_;
690
692 RefinementType refinementType_;
693
695 ClosureType closureType_;
696
704 static int numOfUGGrids;
705
711 bool someElementHasBeenMarkedForRefinement_;
712
718 bool someElementHasBeenMarkedForCoarsening_;
719
721 std::vector<std::shared_ptr<BoundarySegment<dim> > > boundarySegments_;
722
728 unsigned int numBoundarySegments_;
729
730 }; // end Class UGGrid
731
732 namespace Capabilities
733 {
749 template<int dim, int codim>
750 struct hasEntity< UGGrid<dim>, codim>
751 {
752 static const bool v = true;
753 };
754
760 template<int dim, int codim>
761 struct hasEntityIterator<UGGrid<dim>, codim>
762 {
763 static const bool v = false;
764 };
765
770 template<int dim>
771 struct hasEntityIterator<UGGrid<dim>, 0>
772 {
773 static const bool v = true;
774 };
775
780 template<int dim>
781 struct hasEntityIterator<UGGrid<dim>, dim>
782 {
783 static const bool v = true;
784 };
785
789 template<int dim, int codim>
790 struct canCommunicate<UGGrid<dim>, codim>
791 {
792 static const bool v = (codim>=0 && codim<=dim);
793 };
794
798 template<int dim>
800 {
801 static const bool v = true;
802 };
803
807 template<int dim>
809 {
810 static const bool v = false;
811 };
812
816 template<int dim>
817 struct viewThreadSafe< UGGrid<dim> > {
818 static const bool v = true;
819 };
820
821 }
822
823} // namespace Dune
824
825#endif // HAVE_DUNE_UGGRID || DOXYGEN
826#endif // DUNE_UGGRID_HH
Base class for grid boundary segments of arbitrary geometry.
Wrapper class for entities.
Definition: entity.hh:66
vector space out of a tensor product of fields.
Definition: fvector.hh:91
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
Definition: grid.hh:848
Traits::LeafGridView leafGridView() const
View for the leaf grid for All_Partition.
Definition: grid.hh:868
bool loadBalance()
default implementation of load balance does nothing and returns false
Definition: grid.hh:962
Base class for exceptions in Dune grid modules.
Definition: exceptions.hh:20
Provide a generic factory class for unstructured grids.
Definition: gridfactory.hh:275
Grid view abstract base class.
Definition: gridview.hh:66
Id Set Interface.
Definition: indexidset.hh:447
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
Front-end for the grid manager of the finite element toolbox UG3.
Definition: uggrid.hh:208
void postAdapt()
Clean up refinement markers.
size_t numBoundarySegments() const
Return the number of boundary segments.
Definition: uggrid.hh:312
void setRefinementType(RefinementType type)
Sets the type of grid refinement.
Definition: uggrid.hh:630
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.
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:521
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:258
bool loadBalance(const std::vector< Rank > &targetProcessors, unsigned int fromLevel)
Distribute this grid over a distributed machine.
const UGCommunication & comm() const
Definition: uggrid.hh:551
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:331
void setClosureType(ClosureType type)
Sets the type of grid refinement closure.
Definition: uggrid.hh:635
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition: uggrid.hh:325
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition: uggrid.hh:300
RefinementType
The different forms of grid refinement that UG supports.
Definition: uggrid.hh:614
@ COPY
New level consists of the refined elements and the unrefined ones, too.
Definition: uggrid.hh:618
@ LOCAL
New level consists only of the refined elements and the closure.
Definition: uggrid.hh:616
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: uggrid.hh:306
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:264
bool adapt()
Triggers the grid refinement process.
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition: uggrid.hh:319
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:442
Traits::template Codim< Seed::codimension >::Entity entity(const Seed &seed) const
Create an Entity from an EntitySeed.
Definition: uggrid.hh:283
unsigned int Rank
The type used for process ranks.
Definition: uggrid.hh:267
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:622
@ GREEN
Standard red/green refinement.
Definition: uggrid.hh:624
@ NONE
No closure, results in nonconforming meshes.
Definition: uggrid.hh:626
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition: uggrid.hh:339
int size(int codim) const
number of leaf entities per codim in this process
Definition: uggrid.hh:294
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:218
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
Implements an utility class that provides MPI's collective communication methods.
Helpers for dealing with MPI.
Dune namespace.
Definition: alignedallocator.hh:13
specialize with 'true' for all codims that a grid can communicate data on (default=false)
Definition: capabilities.hh:97
specialize with 'true' for all codims that a grid provides an iterator for (default=hasEntity<codim>:...
Definition: capabilities.hh:74
Specialize with 'true' for all codims that a grid implements entities for. (default=false)
Definition: capabilities.hh:58
Specialize with 'true' if implementation guarantees a conforming leaf grid. (default=false)
Definition: capabilities.hh:115
Specialize with 'true' if implementation guarantees conforming level grids. (default=false)
Definition: capabilities.hh:106
Specialize with 'true' if the grid implementation is thread safe, while it is not modified....
Definition: capabilities.hh:169
Static tag representing a codimension.
Definition: dimension.hh:24
A traits struct that collects all associated types of one grid model.
Definition: grid.hh:1013
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 23, 23:29, 2024)