Dune Core Modules (2.5.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
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>
137DataHandle *Dune::UGMessageBufferBase<DataHandle,GridDim,codim>::duneDataHandle_ = 0;
138
139template <class DataHandle, int GridDim, int codim>
140int Dune::UGMessageBufferBase<DataHandle,GridDim,codim>::level = -1;
141#endif // ModelP
142
143namespace Dune {
144
145#ifdef ModelP
146 template <int dim>
147 class CollectiveCommunication<Dune::UGGrid<dim> > :
148 public CollectiveCommunication< Dune::MPIHelper::MPICommunicator >
149 {
150 typedef CollectiveCommunication< Dune::MPIHelper::MPICommunicator > ParentType;
151 public:
153 : ParentType(MPIHelper::getCommunicator())
154 {}
155 };
156#endif // ModelP
157
158 template<int dim>
159 struct UGGridFamily
160 {
161 typedef GridTraits<dim,dim,Dune::UGGrid<dim>,
162 UGGridGeometry,
163 UGGridEntity,
164 UGGridLevelIterator,
165 UGGridLeafIntersection,
166 UGGridLevelIntersection,
167 UGGridLeafIntersectionIterator,
168 UGGridLevelIntersectionIterator,
169 UGGridHierarchicIterator,
170 UGGridLeafIterator,
171 UGGridLevelIndexSet< const UGGrid<dim> >,
172 UGGridLeafIndexSet< const UGGrid<dim> >,
173 UGGridIdSet< const UGGrid<dim> >,
174 typename UG_NS<dim>::UG_ID_TYPE,
175 UGGridIdSet< const UGGrid<dim> >,
176 typename UG_NS<dim>::UG_ID_TYPE,
177 CollectiveCommunication<Dune::UGGrid<dim> >,
178 UGGridLevelGridViewTraits,
179 UGGridLeafGridViewTraits,
180 UGGridEntitySeed,
181 UGGridLocalGeometry>
182 Traits;
183 };
184
185
186 //**********************************************************************
187 //
188 // --UGGrid
189 //
190 //**********************************************************************
191
229 template <int dim>
230 class UGGrid : public GridDefaultImplementation <dim, dim, double, UGGridFamily<dim> >
231 {
233
234 friend class UGGridGeometry<0,dim,const UGGrid<dim> >;
235 friend class UGGridGeometry<dim,dim,const UGGrid<dim> >;
236 friend class UGGridGeometry<1,2,const UGGrid<dim> >;
237 friend class UGGridGeometry<2,3,const UGGrid<dim> >;
238
239 friend class UGGridEntity <0,dim,const UGGrid<dim> >;
240 friend class UGGridEntity <dim,dim,const UGGrid<dim> >;
241 friend class UGGridHierarchicIterator<const UGGrid<dim> >;
242 friend class UGGridLeafIntersection<const UGGrid<dim> >;
243 friend class UGGridLevelIntersection<const UGGrid<dim> >;
244 friend class UGGridLeafIntersectionIterator<const UGGrid<dim> >;
245 friend class UGGridLevelIntersectionIterator<const UGGrid<dim> >;
246
247 friend class UGGridLevelIndexSet<const UGGrid<dim> >;
248 friend class UGGridLeafIndexSet<const UGGrid<dim> >;
249 friend class UGGridIdSet<const UGGrid<dim> >;
250 template <class GridImp_>
251 friend class UGGridLeafGridView;
252 template <class GridImp_>
253 friend class UGGridLevelGridView;
254
255 friend class GridFactory<UGGrid<dim> >;
256
257#ifdef ModelP
258 friend class UGLBGatherScatter;
259#endif
260
261 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
262 friend class UGGridLeafIterator;
263 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
264 friend class UGGridLevelIterator;
265
267 static_assert(dim==2 || dim==3, "Use UGGrid only for 2d and 3d!");
268
269 // The different instantiations are mutual friends so they can access
270 // each others numOfUGGrids field
271 friend class UGGrid<2>;
272 friend class UGGrid<3>;
273
274 //**********************************************************
275 // The Interface Methods
276 //**********************************************************
277 public:
279 typedef UGGridFamily<dim> GridFamily;
280
281 // the Traits
282 typedef typename UGGridFamily<dim>::Traits Traits;
283
285 typedef UG::DOUBLE ctype;
286
288 typedef unsigned int Rank;
289
296
298 ~UGGrid() noexcept(false);
299
302 int maxLevel() const;
303
305 template <typename Seed>
306 typename Traits::template Codim<Seed::codimension>::Entity
307 entity(const Seed& seed) const
308 {
309 const int codim = Seed::codimension;
310 return typename Traits::template Codim<codim>::Entity(UGGridEntity<codim,dim,const UGGrid<dim> >(this->getRealImplementation(seed).target(),this));
311 }
312
315 int size (int level, int codim) const;
316
318 int size (int codim) const
319 {
320 return leafIndexSet().size(codim);
321 }
322
324 int size (int level, GeometryType type) const
325 {
326 return this->levelIndexSet(level).size(type);
327 }
328
330 int size (GeometryType type) const
331 {
332 return this->leafIndexSet().size(type);
333 }
334
336 size_t numBoundarySegments() const {
337 // The number is stored as a member of UGGrid upon grid creation.
338 // The corresponding data structure is not exported by UG. (It is in ug/dom/std/std_internal.h)
339 return numBoundarySegments_;
340 }
341
343 const typename Traits::GlobalIdSet& globalIdSet() const
344 {
345 return idSet_;
346 }
347
349 const typename Traits::LocalIdSet& localIdSet() const
350 {
351 return idSet_;
352 }
353
355 const typename Traits::LevelIndexSet& levelIndexSet(int level) const
356 {
357 if (level<0 || level>maxLevel())
358 DUNE_THROW(GridError, "levelIndexSet of nonexisting level " << level << " requested!");
359 return *levelIndexSets_[level];
360 }
361
363 const typename Traits::LeafIndexSet& leafIndexSet() const
364 {
365 return leafIndexSet_;
366 }
367
370
383 bool mark(int refCount, const typename Traits::template Codim<0>::Entity & e );
384
440 bool mark(const typename Traits::template Codim<0>::Entity & e,
441 typename UG_NS<dim>::RefinementRule rule,
442 int side=0);
443
445 int getMark(const typename Traits::template Codim<0>::Entity& e) const;
446
449 bool preAdapt();
450
452 bool adapt();
453
455 void postAdapt();
459 unsigned int overlapSize(int codim) const {
460 return 0;
461 }
462
464 unsigned int ghostSize(int codim) const {
465 return (codim==0) ? 1 : 0;
466 }
467
469 unsigned int overlapSize(int level, int codim) const {
470 return 0;
471 }
472
474 unsigned int ghostSize(int level, int codim) const {
475 return (codim==0) ? 1 : 0;
476 }
477
485 template<class DataHandle>
486 bool loadBalance (DataHandle& dataHandle)
487 {
488#ifdef ModelP
489 // gather element data
490 // UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
491
492 // gather node data
493 if (dataHandle.contains(dim,dim))
494 UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
495#endif
496
497 // the load balancing step now also attaches
498 // the data to the entities and distributes it
499 loadBalance();
500
501#ifdef ModelP
502 // scatter element data
503 // UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
504
505 // scatter node data
506 if (dataHandle.contains(dim,dim))
507 UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
508#endif
509
510 return true;
511 }
512
519 bool loadBalance(int minlevel=0);
520
551 bool loadBalance(const std::vector<Rank>& targetProcessors, unsigned int fromLevel);
552
562 template<class DataHandle>
563 bool loadBalance (const std::vector<Rank>& targetProcessors, unsigned int fromLevel, DataHandle& dataHandle)
564 {
565#ifdef ModelP
566 // gather element data
567 // UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
568
569 // gather node data
570 UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
571#endif
572
573 // the load balancing step now also attaches
574 // the data to the entities and distributes it
575 loadBalance(targetProcessors,fromLevel);
576
577#ifdef ModelP
578 // scatter element data
579 // UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
580
581 // scatter node data
582 UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
583#endif
584
585 return true;
586 }
587
588
599 template<class DataHandle>
600 void communicate (DataHandle& dataHandle, InterfaceType iftype, CommunicationDirection dir, int level) const
601 {
602#ifdef ModelP
603 typedef typename UGGrid::LevelGridView LevelGridView;
604
605 for (int curCodim = 0; curCodim <= dim; ++curCodim) {
606 if (!dataHandle.contains(dim, curCodim))
607 continue;
608
609 if (curCodim == 0)
610 communicateUG_<LevelGridView, DataHandle, 0>(this->levelGridView(level), level, dataHandle, iftype, dir);
611 else if (curCodim == dim)
612 communicateUG_<LevelGridView, DataHandle, dim>(this->levelGridView(level), level, dataHandle, iftype, dir);
613 else if (curCodim == dim - 1)
614 communicateUG_<LevelGridView, DataHandle, dim-1>(this->levelGridView(level), level, dataHandle, iftype, dir);
615 else if (curCodim == 1)
616 communicateUG_<LevelGridView, DataHandle, 1>(this->levelGridView(level), level, dataHandle, iftype, dir);
617 else
619 className(*this) << "::communicate(): Not "
620 "supported for dim " << dim << " and codim " << curCodim);
621 }
622#endif // ModelP
623 }
624
634 template<class DataHandle>
635 void communicate(DataHandle& dataHandle, InterfaceType iftype, CommunicationDirection dir) const
636 {
637#ifdef ModelP
638 typedef typename UGGrid::LeafGridView LeafGridView;
639
640 for (int curCodim = 0; curCodim <= dim; ++curCodim) {
641 if (!dataHandle.contains(dim, curCodim))
642 continue;
643 int level = -1;
644 if (curCodim == 0)
645 communicateUG_<LeafGridView, DataHandle, 0>(this->leafGridView(), level, dataHandle, iftype, dir);
646 else if (curCodim == dim)
647 communicateUG_<LeafGridView, DataHandle, dim>(this->leafGridView(), level, dataHandle, iftype, dir);
648 else if (curCodim == dim - 1)
649 communicateUG_<LeafGridView, DataHandle, dim-1>(this->leafGridView(), level, dataHandle, iftype, dir);
650 else if (curCodim == 1)
651 communicateUG_<LeafGridView, DataHandle, 1>(this->leafGridView(), level, dataHandle, iftype, dir);
652 else
654 className(*this) << "::communicate(): Not "
655 "supported for dim " << dim << " and codim " << curCodim);
656 }
657#endif // ModelP
658 }
659
662 {
663 return ccobj_;
664 }
665
666 protected:
667#ifdef ModelP
668 template <class GridView, class DataHandle, int codim>
669 void communicateUG_(const GridView& gv, int level,
670 DataHandle &dataHandle,
671 InterfaceType iftype,
672 CommunicationDirection dir) const
673 {
674 typename UG_NS<dim>::DDD_IF_DIR ugIfDir;
675 // Translate the communication direction from Dune-Speak to UG-Speak
676 if (dir==ForwardCommunication)
677 ugIfDir = UG_NS<dim>::IF_FORWARD();
678 else
679 ugIfDir = UG_NS<dim>::IF_BACKWARD();
680
681 typedef UGMessageBuffer<DataHandle,dim,codim> UGMsgBuf;
682 UGMsgBuf::duneDataHandle_ = &dataHandle;
683
684 UGMsgBuf::level = level;
685
686 std::vector<typename UG_NS<dim>::DDD_IF> ugIfs;
687 findDDDInterfaces_(ugIfs, iftype, codim);
688
689 unsigned bufSize = UGMsgBuf::ugBufferSize_(gv);
690 if (!bufSize)
691 return; // we don't need to communicate if we don't have any data!
692 for (unsigned i=0; i < ugIfs.size(); ++i)
693 UG_NS<dim>::DDD_IFOneway(ugIfs[i],
694 ugIfDir,
695 bufSize,
696 &UGMsgBuf::ugGather_,
697 &UGMsgBuf::ugScatter_);
698 }
699
700 void findDDDInterfaces_(std::vector<typename UG_NS<dim>::DDD_IF > &dddIfaces,
701 InterfaceType iftype,
702 int codim) const
703 {
704 dddIfaces.clear();
705 if (codim == 0)
706 {
707 switch (iftype) {
709 // do not communicate anything: Elements can not be in
710 // the interior of two processes at the same time
711 return;
713 // The communicated elements are in the sender's
714 // interior and it does not matter what they are for
715 // the receiver
716 dddIfaces.push_back(UG_NS<dim>::ElementVHIF());
717 return;
718 case All_All_Interface :
719 // It does neither matter what the communicated
720 // elements are for sender nor for the receiver. If
721 // they are seen by these two processes, data is send
722 // and received.
723 dddIfaces.push_back(UG_NS<dim>::ElementSymmVHIF());
724 return;
725 default :
726 DUNE_THROW(GridError,
727 "Element communication not supported for "
728 "interfaces of type "
729 << iftype);
730 }
731 }
732 else if (codim == dim)
733 {
734 switch (iftype)
735 {
737 dddIfaces.push_back(UG_NS<dim>::BorderNodeSymmIF());
738 return;
740 dddIfaces.push_back(UG_NS<dim>::BorderNodeSymmIF());
741 dddIfaces.push_back(UG_NS<dim>::NodeIF());
742 return;
743 case All_All_Interface :
744 dddIfaces.push_back(UG_NS<dim>::NodeAllIF());
745 return;
746 default :
747 DUNE_THROW(GridError,
748 "Node communication not supported for "
749 "interfaces of type "
750 << iftype);
751 }
752 }
753 else if (codim == dim-1)
754 {
755 switch (iftype)
756 {
758 dddIfaces.push_back(UG_NS<dim>::BorderEdgeSymmIF());
759 return;
761 dddIfaces.push_back(UG_NS<dim>::BorderEdgeSymmIF());
762 // Is the following line needed or not?
763 // dddIfaces.push_back(UG_NS<dim>::EdgeIF());
764 return;
765 case All_All_Interface :
766 dddIfaces.push_back(UG_NS<dim>::EdgeSymmVHIF());
767 return;
768 default :
769 DUNE_THROW(GridError,
770 "Edge communication not supported for "
771 "interfaces of type "
772 << iftype);
773 }
774 }
775 else if (codim == 1)
776 {
777 switch (iftype)
778 {
781 dddIfaces.push_back(UG_NS<dim>::BorderVectorSymmIF());
782 return;
783 default :
784 DUNE_THROW(GridError,
785 "Face communication not supported for "
786 "interfaces of type "
787 << iftype);
788 }
789 }
790 else
791 {
792 DUNE_THROW(GridError,
793 "Communication for codim "
794 << codim
795 << " entities is not yet supported "
796 << " by the DUNE UGGrid interface!");
797 }
798 };
799#endif // ModelP
800
801 public:
802 // **********************************************************
803 // End of Interface Methods
804 // **********************************************************
805
813 void getChildrenOfSubface(const typename Traits::template Codim<0>::Entity & e,
814 int elementSide,
815 int maxl,
816 std::vector<typename Traits::template Codim<0>::Entity>& childElements,
817 std::vector<unsigned char>& childElementSides) const;
818
824 COPY
825 };
826
832 NONE
833 };
834
837 refinementType_ = type;
838 }
839
842 closureType_ = type;
843 }
844
851 static void setDefaultHeapSize(unsigned size) {
852 heapSize_ = size;
853 }
854
858 void setPosition(const typename Traits::template Codim<dim>::Entity& e,
859 const FieldVector<double, dim>& pos);
860
865 void globalRefine(int n);
866
871 void saveState(const std::string& filename) const;
872
877 void loadState(const std::string& filename);
878
879 private:
881 typename UG_NS<dim>::MultiGrid* multigrid_;
882
885
891 void setIndices(bool setLevelZero,
892 std::vector<unsigned int>* nodePermutation);
893
894 // Each UGGrid object has a unique name to identify it in the
895 // UG environment structure
896 std::string name_;
897
898 // Our set of level indices
899 std::vector<std::shared_ptr<UGGridLevelIndexSet<const UGGrid<dim> > > > levelIndexSets_;
900
901 UGGridLeafIndexSet<const UGGrid<dim> > leafIndexSet_;
902
903 // One id set implementation
904 // Used for both the local and the global UGGrid id sets
905 UGGridIdSet<const UGGrid<dim> > idSet_;
906
908 RefinementType refinementType_;
909
911 ClosureType closureType_;
912
920 static int numOfUGGrids;
921
927 bool someElementHasBeenMarkedForRefinement_;
928
934 bool someElementHasBeenMarkedForCoarsening_;
935
940 static unsigned int heapSize_;
941
943 std::vector<std::shared_ptr<BoundarySegment<dim> > > boundarySegments_;
944
950 unsigned int numBoundarySegments_;
951
952 }; // end Class UGGrid
953
954 namespace Capabilities
955 {
971 template<int dim>
972 struct hasEntity< UGGrid<dim>, 0>
973 {
974 static const bool v = true;
975 };
976
980 template<int dim>
981 struct hasEntity< UGGrid<dim>, dim>
982 {
983 static const bool v = true;
984 };
985
989 template<int dim>
991 {
992 static const bool v = true;
993 };
994
998 template<int dim>
1000 {
1001 static const bool v = false;
1002 };
1003
1004 }
1005
1006} // namespace Dune
1007
1008#endif // HAVE_UG || DOXYGEN
1009#endif // DUNE_UGGRID_HH
Base class for grid boundary segments of arbitrary geometry.
CollectiveCommunication()
Construct default object.
Definition: collectivecommunication.hh:82
Wrapper class for entities.
Definition: entity.hh:65
vector space out of a tensor product of fields.
Definition: fvector.hh:93
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:268
Definition: grid.hh:920
static std::conditional< std::is_reference< InterfaceType >::value, typenamestd::add_lvalue_reference< typenameReturnImplementationType< typenamestd::remove_reference< InterfaceType >::type >::ImplementationType >::type, typenamestd::remove_const< typenameReturnImplementationType< typenamestd::remove_reference< InterfaceType >::type >::ImplementationType >::type >::type getRealImplementation(InterfaceType &&i)
return real implementation of interface class
Definition: grid.hh:1115
Traits::LeafGridView leafGridView() const
View for the leaf grid for All_Partition.
Definition: grid.hh:940
bool loadBalance()
default implementation of load balance does nothing and returns false
Definition: grid.hh:1030
Traits::LevelGridView levelGridView(int level) const
View for a grid level for All_Partition.
Definition: grid.hh:932
Base class for exceptions in Dune grid modules.
Definition: exceptions.hh:18
Provide a generic factory class for unstructured grids.
Definition: gridfactory.hh:263
Grid view abstract base class.
Definition: gridview.hh:60
GridFamily::Traits::CollectiveCommunication CollectiveCommunication
A type that is a model of Dune::CollectiveCommunication. It provides a portable way for collective co...
Definition: grid.hh:519
GridFamily::Traits::LevelGridView LevelGridView
type of view for level grid
Definition: grid.hh:406
GridFamily::Traits::LeafGridView LeafGridView
type of view for leaf grid
Definition: grid.hh:404
Id Set Interface.
Definition: indexidset.hh:402
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
Default exception for dummy implementations.
Definition: exceptions.hh:261
Front-end for the grid manager of the finite element toolbox UG.
Definition: uggrid.hh:231
void communicate(DataHandle &dataHandle, InterfaceType iftype, CommunicationDirection dir) const
The communication interface for all codims on the leaf level.
Definition: uggrid.hh:635
void postAdapt()
Clean up refinement markers.
size_t numBoundarySegments() const
Return the number of boundary segments.
Definition: uggrid.hh:336
void setRefinementType(RefinementType type)
Sets the type of grid refinement.
Definition: uggrid.hh:836
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.
unsigned int overlapSize(int level, int codim) const
Size of the overlap on a given level.
Definition: uggrid.hh:469
bool loadBalance(int minlevel=0)
Distributes this grid over the available nodes in a distributed machine.
void communicate(DataHandle &dataHandle, InterfaceType iftype, CommunicationDirection dir, int level) const
The communication interface for all codims on a given level.
Definition: uggrid.hh:600
UGGrid()
Default constructor.
~UGGrid() noexcept(false)
Destructor.
const CollectiveCommunication< UGGrid > & comm() const
Definition: uggrid.hh:661
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:563
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:279
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.
int maxLevel() const
const Traits::LevelIndexSet & levelIndexSet(int level) const
Access to the LevelIndexSets.
Definition: uggrid.hh:355
void setClosureType(ClosureType type)
Sets the type of grid refinement closure.
Definition: uggrid.hh:841
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition: uggrid.hh:349
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition: uggrid.hh:324
RefinementType
The different forms of grid refinement that UG supports.
Definition: uggrid.hh:820
@ COPY
New level consists of the refined elements and the unrefined ones, too.
Definition: uggrid.hh:824
@ LOCAL
New level consists only of the refined elements and the closure.
Definition: uggrid.hh:822
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: uggrid.hh:330
void setPosition(const typename Traits::template Codim< dim >::Entity &e, const FieldVector< double, dim > &pos)
Sets a vertex to a new position.
unsigned int ghostSize(int codim) const
Size of the ghost cell layer on the leaf level.
Definition: uggrid.hh:464
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:285
unsigned int ghostSize(int level, int codim) const
Size of the ghost cell layer on a given level.
Definition: uggrid.hh:474
bool adapt()
Triggers the grid refinement process.
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition: uggrid.hh:343
void loadState(const std::string &filename)
Read entire grid hierarchy from disk.
unsigned int overlapSize(int codim) const
Size of the overlap on the leaf level.
Definition: uggrid.hh:459
bool loadBalance(DataHandle &dataHandle)
Distributes the grid and some data over the available nodes in a distributed machine.
Definition: uggrid.hh:486
Traits::template Codim< Seed::codimension >::Entity entity(const Seed &seed) const
Create an Entity from an EntitySeed.
Definition: uggrid.hh:307
unsigned int Rank
The type used for process ranks.
Definition: uggrid.hh:288
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:828
@ GREEN
Standard red/green refinement.
Definition: uggrid.hh:830
@ NONE
No closure, results in nonconforming meshes.
Definition: uggrid.hh:832
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition: uggrid.hh:363
int size(int codim) const
number of leaf entities per codim in this process
Definition: uggrid.hh:318
static void setDefaultHeapSize(unsigned size)
Sets the default heap size.
Definition: uggrid.hh:851
A free function to provide the demangled class name of a given object or type as a string.
Implements an utility class that provides collective communication methods for sequential programs.
A set of traits classes to store static information about grid implementation.
Different resources needed by all grid implementations.
Definition of the DUNE_DEPRECATED macro for the case that config.h is not available.
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
@ 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: alignment.hh:11
std::string className()
Provide the demangled class name of a type T as a string.
Definition: classname.hh:23
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:87
Specialize with 'true' if implementation guarantees conforming level grids. (default=false)
Definition: capabilities.hh:78
Static tag representing a codimension.
Definition: dimension.hh:22
A traits struct that collects all associated types of one grid model.
Definition: grid.hh:1153
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 12, 23:30, 2024)