Dune Core Modules (2.4.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 
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 /* 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 #ifdef ModelP
72 #define __PPIF__
73 #endif
74 
75 #define _3
76 #include "uggrid/ugincludes.hh"
77 
78 // Wrap a few large UG macros by functions before they get undef'ed away.
79 // This time it's the 3d-versions.
80 #define UG_DIM 3
81 #include "uggrid/ugwrapper.hh"
82 #undef UG_DIM
83 
84 // undef all macros defined by UG
85 #include "uggrid/ug_undefs.hh"
86 
87 #undef _3
88 #undef FOR_DUNE
89 
90 // The components of the UGGrid interface
91 #include "uggrid/uggridgeometry.hh"
92 #include "uggrid/uggridlocalgeometry.hh"
93 #include "uggrid/uggridentity.hh"
94 #include "uggrid/uggridentitypointer.hh"
95 #include "uggrid/uggridentityseed.hh"
96 #include "uggrid/uggridintersections.hh"
97 #include "uggrid/uggridintersectioniterators.hh"
98 #include "uggrid/uggridleveliterator.hh"
99 #include "uggrid/uggridleafiterator.hh"
100 #include "uggrid/uggridhieriterator.hh"
101 #include "uggrid/uggridindexsets.hh"
102 #include <dune/grid/uggrid/uggridviews.hh>
103 #ifdef ModelP
104 #include "uggrid/ugmessagebuffer.hh"
105 #include "uggrid/uglbgatherscatter.hh"
106 #endif
107 
108 // Not needed here, but included for user convenience
109 #include "uggrid/uggridfactory.hh"
110 
111 #ifdef ModelP
112 template <class DataHandle, int GridDim, int codim>
113 DataHandle *Dune::UGMessageBufferBase<DataHandle,GridDim,codim>::duneDataHandle_ = 0;
114 
115 template <class DataHandle, int GridDim, int codim>
116 int Dune::UGMessageBufferBase<DataHandle,GridDim,codim>::level = -1;
117 #endif // ModelP
118 
119 namespace Dune {
120 
121 #ifdef ModelP
122  template <int dim>
123  class CollectiveCommunication<Dune::UGGrid<dim> > :
124  public CollectiveCommunication< Dune::MPIHelper::MPICommunicator >
125  {
126  typedef CollectiveCommunication< Dune::MPIHelper::MPICommunicator > ParentType;
127  public:
129  : ParentType(MPIHelper::getCommunicator())
130  {}
131  };
132 #endif // ModelP
133 
134  template<int dim>
135  struct UGGridFamily
136  {
137  typedef GridTraits<dim,dim,Dune::UGGrid<dim>,
138  UGGridGeometry,
139  UGGridEntity,
140  UGGridEntityPointer,
141  UGGridLevelIterator,
142  UGGridLeafIntersection,
143  UGGridLevelIntersection,
144  UGGridLeafIntersectionIterator,
145  UGGridLevelIntersectionIterator,
146  UGGridHierarchicIterator,
147  UGGridLeafIterator,
148  UGGridLevelIndexSet< const UGGrid<dim> >,
149  UGGridLeafIndexSet< const UGGrid<dim> >,
150  UGGridIdSet< const UGGrid<dim> >,
151  typename UG_NS<dim>::UG_ID_TYPE,
152  UGGridIdSet< const UGGrid<dim> >,
153  typename UG_NS<dim>::UG_ID_TYPE,
154  CollectiveCommunication<Dune::UGGrid<dim> >,
155  UGGridLevelGridViewTraits,
156  UGGridLeafGridViewTraits,
157  UGGridEntitySeed,
158  UGGridLocalGeometry>
159  Traits;
160  };
161 
162 
163  //**********************************************************************
164  //
165  // --UGGrid
166  //
167  //**********************************************************************
168 
206  template <int dim>
207  class UGGrid : public GridDefaultImplementation <dim, dim, double, UGGridFamily<dim> >
208  {
210 
211  friend class UGGridGeometry<0,dim,const UGGrid<dim> >;
212  friend class UGGridGeometry<dim,dim,const UGGrid<dim> >;
213  friend class UGGridGeometry<1,2,const UGGrid<dim> >;
214  friend class UGGridGeometry<2,3,const UGGrid<dim> >;
215 
216  friend class UGGridEntity <0,dim,const UGGrid<dim> >;
217  friend class UGGridEntity <dim,dim,const UGGrid<dim> >;
218  friend class UGGridHierarchicIterator<const UGGrid<dim> >;
219  friend class UGGridLeafIntersection<const UGGrid<dim> >;
220  friend class UGGridLevelIntersection<const UGGrid<dim> >;
221  friend class UGGridLeafIntersectionIterator<const UGGrid<dim> >;
222  friend class UGGridLevelIntersectionIterator<const UGGrid<dim> >;
223 
224  friend class UGGridLevelIndexSet<const UGGrid<dim> >;
225  friend class UGGridLeafIndexSet<const UGGrid<dim> >;
226  friend class UGGridIdSet<const UGGrid<dim> >;
227  template <class GridImp_, PartitionIteratorType PiType_>
228  friend class UGGridLeafGridView;
229  template <class GridImp_, PartitionIteratorType PiType_>
230  friend class UGGridLevelGridView;
231 
232  friend class GridFactory<UGGrid<dim> >;
233 
234 #ifdef ModelP
235  friend class UGLBGatherScatter;
236 #endif
237 
238  template <int codim_, PartitionIteratorType PiType_, class GridImp_>
239  friend class UGGridLeafIterator;
240  template <int codim_, PartitionIteratorType PiType_, class GridImp_>
241  friend class UGGridLevelIterator;
242  template <int codim_, class GridImp_>
243  friend class UGGridEntityPointer;
244 
246  static_assert(dim==2 || dim==3, "Use UGGrid only for 2d and 3d!");
247 
248  // The different instantiations are mutual friends so they can access
249  // each others numOfUGGrids field
250  friend class UGGrid<2>;
251  friend class UGGrid<3>;
252 
253  //**********************************************************
254  // The Interface Methods
255  //**********************************************************
256  public:
258  typedef UGGridFamily<dim> GridFamily;
259 
260  // the Traits
261  typedef typename UGGridFamily<dim>::Traits Traits;
262 
264  typedef UG::DOUBLE ctype;
265 
272 
275 #if HAVE_NOEXCEPT_SPECIFIER
276  noexcept(false)
277 #endif
278  ;
279 
282  int maxLevel() const;
283 
285  template<int codim>
286  typename Traits::template Codim<codim>::LevelIterator lbegin (int level) const;
287 
289  template<int codim>
290  typename Traits::template Codim<codim>::LevelIterator lend (int level) const;
291 
293  template<int codim, PartitionIteratorType PiType>
294  typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lbegin (int level) const;
295 
297  template<int codim, PartitionIteratorType PiType>
298  typename Traits::template Codim<codim>::template Partition<PiType>::LevelIterator lend (int level) const;
299 
301  template<int codim>
302  typename Traits::template Codim<codim>::LeafIterator leafbegin() const {
303  return typename Traits::template Codim<codim>::template Partition<All_Partition>::LeafIterator(*this);
304  }
305 
307  template<int codim>
308  typename Traits::template Codim<codim>::LeafIterator leafend() const {
309  return UGGridLeafIterator<codim,All_Partition, const UGGrid<dim> >();
310  }
311 
313  template<int codim, PartitionIteratorType PiType>
314  typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafbegin() const {
315  return typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator(*this);
316  }
317 
319  template<int codim, PartitionIteratorType PiType>
320  typename Traits::template Codim<codim>::template Partition<PiType>::LeafIterator leafend() const {
321  return UGGridLeafIterator<codim,PiType, const UGGrid<dim> >();
322  }
323 
325  template <typename Seed>
326  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.")
327  typename Traits::template Codim<Seed::codimension>::EntityPointer
328  entityPointer(const Seed& seed) const
329  {
330  enum {codim = Seed::codimension};
331  return typename Traits::template Codim<codim>::EntityPointer(UGGridEntityPointer<codim,const UGGrid<dim> >(this->getRealImplementation(seed).target(),this));
332  }
333 
335  template <typename Seed>
336  typename Traits::template Codim<Seed::codimension>::Entity
337  entity(const Seed& seed) const
338  {
339  const int codim = Seed::codimension;
340  return typename Traits::template Codim<codim>::Entity(UGGridEntity<codim,dim,const UGGrid<dim> >(this->getRealImplementation(seed).target(),this));
341  }
342 
345  int size (int level, int codim) const;
346 
348  int size (int codim) const
349  {
350  return leafIndexSet().size(codim);
351  }
352 
354  int size (int level, GeometryType type) const
355  {
356  return this->levelIndexSet(level).size(type);
357  }
358 
360  int size (GeometryType type) const
361  {
362  return this->leafIndexSet().size(type);
363  }
364 
366  size_t numBoundarySegments() const {
367  // The number is stored as a member of UGGrid upon grid creation.
368  // The corresponding data structure is not exported by UG. (It is in ug/dom/std/std_internal.h)
369  return numBoundarySegments_;
370  }
371 
373  const typename Traits::GlobalIdSet& globalIdSet() const
374  {
375  return idSet_;
376  }
377 
379  const typename Traits::LocalIdSet& localIdSet() const
380  {
381  return idSet_;
382  }
383 
385  const typename Traits::LevelIndexSet& levelIndexSet(int level) const
386  {
387  if (level<0 || level>maxLevel())
388  DUNE_THROW(GridError, "levelIndexSet of nonexisting level " << level << " requested!");
389  return *levelIndexSets_[level];
390  }
391 
393  const typename Traits::LeafIndexSet& leafIndexSet() const
394  {
395  return leafIndexSet_;
396  }
397 
400 
413  bool mark(int refCount, const typename Traits::template Codim<0>::Entity & e );
414 
471  bool mark(const typename Traits::template Codim<0>::Entity & e,
472  typename UG_NS<dim>::RefinementRule rule,
473  int side=0);
474 
476  int getMark(const typename Traits::template Codim<0>::Entity& e) const;
477 
480  bool preAdapt();
481 
483  bool adapt();
484 
486  void postAdapt();
490  unsigned int overlapSize(int codim) const {
491  return 0;
492  }
493 
495  unsigned int ghostSize(int codim) const {
496  return (codim==0) ? 1 : 0;
497  }
498 
500  unsigned int overlapSize(int level, int codim) const {
501  return 0;
502  }
503 
505  unsigned int ghostSize(int level, int codim) const {
506  return (codim==0) ? 1 : 0;
507  }
508 
516  template<class DataHandle>
517  bool loadBalance (DataHandle& dataHandle)
518  {
519 #ifdef ModelP
520  // gather element data
521  // UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
522 
523  // gather node data
524  UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
525 #endif
526 
527  // the load balancing step now also attaches
528  // the data to the entities and distributes it
529  loadBalance();
530 
531 #ifdef ModelP
532  // scatter element data
533  // UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
534 
535  // scatter node data
536  UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
537 #endif
538 
539  return true;
540  }
541 
548  bool loadBalance(int minlevel=0);
549 
580  bool loadBalance(const std::vector<unsigned int>& targetProcessors, unsigned int fromLevel);
581 
591  template<class DataHandle>
592  bool loadBalance (const std::vector<unsigned int>& targetProcessors, unsigned int fromLevel, DataHandle& dataHandle)
593  {
594 #ifdef ModelP
595  // gather element data
596  // UGLBGatherScatter::template gather<0>(this->leafGridView(), dataHandle);
597 
598  // gather node data
599  UGLBGatherScatter::template gather<dim>(this->leafGridView(), dataHandle);
600 #endif
601 
602  // the load balancing step now also attaches
603  // the data to the entities and distributes it
604  loadBalance(targetProcessors,fromLevel);
605 
606 #ifdef ModelP
607  // scatter element data
608  // UGLBGatherScatter::template scatter<0>(this->leafGridView(), dataHandle);
609 
610  // scatter node data
611  UGLBGatherScatter::template scatter<dim>(this->leafGridView(), dataHandle);
612 #endif
613 
614  return true;
615  }
616 
617 
628  template<class DataHandle>
629  void communicate (DataHandle& dataHandle, InterfaceType iftype, CommunicationDirection dir, int level) const
630  {
631 #ifdef ModelP
632  typedef typename UGGrid::LevelGridView LevelGridView;
633 
634  for (int curCodim = 0; curCodim <= dim; ++curCodim) {
635  if (!dataHandle.contains(dim, curCodim))
636  continue;
637 
638  if (curCodim == 0)
639  communicateUG_<LevelGridView, DataHandle, 0>(this->levelGridView(level), level, dataHandle, iftype, dir);
640  else if (curCodim == dim)
641  communicateUG_<LevelGridView, DataHandle, dim>(this->levelGridView(level), level, dataHandle, iftype, dir);
642  else if (curCodim == dim - 1)
643  communicateUG_<LevelGridView, DataHandle, dim-1>(this->levelGridView(level), level, dataHandle, iftype, dir);
644  else if (curCodim == 1)
645  communicateUG_<LevelGridView, DataHandle, 1>(this->levelGridView(level), level, dataHandle, iftype, dir);
646  else
648  className(*this) << "::communicate(): Not "
649  "supported for dim " << dim << " and codim " << curCodim);
650  }
651 #endif // ModelP
652  }
653 
663  template<class DataHandle>
664  void communicate(DataHandle& dataHandle, InterfaceType iftype, CommunicationDirection dir) const
665  {
666 #ifdef ModelP
667  typedef typename UGGrid::LeafGridView LeafGridView;
668 
669  for (int curCodim = 0; curCodim <= dim; ++curCodim) {
670  if (!dataHandle.contains(dim, curCodim))
671  continue;
672  int level = -1;
673  if (curCodim == 0)
674  communicateUG_<LeafGridView, DataHandle, 0>(this->leafGridView(), level, dataHandle, iftype, dir);
675  else if (curCodim == dim)
676  communicateUG_<LeafGridView, DataHandle, dim>(this->leafGridView(), level, dataHandle, iftype, dir);
677  else if (curCodim == dim - 1)
678  communicateUG_<LeafGridView, DataHandle, dim-1>(this->leafGridView(), level, dataHandle, iftype, dir);
679  else if (curCodim == 1)
680  communicateUG_<LeafGridView, DataHandle, 1>(this->leafGridView(), level, dataHandle, iftype, dir);
681  else
683  className(*this) << "::communicate(): Not "
684  "supported for dim " << dim << " and codim " << curCodim);
685  }
686 #endif // ModelP
687  }
688 
691  {
692  return ccobj_;
693  }
694 
695  protected:
696 #ifdef ModelP
697  template <class GridView, class DataHandle, int codim>
698  void communicateUG_(const GridView& gv, int level,
699  DataHandle &dataHandle,
700  InterfaceType iftype,
701  CommunicationDirection dir) const
702  {
703  typename UG_NS<dim>::DDD_IF_DIR ugIfDir;
704  // Translate the communication direction from Dune-Speak to UG-Speak
705  if (dir==ForwardCommunication)
706  ugIfDir = UG_NS<dim>::IF_FORWARD();
707  else
708  ugIfDir = UG_NS<dim>::IF_BACKWARD();
709 
710  typedef UGMessageBuffer<DataHandle,dim,codim> UGMsgBuf;
711  UGMsgBuf::duneDataHandle_ = &dataHandle;
712 
713  UGMsgBuf::level = level;
714 
715  std::vector<typename UG_NS<dim>::DDD_IF> ugIfs;
716  findDDDInterfaces_(ugIfs, iftype, codim);
717 
718  unsigned bufSize = UGMsgBuf::ugBufferSize_(gv);
719  if (!bufSize)
720  return; // we don't need to communicate if we don't have any data!
721  for (unsigned i=0; i < ugIfs.size(); ++i)
722  UG_NS<dim>::DDD_IFOneway(ugIfs[i],
723  ugIfDir,
724  bufSize,
725  &UGMsgBuf::ugGather_,
726  &UGMsgBuf::ugScatter_);
727  }
728 
729  void findDDDInterfaces_(std::vector<typename UG_NS<dim>::DDD_IF > &dddIfaces,
730  InterfaceType iftype,
731  int codim) const
732  {
733  dddIfaces.clear();
734  if (codim == 0)
735  {
736  switch (iftype) {
738  // do not communicate anything: Elements can not be in
739  // the interior of two processes at the same time
740  return;
742  // The communicated elements are in the sender's
743  // interior and it does not matter what they are for
744  // the receiver
745  dddIfaces.push_back(UG_NS<dim>::ElementVHIF());
746  return;
747  case All_All_Interface :
748  // It does neither matter what the communicated
749  // elements are for sender nor for the receiver. If
750  // they are seen by these two processes, data is send
751  // and received.
752  dddIfaces.push_back(UG_NS<dim>::ElementSymmVHIF());
753  return;
754  default :
755  DUNE_THROW(GridError,
756  "Element communication not supported for "
757  "interfaces of type "
758  << iftype);
759  }
760  }
761  else if (codim == dim)
762  {
763  switch (iftype)
764  {
766  dddIfaces.push_back(UG_NS<dim>::BorderNodeSymmIF());
767  return;
769  dddIfaces.push_back(UG_NS<dim>::BorderNodeSymmIF());
770  dddIfaces.push_back(UG_NS<dim>::NodeIF());
771  return;
772  case All_All_Interface :
773  dddIfaces.push_back(UG_NS<dim>::NodeAllIF());
774  return;
775  default :
776  DUNE_THROW(GridError,
777  "Node communication not supported for "
778  "interfaces of type "
779  << iftype);
780  }
781  }
782  else if (codim == dim-1)
783  {
784  switch (iftype)
785  {
787  dddIfaces.push_back(UG_NS<dim>::BorderEdgeSymmIF());
788  return;
790  dddIfaces.push_back(UG_NS<dim>::BorderEdgeSymmIF());
791  // Is the following line needed or not?
792  // dddIfaces.push_back(UG_NS<dim>::EdgeIF());
793  return;
794  case All_All_Interface :
795  dddIfaces.push_back(UG_NS<dim>::EdgeSymmVHIF());
796  return;
797  default :
798  DUNE_THROW(GridError,
799  "Edge communication not supported for "
800  "interfaces of type "
801  << iftype);
802  }
803  }
804  else if (codim == 1)
805  {
806  switch (iftype)
807  {
810  dddIfaces.push_back(UG_NS<dim>::BorderVectorSymmIF());
811  return;
812  default :
813  DUNE_THROW(GridError,
814  "Face communication not supported for "
815  "interfaces of type "
816  << iftype);
817  }
818  }
819  else
820  {
821  DUNE_THROW(GridError,
822  "Communication for codim "
823  << codim
824  << " entities is not yet supported "
825  << " by the DUNE UGGrid interface!");
826  }
827  };
828 #endif // ModelP
829 
830  public:
831  // **********************************************************
832  // End of Interface Methods
833  // **********************************************************
834 
835 #ifdef DOXYGEN
836 
848  void getChildrenOfSubface(const typename Traits::template Codim<0>::EntityPointer & e,
849  int elementSide,
850  int maxl,
851  std::vector<typename Traits::template Codim<0>::EntityPointer>& childElements,
852  std::vector<unsigned char>& childElementSides) const;
853 
854 #else
855 
856  // This is the actual implementation of the deprecated method. We need this ugly trick of turning
857  // it into a template to avoid triggering a deprecation warning every time UGGrid is used.
858  template< typename T >
859  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.")
860  typename std::enable_if<
861  std::is_same<
862  T,
863  typename UGGrid<dim>::Traits::template Codim<0>::EntityPointer
864  >::value
865  >::type
866  getChildrenOfSubface(const T & e,
867  int elementSide,
868  int maxl,
869  std::vector<typename Traits::template Codim<0>::EntityPointer>& childElements,
870  std::vector<unsigned char>& childElementSides) const
871  {
872 
873  typedef std::pair<typename UG_NS<dim>::Element*,int> ListEntryType;
874 
876 
877  // //////////////////////////////////////////////////////////////////////
878  // Change the input face number from Dune numbering to UG numbering
879  // //////////////////////////////////////////////////////////////////////
880 
881  elementSide = UGGridRenumberer<dim>::facesDUNEtoUG(elementSide, e->type());
882 
883  // ///////////////
884  // init list
885  // ///////////////
886  if (!e->isLeaf() // Get_Sons_of_ElementSide returns GM_FATAL when called for a leaf !?!
887  && e->level() < maxl) {
888 
889  typename UG_NS<dim>::Element* theElement = this->getRealImplementation(*e).target_;
890 
891  int Sons_of_Side = 0;
892  typename UG_NS<dim>::Element* SonList[UG_NS<dim>::MAX_SONS];
893  int SonSides[UG_NS<dim>::MAX_SONS];
894 
895  int rv = Get_Sons_of_ElementSide(theElement,
896  elementSide,
897  &Sons_of_Side,
898  SonList, // the output elements
899  SonSides, // Output element side numbers
900  true, // Element sons are not precomputed
901  true); // ioflag: I have no idea what this is supposed to do
902  if (rv!=0)
903  DUNE_THROW(GridError, "Get_Sons_of_ElementSide returned with error value " << rv);
904 
905  for (int i=0; i<Sons_of_Side; i++)
906  list.push_back(ListEntryType(SonList[i],SonSides[i]));
907 
908  }
909 
910  // //////////////////////////////////////////////////
911  // Traverse and collect all children of the side
912  // //////////////////////////////////////////////////
913 
914  typename SLList<ListEntryType>::iterator f = list.begin();
915  for (; f!=list.end(); ++f) {
916 
917  typename UG_NS<dim>::Element* theElement = f->first;
918  int side = f->second;
919 
920  int Sons_of_Side = 0;
921  typename UG_NS<dim>::Element* SonList[UG_NS<dim>::MAX_SONS];
922  int SonSides[UG_NS<dim>::MAX_SONS];
923 
924  if (UG_NS<dim>::myLevel(theElement) < maxl) {
925 
926  Get_Sons_of_ElementSide(theElement,
927  side, // Input element side number
928  &Sons_of_Side, // Number of topological sons of the element side
929  SonList, // Output elements
930  SonSides, // Output element side numbers
931  true,
932  true);
933 
934  for (int i=0; i<Sons_of_Side; i++)
935  list.push_back(ListEntryType(SonList[i],SonSides[i]));
936 
937  }
938 
939  }
940 
941  // //////////////////////////////
942  // Extract result from stack
943  // //////////////////////////////
944 
945  // Use reserve / push_back since EntityPointer is not default constructable
946  childElements.clear();
947  childElements.reserve( list.size() );
948  childElementSides.resize(list.size());
949 
950  int i=0;
951  for (f = list.begin(); f!=list.end(); ++f, ++i)
952  {
953 
954  // Set element
955  typedef typename Traits::template Codim< 0 >::EntityPointer EntityPointer;
956  childElements.push_back( EntityPointer( UGGridEntityPointer< 0, const UGGrid< dim > >( f->first, this ) ) );
957 
958  int side = f->second;
959 
960  // Dune numbers the faces of several elements differently than UG.
961  // The following switch does the transformation
962  childElementSides[i] = UGGridRenumberer<dim>::facesUGtoDUNE(side, childElements[i]->type());
963 
964  }
965 
966  }
967 #endif
968 
976  void getChildrenOfSubface(const typename Traits::template Codim<0>::Entity & e,
977  int elementSide,
978  int maxl,
979  std::vector<typename Traits::template Codim<0>::Entity>& childElements,
980  std::vector<unsigned char>& childElementSides) const;
981 
987  COPY
988  };
989 
991  enum ClosureType {
995  NONE
996  };
997 
1000  refinementType_ = type;
1001  }
1002 
1005  closureType_ = type;
1006  }
1007 
1014  static void setDefaultHeapSize(unsigned size) {
1015  heapSize_ = size;
1016  }
1017 
1018 #ifdef DOXYGEN
1019 
1028  void setPosition(const typename Traits::template Codim<dim>::EntityPointer& e,
1029  const FieldVector<double, dim>& pos);
1030 
1031 #else
1032 
1033  // This is the actual implementation of the deprecated method. We need this ugly trick of turning
1034  // it into a template to avoid triggering a deprecation warning every time UGGrid is used.
1035  template< typename T >
1036  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.")
1037  typename std::enable_if<
1038  std::is_same<
1039  T,
1040  typename UGGrid<dim>::Traits::template Codim<dim>::EntityPointer
1041  >::value
1042  >::type
1043  setPosition(const T& e,
1044  const FieldVector<double, dim>& pos)
1045 {
1046  typename UG_NS<dim>::Node* target = this->getRealImplementation(*e).target_;
1047 
1048  for (int i=0; i<dim; i++)
1049  target->myvertex->iv.x[i] = pos[i];
1050 }
1051 
1052 #endif
1053 
1057  void setPosition(const typename Traits::template Codim<dim>::Entity& e,
1058  const FieldVector<double, dim>& pos);
1059 
1064  void globalRefine(int n);
1065 
1070  void saveState(const std::string& filename) const;
1071 
1076  void loadState(const std::string& filename);
1077 
1078  private:
1080  typename UG_NS<dim>::MultiGrid* multigrid_;
1081 
1084 
1090  void setIndices(bool setLevelZero,
1091  std::vector<unsigned int>* nodePermutation);
1092 
1093  // Each UGGrid object has a unique name to identify it in the
1094  // UG environment structure
1095  std::string name_;
1096 
1097  // Our set of level indices
1098  std::vector<shared_ptr<UGGridLevelIndexSet<const UGGrid<dim> > > > levelIndexSets_;
1099 
1100  UGGridLeafIndexSet<const UGGrid<dim> > leafIndexSet_;
1101 
1102  // One id set implementation
1103  // Used for both the local and the global UGGrid id sets
1104  UGGridIdSet<const UGGrid<dim> > idSet_;
1105 
1107  RefinementType refinementType_;
1108 
1110  ClosureType closureType_;
1111 
1119  static int numOfUGGrids;
1120 
1126  bool someElementHasBeenMarkedForRefinement_;
1127 
1133  bool someElementHasBeenMarkedForCoarsening_;
1134 
1139  static unsigned int heapSize_;
1140 
1142  std::vector<shared_ptr<BoundarySegment<dim> > > boundarySegments_;
1143 
1149  unsigned int numBoundarySegments_;
1150 
1151  }; // end Class UGGrid
1152 
1153  namespace Capabilities
1154  {
1170  template<int dim>
1171  struct hasEntity< UGGrid<dim>, 0>
1172  {
1173  static const bool v = true;
1174  };
1175 
1179  template<int dim>
1180  struct hasEntity< UGGrid<dim>, dim>
1181  {
1182  static const bool v = true;
1183  };
1184 
1188  template<int dim>
1189  struct DUNE_DEPRECATED_MSG("Capabilities::isParallel will be removed after dune-grid-2.4.") isParallel< UGGrid<dim> >
1190  {
1191 #ifdef ModelP
1192  static const bool DUNE_DEPRECATED_MSG("Capabilities::isParallel will be removed after dune-grid-2.4.") v = true;
1193 #else // !ModelP
1194  static const bool DUNE_DEPRECATED_MSG("Capabilities::isParallel will be removed after dune-grid-2.4.") v = false;
1195 #endif // !ModelP
1196  };
1197 
1201  template<int dim>
1203  {
1204  static const bool v = true;
1205  };
1206 
1210  template<int dim>
1212  {
1213  static const bool v = false;
1214  };
1215 
1216  }
1217 
1218 } // namespace Dune
1219 
1220 #endif // HAVE_UG || DOXYGEN
1221 #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:208
void communicate(DataHandle &dataHandle, InterfaceType iftype, CommunicationDirection dir) const
The communication interface for all codims on the leaf level.
Definition: uggrid.hh:664
~UGGrid()
Destructor.
void postAdapt()
Clean up refinement markers.
size_t numBoundarySegments() const
Return the number of boundary segments.
Definition: uggrid.hh:366
void setRefinementType(RefinementType type)
Sets the type of grid refinement.
Definition: uggrid.hh:999
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:500
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:629
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:690
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:258
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:592
void setClosureType(ClosureType type)
Sets the type of grid refinement closure.
Definition: uggrid.hh:1004
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:348
const Traits::GlobalIdSet & globalIdSet() const
Access to the GlobalIdSet.
Definition: uggrid.hh:373
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition: uggrid.hh:354
RefinementType
The different forms of grid refinement that UG supports.
Definition: uggrid.hh:983
@ COPY
New level consists of the refined elements and the unrefined ones, too.
Definition: uggrid.hh:987
@ LOCAL
New level consists only of the refined elements and the closure.
Definition: uggrid.hh:985
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: uggrid.hh:360
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:495
UG::DOUBLE ctype
The type used to store coordinates.
Definition: uggrid.hh:264
const Traits::LocalIdSet & localIdSet() const
Access to the LocalIdSet.
Definition: uggrid.hh:379
unsigned int ghostSize(int level, int codim) const
Size of the ghost cell layer on a given level.
Definition: uggrid.hh:505
bool adapt()
Triggers the grid refinement process.
const Traits::LevelIndexSet & levelIndexSet(int level) const
Access to the LevelIndexSets.
Definition: uggrid.hh:385
Traits::template Codim< codim >::LeafIterator leafbegin() const
Iterator to first leaf entity of given codim.
Definition: uggrid.hh:302
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:490
bool loadBalance(DataHandle &dataHandle)
Distributes the grid and some data over the available nodes in a distributed machine.
Definition: uggrid.hh:517
Traits::template Codim< codim >::LeafIterator leafend() const
one past the end of the sequence of leaf entities
Definition: uggrid.hh:308
const Traits::LeafIndexSet & leafIndexSet() const
Access to the LeafIndexSet.
Definition: uggrid.hh:393
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:314
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:991
@ GREEN
Standard red/green refinement.
Definition: uggrid.hh:993
@ NONE
No closure, results in nonconforming meshes.
Definition: uggrid.hh:995
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:320
static void setDefaultHeapSize(unsigned size)
Sets the default heap size.
Definition: uggrid.hh:1014
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 15, 22:30, 2024)