Dune Core Modules (2.4.2)

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 <dune/common/classname.hh>
16 
19 #include <dune/grid/common/grid.hh>
20 
21 #if HAVE_UG || DOXYGEN
22 
23 #ifdef ModelP
25 #endif
26 
27 /* [Before reading the following: the macros UG_DIM_2 and UG_DIM_3 where named
28  * _2 and _3, respectively, up until ug-3.12.0.]
29  *
30  * The following lines including the necessary UG headers are somewhat
31  tricky. Here's what's happening:
32  UG can support two- and three-dimensional grids. You choose by setting
33  either UG_DIM_2 or UG_DIM_3 while compiling. This changes all sorts of stuff, in
34  particular data structures in the headers.
35  UG was never supposed to provide 2d and 3d grids at the same time.
36  However, when compiling it as c++, the dimension-dependent parts are
37  wrapped up cleanly in the namespaces UG::D2 and UG::D3, respectively. That
38  way it is possible to link together the UG lib for 2d and the one for 3d.
39  But we also need the headers twice! Once with UG_DIM_2 set and once with UG_DIM_3!
40  So here we go:*/
41 
42 /* The following define tells the UG headers that we want access to a few
43  special fields, for example the extra index fields in the element data structures. */
44 #define FOR_DUNE
45 
46 // Set UG's space-dimension flag to 2d
47 #ifdef UG_USE_NEW_DIMENSION_DEFINES
48 #define UG_DIM_2
49 #else
50 #define _2
51 #endif
52 // And include all necessary UG headers
53 #include "uggrid/ugincludes.hh"
54 
55 // Wrap a few large UG macros by functions before they get undef'ed away.
56 // Here: The 2d-version of the macros
57 #define UG_DIM 2
58 #include "uggrid/ugwrapper.hh"
59 #undef UG_DIM
60 
61 // UG defines a whole load of preprocessor macros. ug_undefs.hh undefines
62 // them all, so we don't get name clashes.
63 #include "uggrid/ug_undefs.hh"
64 #ifdef UG_USE_NEW_DIMENSION_DEFINES
65 #undef UG_DIM_2
66 #else
67 #undef _2
68 #endif
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 #ifdef UG_USE_NEW_DIMENSION_DEFINES
87 #define UG_DIM_3
88 #else
89 #define _3
90 #endif
91 
92 #include "uggrid/ugincludes.hh"
93 
94 // Wrap a few large UG macros by functions before they get undef'ed away.
95 // This time it's the 3d-versions.
96 #define UG_DIM 3
97 #include "uggrid/ugwrapper.hh"
98 #undef UG_DIM
99 
100 // undef all macros defined by UG
101 #include "uggrid/ug_undefs.hh"
102 
103 #ifdef UG_USE_NEW_DIMENSION_DEFINES
104 #undef UG_DIM_3
105 #else
106 #undef _3
107 #endif
108 #undef FOR_DUNE
109 
110 // The components of the UGGrid interface
111 #include "uggrid/uggridgeometry.hh"
112 #include "uggrid/uggridlocalgeometry.hh"
113 #include "uggrid/uggridentity.hh"
114 #include "uggrid/uggridentitypointer.hh"
115 #include "uggrid/uggridentityseed.hh"
116 #include "uggrid/uggridintersections.hh"
117 #include "uggrid/uggridintersectioniterators.hh"
118 #include "uggrid/uggridleveliterator.hh"
119 #include "uggrid/uggridleafiterator.hh"
120 #include "uggrid/uggridhieriterator.hh"
121 #include "uggrid/uggridindexsets.hh"
122 #include <dune/grid/uggrid/uggridviews.hh>
123 #ifdef ModelP
124 #include "uggrid/ugmessagebuffer.hh"
125 #include "uggrid/uglbgatherscatter.hh"
126 #endif
127 
128 // Not needed here, but included for user convenience
129 #include "uggrid/uggridfactory.hh"
130 
131 #ifdef ModelP
132 template <class DataHandle, int GridDim, int codim>
133 DataHandle *Dune::UGMessageBufferBase<DataHandle,GridDim,codim>::duneDataHandle_ = 0;
134 
135 template <class DataHandle, int GridDim, int codim>
136 int Dune::UGMessageBufferBase<DataHandle,GridDim,codim>::level = -1;
137 #endif // ModelP
138 
139 namespace Dune {
140 
141 #ifdef ModelP
142  template <int dim>
143  class CollectiveCommunication<Dune::UGGrid<dim> > :
144  public CollectiveCommunication< Dune::MPIHelper::MPICommunicator >
145  {
146  typedef CollectiveCommunication< Dune::MPIHelper::MPICommunicator > ParentType;
147  public:
149  : ParentType(MPIHelper::getCommunicator())
150  {}
151  };
152 #endif // ModelP
153 
154  template<int dim>
155  struct UGGridFamily
156  {
157  typedef GridTraits<dim,dim,Dune::UGGrid<dim>,
158  UGGridGeometry,
159  UGGridEntity,
160  UGGridEntityPointer,
161  UGGridLevelIterator,
162  UGGridLeafIntersection,
163  UGGridLevelIntersection,
164  UGGridLeafIntersectionIterator,
165  UGGridLevelIntersectionIterator,
166  UGGridHierarchicIterator,
167  UGGridLeafIterator,
168  UGGridLevelIndexSet< const UGGrid<dim> >,
169  UGGridLeafIndexSet< const UGGrid<dim> >,
170  UGGridIdSet< const UGGrid<dim> >,
171  typename UG_NS<dim>::UG_ID_TYPE,
172  UGGridIdSet< const UGGrid<dim> >,
173  typename UG_NS<dim>::UG_ID_TYPE,
174  CollectiveCommunication<Dune::UGGrid<dim> >,
175  UGGridLevelGridViewTraits,
176  UGGridLeafGridViewTraits,
177  UGGridEntitySeed,
178  UGGridLocalGeometry>
179  Traits;
180  };
181 
182 
183  //**********************************************************************
184  //
185  // --UGGrid
186  //
187  //**********************************************************************
188 
226  template <int dim>
227  class UGGrid : public GridDefaultImplementation <dim, dim, double, UGGridFamily<dim> >
228  {
230 
231  friend class UGGridGeometry<0,dim,const UGGrid<dim> >;
232  friend class UGGridGeometry<dim,dim,const UGGrid<dim> >;
233  friend class UGGridGeometry<1,2,const UGGrid<dim> >;
234  friend class UGGridGeometry<2,3,const UGGrid<dim> >;
235 
236  friend class UGGridEntity <0,dim,const UGGrid<dim> >;
237  friend class UGGridEntity <dim,dim,const UGGrid<dim> >;
238  friend class UGGridHierarchicIterator<const UGGrid<dim> >;
239  friend class UGGridLeafIntersection<const UGGrid<dim> >;
240  friend class UGGridLevelIntersection<const UGGrid<dim> >;
241  friend class UGGridLeafIntersectionIterator<const UGGrid<dim> >;
242  friend class UGGridLevelIntersectionIterator<const UGGrid<dim> >;
243 
244  friend class UGGridLevelIndexSet<const UGGrid<dim> >;
245  friend class UGGridLeafIndexSet<const UGGrid<dim> >;
246  friend class UGGridIdSet<const UGGrid<dim> >;
247  template <class GridImp_, PartitionIteratorType PiType_>
248  friend class UGGridLeafGridView;
249  template <class GridImp_, PartitionIteratorType PiType_>
250  friend class UGGridLevelGridView;
251 
252  friend class GridFactory<UGGrid<dim> >;
253 
254 #ifdef ModelP
255  friend class UGLBGatherScatter;
256 #endif
257 
258  template <int codim_, PartitionIteratorType PiType_, class GridImp_>
259  friend class UGGridLeafIterator;
260  template <int codim_, PartitionIteratorType PiType_, class GridImp_>
261  friend class UGGridLevelIterator;
262  template <int codim_, class GridImp_>
263  friend class UGGridEntityPointer;
264 
266  static_assert(dim==2 || dim==3, "Use UGGrid only for 2d and 3d!");
267 
268  // The different instantiations are mutual friends so they can access
269  // each others numOfUGGrids field
270  friend class UGGrid<2>;
271  friend class UGGrid<3>;
272 
273  //**********************************************************
274  // The Interface Methods
275  //**********************************************************
276  public:
278  typedef UGGridFamily<dim> GridFamily;
279 
280  // the Traits
281  typedef typename UGGridFamily<dim>::Traits Traits;
282 
284  typedef UG::DOUBLE ctype;
285 
292 
295 #if HAVE_NOEXCEPT_SPECIFIER
296  noexcept(false)
297 #endif
298  ;
299 
302  int maxLevel() const;
303 
305  template<int codim>
306  typename Traits::template Codim<codim>::LevelIterator lbegin (int level) const;
307 
309  template<int codim>
310  typename Traits::template Codim<codim>::LevelIterator lend (int level) const;
311 
313  template<int codim, PartitionIteratorType PiType>
314  typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lbegin (int level) const;
315 
317  template<int codim, PartitionIteratorType PiType>
318  typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lend (int level) const;
319 
321  template<int codim>
322  typename Traits::template Codim<codim>::LeafIterator leafbegin() const {
323  return typename Traits::template Codim<codim>::template Partition<All_Partition>::LeafIterator(*this);
324  }
325 
327  template<int codim>
328  typename Traits::template Codim<codim>::LeafIterator leafend() const {
329  return UGGridLeafIterator<codim,All_Partition, const UGGrid<dim> >();
330  }
331 
333  template<int codim, PartitionIteratorType PiType>
334  typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafbegin() const {
335  return typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator(*this);
336  }
337 
339  template<int codim, PartitionIteratorType PiType>
340  typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafend() const {
341  return UGGridLeafIterator<codim,PiType, const UGGrid<dim> >();
342  }
343 
345  template <typename Seed>
346  DUNE_DEPRECATED_MSG("entityPointer() is deprecated and will be removed after the release of dune-grid 2.4. Use entity() instead to directly obtain an Entity object.")
347  typename Traits::template Codim<Seed::codimension>::EntityPointer
348  entityPointer(const Seed& seed) const
349  {
350  enum {codim = Seed::codimension};
351  return typename Traits::template Codim<codim>::EntityPointer(UGGridEntityPointer<codim,const UGGrid<dim> >(this->getRealImplementation(seed).target(),this));
352  }
353 
355  template <typename Seed>
356  typename Traits::template Codim<Seed::codimension>::Entity
357  entity(const Seed& seed) const
358  {
359  const int codim = Seed::codimension;
360  return typename Traits::template Codim<codim>::Entity(UGGridEntity<codim,dim,const UGGrid<dim> >(this->getRealImplementation(seed).target(),this));
361  }
362 
365  int size (int level, int codim) const;
366 
368  int size (int codim) const
369  {
370  return leafIndexSet().size(codim);
371  }
372 
374  int size (int level, GeometryType type) const
375  {
376  return this->levelIndexSet(level).size(type);
377  }
378 
380  int size (GeometryType type) const
381  {
382  return this->leafIndexSet().size(type);
383  }
384 
386  size_t numBoundarySegments() const {
387  // The number is stored as a member of UGGrid upon grid creation.
388  // The corresponding data structure is not exported by UG. (It is in ug/dom/std/std_internal.h)
389  return numBoundarySegments_;
390  }
391 
393  const typename Traits::GlobalIdSet& globalIdSet() const
394  {
395  return idSet_;
396  }
397 
399  const typename Traits::LocalIdSet& localIdSet() const
400  {
401  return idSet_;
402  }
403 
405  const typename Traits::LevelIndexSet& levelIndexSet(int level) const
406  {
407  if (level<0 || level>maxLevel())
408  DUNE_THROW(GridError, "levelIndexSet of nonexisting level " << level << " requested!");
409  return *levelIndexSets_[level];
410  }
411 
413  const typename Traits::LeafIndexSet& leafIndexSet() const
414  {
415  return leafIndexSet_;
416  }
417 
420 
433  bool mark(int refCount, const typename Traits::template Codim<0>::Entity & e );
434 
491  bool mark(const typename Traits::template Codim<0>::Entity & e,
492  typename UG_NS<dim>::RefinementRule rule,
493  int side=0);
494 
496  int getMark(const typename Traits::template Codim<0>::Entity& e) const;
497 
500  bool preAdapt();
501 
503  bool adapt();
504 
506  void postAdapt();
510  unsigned int overlapSize(int codim) const {
511  return 0;
512  }
513 
515  unsigned int ghostSize(int codim) const {
516  return (codim==0) ? 1 : 0;
517  }
518 
520  unsigned int overlapSize(int level, int codim) const {
521  return 0;
522  }
523 
525  unsigned int ghostSize(int level, int codim) const {
526  return (codim==0) ? 1 : 0;
527  }
528 
536  template<class DataHandle>
537  bool loadBalance (DataHandle& dataHandle)
538  {
539 #ifdef ModelP
540  // gather element data
541  // UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
542 
543  // gather node data
544  UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
545 #endif
546 
547  // the load balancing step now also attaches
548  // the data to the entities and distributes it
549  loadBalance();
550 
551 #ifdef ModelP
552  // scatter element data
553  // UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
554 
555  // scatter node data
556  UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
557 #endif
558 
559  return true;
560  }
561 
568  bool loadBalance(int minlevel=0);
569 
600  bool loadBalance(const std::vector<unsigned int>& targetProcessors, unsigned int fromLevel);
601 
611  template<class DataHandle>
612  bool loadBalance (const std::vector<unsigned int>& targetProcessors, unsigned int fromLevel, DataHandle& dataHandle)
613  {
614 #ifdef ModelP
615  // gather element data
616  // UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
617 
618  // gather node data
619  UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
620 #endif
621 
622  // the load balancing step now also attaches
623  // the data to the entities and distributes it
624  loadBalance(targetProcessors,fromLevel);
625 
626 #ifdef ModelP
627  // scatter element data
628  // UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
629 
630  // scatter node data
631  UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
632 #endif
633 
634  return true;
635  }
636 
637 
648  template<class DataHandle>
649  void communicate (DataHandle& dataHandle, InterfaceType iftype, CommunicationDirection dir, int level) const
650  {
651 #ifdef ModelP
652  typedef typename UGGrid::LevelGridView LevelGridView;
653 
654  for (int curCodim = 0; curCodim <= dim; ++curCodim) {
655  if (!dataHandle.contains(dim, curCodim))
656  continue;
657 
658  if (curCodim == 0)
659  communicateUG_<LevelGridView, DataHandle, 0>(this->levelGridView(level), level, dataHandle, iftype, dir);
660  else if (curCodim == dim)
661  communicateUG_<LevelGridView, DataHandle, dim>(this->levelGridView(level), level, dataHandle, iftype, dir);
662  else if (curCodim == dim - 1)
663  communicateUG_<LevelGridView, DataHandle, dim-1>(this->levelGridView(level), level, dataHandle, iftype, dir);
664  else if (curCodim == 1)
665  communicateUG_<LevelGridView, DataHandle, 1>(this->levelGridView(level), level, dataHandle, iftype, dir);
666  else
668  className(*this) << "::communicate(): Not "
669  "supported for dim " << dim << " and codim " << curCodim);
670  }
671 #endif // ModelP
672  }
673 
683  template<class DataHandle>
684  void communicate(DataHandle& dataHandle, InterfaceType iftype, CommunicationDirection dir) const
685  {
686 #ifdef ModelP
687  typedef typename UGGrid::LeafGridView LeafGridView;
688 
689  for (int curCodim = 0; curCodim <= dim; ++curCodim) {
690  if (!dataHandle.contains(dim, curCodim))
691  continue;
692  int level = -1;
693  if (curCodim == 0)
694  communicateUG_<LeafGridView, DataHandle, 0>(this->leafGridView(), level, dataHandle, iftype, dir);
695  else if (curCodim == dim)
696  communicateUG_<LeafGridView, DataHandle, dim>(this->leafGridView(), level, dataHandle, iftype, dir);
697  else if (curCodim == dim - 1)
698  communicateUG_<LeafGridView, DataHandle, dim-1>(this->leafGridView(), level, dataHandle, iftype, dir);
699  else if (curCodim == 1)
700  communicateUG_<LeafGridView, DataHandle, 1>(this->leafGridView(), level, dataHandle, iftype, dir);
701  else
703  className(*this) << "::communicate(): Not "
704  "supported for dim " << dim << " and codim " << curCodim);
705  }
706 #endif // ModelP
707  }
708 
711  {
712  return ccobj_;
713  }
714 
715  protected:
716 #ifdef ModelP
717  template <class GridView, class DataHandle, int codim>
718  void communicateUG_(const GridView& gv, int level,
719  DataHandle &dataHandle,
720  InterfaceType iftype,
721  CommunicationDirection dir) const
722  {
723  typename UG_NS<dim>::DDD_IF_DIR ugIfDir;
724  // Translate the communication direction from Dune-Speak to UG-Speak
725  if (dir==ForwardCommunication)
726  ugIfDir = UG_NS<dim>::IF_FORWARD();
727  else
728  ugIfDir = UG_NS<dim>::IF_BACKWARD();
729 
730  typedef UGMessageBuffer<DataHandle,dim,codim> UGMsgBuf;
731  UGMsgBuf::duneDataHandle_ = &dataHandle;
732 
733  UGMsgBuf::level = level;
734 
735  std::vector<typename UG_NS<dim>::DDD_IF> ugIfs;
736  findDDDInterfaces_(ugIfs, iftype, codim);
737 
738  unsigned bufSize = UGMsgBuf::ugBufferSize_(gv);
739  if (!bufSize)
740  return; // we don't need to communicate if we don't have any data!
741  for (unsigned i=0; i < ugIfs.size(); ++i)
742  UG_NS<dim>::DDD_IFOneway(ugIfs[i],
743  ugIfDir,
744  bufSize,
745  &UGMsgBuf::ugGather_,
746  &UGMsgBuf::ugScatter_);
747  }
748 
749  void findDDDInterfaces_(std::vector<typename UG_NS<dim>::DDD_IF > &dddIfaces,
750  InterfaceType iftype,
751  int codim) const
752  {
753  dddIfaces.clear();
754  if (codim == 0)
755  {
756  switch (iftype) {
758  // do not communicate anything: Elements can not be in
759  // the interior of two processes at the same time
760  return;
762  // The communicated elements are in the sender's
763  // interior and it does not matter what they are for
764  // the receiver
765  dddIfaces.push_back(UG_NS<dim>::ElementVHIF());
766  return;
767  case All_All_Interface :
768  // It does neither matter what the communicated
769  // elements are for sender nor for the receiver. If
770  // they are seen by these two processes, data is send
771  // and received.
772  dddIfaces.push_back(UG_NS<dim>::ElementSymmVHIF());
773  return;
774  default :
775  DUNE_THROW(GridError,
776  "Element communication not supported for "
777  "interfaces of type "
778  << iftype);
779  }
780  }
781  else if (codim == dim)
782  {
783  switch (iftype)
784  {
786  dddIfaces.push_back(UG_NS<dim>::BorderNodeSymmIF());
787  return;
789  dddIfaces.push_back(UG_NS<dim>::BorderNodeSymmIF());
790  dddIfaces.push_back(UG_NS<dim>::NodeIF());
791  return;
792  case All_All_Interface :
793  dddIfaces.push_back(UG_NS<dim>::NodeAllIF());
794  return;
795  default :
796  DUNE_THROW(GridError,
797  "Node communication not supported for "
798  "interfaces of type "
799  << iftype);
800  }
801  }
802  else if (codim == dim-1)
803  {
804  switch (iftype)
805  {
807  dddIfaces.push_back(UG_NS<dim>::BorderEdgeSymmIF());
808  return;
810  dddIfaces.push_back(UG_NS<dim>::BorderEdgeSymmIF());
811  // Is the following line needed or not?
812  // dddIfaces.push_back(UG_NS<dim>::EdgeIF());
813  return;
814  case All_All_Interface :
815  dddIfaces.push_back(UG_NS<dim>::EdgeSymmVHIF());
816  return;
817  default :
818  DUNE_THROW(GridError,
819  "Edge communication not supported for "
820  "interfaces of type "
821  << iftype);
822  }
823  }
824  else if (codim == 1)
825  {
826  switch (iftype)
827  {
830  dddIfaces.push_back(UG_NS<dim>::BorderVectorSymmIF());
831  return;
832  default :
833  DUNE_THROW(GridError,
834  "Face communication not supported for "
835  "interfaces of type "
836  << iftype);
837  }
838  }
839  else
840  {
841  DUNE_THROW(GridError,
842  "Communication for codim "
843  << codim
844  << " entities is not yet supported "
845  << " by the DUNE UGGrid interface!");
846  }
847  };
848 #endif // ModelP
849 
850  public:
851  // **********************************************************
852  // End of Interface Methods
853  // **********************************************************
854 
855 #ifdef DOXYGEN
856 
868  void getChildrenOfSubface(const typename Traits::template Codim<0>::EntityPointer & e,
869  int elementSide,
870  int maxl,
871  std::vector<typename Traits::template Codim<0>::EntityPointer>& childElements,
872  std::vector<unsigned char>& childElementSides) const;
873 
874 #else
875 
876  // This is the actual implementation of the deprecated method. We need this ugly trick of turning
877  // it into a template to avoid triggering a deprecation warning every time UGGrid is used.
878  template< typename T >
879  DUNE_DEPRECATED_MSG("This version of getChildrenOfSubface() uses EntityPointer and is deprecated. It will be removed after the release of Dune 2.4. Please use the new version with entities instead.")
880  typename std::enable_if<
881  std::is_same<
882  T,
883  typename UGGrid<dim>::Traits::template Codim<0>::EntityPointer
884  >::value
885  >::type
886  getChildrenOfSubface(const T & e,
887  int elementSide,
888  int maxl,
889  std::vector<typename Traits::template Codim<0>::EntityPointer>& childElements,
890  std::vector<unsigned char>& childElementSides) const
891  {
892 
893  typedef std::pair<typename UG_NS<dim>::Element*,int> ListEntryType;
894 
896 
897  // //////////////////////////////////////////////////////////////////////
898  // Change the input face number from Dune numbering to UG numbering
899  // //////////////////////////////////////////////////////////////////////
900 
901  elementSide = UGGridRenumberer<dim>::facesDUNEtoUG(elementSide, e->type());
902 
903  // ///////////////
904  // init list
905  // ///////////////
906  if (!e->isLeaf() // Get_Sons_of_ElementSide returns GM_FATAL when called for a leaf !?!
907  && e->level() < maxl) {
908 
909  typename UG_NS<dim>::Element* theElement = this->getRealImplementation(*e).target_;
910 
911  int Sons_of_Side = 0;
912  typename UG_NS<dim>::Element* SonList[UG_NS<dim>::MAX_SONS];
913  int SonSides[UG_NS<dim>::MAX_SONS];
914 
915  int rv = Get_Sons_of_ElementSide(theElement,
916  elementSide,
917  &Sons_of_Side,
918  SonList, // the output elements
919  SonSides, // Output element side numbers
920  true, // Element sons are not precomputed
921  true); // ioflag: I have no idea what this is supposed to do
922  if (rv!=0)
923  DUNE_THROW(GridError, "Get_Sons_of_ElementSide returned with error value " << rv);
924 
925  for (int i=0; i<Sons_of_Side; i++)
926  list.push_back(ListEntryType(SonList[i],SonSides[i]));
927 
928  }
929 
930  // //////////////////////////////////////////////////
931  // Traverse and collect all children of the side
932  // //////////////////////////////////////////////////
933 
934  typename SLList<ListEntryType>::iterator f = list.begin();
935  for (; f!=list.end(); ++f) {
936 
937  typename UG_NS<dim>::Element* theElement = f->first;
938  int side = f->second;
939 
940  int Sons_of_Side = 0;
941  typename UG_NS<dim>::Element* SonList[UG_NS<dim>::MAX_SONS];
942  int SonSides[UG_NS<dim>::MAX_SONS];
943 
944  if (UG_NS<dim>::myLevel(theElement) < maxl) {
945 
946  Get_Sons_of_ElementSide(theElement,
947  side, // Input element side number
948  &Sons_of_Side, // Number of topological sons of the element side
949  SonList, // Output elements
950  SonSides, // Output element side numbers
951  true,
952  true);
953 
954  for (int i=0; i<Sons_of_Side; i++)
955  list.push_back(ListEntryType(SonList[i],SonSides[i]));
956 
957  }
958 
959  }
960 
961  // //////////////////////////////
962  // Extract result from stack
963  // //////////////////////////////
964 
965  // Use reserve / push_back since EntityPointer is not default constructable
966  childElements.clear();
967  childElements.reserve( list.size() );
968  childElementSides.resize(list.size());
969 
970  int i=0;
971  for (f = list.begin(); f!=list.end(); ++f, ++i)
972  {
973 
974  // Set element
975  typedef typename Traits::template Codim< 0 >::EntityPointer EntityPointer;
976  childElements.push_back( EntityPointer( UGGridEntityPointer< 0, const UGGrid< dim > >( f->first, this ) ) );
977 
978  int side = f->second;
979 
980  // Dune numbers the faces of several elements differently than UG.
981  // The following switch does the transformation
982  childElementSides[i] = UGGridRenumberer<dim>::facesUGtoDUNE(side, childElements[i]->type());
983 
984  }
985 
986  }
987 #endif
988 
996  void getChildrenOfSubface(const typename Traits::template Codim<0>::Entity & e,
997  int elementSide,
998  int maxl,
999  std::vector<typename Traits::template Codim<0>::Entity>& childElements,
1000  std::vector<unsigned char>& childElementSides) const;
1001 
1007  COPY
1008  };
1009 
1015  NONE
1016  };
1017 
1020  refinementType_ = type;
1021  }
1022 
1025  closureType_ = type;
1026  }
1027 
1034  static void setDefaultHeapSize(unsigned size) {
1035  heapSize_ = size;
1036  }
1037 
1038 #ifdef DOXYGEN
1039 
1048  void setPosition(const typename Traits::template Codim<dim>::EntityPointer& e,
1049  const FieldVector<double, dim>& pos);
1050 
1051 #else
1052 
1053  // This is the actual implementation of the deprecated method. We need this ugly trick of turning
1054  // it into a template to avoid triggering a deprecation warning every time UGGrid is used.
1055  template< typename T >
1056  DUNE_DEPRECATED_MSG("This version of setPosition() uses EntityPointer and is deprecated. It will be removed after the release of Dune 2.4. Please use the new version with entities instead.")
1057  typename std::enable_if<
1058  std::is_same<
1059  T,
1060  typename UGGrid<dim>::Traits::template Codim<dim>::EntityPointer
1061  >::value
1062  >::type
1063  setPosition(const T& e,
1064  const FieldVector<double, dim>& pos)
1065 {
1066  typename UG_NS<dim>::Node* target = this->getRealImplementation(*e).target_;
1067 
1068  for (int i=0; i<dim; i++)
1069  target->myvertex->iv.x[i] = pos[i];
1070 }
1071 
1072 #endif
1073 
1077  void setPosition(const typename Traits::template Codim<dim>::Entity& e,
1078  const FieldVector<double, dim>& pos);
1079 
1084  void globalRefine(int n);
1085 
1090  void saveState(const std::string& filename) const;
1091 
1096  void loadState(const std::string& filename);
1097 
1098  private:
1100  typename UG_NS<dim>::MultiGrid* multigrid_;
1101 
1104 
1110  void setIndices(bool setLevelZero,
1111  std::vector<unsigned int>* nodePermutation);
1112 
1113  // Each UGGrid object has a unique name to identify it in the
1114  // UG environment structure
1115  std::string name_;
1116 
1117  // Our set of level indices
1118  std::vector<shared_ptr<UGGridLevelIndexSet<const UGGrid<dim> > > > levelIndexSets_;
1119 
1120  UGGridLeafIndexSet<const UGGrid<dim> > leafIndexSet_;
1121 
1122  // One id set implementation
1123  // Used for both the local and the global UGGrid id sets
1124  UGGridIdSet<const UGGrid<dim> > idSet_;
1125 
1127  RefinementType refinementType_;
1128 
1130  ClosureType closureType_;
1131 
1139  static int numOfUGGrids;
1140 
1146  bool someElementHasBeenMarkedForRefinement_;
1147 
1153  bool someElementHasBeenMarkedForCoarsening_;
1154 
1159  static unsigned int heapSize_;
1160 
1162  std::vector<shared_ptr<BoundarySegment<dim> > > boundarySegments_;
1163 
1169  unsigned int numBoundarySegments_;
1170 
1171  }; // end Class UGGrid
1172 
1173  namespace Capabilities
1174  {
1190  template<int dim>
1191  struct hasEntity< UGGrid<dim>, 0>
1192  {
1193  static const bool v = true;
1194  };
1195 
1199  template<int dim>
1200  struct hasEntity< UGGrid<dim>, dim>
1201  {
1202  static const bool v = true;
1203  };
1204 
1208  template<int dim>
1209  struct DUNE_DEPRECATED_MSG("Capabilities::isParallel will be removed after dune-grid-2.4.") isParallel< UGGrid<dim> >
1210  {
1211 #ifdef ModelP
1212  static const bool DUNE_DEPRECATED_MSG("Capabilities::isParallel will be removed after dune-grid-2.4.") v = true;
1213 #else // !ModelP
1214  static const bool DUNE_DEPRECATED_MSG("Capabilities::isParallel will be removed after dune-grid-2.4.") v = false;
1215 #endif // !ModelP
1216  };
1217 
1221  template<int dim>
1223  {
1224  static const bool v = true;
1225  };
1226 
1230  template<int dim>
1232  {
1233  static const bool v = false;
1234  };
1235 
1236  }
1237 
1238 } // namespace Dune
1239 
1240 #endif // HAVE_UG || DOXYGEN
1241 #endif // DUNE_UGGRID_HH
Base class for grid boundary segments of arbitrary geometry.
CollectiveCommunication()
Construct default object.
Definition: collectivecommunication.hh:76
Wrapper class for pointers to entities.
Definition: entitypointer.hh:113
vector space out of a tensor product of fields.
Definition: fvector.hh:94
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:25
Definition: grid.hh:1030
static std::conditional< std::is_reference< InterfaceType >::value, typename std::add_lvalue_reference< typename ReturnImplementationType< typename std::remove_reference< InterfaceType >::type >::ImplementationType >::type, typename std::remove_const< typename ReturnImplementationType< typename std::remove_reference< InterfaceType >::type >::ImplementationType >::type >::type getRealImplementation(InterfaceType &&i)
return real implementation of interface class
Definition: grid.hh:1305
Traits::template Codim< codim >::LevelIterator DUNE_DEPRECATED_MSG("The method lbegin( level ) is superseded by levelGridView( level ).begin.") lbegin(int level) const
Iterator to first entity of given codim on level for PartitionType All_Partition.
Definition: grid.hh:1043
Traits::template Partition< pitype >::LevelGridView levelGridView(int level) const
View for a grid level.
Definition: grid.hh:1105
bool loadBalance()
default implementation of load balance does nothing and returns false
Definition: grid.hh:1220
Traits::template Partition< pitype >::LeafGridView leafGridView() const
View for the leaf grid.
Definition: grid.hh:1115
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:59
GridFamily::Traits::CollectiveCommunication CollectiveCommunication
A type that is a model of Dune::CollectiveCommunication. It provides a portable way for collective co...
Definition: grid.hh:545
Partition< All_Partition >::LevelGridView LevelGridView
View types for All_Partition.
Definition: grid.hh:428
Id Set Interface.
Definition: indexidset.hh:414
Index Set Interface base class.
Definition: indexidset.hh:76
Default exception for dummy implementations.
Definition: exceptions.hh:288
A mutable iterator for the SLList.
Definition: sllist.hh:269
A single linked list.
Definition: sllist.hh:42
Front-end for the grid manager of the finite element toolbox UG.
Definition: uggrid.hh:228
void communicate(DataHandle &dataHandle, InterfaceType iftype, CommunicationDirection dir) const
The communication interface for all codims on the leaf level.
Definition: uggrid.hh:684
~UGGrid()
Destructor.
void postAdapt()
Clean up refinement markers.
size_t numBoundarySegments() const
Return the number of boundary segments.
Definition: uggrid.hh:386
void setRefinementType(RefinementType type)
Sets the type of grid refinement.
Definition: uggrid.hh:1019
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:520
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:649
UGGrid()
Default constructor.
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator lbegin(int level) const
Iterator to first entity of given codim on level.
const CollectiveCommunication< UGGrid > & comm() const
Definition: uggrid.hh:710
void setPosition(const typename Traits::template Codim< dim >::EntityPointer &e, const FieldVector< double, dim > &pos)
Sets a vertex to a new position.
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:278
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
bool loadBalance(const std::vector< unsigned int > &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:612
void setClosureType(ClosureType type)
Sets the type of grid refinement closure.
Definition: uggrid.hh:1024
DUNE_DEPRECATED_MSG("entityPointer() is deprecated and will be removed after the release of dune-grid 2.4. Use entity() instead to directly obtain an Entity object.") typename Traits int size(int codim) const
Create an EntityPointer from an EntitySeed.
Definition: uggrid.hh:368
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition: uggrid.hh:393
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition: uggrid.hh:374
RefinementType
The different forms of grid refinement that UG supports.
Definition: uggrid.hh:1003
@ COPY
New level consists of the refined elements and the unrefined ones, too.
Definition: uggrid.hh:1007
@ LOCAL
New level consists only of the refined elements and the closure.
Definition: uggrid.hh:1005
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: uggrid.hh:380
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:515
UG::DOUBLE ctype
The type used to store coordinates.
Definition: uggrid.hh:284
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition: uggrid.hh:399
unsigned int ghostSize(int level, int codim) const
Size of the ghost cell layer on a given level.
Definition: uggrid.hh:525
bool adapt()
Triggers the grid refinement process.
const Traits::LevelIndexSet & levelIndexSet(int level) const
Access to the LevelIndexSets.
Definition: uggrid.hh:405
Traits::template Codim< codim >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim.
Definition: uggrid.hh:322
Traits::template Codim< codim >::template Partition< PiType >::LevelIterator lend(int level) const
one past the end on this level
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:510
bool loadBalance(DataHandle &dataHandle)
Distributes the grid and some data over the available nodes in a distributed machine.
Definition: uggrid.hh:537
Traits::template Codim< codim >::LeafIterator leafend() const
one past the end of the sequence of leaf entities
Definition: uggrid.hh:328
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition: uggrid.hh:413
Traits::template Codim< codim >::LevelIterator lend(int level) const
one past the end on this 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.
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim.
Definition: uggrid.hh:334
Traits::template Codim< codim >::LevelIterator lbegin(int level) const
Iterator to first entity of given codim on level.
ClosureType
Decide whether to add a green closure to locally refined grid sections or not.
Definition: uggrid.hh:1011
@ GREEN
Standard red/green refinement.
Definition: uggrid.hh:1013
@ NONE
No closure, results in nonconforming meshes.
Definition: uggrid.hh:1015
bool loadBalance(const std::vector< unsigned int > &targetProcessors, unsigned int fromLevel)
Distribute this grid over a distributed machine.
Traits::template Codim< codim >::template Partition< PiType >::LeafIterator leafend() const
one past the end of the sequence of leaf entities
Definition: uggrid.hh:340
static void setDefaultHeapSize(unsigned size)
Sets the default heap size.
Definition: uggrid.hh:1034
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.
void push_back(const MemberType &item)
Add a new entry to the end of the list.
Definition: sllist.hh:657
int size() const
Get the number of elements the list contains.
Definition: sllist.hh:770
iterator end()
Get an iterator pointing to the end of the list.
Definition: sllist.hh:788
iterator begin()
Get an iterator pointing to the first element in the list.
Definition: sllist.hh:776
#define DUNE_THROW(E, m)
Definition: exceptions.hh:243
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
struct DUNE_DEPRECATED_MSG("Capabilities::isParallel will be removed after dune-grid-2.4.") isParallel< ALU3dGrid< elType
YaspGrid is parallel.
Implements an utility class that provides MPI's collective communication methods.
Helpers for dealing with MPI.
Dune namespace.
Definition: alignment.hh:10
std::string className(T &t)
Provide the demangled class name of a given object as a string.
Definition: classname.hh:23
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:108
Specialize with 'true' if implementation guarantees conforming level grids. (default=false)
Definition: capabilities.hh:99
Specialize with 'true' if implementation supports parallelism. (default=false)
Definition: capabilities.hh:73
Static tag representing a codimension.
Definition: dimension.hh:22
A traits struct that collects all associated types of one grid model.
Definition: grid.hh:1344
The specialization of the generic GridFactory for UGGrid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 10, 22:30, 2024)