Dune Core Modules (2.4.2)

sgrid.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 #ifndef DUNE_SGRID_HH
4 #define DUNE_SGRID_HH
5 
6 #ifndef DUNE_AVOID_SGRID_DEPRE_WARNING_BECAUSE_I_KNOW_WHAT_IM_DOING
7 #warning The SGrid grid manager has been deprecated, and will be removed after the dune-grid-2.4 release. \
8  Please use YaspGrid instead.
9 #endif
10 
11 #include <limits>
12 #include <vector>
13 #include <stack>
14 
15 #include <dune/common/fvector.hh>
16 #include <dune/common/fmatrix.hh>
20 #include <dune/geometry/genericgeometry/topologytypes.hh>
23 #include <dune/grid/common/grid.hh>
24 #include <dune/grid/sgrid/numbering.hh>
27 
33 namespace Dune {
34 
35  //************************************************************************
39  typedef double sgrid_ctype;
40 
41  // globally define the persistent index type
42  const int sgrid_dim_bits = 24; // bits for encoding each dimension
43  const int sgrid_level_bits = 6; // bits for encoding level number
44  const int sgrid_codim_bits = 4; // bits for encoding codimension
45 
46  //************************************************************************
47  // forward declaration of templates
48 
49  template<int dim, int dimworld, class GridImp> class SGeometry;
50  template<int codim, int dim, class GridImp> class SEntity;
51  template<int codim, class GridImp> class SEntityPointer;
52  template<int codim, class GridImp> class SEntitySeed;
53  template<int codim, PartitionIteratorType, class GridImp> class SLevelIterator;
54  template<int dim, int dimworld, class ctype> class SGrid;
55  template<class GridImp> class SIntersection;
56  template<class GridImp> class SIntersectionIterator;
57  template<class GridImp> class SHierarchicIterator;
58 
59  //************************************************************************
90  template<int mydim, int cdim, class GridImp>
91  class SGeometry
92  : public AxisAlignedCubeGeometry<typename GridImp::ctype,mydim,cdim>
93  {
94  public:
96  typedef typename GridImp::ctype ctype;
97 
105  void make (const FieldVector<ctype,cdim>& lower,
107  {
108  if (mydim==0) {
109  // set up base class
111  return;
112  }
113 
114  // construct the upper right corner of the cube geometry
115  FieldVector<ctype, cdim> upper = lower;
116  for (int i=0; i<mydim; i++)
117  upper += A[i];
118 
119  // look for the directions where the cube is actually extended
120  std::bitset<cdim> axes;
121 
122  for (size_t i=0; i<cdim; i++)
123  if ((upper[i] - lower[i]) > 1e-10)
124  axes[i] = true;
125 
126  // set up base class
127  static_cast< AxisAlignedCubeGeometry<ctype,mydim,cdim> & >( *this ) = AxisAlignedCubeGeometry<ctype,mydim,cdim>(lower, upper, axes);
128  }
129 
132  : AxisAlignedCubeGeometry<ctype,mydim,cdim>(FieldVector<ctype,cdim>(0),FieldVector<ctype,cdim>(0)) // anything
133  {}
134  };
135 
136 
137  //************************************************************************
142  template<int codim, int dim, class GridImp, template<int,int,class> class EntityImp>
143  class SEntityBase :
144  public EntityDefaultImplementation<codim,dim,GridImp,EntityImp>
145  {
146  friend class SEntityPointer<codim,GridImp>;
147  friend class SIntersection<GridImp>;
148  enum { dimworld = GridImp::dimensionworld };
149 
150  typedef typename GridImp::Traits::template Codim< codim >::GeometryImpl GeometryImpl;
151 
152  public:
153  typedef typename GridImp::ctype ctype;
154  typedef typename GridImp::template Codim<codim>::Geometry Geometry;
155  typedef typename GridImp::PersistentIndexType PersistentIndexType;
156 
158  int level () const
159  {
160  return l;
161  }
162 
163  bool equals(const SEntityBase& other) const
164  {
165  return (index==other.index)&&(l==other.l)&&(grid==other.grid);
166  }
167 
169  int globalIndex() const;
170 
176  }
177 
180  {
181  static const GeometryType cubeType(GeometryType::cube,dim-codim);
182  return cubeType;
183  }
184 
186  Geometry geometry () const
187  {
188  if (!builtgeometry) makegeometry();
189 
190  // return result
191  return Geometry( geo );
192  }
193 
194  PartitionType partitionType () const { return InteriorEntity; }
195 
197  SEntityBase (GridImp* _grid, int _l, int _index) :
198  grid(_grid),
199  l(_l),
200  index(_index),
201  z(grid->z(l,index,codim)),
202  builtgeometry(false) {}
203 
206  grid(nullptr),
207  l(-1), // marker for invalid entity
208  index(-1), // marker for invalid entity
209  builtgeometry(false) // mark geometry as not built
210  {}
211 
213  SEntityBase ( const SEntityBase& other ) :
214  grid(other.grid),
215  l(other.l),
216  index(other.index),
217  z(other.z),
218  geo(), // do not copy geometry
219  builtgeometry(false) // mark geometry as not built
220  {}
221 
223  void make (GridImp* _grid, int _l, int _id);
224 
226  void make (int _l, int _id);
227 
229  void makegeometry () const;
230 
232  PersistentIndexType persistentIndex () const
233  {
234  return grid->persistentIndex(l, codim, z);
235  }
236 
238  int compressedIndex () const
239  {
240  return index;
241  }
242 
244  int compressedLeafIndex () const
245  {
246  // codim != dim -> there are no copies of entities
247  // maxlevel -> ids are fine
248  if (codim<dim || l==grid->maxLevel())
249  return compressedIndex();
250 
251  // this is a vertex which is not on the finest level
252  // move coordinates up to maxlevel (multiply by 2 for each level
253  array<int,dim> coord;
254  for (int k=0; k<dim; k++)
255  coord[k] = z[k]*(1<<(grid->maxLevel()-l));
256 
257  // compute number with respect to maxLevel
258  return grid->n(grid->maxLevel(),coord);
259  }
260 
262  int subCompressedIndex (int cd, int i) const
263  {
264  DUNE_THROW(NotImplemented,"subIndex for entities with codimension > 0 is not implemented");
265  return -1;
266  }
267 
269  int subCompressedLeafIndex (int cd, int i) const
270  {
271  DUNE_THROW(NotImplemented,"subIndex for entities with codimension > 0 is not implemented");
272  return -1;
273  }
274 
275  protected:
276  // this is how we implement our elements
277  GridImp* grid;
278  int l;
279  int index;
280  array<int,dim> z;
281  mutable GeometryImpl geo;
282  mutable bool builtgeometry;
283  };
284 
285 
291  template<int codim, int dim, class GridImp>
292  class SEntity : public SEntityBase<codim,dim,GridImp,SEntity>
293  {
295  friend class SEntityPointer<codim,GridImp>;
296  friend class SIntersectionIterator<GridImp>;
297  public:
299  SEntity (GridImp* _grid, int _l, int _id) :
300  SEntityBase(_grid,_l,_id) {}
301 
302  SEntity()
303  {}
304  };
305 
332  template<int dim, class GridImp>
333  class SEntity<0,dim,GridImp> : public SEntityBase<0,dim,GridImp,SEntity>
334  {
335  enum { dimworld = GridImp::dimensionworld };
337  using SEntityBase::grid;
338  using SEntityBase::l;
339  using SEntityBase::index;
340  using SEntityBase::z;
341 
342  typedef typename GridImp::Traits::template Codim< 0 >::GeometryImpl GeometryImpl;
343  typedef typename GridImp::Traits::template Codim< 0 >::LocalGeometryImpl LocalGeometryImpl;
344 
345  friend class SEntityPointer<0,GridImp>;
346  friend class SIntersectionIterator<GridImp>;
347  friend class SIntersection<GridImp>;
348 
349  public:
350  typedef typename GridImp::ctype ctype;
351  typedef typename GridImp::template Codim<0>::Geometry Geometry;
352  typedef typename GridImp::template Codim<0>::LocalGeometry LocalGeometry;
353  template <int cd>
354  struct Codim
355  {
356  typedef typename GridImp::template Codim<cd>::EntityPointer EntityPointer;
357  typedef typename GridImp::template Codim<cd>::Entity Entity;
358  };
359  typedef typename GridImp::template Codim<0>::EntityPointer EntityPointer;
360  typedef typename GridImp::template Codim<0>::Entity Entity;
361  typedef typename GridImp::LeafIntersectionIterator IntersectionIterator;
362  typedef typename GridImp::HierarchicIterator HierarchicIterator;
363  typedef typename GridImp::PersistentIndexType PersistentIndexType;
364 
366  friend class SHierarchicIterator<GridImp>;
367 
372  template<int cc> int count () const;
373 
376  unsigned int subEntities (unsigned int codim) const;
377 
382  template<int cc> typename Codim<cc>::Entity subEntity (int i) const;
383 
385  int subCompressedIndex (int codim, int i) const
386  {
387  if (codim==0) return this->compressedIndex();
388  // compute subIndex
389  return (this->grid)->n(this->l, this->grid->subz(this->z,i,codim));
390  }
391 
395  int subCompressedLeafIndex (int codim, int i) const
396  {
397  if (codim==0) return this->compressedLeafIndex();
398 
399  assert(this->l == this->grid->maxLevel());
400  // compute subIndex
401  return (this->grid)->n(this->l, this->grid->subz(this->z,i,codim));
402  }
403 
405  PersistentIndexType subPersistentIndex (int codim, int i) const
406  {
407  if (codim==0) return this->persistentIndex();
408  // compute subId
409  return this->grid->persistentIndex(this->l, codim, this->grid->subz(this->z,i,codim));
410  }
411 
419  IntersectionIterator ibegin () const;
420  IntersectionIterator ileafbegin () const;
421  IntersectionIterator ilevelbegin () const;
423  IntersectionIterator iend () const;
424  IntersectionIterator ileafend () const;
425  IntersectionIterator ilevelend () const;
426 
432  Entity father () const;
433 
435  bool hasFather () const
436  {
437  return (this->level()>0);
438  }
439 
441  bool isLeaf () const
442  {
443  return ( this->grid->maxLevel() == this->level() );
444  }
445 
457  LocalGeometry geometryInFather () const;
458 
465  HierarchicIterator hbegin (int maxLevel) const;
466 
468  HierarchicIterator hend (int maxLevel) const;
469 
470  SEntity()
471  {}
472 
473  // members specific to SEntity
475  SEntity (GridImp* _grid, int _l, int _index) :
476  SEntityBase(_grid,_l,_index),
477  built_father(false)
478  {}
479 
480  SEntity (const SEntity& other ) :
481  SEntityBase(other.grid, other.l, other.index ),
482  built_father(false)
483  {}
484 
486  void make (GridImp* _grid, int _l, int _id)
487  {
488  SEntityBase::make(_grid,_l,_id);
489  built_father = false;
490  }
491 
493  void make (int _l, int _id)
494  {
495  SEntityBase::make(_l,_id);
496  built_father = false;
497  }
498 
499 
500  private:
501 
502  mutable bool built_father;
503  mutable int father_index;
504  mutable LocalGeometryImpl in_father_local;
505  void make_father() const;
506  };
507 
508 
509  //************************************************************************
518  int l;
519  int index;
520  SHierarchicStackElem () : l(-1), index(-1) {}
521  SHierarchicStackElem (int _l, int _index) {l=_l; index=_index;}
522  bool operator== (const SHierarchicStackElem& s) const {return !operator!=(s);}
523  bool operator!= (const SHierarchicStackElem& s) const {return l!=s.l || index!=s.index;}
524  };
525 
526  template<class GridImp>
527  class SHierarchicIterator :
528  public Dune::SEntityPointer <0,GridImp>
529  {
530  friend class SHierarchicIterator<const GridImp>;
531  enum { dim = GridImp::dimension };
532  enum { dimworld = GridImp::dimensionworld };
534  using SEntityPointer::realEntity;
535  using SEntityPointer::grid;
536  using SEntityPointer::l;
537  using SEntityPointer::index;
538  public:
539  typedef typename GridImp::template Codim<0>::Entity Entity;
540  typedef typename GridImp::ctype ctype;
541 
543  void increment();
544 
551  SHierarchicIterator (GridImp* _grid,
553  int _maxLevel, bool makeend) :
554  SEntityPointer(_grid,_e.level(),_e.compressedIndex())
555  {
556  // without sons, we are done
557  // (the end iterator is equal to the calling iterator)
558  if (makeend) return;
559 
560  // remember element where begin has been called
561  orig_l = this->entity().level();
562  orig_index = _grid->getRealImplementation(this->entity()).compressedIndex();
563 
564  // push original element on stack
565  SHierarchicStackElem originalElement(orig_l, orig_index);
566  stack.push(originalElement);
567 
568  // compute maxLevel
569  maxLevel = std::min(_maxLevel,this->grid->maxLevel());
570 
571  // ok, push all the sons as well
572  push_sons(orig_l,orig_index);
573 
574  // and pop the first son
575  increment();
576  }
577 
578  private:
579  int maxLevel;
580  int orig_l, orig_index;
581 
583  std::stack<SHierarchicStackElem, Dune::ReservedVector<SHierarchicStackElem,GridImp::MAXL> > stack;
584 
585  void push_sons (int level, int fatherid);
586  };
587 
588  //************************************************************************
595  template<class GridImp>
597  {
598  enum { dim=GridImp::dimension };
599  enum { dimworld=GridImp::dimensionworld };
600 
601 
602  friend class SIntersection<GridImp>;
603 
604  public:
605  typedef typename GridImp::template Codim<0>::Entity Entity;
606  typedef typename GridImp::template Codim<0>::EntityPointer EntityPointer;
607  typedef typename GridImp::template Codim<1>::Geometry Geometry;
608  typedef typename GridImp::template Codim<1>::LocalGeometry LocalGeometry;
609  typedef Dune::SIntersection<GridImp> IntersectionImp;
612  enum { dimension=dim };
614  enum { dimensionworld=dimworld };
616  typedef typename GridImp::ctype ctype;
617 
620  {
621  return intersection == i.intersection;
622  }
623 
625  void increment();
626 
628  const Intersection & dereference() const
629  {
630  return intersection;
631  }
632 
634  {}
635 
636  SIntersectionIterator (GridImp* _grid, const SEntity<0,dim,GridImp >* _self, int _count) :
637  intersection(IntersectionImp(_grid,_self,_count))
638  {}
639 
640  private:
641  Intersection intersection;
642  };
643 
644  template<class GridImp>
645  class SIntersection
646  {
647  enum { dim=GridImp::dimension };
648  enum { dimworld=GridImp::dimensionworld };
649 
650  friend class SIntersectionIterator<GridImp>;
651 
652  typedef typename GridImp::Traits::template Codim< 1 >::GeometryImpl GeometryImpl;
653  typedef typename GridImp::Traits::template Codim< 1 >::LocalGeometryImpl LocalGeometryImpl;
654 
655  public:
656  typedef typename GridImp::template Codim<0>::Entity Entity;
657  typedef typename GridImp::template Codim<0>::EntityPointer EntityPointer;
658  typedef typename GridImp::template Codim<1>::Geometry Geometry;
659  typedef typename Geometry::LocalCoordinate LocalCoordinate;
660  typedef typename Geometry::GlobalCoordinate GlobalCoordinate;
661  typedef typename GridImp::template Codim<1>::LocalGeometry LocalGeometry;
664  enum { dimension=dim };
666  enum { dimensionworld=dimworld };
668  typedef typename GridImp::ctype ctype;
669 
670  bool equals(const SIntersection& other) const;
671 
674  Entity inside() const;
675 
678  Entity outside() const;
679 
681  bool boundary () const;
682 
683  int boundaryId () const {
684  if (boundary()) return count + 1;
685  return 0;
686  }
687 
688  int boundarySegmentIndex () const {
689  if (boundary())
690  return grid->boundarySegmentIndex(self.level(), count, zred);
691  return -1;
692  }
693 
695  bool neighbor () const;
696 
700  LocalGeometry geometryInInside () const;
704  LocalGeometry geometryInOutside () const;
708  Geometry geometry () const;
709 
711  GeometryType type () const
712  {
713  return GeometryType(GeometryType::cube,dim-1);
714  }
715 
717  int indexInInside () const;
719  int indexInOutside () const;
720 
721  SIntersection()
722  : grid(nullptr)
723  , partition(-1) // marker for invalid intersection
724  , count(-1) // marker for invalid intersection
725  , valid_count(false)
726  , valid_nb(false)
727  , is_on_boundary(false)
728  , built_intersections(false)
729  {}
730 
732  SIntersection (GridImp* _grid, const SEntity<0,dim,GridImp >* _self, int _count) :
733  self(*_self), ne(self), grid(_grid),
734  partition(_grid->partition(grid->getRealImplementation(ne).level(),_self->z)),
735  zred(_grid->compress(grid->getRealImplementation(ne).level(),_self->z))
736  {
737  // make neighbor
738  make(_count);
739  }
740 
741  SIntersection (const SIntersection & other) :
742  self(other.self), ne(other.ne), grid(other.grid),
743  partition(other.partition), zred(other.zred),
744  count(other.count), valid_count(other.valid_count),
745  valid_nb(other.valid_nb), is_on_boundary(other.is_on_boundary),
746  built_intersections(false)
747  {}
748 
750  SIntersection& operator = (const SIntersection& other)
751  {
752  /* Assign data from other */
753  self = other.self;
754  ne = other.ne;
755  grid = other.grid;
756  partition = other.partition;
757  zred = other.zred;
758  count = other.count;
759  valid_count = other.valid_count;
760  valid_nb = other.valid_nb;
761  is_on_boundary = other.is_on_boundary;
762 
763  /* mark cached data as invalid */
764  built_intersections = false;
765 
766  return *this;
767  }
768 
770  bool conforming () const
771  {
772  return true;
773  }
774 
776  GlobalCoordinate outerNormal (const LocalCoordinate& local) const
777  {
778  return centerUnitOuterNormal();
779  }
780 
782  GlobalCoordinate integrationOuterNormal (const LocalCoordinate& local) const
783  {
784  FieldVector<ctype, dimworld> n = centerUnitOuterNormal();
785  n *= geometry().integrationElement(local);
786  return n;
787  }
788 
790  GlobalCoordinate unitOuterNormal (const LocalCoordinate& local) const
791  {
792  return centerUnitOuterNormal();
793  }
794 
796  GlobalCoordinate centerUnitOuterNormal () const
797  {
798  FieldVector<ctype, dimworld> normal(0.0);
799  normal[count/2] = (count%2) ? 1.0 : -1.0;
800  return normal;
801  }
802 
803  private:
804 
805  void make (int _count) const;
806  void makeintersections () const;
807  Entity self;
808  mutable Entity ne;
809  const GridImp * grid;
810  int partition;
811  array<int,dim> zred;
812  mutable int count;
813  mutable bool valid_count;
814  mutable bool valid_nb;
815  mutable bool is_on_boundary;
816  mutable bool built_intersections;
817  mutable LocalGeometryImpl is_self_local;
818  mutable GeometryImpl is_global;
819  mutable LocalGeometryImpl is_nb_local;
820 
821  };
822 
823  //************************************************************************
824 
828  template <class T>
829  class AutoPtrStack : public std::stack<T*>
830  {
831  public:
832  ~AutoPtrStack()
833  {
834  while(! this->empty())
835  {
836  T* e = this->top();
837  delete e;
838  this->pop();
839  }
840  }
841  };
842 
845  template<int codim, class GridImp>
847  {
848  enum { dim = GridImp::dimension };
849  friend class SIntersectionIterator<GridImp>;
850  public:
852  typedef typename GridImp::template Codim<codim>::Entity Entity;
854  enum { codimension = codim };
855 
857  bool equals(const SEntityPointer<codim,GridImp>& i) const;
859  Entity& dereference() const;
861  int level () const;
862 
864  SEntityPointer (GridImp * _grid, int _l, int _index) :
865  grid(_grid), l(_l), index(_index),
866  e(0)
867  {}
868 
871  grid(_e.grid), l(_e.l), index(_e.index),
872  e(0)
873  {}
874 
877  grid(other.grid), l(other.l), index(other.index),
878  e( 0 )
879  {}
880 
883  {
884  if( e )
885  enStack().push( e );
886 #ifndef NDEBUG
887  index = -1;
888 #endif
889  }
890 
893  {
894  grid = other.grid;
895  l = other.l;
896  index = other.index;
897 
898  // free current entity
899  if( e )
900  enStack().push( e );
901  e = 0;
902 
903  return *this;
904  }
905 
906  protected:
907  inline SEntity<codim,dim,GridImp>& realEntity() const
908  {
909  return grid->getRealImplementation(entity());
910  }
911 
912  inline Entity& entity() const
913  {
914  if( ! e )
915  {
916  e = getEntity( grid, l, index );
917  }
918  return *e;
919  }
920 
921  typedef AutoPtrStack< Entity > EntityStackType;
922  static inline EntityStackType& enStack()
923  {
924  static EntityStackType eStack;
925  return eStack;
926  }
927 
928  inline Entity* getEntity(GridImp* _grid, int _l, int _id ) const
929  {
930  // get stack reference
931  EntityStackType& enSt = enStack();
932 
933  if( enSt.empty() )
934  {
935  return (new Entity(SEntity<codim,dim,GridImp>(_grid, _l, _id)));
936  }
937  else
938  {
939  Entity* e = enSt.top();
940  enSt.pop();
941  grid->getRealImplementation(*e).make(_grid, _l,_id);
942  return e;
943  }
944  }
945 
946  GridImp* grid;
947  int l;
948  mutable int index;
949  mutable Entity* e;
950  };
951 
954  template<int codim, class GridImp>
956  {
957  enum { dim = GridImp::dimension };
958  public:
959  enum { codimension = codim };
960 
963  _l(-1), _index(0)
964  {}
965 
967  SEntitySeed (int l, int index) :
968  _l(l), _index(index)
969  {}
970 
972  bool isValid() const
973  {
974  return _l != -1;
975  }
976 
977  int level () const { return this->_l; }
978  int index () const { return this->_index; }
979 
980  private:
981  int _l;
982  int _index;
983  };
984 
985  //************************************************************************
986 
987 
990  template<int codim, PartitionIteratorType pitype, class GridImp>
992  public Dune::SEntityPointer <codim,GridImp>
993  {
994  friend class SLevelIterator<codim, pitype,const GridImp>;
995  enum { dim = GridImp::dimension };
997  using SEntityPointer::realEntity;
998  using SEntityPointer::l;
999  using SEntityPointer::index;
1000  public:
1001  typedef typename GridImp::template Codim<codim>::Entity Entity;
1002 
1004  void increment();
1005 
1007  SLevelIterator (GridImp * _grid, int _l, int _id) :
1008  SEntityPointer(_grid,_l,_id) {}
1009  };
1010 
1011 
1012  //========================================================================
1017  //========================================================================
1018 
1019  template<class GridImp>
1020  class SGridLevelIndexSet : public IndexSet<GridImp,SGridLevelIndexSet<GridImp> >
1021  {
1024 
1025  enum { dim = GridImp::dimension };
1026 
1027  public:
1028 
1030  SGridLevelIndexSet ( const GridImp &g, int l )
1031  : grid( g ),
1032  level( l )
1033  {
1034  // TODO move list of geometrytypes to grid, can be computed static (singleton)
1035  // contains a single element type;
1036  for (int codim=0; codim<=GridImp::dimension; codim++)
1037  mytypes[codim].push_back(GeometryType(GeometryType::cube,GridImp::dimension-codim));
1038  }
1039 
1041  template<int cd>
1042  int index (const typename GridImp::Traits::template Codim<cd>::Entity& e) const
1043  {
1044  return grid.getRealImplementation(e).compressedIndex();
1045  }
1046 
1047  template< int cc >
1048  int subIndex ( const typename GridImp::Traits::template Codim< cc >::Entity &e,
1049  int i, unsigned int codim ) const
1050  {
1051  if( cc == 0 )
1052  return grid.getRealImplementation(e).subCompressedIndex(codim, i);
1053  else
1054  DUNE_THROW( NotImplemented, "subIndex for higher codimension entity not implemented for SGrid." );
1055  }
1056 
1057  // return true if the given entity is contained in \f$E\f$.
1058  template< class EntityType >
1059  bool contains ( const EntityType &e ) const
1060  {
1061  return (e.level() == level);
1062  }
1063 
1065  int size (GeometryType type) const
1066  {
1067  return grid.size( level, type );
1068  }
1069 
1071  int size (int codim) const
1072  {
1073  return grid.size( level, codim );
1074  }
1075 
1077  const std::vector<GeometryType>& geomTypes (int codim) const
1078  {
1079  return mytypes[codim];
1080  }
1081 
1083  const std::vector<GeometryType>& types (int codim) const
1084  {
1085  return geomTypes(codim);
1086  }
1087 
1088  private:
1089  const GridImp& grid;
1090  int level;
1091  std::vector<GeometryType> mytypes[GridImp::dimension+1];
1092  };
1093 
1094 
1095 
1096  //========================================================================
1101  //========================================================================
1102 
1103  template<class GridImp>
1105  public IdSet<GridImp,SGridGlobalIdSet<GridImp>, typename remove_const<GridImp>::type::PersistentIndexType>
1106  /*
1107  We used the remove_const to extract the Type from the mutable class,
1108  because the const class is not instantiated yet.
1109  */
1110  {
1112 
1113  public:
1114 
1116  /*
1117  We use the remove_const to extract the Type from the mutable class,
1118  because the const class is not instantiated yet.
1119  */
1120  typedef typename remove_const<GridImp>::type::PersistentIndexType IdType;
1121 
1123  /*
1124  We use the remove_const to extract the Type from the mutable class,
1125  because the const class is not instantiated yet.
1126  */
1127  template<int cd>
1128  IdType id (const typename remove_const<GridImp>::type::Traits::template Codim<cd>::Entity& e) const
1129  {
1130  return GridImp::getRealImplementation(e).persistentIndex();
1131  }
1132 
1134  /*
1135  We use the remove_const to extract the Type from the mutable class,
1136  because the const class is not instantiated yet.
1137  */
1138  IdType subId ( const typename remove_const< GridImp >::type::Traits::template Codim< 0 >::Entity &e,
1139  int i, unsigned int codim ) const
1140  {
1141  return GridImp::getRealImplementation(e).subPersistentIndex(codim, i);
1142  }
1143  };
1144 
1145 
1146  template<int dim, int dimworld, class ctype>
1147  struct SGridFamily
1148  {
1149  typedef GridTraits<dim,dimworld,Dune::SGrid<dim,dimworld,ctype>,
1150  SGeometry,SEntity,
1151  SEntityPointer,SLevelIterator,
1152  SIntersection, // leaf intersection
1153  SIntersection, // level intersection
1154  SIntersectionIterator, // leaf intersection iter
1155  SIntersectionIterator, // level intersection iter
1156  SHierarchicIterator,
1157  SLevelIterator,
1158  SGridLevelIndexSet<const SGrid<dim,dimworld,ctype> >,
1159  SGridLevelIndexSet<const SGrid<dim,dimworld,ctype> >,
1160  SGridGlobalIdSet<const SGrid<dim,dimworld,ctype> >,
1161  bigunsignedint<dim*sgrid_dim_bits+sgrid_level_bits+sgrid_codim_bits>,
1162  SGridGlobalIdSet<const SGrid<dim,dimworld,ctype> >,
1163  bigunsignedint<dim*sgrid_dim_bits+sgrid_level_bits+sgrid_codim_bits>,
1164  CollectiveCommunication<Dune::SGrid<dim,dimworld,ctype> >,
1165  DefaultLevelGridViewTraits, DefaultLeafGridViewTraits,
1166  SEntitySeed>
1167  Traits;
1168  };
1169 
1170 
1171  //************************************************************************
1226  template<int dim, int dimworld, typename _ctype = sgrid_ctype>
1227  class SGrid : public GridDefaultImplementation <dim,dimworld,_ctype,SGridFamily<dim,dimworld,_ctype> >
1228  {
1229  public:
1230  typedef SGridFamily<dim,dimworld,_ctype> GridFamily;
1232 
1233  // need for friend declarations in entity
1237 
1239 
1241  enum { MAXL=32 };
1242 
1244  typedef _ctype ctype;
1245 
1246  // constructors
1247 
1255  SGrid (const int * const N_, const ctype * const H_);
1256 
1264  SGrid (const int * const N_, const ctype * const L_, const ctype * const H_);
1265 
1276 
1278  SGrid ();
1279 
1282 
1285  int maxLevel() const;
1286 
1288  template<int cd, PartitionIteratorType pitype>
1289  typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lbegin (int level) const;
1290 
1292  template<int cd, PartitionIteratorType pitype>
1293  typename Traits::template Codim<cd>::template Partition<pitype>::LevelIterator lend (int level) const;
1294 
1296  template<int cd>
1297  typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lbegin (int level) const
1298  {
1299  return lbegin<cd,All_Partition>(level);
1300  }
1301 
1303  template<int cd>
1304  typename Traits::template Codim<cd>::template Partition<All_Partition>::LevelIterator lend (int level) const
1305  {
1306  return lend<cd,All_Partition>(level);
1307  }
1308 
1310  template<int cd, PartitionIteratorType pitype>
1311  typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator leafbegin () const;
1312 
1314  template<int cd, PartitionIteratorType pitype>
1315  typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator leafend () const;
1316 
1318  template<int cd>
1319  typename Traits::template Codim<cd>::template Partition<All_Partition>::LeafIterator leafbegin () const
1320  {
1321  return leafbegin<cd,All_Partition>();
1322  }
1323 
1325  template<int cd>
1326  typename Traits::template Codim<cd>::template Partition<All_Partition>::LeafIterator leafend () const
1327  {
1328  return leafend<cd,All_Partition>();
1329  }
1330 
1331  // \brief obtain EntityPointer from EntitySeed. */
1332  template <typename Seed>
1333  typename Traits::template Codim<Seed::codimension>::EntityPointer
1334  entityPointer(const Seed& seed) const
1335  {
1336  enum { codim = Seed::codimension };
1337  return SEntityPointer<codim,const SGrid<dim,dimworld> >(this,
1338  this->getRealImplementation(seed).level(),
1339  this->getRealImplementation(seed).index());
1340  }
1341 
1342  // \brief obtain Entity from EntitySeed. */
1343  template <typename Seed>
1344  typename Traits::template Codim<Seed::codimension>::Entity
1345  entity(const Seed& seed) const
1346  {
1347  enum { codim = Seed::codimension };
1348  return typename Traits::template Codim<Seed::codimension>::Entity(
1349  SEntity<codim,dim,const SGrid<dim,dimworld> >(this,
1350  this->getRealImplementation(seed).level(),
1351  this->getRealImplementation(seed).index())
1352  );
1353  }
1354 
1368  template<class T, template<class> class P, int codim>
1369  void communicate (T& t, InterfaceType iftype, CommunicationDirection dir, int level)
1370  {
1371  // SGrid is sequential and has no periodic boundaries, so do nothing ...
1372  return;
1373  }
1374 
1376  int size (int level, int codim) const;
1377 
1379  int size (int codim) const
1380  {
1381  return size(maxLevel(),codim);
1382  }
1383 
1385  int size (int level, GeometryType type) const
1386  {
1387  return (type.isCube()) ? size(level,dim-type.dim()) : 0;
1388  }
1389 
1391  int size (GeometryType type) const
1392  {
1393  return size(maxLevel(),type);
1394  }
1395 
1397  size_t numBoundarySegments () const
1398  {
1399  return boundarysize;
1400  }
1401 
1403  int global_size (int codim) const;
1404 
1406  int overlapSize (int level, int codim)
1407  {
1408  return 0;
1409  }
1410 
1412  int ghostSize (int level, int codim)
1413  {
1414  return 0;
1415  }
1416 
1417  // these are all members specific to sgrid
1418 
1420  void globalRefine (int refCount);
1421 
1423  const array<int, dim>& dims(int level) const {
1424  return N[level];
1425  }
1426 
1429  return low;
1430  }
1431 
1434  return H;
1435  }
1436 
1438  bool adapt ()
1439  {
1440  globalRefine(1);
1441  return true;
1442  }
1443 
1444  // The new index sets from DDM 11.07.2005
1445  const typename Traits::GlobalIdSet& globalIdSet() const
1446  {
1447  return theglobalidset;
1448  }
1449 
1450  const typename Traits::LocalIdSet& localIdSet() const
1451  {
1452  return theglobalidset;
1453  }
1454 
1455  const typename Traits::LevelIndexSet& levelIndexSet(int level) const
1456  {
1457  assert(level>=0 && level<=maxLevel());
1458  return *(indexsets[level]);
1459  }
1460 
1461  const typename Traits::LeafIndexSet& leafIndexSet() const
1462  {
1463  return *indexsets.back();
1464  }
1465 
1470  template<class DataHandle>
1471  void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir, int level) const
1472  {}
1473 
1474  template<class DataHandle>
1475  void communicate (DataHandle& data, InterfaceType iftype, CommunicationDirection dir) const
1476  {}
1477 
1478  const CollectiveCommunication<SGrid>& comm () const
1479  {
1480  return ccobj;
1481  }
1482 
1484  int overlapSize (int level, int codim) const
1485  {
1486  return 0;
1487  }
1488 
1490  int overlapSize (int codim) const
1491  {
1492  return 0;
1493  }
1494 
1496  int ghostSize (int level, int codim) const
1497  {
1498  return 0;
1499  }
1500 
1502  int ghostSize (int codim) const
1503  {
1504  return 0;
1505  }
1506 
1507  /*
1508  @}
1509  */
1510 
1511  private:
1512  /*
1513  Make associated classes friends to grant access to the real entity
1514  */
1515  friend class Dune::SGridLevelIndexSet<Dune::SGrid<dim,dimworld> >;
1516  friend class Dune::SGridGlobalIdSet<Dune::SGrid<dim,dimworld> >;
1517  friend class Dune::SIntersection<const Dune::SGrid<dim,dimworld> >;
1518  friend class Dune::SHierarchicIterator<Dune::SGrid<dim,dimworld> >;
1519  friend class Dune::SEntity<0,dim,Dune::SGrid<dim,dimworld> >;
1520 
1521  friend class Dune::SGridLevelIndexSet<const Dune::SGrid<dim,dimworld> >;
1522  friend class Dune::SGridGlobalIdSet<const Dune::SGrid<dim,dimworld> >;
1523  friend class Dune::SIntersectionIterator<const Dune::SGrid<dim,dimworld> >;
1524  friend class Dune::SHierarchicIterator<const Dune::SGrid<dim,dimworld> >;
1525  friend class Dune::SEntity<0,dim,const Dune::SGrid<dim,dimworld> >;
1526 
1527  template<int codim_, int dim_, class GridImp_, template<int,int,class> class EntityImp_>
1528  friend class Dune::SEntityBase;
1529 
1530  template<int codim_, class GridImp_>
1531  friend class Dune::SEntityPointer;
1532 
1533  template<int codim_, int dim_, class GridImp_, template<int,int,class> class EntityImp_>
1534  friend class Entity;
1535 
1537  FieldVector<ctype, dimworld> pos (int level, array<int,dim>& z) const;
1538 
1540  int calc_codim (int level, const array<int,dim>& z) const;
1541 
1543  int n (int level, const array<int,dim>& z) const;
1544 
1546  array<int,dim> z (int level, int i, int codim) const;
1547 
1549  array<int,dim> subz (const array<int,dim> & z, int i, int codim) const;
1550 
1552  array<int,dim> compress (int level, const array<int,dim>& z) const;
1553 
1555  array<int,dim> expand (int level, const array<int,dim>& r, int b) const;
1556 
1560  int partition (int level, const array<int,dim>& z) const;
1561 
1563  bool exists (int level, const array<int,dim>& zred) const;
1564 
1565  // compute boundary segment index for a given zentity and a face
1566  int boundarySegmentIndex (int l, int face, const array<int,dim> & zentity) const
1567  {
1568  array<int,dim-1> zface;
1569  int dir = face/2;
1570  int side = face%2;
1571  // compute z inside the global face
1572  for (int i=0; i<dir; i++) zface[i] = zentity[i]/(1<<l);
1573  for (int i=dir+1; i<dim; i++) zface[i-1] = zentity[i]/(1<<l);
1574  zface = boundarymapper[dir].expand(zface, 0);
1575  // compute index in the face
1576  int index = boundarymapper[dir].n(zface);
1577  // compute offset
1578  for (int i=0; i<dir; i++)
1579  index += 2*boundarymapper[i].elements(0);
1580  index += side*boundarymapper[dir].elements(0);
1581  return index;
1582  }
1583 
1584  // compute persistent index for a given zentity
1585  PersistentIndexType persistentIndex (int l, int codim, const array<int,dim> & zentity) const
1586  {
1587  if (codim!=dim)
1588  {
1589  // encode codim, this would actually not be necessary
1590  // because z is unique in codim
1591  PersistentIndexType id(codim);
1592 
1593  // encode level
1594  id = id << sgrid_level_bits;
1595  id = id+PersistentIndexType(l);
1596 
1597  // encode coordinates
1598  for (int i=dim-1; i>=0; i--)
1599  {
1600  id = id << sgrid_dim_bits;
1601  id = id+PersistentIndexType(zentity[i]);
1602  }
1603 
1604  return id;
1605  }
1606  else
1607  {
1608  // determine min number of trailing zeroes
1609  // consider that z is on the doubled grid !
1610  int trailing = 1000;
1611  for (int i=0; i<dim; i++)
1612  {
1613  // count trailing zeros
1614  int zeros = 0;
1615  for (int j=0; j<l; j++)
1616  if (zentity[i]&(1<<(j+1)))
1617  break;
1618  else
1619  zeros++;
1620  trailing = std::min(trailing,zeros);
1621  }
1622 
1623  // determine the level of this vertex
1624  int level = l-trailing;
1625 
1626  // encode codim
1627  PersistentIndexType id(dim);
1628 
1629  // encode level
1630  id = id << sgrid_level_bits;
1631  id = id+PersistentIndexType(level);
1632 
1633  // encode coordinates
1634  for (int i=dim-1; i>=0; i--)
1635  {
1636  id = id << sgrid_dim_bits;
1637  id = id+PersistentIndexType(zentity[i]>>trailing);
1638  }
1639 
1640  return id;
1641  }
1642  }
1643 
1644  // disable copy and assign
1645  SGrid(const SGrid &) {}
1646  SGrid & operator = (const SGrid &) { return *this; }
1647  // generate SGrid
1648  void makeSGrid (const array<int,dim>& N_, const FieldVector<ctype, dim>& L_, const FieldVector<ctype, dim>& H_);
1649 
1650  /*
1651  internal data
1652  */
1653  CollectiveCommunication<SGrid> ccobj;
1654 
1655  ReservedVector<SGridLevelIndexSet<const SGrid<dim,dimworld> >*, MAXL> indexsets;
1656  SGridGlobalIdSet<const SGrid<dim,dimworld> > theglobalidset;
1657 
1658  int L; // number of levels in hierarchic mesh 0<=level<L
1659  FieldVector<ctype, dim> low; // lower left corner of the grid
1660  FieldVector<ctype, dim> H; // length of cube per direction
1661  std::vector<array<int,dim> > N; // number of elements per direction for each level
1662  std::vector<FieldVector<ctype, dim> > h; // mesh size per direction for each level
1663  mutable CubeMapper<dim> *mapper; // a mapper for each level
1664 
1665  // boundary segement index set
1666  array<CubeMapper<dim-1>, dim> boundarymapper; // a mapper for each coarse grid face
1667  int boundarysize;
1668  };
1669 
1670  namespace Capabilities
1671  {
1672 
1684  template<int dim, int dimw>
1685  struct hasSingleGeometryType< SGrid<dim,dimw> >
1686  {
1687  static const bool v = true;
1688  static const unsigned int topologyId = GenericGeometry :: CubeTopology< dim > :: type :: id ;
1689  };
1690 
1694  template<int dim, int dimw>
1695  struct isCartesian< SGrid<dim,dimw> >
1696  {
1697  static const bool v = true;
1698  };
1699 
1703  template<int dim, int dimw, int cdim>
1704  struct hasEntity< SGrid<dim,dimw>, cdim>
1705  {
1706  static const bool v = true;
1707  };
1708 
1712  template<int dim, int dimw>
1713  struct isLevelwiseConforming< SGrid<dim,dimw> >
1714  {
1715  static const bool v = true;
1716  };
1717 
1721  template<int dim, int dimw>
1722  struct isLeafwiseConforming< SGrid<dim,dimw> >
1723  {
1724  static const bool v = true;
1725  };
1726 
1727  } // end namespace Capabilities
1728 
1735  template<int dim>
1736  class StructuredGridFactory<SGrid<dim, dim> > {
1737  typedef SGrid<dim, dim> GridType;
1738  typedef typename GridType::ctype ctype;
1739  static const int dimworld = GridType::dimensionworld;
1740 
1741  public:
1748  static shared_ptr<GridType>
1750  const FieldVector<ctype,dimworld>& upperRight,
1751  const array<unsigned int,dim>& elements)
1752  {
1753  FieldVector<int, dim> elements_;
1754  std::copy(elements.begin(), elements.end(), elements_.begin());
1755 
1756  return shared_ptr<GridType>
1757  (new GridType(elements_, lowerLeft, upperRight));
1758  }
1759 
1769  static shared_ptr<GridType>
1771  const FieldVector<ctype,dimworld>& upperRight,
1772  const array<unsigned int,dim>& elements)
1773  {
1774  DUNE_THROW(GridError, className<StructuredGridFactory>()
1775  << "::createSimplexGrid(): Simplices are not supported "
1776  "by SGrid.");
1777  }
1778  };
1779 
1780 } // end namespace Dune
1781 
1782 #include "sgrid/sgrid.cc"
1783 
1784 #endif
A geometry implementation for axis-aligned hypercubes.
Portable very large unsigned integers.
a stack of pointers with auto destruction if the stack is destructed
Definition: sgrid.hh:830
A geometry implementation for axis-aligned hypercubes.
Definition: axisalignedcubegeometry.hh:50
Iterator begin()
begin iterator
Definition: densevector.hh:307
Default Implementations for EntityImp.
Definition: entity.hh:740
GridImp::template Codim< cd >::EntityPointer EntityPointer
The corresponding entity seed (for storage of entities)
Definition: entity.hh:755
Wrapper class for entities.
Definition: entity.hh:62
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:25
unsigned int dim() const
Return dimension of the type.
Definition: type.hh:321
bool isCube() const
Return true if entity is a cube of any dimension.
Definition: type.hh:311
@ cube
Cube element in any nonnegative dimension.
Definition: type.hh:31
FieldVector< ctype, cdim > GlobalCoordinate
type of the global coordinates
Definition: geometry.hh:104
FieldVector< ctype, mydim > LocalCoordinate
type of local coordinates
Definition: geometry.hh:101
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
Base class for exceptions in Dune grid modules.
Definition: exceptions.hh:18
Id Set Interface.
Definition: indexidset.hh:414
Index Set Interface base class.
Definition: indexidset.hh:76
Intersection of a mesh entities of codimension 0 ("elements") with a "neighboring" element or with th...
Definition: intersection.hh:161
Default exception for dummy implementations.
Definition: exceptions.hh:288
Definition: sgrid.hh:145
int level() const
level of this element
Definition: sgrid.hh:158
void make(int _l, int _id)
Reinitialization.
Geometry geometry() const
geometry of this entity
Definition: sgrid.hh:186
GridImp * grid
grid containes mapper, geometry, etc.
Definition: sgrid.hh:277
int compressedIndex() const
consecutive, codim-wise, level-wise index
Definition: sgrid.hh:238
void make(GridImp *_grid, int _l, int _id)
Reinitialization.
void makegeometry() const
geometry of this entity
int index
my consecutive index
Definition: sgrid.hh:279
int subCompressedLeafIndex(int cd, int i) const
subentity compressed leaf index (not available here)
Definition: sgrid.hh:269
array< int, dim > z
my coordinate, number of even components = codim
Definition: sgrid.hh:280
GeometryImpl geo
geometry, is only built on demand
Definition: sgrid.hh:281
int l
level where element is on
Definition: sgrid.hh:278
PersistentIndexType persistentIndex() const
globally unique, persistent index
Definition: sgrid.hh:232
SEntitySeed< codim, GridImp > seed() const
Return the entity seed which contains sufficient information to generate the entity again and uses as...
Definition: sgrid.hh:174
int globalIndex() const
global index is calculated from the index and grid size
int compressedLeafIndex() const
consecutive, codim-wise, level-wise index
Definition: sgrid.hh:244
GeometryType type() const
return the element type identifier
Definition: sgrid.hh:179
SEntityBase(const SEntityBase &other)
copy constructor
Definition: sgrid.hh:213
SEntityBase(GridImp *_grid, int _l, int _index)
constructor
Definition: sgrid.hh:197
bool builtgeometry
true if geometry has been constructed
Definition: sgrid.hh:282
SEntityBase()
empty constructor
Definition: sgrid.hh:205
int subCompressedIndex(int cd, int i) const
subentity compressed index (not available here)
Definition: sgrid.hh:262
Definition: sgrid.hh:847
int index
my consecutive index
Definition: sgrid.hh:948
SEntityPointer(const SEntityPointer< codim, GridImp > &other)
constructor
Definition: sgrid.hh:876
bool equals(const SEntityPointer< codim, GridImp > &i) const
equality
SEntityPointer(GridImp *_grid, int _l, int _index)
constructor
Definition: sgrid.hh:864
int l
level where element is on
Definition: sgrid.hh:947
SEntityPointer(const SEntity< codim, dim, GridImp > &_e)
constructor
Definition: sgrid.hh:870
Entity & dereference() const
dereferencing
SEntityPointer & operator=(const SEntityPointer &other)
assignment operator
Definition: sgrid.hh:892
Entity * e
virtual entity
Definition: sgrid.hh:949
GridImp * grid
my grid
Definition: sgrid.hh:946
~SEntityPointer()
destructor pointer
Definition: sgrid.hh:882
int level() const
ask for level of entity
Definition: sgrid.hh:956
bool isValid() const
check whether the EntitySeed refers to a valid Entity
Definition: sgrid.hh:972
SEntitySeed()
default constructor (invalid)
Definition: sgrid.hh:962
SEntitySeed(int l, int index)
constructor
Definition: sgrid.hh:967
bool hasFather() const
returns true if father entity exists
Definition: sgrid.hh:435
Entity father() const
Inter-level access to father element on coarser grid.
Codim< cc >::Entity subEntity(int i) const
bool isLeaf() const
return true if the entity is leaf
Definition: sgrid.hh:441
IntersectionIterator ibegin() const
int subCompressedLeafIndex(int codim, int i) const
Definition: sgrid.hh:395
PersistentIndexType subPersistentIndex(int codim, int i) const
subentity persistent index
Definition: sgrid.hh:405
int subCompressedIndex(int codim, int i) const
subentity compressed index
Definition: sgrid.hh:385
HierarchicIterator hbegin(int maxLevel) const
Inter-level access to son elements on higher levels<=maxLevel.
IntersectionIterator iend() const
Reference to one past the last intersection.
LocalGeometry geometryInFather() const
Location of this element relative to the reference element element of the father.
SEntity(GridImp *_grid, int _l, int _index)
constructor
Definition: sgrid.hh:475
void make(int _l, int _id)
Reinitialization.
Definition: sgrid.hh:493
void make(GridImp *_grid, int _l, int _id)
Reinitialization.
Definition: sgrid.hh:486
HierarchicIterator hend(int maxLevel) const
Returns iterator to one past the last son.
unsigned int subEntities(unsigned int codim) const
Return number of subentities with codimension codim.
Definition: sgrid.hh:293
SEntity(GridImp *_grid, int _l, int _id)
constructor
Definition: sgrid.hh:299
Definition: sgrid.hh:93
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: sgrid.hh:96
SGeometry()
constructor
Definition: sgrid.hh:131
void make(const FieldVector< ctype, cdim > &lower, const FieldMatrix< ctype, mydim, cdim > &A)
Set up the geometry.
Definition: sgrid.hh:105
persistent, globally unique Ids
Definition: sgrid.hh:1110
IdType subId(const typename remove_const< GridImp >::type::Traits::template Codim< 0 >::Entity &e, int i, unsigned int codim) const
get id of subentity
Definition: sgrid.hh:1138
IdType id(const typename remove_const< GridImp >::type::Traits::template Codim< cd >::Entity &e) const
get id of an entity
Definition: sgrid.hh:1128
remove_const< GridImp >::type::PersistentIndexType IdType
define the type used for persistent indices
Definition: sgrid.hh:1120
implementation of index set
Definition: sgrid.hh:1021
const std::vector< GeometryType > & geomTypes(int codim) const
deliver all geometry types used in this grid
Definition: sgrid.hh:1077
int index(const typename GridImp::Traits::template Codim< cd >::Entity &e) const
get index of an entity
Definition: sgrid.hh:1042
int size(int codim) const
return size of set for a given codim
Definition: sgrid.hh:1071
const std::vector< GeometryType > & types(int codim) const
deliver all geometry types used in this grid
Definition: sgrid.hh:1083
SGridLevelIndexSet(const GridImp &g, int l)
constructor stores reference to a grid and level
Definition: sgrid.hh:1030
int size(GeometryType type) const
get number of entities of given type and level (the level is known to the object)
Definition: sgrid.hh:1065
[ provides Dune::Grid ]
Definition: sgrid.hh:1228
Traits::template Codim< cd >::template Partition< All_Partition >::LevelIterator lend(int level) const
one past the end on this level
Definition: sgrid.hh:1304
Traits::template Codim< cd >::template Partition< All_Partition >::LevelIterator lbegin(int level) const
Iterator to first entity of given codim on level.
Definition: sgrid.hh:1297
Traits::template Codim< cd >::template Partition< All_Partition >::LeafIterator leafend() const
return LeafIterator which points behind the last entity
Definition: sgrid.hh:1326
int ghostSize(int level, int codim) const
return size (= distance in graph) of ghost region
Definition: sgrid.hh:1496
int size(int level, int codim) const
number of grid entities per level and codim
int overlapSize(int level, int codim)
return size (= distance in graph) of overlap region
Definition: sgrid.hh:1406
SGrid()
empty constructor making grid of unit square discretized with one cell
bool adapt()
map adapt to global refine
Definition: sgrid.hh:1438
int ghostSize(int level, int codim)
return size (= distance in graph) of ghost region
Definition: sgrid.hh:1412
Traits::template Codim< cd >::template Partition< pitype >::LeafIterator leafbegin() const
return LeafIterator which points to the first entity
int size(GeometryType type) const
number of leaf entities per codim and geometry type in this process
Definition: sgrid.hh:1391
int size(int level, GeometryType type) const
number of entities per level and geometry type in this process
Definition: sgrid.hh:1385
FieldVector< ctype, dimworld > upperRight() const
Get upper right corner.
Definition: sgrid.hh:1433
int overlapSize(int level, int codim) const
return size (= distance in graph) of overlap region
Definition: sgrid.hh:1484
int global_size(int codim) const
number of grid entities of all level for given codim
const FieldVector< ctype, dimworld > & lowerLeft() const
Get lower left corner.
Definition: sgrid.hh:1428
int overlapSize(int codim) const
return size (= distance in graph) of overlap region
Definition: sgrid.hh:1490
_ctype ctype
define type used for coordinates in grid module
Definition: sgrid.hh:1244
SGrid(const int *const N_, const ctype *const L_, const ctype *const H_)
Make an SGrid from position, extend and number of cells per direction.
const array< int, dim > & dims(int level) const
Get number of elements in each coordinate direction.
Definition: sgrid.hh:1423
SGrid(const int *const N_, const ctype *const H_)
Make an SGrid from extend and number of cells per direction.
int size(int codim) const
number of leaf entities per codim in this process
Definition: sgrid.hh:1379
~SGrid()
SGrid destructor.
Traits::template Codim< cd >::template Partition< pitype >::LevelIterator lend(int level) const
one past the end on this level
Traits::template Codim< cd >::template Partition< All_Partition >::LeafIterator leafbegin() const
return LeafIterator which points to the first entity
Definition: sgrid.hh:1319
void globalRefine(int refCount)
Refine mesh globally by one refCount levels.
int ghostSize(int codim) const
return size (= distance in graph) of ghost region
Definition: sgrid.hh:1502
SGrid(FieldVector< int, dim > N_, FieldVector< ctype, dim > L_, FieldVector< ctype, dim > H_)
Make an SGrid from position, extend and number of cells per direction.
Traits::template Codim< cd >::template Partition< pitype >::LevelIterator lbegin(int level) const
Iterator to first entity of given codim on level.
int maxLevel() const
Traits::template Codim< cd >::template Partition< pitype >::LeafIterator leafend() const
one past the end on the leaf level
void communicate(T &t, InterfaceType iftype, CommunicationDirection dir, int level)
Definition: sgrid.hh:1369
size_t numBoundarySegments() const
returns the number of boundary segments within the macro grid
Definition: sgrid.hh:1397
Definition: sgrid.hh:597
bool equals(const SIntersectionIterator< GridImp > &i) const
equality
Definition: sgrid.hh:619
const Intersection & dereference() const
dereferencing
Definition: sgrid.hh:628
GridImp::ctype ctype
define type used for coordinates in grid module
Definition: sgrid.hh:616
Definition: sgrid.hh:993
SLevelIterator(GridImp *_grid, int _l, int _id)
constructor
Definition: sgrid.hh:1007
void increment()
increment
static shared_ptr< GridType > createSimplexGrid(const FieldVector< ctype, dimworld > &lowerLeft, const FieldVector< ctype, dimworld > &upperRight, const array< unsigned int, dim > &elements)
Create a structured simplex grid.
Definition: sgrid.hh:1770
static shared_ptr< GridType > createCubeGrid(const FieldVector< ctype, dimworld > &lowerLeft, const FieldVector< ctype, dimworld > &upperRight, const array< unsigned int, dim > &elements)
Create a structured cube grid.
Definition: sgrid.hh:1749
Construct structured cube and simplex grids in unstructured grid managers.
Definition: structuredgridfactory.hh:31
Portable very large unsigned integers.
Definition: bigunsignedint.hh:70
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.
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:178
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:243
bool ne(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test for inequality using epsilon
Definition: float_cmp.cc:125
PartitionType
Attributes used in the generic overlap model.
Definition: gridenums.hh:28
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
@ InteriorEntity
all interior entities
Definition: gridenums.hh:29
Provides base classes for index and id sets.
Dune namespace.
Definition: alignment.hh:10
double sgrid_ctype
Definition: sgrid.hh:39
An stl-compliant random-access container which stores everything on the stack.
Specialize with 'true' for all codims that a grid implements entities for. (default=false)
Definition: capabilities.hh:58
Specialize with 'true' for if the codimension 0 entity of the grid has only one possible geometry typ...
Definition: capabilities.hh:27
Specialize with 'true' if the grid is a Cartesian grid. Cartesian grids satisfy the following propert...
Definition: capabilities.hh:48
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
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
IdSet< const GridImp, GlobalIdSetImp, GIDType > GlobalIdSet
The type of the global id set.
Definition: grid.hh:1430
IndexSet< const GridImp, LevelIndexSetImp > LevelIndexSet
The type of the level index set.
Definition: grid.hh:1426
IndexSet< const GridImp, LeafIndexSetImp > LeafIndexSet
The type of the leaf index set.
Definition: grid.hh:1428
IdSet< const GridImp, LocalIdSetImp, LIDType > LocalIdSet
The type of the local id set.
Definition: grid.hh:1432
Definition: sgrid.hh:517
A class to construct structured cube and simplex grids using the grid factory.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 10, 22:30, 2024)