Dune Core Modules (2.3.1)

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
16
20
21#if HAVE_UG || DOXYGEN
22
23#ifdef ModelP
25#endif
26
27/* The following lines including the necessary UG headers are somewhat
28 tricky. Here's what's happening:
29 UG can support two- and three-dimensional grids. You choose be setting
30 either _2 oder _3 while compiling. This changes all sorts of stuff, in
31 particular data structures in the headers.
32 UG was never supposed to provide 2d and 3d grids at the same time.
33 However, when compiling it as c++, the dimension-dependent parts are
34 wrapped up cleanly in the namespaces UG::D2 and UG::D3, respectively. That
35 way it is possible to link together the UG lib for 2d and the one for 3d.
36 But we also need the headers twice! Once with _2 set and once with _3!
37 So here we go:*/
38
39/* The following define tells the UG headers that we want access to a few
40 special fields, for example the extra index fields in the element data structures. */
41#define FOR_DUNE
42
43// Set UG's space-dimension flag to 2d
44#define _2
45// And include all necessary UG headers
46#include "uggrid/ugincludes.hh"
47
48// Wrap a few large UG macros by functions before they get undef'ed away.
49// Here: The 2d-version of the macros
50#define UG_DIM 2
51#include "uggrid/ugwrapper.hh"
52#undef UG_DIM
53
54// UG defines a whole load of preprocessor macros. ug_undefs.hh undefines
55// them all, so we don't get name clashes.
56#include "uggrid/ug_undefs.hh"
57#undef _2
58
59/* Now we're done with 2d, and we can do the whole thing over again for 3d */
60
61/* All macros set by UG have been unset. This includes the macros that ensure
62 single inclusion of headers. We can thus include them again. However, we
63 only want to include those headers again that contain dimension-dependent stuff.
64 Therefore, we set a few single-inclusion defines manually before including
65 ugincludes.hh again.
66 */
67#define UGTYPES_H
68#define __HEAPS__
69#define __UGENV__
70#define __DEVICESH__
71
72#define _3
73#include "uggrid/ugincludes.hh"
74
75// Wrap a few large UG macros by functions before they get undef'ed away.
76// This time it's the 3d-versions.
77#define UG_DIM 3
78#include "uggrid/ugwrapper.hh"
79#undef UG_DIM
80
81// undef all macros defined by UG
82#include "uggrid/ug_undefs.hh"
83
84#undef _3
85#undef FOR_DUNE
86
87// The components of the UGGrid interface
88#include "uggrid/uggridgeometry.hh"
89#include "uggrid/uggridlocalgeometry.hh"
90#include "uggrid/uggridentity.hh"
91#include "uggrid/uggridentitypointer.hh"
92#include "uggrid/uggridentityseed.hh"
93#include "uggrid/uggridintersections.hh"
94#include "uggrid/uggridintersectioniterators.hh"
95#include "uggrid/uggridleveliterator.hh"
96#include "uggrid/uggridleafiterator.hh"
97#include "uggrid/uggridhieriterator.hh"
98#include "uggrid/uggridindexsets.hh"
99#ifdef ModelP
100#include "uggrid/ugmessagebuffer.hh"
101#include "uggrid/uglbgatherscatter.hh"
102#endif
103
104// Not needed here, but included for user convenience
106
107#ifdef ModelP
108template <class DataHandle, int GridDim, int codim>
109DataHandle *Dune::UGMessageBufferBase<DataHandle,GridDim,codim>::duneDataHandle_ = 0;
110
111template <class DataHandle, int GridDim, int codim>
112int Dune::UGMessageBufferBase<DataHandle,GridDim,codim>::level = -1;
113#endif // ModelP
114
115namespace Dune {
116
117#ifdef ModelP
118 template <int dim>
119 class CollectiveCommunication<Dune::UGGrid<dim> > :
120 public CollectiveCommunication< Dune::MPIHelper::MPICommunicator >
121 {
122 typedef CollectiveCommunication< Dune::MPIHelper::MPICommunicator > ParentType;
123 public:
125 : ParentType(MPIHelper::getCommunicator())
126 {}
127 };
128#endif // ModelP
129
130 template<int dim>
131 struct UGGridFamily
132 {
133 typedef GridTraits<dim,dim,Dune::UGGrid<dim>,
134 UGGridGeometry,
135 UGGridEntity,
136 UGGridEntityPointer,
137 UGGridLevelIterator,
138 UGGridLeafIntersection,
139 UGGridLevelIntersection,
140 UGGridLeafIntersectionIterator,
141 UGGridLevelIntersectionIterator,
142 UGGridHierarchicIterator,
143 UGGridLeafIterator,
144 UGGridLevelIndexSet< const UGGrid<dim> >,
145 UGGridLeafIndexSet< const UGGrid<dim> >,
146 UGGridIdSet< const UGGrid<dim> >,
147 typename UG_NS<dim>::UG_ID_TYPE,
148 UGGridIdSet< const UGGrid<dim> >,
149 typename UG_NS<dim>::UG_ID_TYPE,
150 CollectiveCommunication<Dune::UGGrid<dim> >,
151 DefaultLevelGridViewTraits,
152 DefaultLeafGridViewTraits,
153 UGGridEntitySeed,
154 UGGridLocalGeometry>
155 Traits;
156 };
157
158
159 //**********************************************************************
160 //
161 // --UGGrid
162 //
163 //**********************************************************************
164
202 template <int dim>
203 class UGGrid : public GridDefaultImplementation <dim, dim, double, UGGridFamily<dim> >
204 {
206
207 friend class UGGridGeometry<0,dim,const UGGrid<dim> >;
208 friend class UGGridGeometry<dim,dim,const UGGrid<dim> >;
209 friend class UGGridGeometry<1,2,const UGGrid<dim> >;
210 friend class UGGridGeometry<2,3,const UGGrid<dim> >;
211
212 friend class UGGridEntity <0,dim,const UGGrid<dim> >;
213 friend class UGGridEntity <dim,dim,const UGGrid<dim> >;
214 friend class UGGridHierarchicIterator<const UGGrid<dim> >;
215 friend class UGGridLeafIntersection<const UGGrid<dim> >;
216 friend class UGGridLevelIntersection<const UGGrid<dim> >;
217 friend class UGGridLeafIntersectionIterator<const UGGrid<dim> >;
218 friend class UGGridLevelIntersectionIterator<const UGGrid<dim> >;
219
220 friend class UGGridLevelIndexSet<const UGGrid<dim> >;
221 friend class UGGridLeafIndexSet<const UGGrid<dim> >;
222 friend class UGGridIdSet<const UGGrid<dim> >;
223
224 friend class GridFactory<UGGrid<dim> >;
225
226#ifdef ModelP
227 friend class UGLBGatherScatter;
228#endif
229
230 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
231 friend class UGGridLeafIterator;
232 template <int codim_, PartitionIteratorType PiType_, class GridImp_>
233 friend class UGGridLevelIterator;
234 template <int codim_, class GridImp_>
235 friend class UGGridEntityPointer;
236
238 dune_static_assert(dim==2 || dim==3, "Use UGGrid only for 2d and 3d!");
239
240 // The different instantiations are mutual friends so they can access
241 // each others numOfUGGrids field
242 friend class UGGrid<2>;
243 friend class UGGrid<3>;
244
245 //**********************************************************
246 // The Interface Methods
247 //**********************************************************
248 public:
250 typedef UGGridFamily<dim> GridFamily;
251
252 // the Traits
253 typedef typename UGGridFamily<dim>::Traits Traits;
254
256 typedef UG::DOUBLE ctype;
257
264
267
270 int maxLevel() const;
271
273 template<int codim>
274 typename Traits::template Codim<codim>::LevelIterator lbegin (int level) const;
275
277 template<int codim>
278 typename Traits::template Codim<codim>::LevelIterator lend (int level) const;
279
281 template<int codim, PartitionIteratorType PiType>
282 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lbegin (int level) const;
283
285 template<int codim, PartitionIteratorType PiType>
286 typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lend (int level) const;
287
289 template<int codim>
290 typename Traits::template Codim<codim>::LeafIterator leafbegin() const {
291 return typename Traits::template Codim<codim>::template Partition<All_Partition>::LeafIterator(*this);
292 }
293
295 template<int codim>
296 typename Traits::template Codim<codim>::LeafIterator leafend() const {
297 return UGGridLeafIterator<codim,All_Partition, const UGGrid<dim> >();
298 }
299
301 template<int codim, PartitionIteratorType PiType>
302 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafbegin() const {
303 return typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator(*this);
304 }
305
307 template<int codim, PartitionIteratorType PiType>
308 typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafend() const {
309 return UGGridLeafIterator<codim,PiType, const UGGrid<dim> >();
310 }
311
313 template <typename Seed>
314 typename Traits::template Codim<Seed::codimension>::EntityPointer
315 entityPointer(const Seed& seed) const
316 {
317 enum {codim = Seed::codimension};
318 return typename Traits::template Codim<codim>::EntityPointer(UGGridEntityPointer<codim,const UGGrid<dim> >(this->getRealImplementation(seed).target(),this));
319 }
320
323 int size (int level, int codim) const;
324
326 int size (int codim) const
327 {
328 return leafIndexSet().size(codim);
329 }
330
332 int size (int level, GeometryType type) const
333 {
334 return this->levelIndexSet(level).size(type);
335 }
336
338 int size (GeometryType type) const
339 {
340 return this->leafIndexSet().size(type);
341 }
342
344 size_t numBoundarySegments() const {
345 // The number is stored as a member of UGGrid upon grid creation.
346 // The corresponding data structure is not exported by UG. (It is in ug/dom/std/std_internal.h)
347 return numBoundarySegments_;
348 }
349
351 const typename Traits::GlobalIdSet& globalIdSet() const
352 {
353 return idSet_;
354 }
355
357 const typename Traits::LocalIdSet& localIdSet() const
358 {
359 return idSet_;
360 }
361
363 const typename Traits::LevelIndexSet& levelIndexSet(int level) const
364 {
365 if (level<0 || level>maxLevel())
366 DUNE_THROW(GridError, "levelIndexSet of nonexisting level " << level << " requested!");
367 return *levelIndexSets_[level];
368 }
369
371 const typename Traits::LeafIndexSet& leafIndexSet() const
372 {
373 return leafIndexSet_;
374 }
375
378
391 bool mark(int refCount, const typename Traits::template Codim<0>::Entity & e );
392
449 bool mark(const typename Traits::template Codim<0>::Entity & e,
450 typename UG_NS<dim>::RefinementRule rule,
451 int side=0);
452
454 int getMark(const typename Traits::template Codim<0>::Entity& e) const;
455
458 bool preAdapt();
459
461 bool adapt();
462
464 void postAdapt();
468 unsigned int overlapSize(int codim) const {
469 return 0;
470 }
471
473 unsigned int ghostSize(int codim) const {
474 return (codim==0) ? 1 : 0;
475 }
476
478 unsigned int overlapSize(int level, int codim) const {
479 return 0;
480 }
481
483 unsigned int ghostSize(int level, int codim) const {
484 return (codim==0) ? 1 : 0;
485 }
486
491 bool loadBalance() {
492 return loadBalance(0,0,2,32,1);
493 }
494
502 template<class DataHandle>
503 bool loadBalance (DataHandle& dataHandle)
504 {
505#if !HAVE_UG_PATCH10
506 DUNE_THROW(NotImplemented, "load balancing with data attached");
507#else
508#ifdef ModelP
509 // gather element data
510 // UGLBGatherScatter::template gather<0>(this->leafView(), dataHandle);
511
512 // gather node data
513 UGLBGatherScatter::template gather<dim>(this->leafView(), dataHandle);
514#endif
515
516 // the load balancing step now also attaches
517 // the data to the entities and distributes it
518 loadBalance();
519
520#ifdef ModelP
521 // scatter element data
522 // UGLBGatherScatter::template scatter<0>(this->leafView(), dataHandle);
523
524 // scatter node data
525 UGLBGatherScatter::template scatter<dim>(this->leafView(), dataHandle);
526#endif
527
528 return true;
529#endif // HAVE_UG_PATCH10
530 }
531
548 bool loadBalance(int strategy, int minlevel, int depth, int maxlevel, int minelement);
549
560 template<class DataHandle>
561 void communicate (DataHandle& dataHandle, InterfaceType iftype, CommunicationDirection dir, int level) const
562 {
563#ifdef ModelP
564 typedef typename UGGrid::LevelGridView LevelGridView;
565
566 for (int curCodim = 0; curCodim <= dim; ++curCodim) {
567 if (!dataHandle.contains(dim, curCodim))
568 continue;
569
570 if (curCodim == 0)
571 communicateUG_<LevelGridView, DataHandle, 0>(this->levelGridView(level), level, dataHandle, iftype, dir);
572 else if (curCodim == dim)
573 communicateUG_<LevelGridView, DataHandle, dim>(this->levelGridView(level), level, dataHandle, iftype, dir);
574 else if (curCodim == dim - 1)
575 communicateUG_<LevelGridView, DataHandle, dim-1>(this->levelGridView(level), level, dataHandle, iftype, dir);
576 else if (curCodim == 1)
577 communicateUG_<LevelGridView, DataHandle, 1>(this->levelGridView(level), level, dataHandle, iftype, dir);
578 else
580 className(*this) << "::communicate(): Not "
581 "supported for dim " << dim << " and codim " << curCodim);
582 }
583#endif // ModelP
584 }
585
595 template<class DataHandle>
596 void communicate(DataHandle& dataHandle, InterfaceType iftype, CommunicationDirection dir) const
597 {
598#ifdef ModelP
599 typedef typename UGGrid::LeafGridView LeafGridView;
600
601 for (int curCodim = 0; curCodim <= dim; ++curCodim) {
602 if (!dataHandle.contains(dim, curCodim))
603 continue;
604 int level = -1;
605 if (curCodim == 0)
606 communicateUG_<LeafGridView, DataHandle, 0>(this->leafView(), level, dataHandle, iftype, dir);
607 else if (curCodim == dim)
608 communicateUG_<LeafGridView, DataHandle, dim>(this->leafView(), level, dataHandle, iftype, dir);
609 else if (curCodim == dim - 1)
610 communicateUG_<LeafGridView, DataHandle, dim-1>(this->leafView(), level, dataHandle, iftype, dir);
611 else if (curCodim == 1)
612 communicateUG_<LeafGridView, DataHandle, 1>(this->leafView(), level, dataHandle, iftype, dir);
613 else
615 className(*this) << "::communicate(): Not "
616 "supported for dim " << dim << " and codim " << curCodim);
617 }
618#endif // ModelP
619 }
620
623 {
624 return ccobj_;
625 }
626
627 protected:
628#ifdef ModelP
629 template <class GridView, class DataHandle, int codim>
630 void communicateUG_(const GridView& gv, int level,
631 DataHandle &dataHandle,
632 InterfaceType iftype,
633 CommunicationDirection dir) const
634 {
635 typename UG_NS<dim>::DDD_IF_DIR ugIfDir;
636 // Translate the communication direction from Dune-Speak to UG-Speak
637 if (dir==ForwardCommunication)
638 ugIfDir = UG_NS<dim>::IF_FORWARD();
639 else
640 ugIfDir = UG_NS<dim>::IF_BACKWARD();
641
642 typedef UGMessageBuffer<DataHandle,dim,codim> UGMsgBuf;
643 UGMsgBuf::duneDataHandle_ = &dataHandle;
644
645 UGMsgBuf::level = level;
646
647 std::vector<typename UG_NS<dim>::DDD_IF> ugIfs;
648 findDDDInterfaces_(ugIfs, iftype, codim);
649
650 unsigned bufSize = UGMsgBuf::ugBufferSize_(gv);
651 if (!bufSize)
652 return; // we don't need to communicate if we don't have any data!
653 for (unsigned i=0; i < ugIfs.size(); ++i)
654 UG_NS<dim>::DDD_IFOneway(ugIfs[i],
655 ugIfDir,
656 bufSize,
657 &UGMsgBuf::ugGather_,
658 &UGMsgBuf::ugScatter_);
659 }
660
661 void findDDDInterfaces_(std::vector<typename UG_NS<dim>::DDD_IF > &dddIfaces,
662 InterfaceType iftype,
663 int codim) const
664 {
665 dddIfaces.clear();
666 if (codim == 0)
667 {
668 switch (iftype) {
670 // do not communicate anything: Elements can not be in
671 // the interior of two processes at the same time
672 return;
674 // The communicated elements are in the sender's
675 // interior and it does not matter what they are for
676 // the receiver
677 dddIfaces.push_back(UG_NS<dim>::ElementVHIF());
678 return;
679 case All_All_Interface :
680 // It does neither matter what the communicated
681 // elements are for sender nor for the receiver. If
682 // they are seen by these two processes, data is send
683 // and received.
684 dddIfaces.push_back(UG_NS<dim>::ElementSymmVHIF());
685 return;
686 default :
687 DUNE_THROW(GridError,
688 "Element communication not supported for "
689 "interfaces of type "
690 << iftype);
691 }
692 }
693 else if (codim == dim)
694 {
695 switch (iftype)
696 {
698 dddIfaces.push_back(UG_NS<dim>::BorderNodeSymmIF());
699 return;
701 dddIfaces.push_back(UG_NS<dim>::BorderNodeSymmIF());
702 dddIfaces.push_back(UG_NS<dim>::NodeIF());
703 return;
704 case All_All_Interface :
705 dddIfaces.push_back(UG_NS<dim>::NodeAllIF());
706 return;
707 default :
708 DUNE_THROW(GridError,
709 "Node communication not supported for "
710 "interfaces of type "
711 << iftype);
712 }
713 }
714 else if (codim == dim-1)
715 {
716 switch (iftype)
717 {
719 dddIfaces.push_back(UG_NS<dim>::BorderEdgeSymmIF());
720 return;
722 dddIfaces.push_back(UG_NS<dim>::BorderEdgeSymmIF());
723 // Is the following line needed or not?
724 // dddIfaces.push_back(UG_NS<dim>::EdgeIF());
725 return;
726 default :
727 DUNE_THROW(GridError,
728 "Edge communication not supported for "
729 "interfaces of type "
730 << iftype);
731 }
732 }
733 else if (codim == 1)
734 {
735 switch (iftype)
736 {
739 dddIfaces.push_back(UG_NS<dim>::BorderVectorSymmIF());
740 return;
741 default :
742 DUNE_THROW(GridError,
743 "Face communication not supported for "
744 "interfaces of type "
745 << iftype);
746 }
747 }
748 else
749 {
750 DUNE_THROW(GridError,
751 "Communication for codim "
752 << codim
753 << " entities is not yet supported "
754 << " by the DUNE UGGrid interface!");
755 }
756 };
757#endif // ModelP
758
759 public:
760 // **********************************************************
761 // End of Interface Methods
762 // **********************************************************
763
771 void getChildrenOfSubface(const typename Traits::template Codim<0>::EntityPointer & e,
772 int elementSide,
773 int maxl,
774 std::vector<typename Traits::template Codim<0>::EntityPointer>& childElements,
775 std::vector<unsigned char>& childElementSides) const;
776
782 COPY
783 };
784
790 NONE
791 };
792
795 refinementType_ = type;
796 }
797
800 closureType_ = type;
801 }
802
809 static void setDefaultHeapSize(unsigned size) {
810 heapSize_ = size;
811 }
812
816 void setPosition(const typename Traits::template Codim<dim>::EntityPointer& e,
817 const FieldVector<double, dim>& pos);
818
823 void globalRefine(int n);
824
829 void saveState(const std::string& filename) const;
830
835 void loadState(const std::string& filename);
836
837 private:
839 typename UG_NS<dim>::MultiGrid* multigrid_;
840
843
849 void setIndices(bool setLevelZero,
850 std::vector<unsigned int>* nodePermutation);
851
852 // Each UGGrid object has a unique name to identify it in the
853 // UG environment structure
854 std::string name_;
855
856 // Our set of level indices
857 std::vector<shared_ptr<UGGridLevelIndexSet<const UGGrid<dim> > > > levelIndexSets_;
858
859 UGGridLeafIndexSet<const UGGrid<dim> > leafIndexSet_;
860
861 // One id set implementation
862 // Used for both the local and the global UGGrid id sets
863 UGGridIdSet<const UGGrid<dim> > idSet_;
864
866 RefinementType refinementType_;
867
869 ClosureType closureType_;
870
878 static int numOfUGGrids;
879
885 bool someElementHasBeenMarkedForRefinement_;
886
892 bool someElementHasBeenMarkedForCoarsening_;
893
898 static unsigned int heapSize_;
899
901 std::vector<shared_ptr<BoundarySegment<dim> > > boundarySegments_;
902
908 unsigned int numBoundarySegments_;
909
910 }; // end Class UGGrid
911
912 namespace Capabilities
913 {
929 template<int dim>
930 struct hasEntity< UGGrid<dim>, 0>
931 {
932 static const bool v = true;
933 };
934
938 template<int dim>
939 struct hasEntity< UGGrid<dim>, dim>
940 {
941 static const bool v = true;
942 };
943
947 template<int dim>
948 struct isParallel< UGGrid<dim> >
949 {
950#ifdef ModelP
951 static const bool v = true;
952#else // !ModelP
953 static const bool v = false;
954#endif // !ModelP
955 };
956
960 template<int dim>
962 {
963 static const bool v = true;
964 };
965
969 template<int dim>
971 {
972 static const bool v = false;
973 };
974
975 }
976
977} // namespace Dune
978
979#endif // HAVE_UG || DOXYGEN
980#endif // DUNE_UGGRID_HH
Base class for grid boundary segments of arbitrary geometry.
CollectiveCommunication()
Construct default object.
Definition: collectivecommunication.hh:75
vector space out of a tensor product of fields.
Definition: fvector.hh:92
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:25
Definition: grid.hh:1017
Traits::template Partition< All_Partition >::LeafGridView leafView() const DUNE_DEPRECATED_MSG("The method leafView has been renamed to leafGridView.")
View for the leaf grid for All_Partition.
Definition: grid.hh:1055
static ReturnImplementationType< InterfaceType >::ImplementationType & getRealImplementation(InterfaceType &i)
return real implementation of interface class
Definition: grid.hh:1223
Traits::template Partition< pitype >::LevelGridView levelGridView(int level) const
View for a grid level.
Definition: grid.hh:1064
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:57
GridFamily::Traits::CollectiveCommunication CollectiveCommunication
A type that is a model of Dune::CollectiveCommunication. It provides a portable way for collective co...
Definition: grid.hh:543
Partition< All_Partition >::LevelGridView LevelGridView
View types for All_Partition.
Definition: grid.hh:426
Id Set Interface.
Definition: indexidset.hh:403
Index Set Interface base class.
Definition: indexidset.hh:78
IndexType size(GeometryType type) const
Return total number of entities of given geometry type in entity set .
Definition: indexidset.hh:204
Default exception for dummy implementations.
Definition: exceptions.hh:289
Front-end for the grid manager of the finite element toolbox UG.
Definition: uggrid.hh:204
Traits::template Codim< codim >::LevelIterator lend(int level) const
one past the end on this level
void communicate(DataHandle &dataHandle, InterfaceType iftype, CommunicationDirection dir) const
The communication interface for all codims on the leaf level.
Definition: uggrid.hh:596
~UGGrid()
Destructor.
void postAdapt()
Clean up refinement markers.
size_t numBoundarySegments() const
Return the number of boundary segments.
Definition: uggrid.hh:344
void setRefinementType(RefinementType type)
Sets the type of grid refinement.
Definition: uggrid.hh:794
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:478
void communicate(DataHandle &dataHandle, InterfaceType iftype, CommunicationDirection dir, int level) const
The communication interface for all codims on a given level.
Definition: uggrid.hh:561
UGGrid()
Default constructor.
void setPosition(const typename Traits::template Codim< dim >::EntityPointer &e, const FieldVector< double, dim > &pos)
Sets a vertex to a new position.
Traits::template Codim< Seed::codimension >::EntityPointer entityPointer(const Seed &seed) const
Create an EntityPointer from an EntitySeed.
Definition: uggrid.hh:315
const CollectiveCommunication< UGGrid > & comm() const
Definition: uggrid.hh:622
UGGridFamily< dim > GridFamily
type of the used GridFamily for this grid
Definition: uggrid.hh:250
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.
void getChildrenOfSubface(const typename Traits::template Codim< 0 >::EntityPointer &e, int elementSide, int maxl, std::vector< typename Traits::template Codim< 0 >::EntityPointer > &childElements, std::vector< unsigned char > &childElementSides) const
Rudimentary substitute for a hierarchic iterator on faces.
int maxLevel() const
const Traits::LevelIndexSet & levelIndexSet(int level) const
Access to the LevelIndexSets.
Definition: uggrid.hh:363
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafend() const
one past the end of the sequence of leaf entities
Definition: uggrid.hh:308
void setClosureType(ClosureType type)
Sets the type of grid refinement closure.
Definition: uggrid.hh:799
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator lend(int level) const
one past the end on this level
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition: uggrid.hh:357
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition: uggrid.hh:332
bool loadBalance(int strategy, int minlevel, int depth, int maxlevel, int minelement)
Distributes this grid over the available nodes in a distributed machine.
RefinementType
The different forms of grid refinement that UG supports.
Definition: uggrid.hh:778
@ COPY
New level consists of the refined elements and the unrefined ones, too.
Definition: uggrid.hh:782
@ LOCAL
New level consists only of the refined elements and the closure.
Definition: uggrid.hh:780
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: uggrid.hh:338
unsigned int ghostSize(int codim) const
Size of the ghost cell layer on the leaf level.
Definition: uggrid.hh:473
Traits::template Codim< codim >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim.
Definition: uggrid.hh:290
Traits::template Codim< codim >::LeafIterator leafend() const
one past the end of the sequence of leaf entities
Definition: uggrid.hh:296
int size(int level, int codim) const
Number of grid entities per level and codim.
bool loadBalance()
Default load balancing.
Definition: uggrid.hh:491
UG::DOUBLE ctype
The type used to store coordinates.
Definition: uggrid.hh:256
unsigned int ghostSize(int level, int codim) const
Size of the ghost cell layer on a given level.
Definition: uggrid.hh:483
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator lbegin(int level) const
Iterator to first entity of given codim on level.
bool adapt()
Triggers the grid refinement process.
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition: uggrid.hh:351
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:468
bool loadBalance(DataHandle &dataHandle)
Distributes the grid and some data over the available nodes in a distributed machine.
Definition: uggrid.hh:503
Traits::template Codim< codim >::LevelIterator lbegin(int level) const
Iterator to first entity of given codim on level.
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:786
@ GREEN
Standard red/green refinement.
Definition: uggrid.hh:788
@ NONE
No closure, results in nonconforming meshes.
Definition: uggrid.hh:790
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim.
Definition: uggrid.hh:302
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition: uggrid.hh:371
int size(int codim) const
number of leaf entities per codim in this process
Definition: uggrid.hh:326
static void setDefaultHeapSize(unsigned size)
Sets the default heap size.
Definition: uggrid.hh:809
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.
A few common exception classes.
#define dune_static_assert(COND, MSG)
Helper template so that compilation fails if condition is not true.
Definition: static_assert.hh:79
#define DUNE_THROW(E, m)
Definition: exceptions.hh:244
Dune namespace.
Definition: alignment.hh:14
std::string className(T &t)
Provide the demangled class name of a given object as a string.
Definition: classname.hh:23
CommunicationDirection
Define a type for communication direction parameter.
Definition: gridenums.hh:164
@ ForwardCommunication
communicate as given in InterfaceType
Definition: gridenums.hh:165
InterfaceType
Parameter to be used for the communication functions.
Definition: gridenums.hh:80
@ InteriorBorder_All_Interface
send interior and border, receive all entities
Definition: gridenums.hh:82
@ All_All_Interface
send all and receive all entities
Definition: gridenums.hh:85
@ InteriorBorder_InteriorBorder_Interface
send/receive interior and border entities
Definition: gridenums.hh:81
Implements an utility class that provides collective communication methods for sequential programs.
Implements an utility class that provides MPI's collective communication methods.
Helpers for dealing with MPI.
Fallback implementation of the C++0x static_assert feature.
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:96
Specialize with 'true' if implementation guarantees conforming level grids. (default=false)
Definition: capabilities.hh:87
Specialize with 'true' if implementation supports parallelism. (default=false)
Definition: capabilities.hh:65
A traits struct that collects all associated types of one grid model.
Definition: grid.hh:1261
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)