1 #ifndef DUNE_YASPGRID_HH
2 #define DUNE_YASPGRID_HH
17 #include <dune/grid/yaspgrid/grids.hh>
19 #include <dune/common/misc.hh>
20 #include <dune/common/bigunsignedint.hh>
21 #include <dune/common/typetraits.hh>
22 #include <dune/common/collectivecommunication.hh>
23 #include <dune/common/mpihelper.hh>
24 #include <dune/geometry/genericgeometry/topologytypes.hh>
30 #include <dune/common/mpicollectivecommunication.hh>
60 template<
int mydim,
int cdim,
class Gr
idImp>
class YaspGeometry;
61 template<
int codim,
int dim,
class Gr
idImp>
class YaspEntity;
64 template<
int codim, PartitionIteratorType pitype,
class Gr
idImp>
class YaspLevelIterator;
72 namespace FacadeOptions
75 template<
int dim,
int mydim,
int cdim>
78 static const bool v =
false;
81 template<
int dim,
int mydim,
int cdim>
84 static const bool v =
false;
92 template<
int dim,
class Gr
idImp>
94 static const array<YaspGeometry<dim,dim,GridImp>, (1<<dim) >
_geo;
95 static array<YaspGeometry<dim,dim,GridImp>, (1<<dim) >
initSons()
97 array<YaspGeometry<dim,dim,GridImp>, (1<<dim) > geo;
98 FieldVector<yaspgrid_ctype,dim> midpoint(0.25);
99 FieldVector<yaspgrid_ctype,dim> extension(0.5);
100 for (
int i=0; i<(1<<dim); i++)
103 for (
int k=0; k<dim; k++)
114 template<
int dim,
class Gr
idImp>
115 const array<YaspGeometry<dim,dim,GridImp>, (1<<dim)>
116 YaspFatherRelativeLocalElement<dim, GridImp>::_geo =
117 YaspFatherRelativeLocalElement<dim, GridImp>::initSons();
130 template<
int mydim,
int cdim,
class Gr
idImp>
135 typedef typename GridImp::ctype
ctype;
153 FieldVector< ctype, cdim >
corner (
const int i )
const
155 assert( i >= 0 && i < (
int) coord_.N() );
156 FieldVector<ctype, cdim>& c = coord_[i];
158 for (
int k=0; k<cdim; k++)
167 c[k] = midpoint[k]+0.5*extension[k];
169 c[k] = midpoint[k]-0.5*extension[k];
177 FieldVector< ctype, cdim >
center ( )
const
183 FieldVector<ctype, cdim>
global (
const FieldVector<ctype, mydim>&
local)
const
185 FieldVector<ctype, cdim> g;
187 for (
int k=0; k<cdim; k++)
192 g[k] = midpoint[k] + (local[bit]-0.5)*extension[k];
199 FieldVector<ctype, mydim>
local (
const FieldVector<ctype, cdim>&
global)
const
201 FieldVector<ctype, mydim> l;
203 for (
int k=0; k<cdim; k++)
206 l[bit] = (global[k]-midpoint[k])/extension[k] + 0.5;
216 for (
int k=0; k<cdim; k++)
217 if (k!=missing) volume *= extension[k];
233 for (
int i=0; i<cdim; ++i)
237 JT[k][i] = extension[i];
248 for (
int i=0; i<cdim; ++i)
252 Jinv[i][k] = 1.0/extension[i];
264 : midpoint(p), extension(h), missing(m)
267 DUNE_THROW(
GridError,
"general YaspGeometry assumes cdim=mydim+1");
272 : midpoint(other.midpoint),
273 extension(other.extension),
274 missing(other.missing)
281 s <<
"YaspGeometry<"<<mydim<<
","<<cdim<<
"> ";
283 for (
int i=0; i<cdim; i++)
284 s <<
" " << midpoint[i];
286 for (
int i=0; i<cdim; i++)
287 s <<
" " << extension[i];
288 s <<
" missing is " << missing;
299 FieldVector<ctype, cdim> midpoint;
300 FieldVector<ctype, cdim> extension;
305 mutable FieldMatrix<ctype, mydim, cdim> JT;
306 mutable FieldMatrix<ctype, cdim, mydim> Jinv;
307 mutable FieldMatrix<ctype, Power_m_p<2,mydim>::power, cdim> coord_;
314 template<
int mydim,
class Gr
idImp>
318 typedef typename GridImp::ctype
ctype;
336 const FieldVector<ctype, mydim>& operator[] (
int i)
const
342 FieldVector< ctype, mydim >
corner (
const int i )
const
344 assert( i >= 0 && i < (
int) coord_.N() );
345 FieldVector<ctype, mydim>& c = coord_[i];
346 for (
int k=0; k<mydim; k++)
348 c[k] = midpoint[k]+0.5*extension[k];
350 c[k] = midpoint[k]-0.5*extension[k];
355 FieldVector< ctype, mydim >
center ( )
const
361 FieldVector<ctype, mydim>
global (
const FieldVector<ctype, mydim>&
local)
const
363 FieldVector<ctype,mydim> g;
364 for (
int k=0; k<mydim; k++)
365 g[k] = midpoint[k] + (local[k]-0.5)*extension[k];
370 FieldVector<ctype, mydim>
local (
const FieldVector<ctype,mydim>&
global)
const
372 FieldVector<ctype, mydim> l;
373 for (
int k=0; k<mydim; k++)
374 l[k] = (global[k]-midpoint[k])/extension[k] + 0.5;
389 for (
int k=0; k<mydim; k++) vol *= extension[k];
396 for (
int i=0; i<mydim; ++i)
399 JT[i][i] = extension[i];
406 for (
int i=0; i<mydim; ++i)
409 Jinv[i][i] = 1.0/extension[i];
418 YaspGeometry (
const FieldVector<ctype, mydim>& p,
const FieldVector<ctype, mydim>& h)
419 : midpoint(p), extension(h)
424 : midpoint(other.midpoint),
425 extension(other.extension)
432 s <<
"YaspGeometry<"<<mydim<<
","<<mydim<<
"> ";
434 for (
int i=0; i<mydim; i++)
435 s <<
" " << midpoint[i];
437 for (
int i=0; i<mydim; i++)
438 s <<
" " << extension[i];
450 FieldVector<ctype, mydim> midpoint;
451 FieldVector<ctype, mydim> extension;
455 mutable FieldMatrix<ctype, mydim, mydim> Jinv,JT;
456 mutable FieldMatrix<ctype, Power_m_p<2,mydim>::power, mydim> coord_;
460 template<
int cdim,
class Gr
idImp>
464 typedef typename GridImp::ctype
ctype;
482 const FieldVector<ctype, cdim>& operator[] (
int i)
const
488 FieldVector< ctype, cdim >
corner (
const int i )
const
494 FieldVector< ctype, cdim >
center ( )
const
509 static FieldMatrix<ctype,0,cdim> JT(0.0);
515 static FieldMatrix<ctype,cdim,0> Jinv(0.0);
535 s <<
"YaspGeometry<"<<0<<
","<<cdim<<
"> ";
536 s <<
"position " << position;
543 FieldVector<ctype, cdim> position;
547 template <
int mydim,
int cdim,
class Gr
idImp>
549 std::ostream& operator<< (std::ostream& s, YaspGeometry<mydim,cdim,GridImp>& e)
565 template<
int codim,
int dim,
class Gr
idImp>
567 public GridImp::template Codim<codim>
::Entity
570 typedef typename GridImp::ctype
ctype;
572 typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator
YGLI;
573 typedef typename SubYGrid<dim,ctype>::TransformingSubIterator
TSI;
576 GridImp::template Codim<codim>::
Entity (
YaspEntity<codim, dim, GridImp>(yg,g,it))
579 GridImp::template Codim<codim>::
Entity (e)
583 return this->realEntity.transformingsubiterator();
587 return this->realEntity.gridlevel();
591 return this->realEntity.yaspgrid();
595 template<
int codim,
int dim,
class Gr
idImp>
600 typedef typename GridImp::ctype
ctype;
602 typedef typename GridImp::template Codim<codim>::Geometry
Geometry;
606 DUNE_THROW(
GridError,
"YaspEntity not implemented");
612 DUNE_THROW(
GridError,
"YaspEntity not implemented");
618 DUNE_THROW(
GridError,
"YaspEntity not implemented");
624 DUNE_THROW(
GridError,
"YaspEntity not implemented");
629 DUNE_THROW(
GridError,
"YaspEntity not implemented");
632 typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator
YGLI;
633 typedef typename SubYGrid<dim,ctype>::TransformingSubIterator
TSI;
636 DUNE_THROW(
GridError,
"YaspEntity not implemented");
648 DUNE_THROW(
GridError,
"YaspEntity not implemented");
654 DUNE_THROW(
GridError,
"YaspEntity not implemented");
660 DUNE_THROW(
GridError,
"YaspEntity not implemented");
666 template<
int dim,
class Gr
idImp>
670 enum { dimworld = GridImp::dimensionworld };
672 typedef typename GridImp::Traits::template Codim< 0 >::GeometryImpl GeometryImpl;
675 typedef typename GridImp::ctype
ctype;
677 typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator
YGLI;
678 typedef typename SubYGrid<dim,ctype>::TransformingSubIterator
TSI;
680 typedef typename GridImp::template Codim< 0 >::Geometry
Geometry;
690 typedef typename GridImp::template Codim<0>::EntitySeed
EntitySeed;
700 typedef typename YGrid<dim,ctype>::iTupel
iTupel;
704 : _yg(yg), _it(it), _g(g)
709 int level ()
const {
return _g.level(); }
712 int index ()
const {
return _it.superindex(); }
716 return _g.cell_global().index(_it.coord());
729 if (_g.cell_interior().inside(_it.coord()))
return InteriorEntity;
730 if (_g.cell_overlap().inside(_it.coord()))
return OverlapEntity;
731 DUNE_THROW(
GridError,
"Impossible GhostEntity " << _it.coord() <<
"\t"
732 << _g.cell_interior().origin() <<
"/" << _g.cell_interior().size());
739 GeometryImpl _geometry(_it.position(),_it.meshsize());
745 template<
int cc>
int count ()
const
747 if (cc==dim)
return 1<<dim;
748 if (cc==1)
return 2*dim;
749 if (cc==dim-1)
return dim*(1<<(dim-1));
751 DUNE_THROW(
GridError,
"codim " << cc <<
" (dim=" << dim <<
") not (yet) implemented");
759 dune_static_assert( cc == dim || cc == 0 ,
760 "YaspGrid only supports Entities with codim=dim and codim=0");
764 iTupel coord = _it.coord();
767 for (
int k=0; k<dim; k++)
768 if (i&(1<<k)) (coord[k])++;
776 DUNE_THROW(
GridError,
"codim " << cc <<
" (dim=" << dim <<
") not (yet) implemented");
784 DUNE_THROW(
GridError,
"tried to call father on level 0");
787 YGLI cg = _g.coarser();
790 iTupel coord = _it.coord();
793 for (
int k=0; k<dim; k++) coord[k] = coord[k]/2;
799 bool hasFather ()
const
801 return (_g.level()>0);
817 for (
int k=0; k<dim; k++)
825 const TSI& transformingsubiterator ()
const
842 return (_g.level() == _g.mg()->maxlevel());
847 bool isNew ()
const {
return _yg->adaptRefCount > 0 && _g.mg()->maxlevel() < _g.level() + _yg->adaptRefCount; }
916 const iTupel& size = _g.cell_global().size();
920 for (
int i=0; i<dim; i++)
922 coord[i] = _it.coord(i);
923 if (coord[i]<0) coord[i] += size[i];
924 if (coord[i]>=size[i]) coord[i] -= size[i];
936 for (
int i=dim-1; i>=0; i--)
948 return _it.superindex();
954 return _it.superindex();
966 for (
int k=0; k<dim; k++)
968 coord[k] = _it.coord(k);
969 if (coord[k]<0) coord[k] += _g.cell_global().size(k);
970 if (coord[k]>=_g.cell_global().size(k)) coord[k] -= _g.cell_global().size(k);
976 for (
int k=0; k<dim; k++)
977 if (i&(1<<k)) (coord[k])++;
981 for (
int i=0; i<dim; i++)
985 for (
int j=0; j<_g.level(); j++)
990 trailing =
std::min(trailing,zeros);
994 int level = _g.level()-trailing;
1004 for (
int i=dim-1; i>=0; i--)
1021 for (
int k=0; k<dim; k++)
1022 coord[k] = coord[k]*2 + 1;
1036 for (
int i=dim-1; i>=0; i--)
1046 static unsigned int edge[ 12 ] = { 0, 1, 2, 3, 4, 5, 8, 9, 6, 7, 10, 11 };
1057 int ifix=(dim-1)-(i/m);
1061 for (
int k=0; k<dim; k++)
1063 coord[k] = coord[k]*2+1;
1064 if (k==ifix)
continue;
1065 if ((i%m)&bit) coord[k] += 1;
else coord[k] -= 1;
1077 for (
int i=dim-1; i>=0; i--)
1086 DUNE_THROW(GridError,
"codim " << cc <<
" (dim=" << dim <<
") not (yet) implemented");
1090 int subCompressedIndex (
int i,
int cc)
const
1097 for (
int k=0; k<dim; ++k)
1098 coord[k] = _it.coord(k)-_g.cell_overlap().origin(k);
1103 for (
int k=0; k<dim; k++)
1104 if (i&(1<<k)) (coord[k])++;
1107 int index = coord[dim-1];
1108 for (
int k=dim-2; k>=0; --k)
1109 index = (index*(_g.cell_overlap().size(k)+1))+coord[k];
1121 if (i%2) coord[ivar] += 1;
1124 int index = coord[dim-1];
1125 for (
int k=dim-2; k>=0; --k)
1127 index = (index*(_g.cell_overlap().size(k)+1))+coord[k];
1129 index = (index*(_g.cell_overlap().size(k)))+coord[k];
1132 for (
int j=0; j<ivar; j++)
1134 int n=_g.cell_overlap().size(j)+1;
1135 for (
int l=0; l<dim; l++)
1136 if (l!=j) n *= _g.cell_overlap().size(l);
1144 static unsigned int edge[ 12 ] = { 0, 1, 2, 3, 4, 5, 8, 9, 6, 7, 10, 11 };
1155 int ifix=(dim-1)-(i/m);
1159 for (
int k=0; k<dim; k++)
1161 if (k==ifix)
continue;
1162 if ((i%m)&bit) coord[k] += 1;
1167 int index = coord[dim-1];
1168 for (
int k=dim-2; k>=0; --k)
1170 index = (index*(_g.cell_overlap().size(k)+1))+coord[k];
1172 index = (index*(_g.cell_overlap().size(k)))+coord[k];
1175 for (
int j=dim-1; j>ifix; j--)
1177 int n=_g.cell_overlap().size(j);
1178 for (
int l=0; l<dim; l++)
1179 if (l!=j) n *= _g.cell_overlap().size(l)+1;
1186 DUNE_THROW(GridError,
"codim " << cc <<
" (dim=" << dim <<
") not (yet) implemented");
1190 int subCompressedLeafIndex (
int i,
int cc)
const
1197 for (
int k=0; k<dim; ++k)
1198 coord[k] = _it.coord(k)-_g.cell_overlap().origin(k);
1203 for (
int k=0; k<dim; k++)
1204 if (i&(1<<k)) (coord[k])++;
1207 for (
int k=0; k<dim; k++)
1208 coord[k] = coord[k]<<(_g.mg()->maxlevel()-_g.level());
1211 int index = coord[dim-1];
1212 for (
int k=dim-2; k>=0; --k)
1213 index = (index*(_g.mg()->rbegin().cell_overlap().size(k)+1))+coord[k];
1225 if (i%2) coord[ivar] += 1;
1228 int index = coord[dim-1];
1229 for (
int k=dim-2; k>=0; --k)
1231 index = (index*(_g.cell_overlap().size(k)+1))+coord[k];
1233 index = (index*(_g.cell_overlap().size(k)))+coord[k];
1236 for (
int j=0; j<ivar; j++)
1238 int n=_g.cell_overlap().size(j)+1;
1239 for (
int l=0; l<dim; l++)
1240 if (l!=j) n *= _g.cell_overlap().size(l);
1248 static unsigned int edge[ 12 ] = { 0, 1, 2, 3, 4, 5, 8, 9, 6, 7, 10, 11 };
1259 int ifix=(dim-1)-(i/m);
1263 for (
int k=0; k<dim; k++)
1265 if (k==ifix)
continue;
1266 if ((i%m)&bit) coord[k] += 1;
1271 int index = coord[dim-1];
1272 for (
int k=dim-2; k>=0; --k)
1274 index = (index*(_g.cell_overlap().size(k)+1))+coord[k];
1276 index = (index*(_g.cell_overlap().size(k)))+coord[k];
1279 for (
int j=dim-1; j>ifix; j--)
1281 int n=_g.cell_overlap().size(j);
1282 for (
int l=0; l<dim; l++)
1283 if (l!=j) n *= _g.cell_overlap().size(l)+1;
1290 DUNE_THROW(GridError,
"codim " << cc <<
" (dim=" << dim <<
") not (yet) implemented");
1293 const GridImp * _yg;
1300 template<
int dim,
class Gr
idImp>
1304 enum { dimworld = GridImp::dimensionworld };
1306 typedef typename GridImp::Traits::template Codim<dim>::GeometryImpl GeometryImpl;
1311 typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator
YGLI;
1312 typedef typename SubYGrid<dim,ctype>::TransformingSubIterator
TSI;
1314 typedef typename GridImp::template Codim<dim>::Geometry
Geometry;
1323 typedef typename GridImp::template Codim<dim>::EntitySeed
EntitySeed;
1329 typedef typename YGrid<dim,ctype>::iTupel
iTupel;
1333 : _yg(yg), _it(it), _g(g)
1340 int index ()
const {
return _it.superindex();}
1343 int globalIndex ()
const {
return _g.cell_global().index(_it.coord()); }
1354 GeometryImpl _geometry(_it.position());
1361 if (_g.vertex_interior().inside(_it.coord()))
return InteriorEntity;
1362 if (_g.vertex_interiorborder().inside(_it.coord()))
return BorderEntity;
1363 if (_g.vertex_overlap().inside(_it.coord()))
return OverlapEntity;
1364 if (_g.vertex_overlapfront().inside(_it.coord()))
return FrontEntity;
1378 const iTupel& size = _g.vertex_global().size();
1382 for (
int i=0; i<dim; i++)
1384 coord[i] = _it.coord(i);
1385 if (coord[i]<0) coord[i] += size[i];
1386 if (coord[i]>=size[i]) coord[i] -= size[i];
1390 int trailing = 1000;
1391 for (
int i=0; i<dim; i++)
1395 for (
int j=0; j<_g.level(); j++)
1396 if (coord[i]&(1<<j))
1400 trailing =
std::min(trailing,zeros);
1404 int level = _g.level()-trailing;
1414 for (
int i=dim-1; i>=0; i--)
1429 if (_g.level()==_g.mg()->maxlevel())
1430 return _it.superindex();
1434 for (
int i=0; i<dim; i++) coord[i] = _it.coord(i)-(_g).vertex_overlap().origin(i);
1437 for (
int k=0; k<dim; k++)
1438 coord[k] = coord[k]*(1<<(_g.mg()->maxlevel()-_g.level()));
1441 int index = coord[dim-1];
1442 for (
int k=dim-2; k>=0; --k)
1443 index = (index*(_g.mg()->rbegin().cell_overlap().size(k)+1))+coord[k];
1456 mutable FieldVector<ctype, dim>
loc;
1466 template<
class Gr
idImp>
1469 enum { dim=GridImp::dimension };
1470 enum { dimworld=GridImp::dimensionworld };
1471 typedef typename GridImp::ctype ctype;
1475 typedef typename GridImp::Traits::template Codim< 1 >::GeometryImpl GeometryImpl;
1476 typedef typename GridImp::Traits::template Codim< 1 >::LocalGeometryImpl LocalGeometryImpl;
1480 typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator
YGLI;
1481 typedef typename SubYGrid<dim,ctype>::TransformingSubIterator
TSI;
1482 typedef typename GridImp::template Codim<0>::Entity
Entity;
1484 typedef typename GridImp::template Codim<1>::Geometry
Geometry;
1493 if (_count == 2*_dir + _face || _count >= 2*dim)
1513 _count += (_count < 2*dim);
1519 return _inside.
equals(other._inside) && _count == other._count;
1591 DUNE_THROW(
GridError,
"called boundarySegmentIndex while boundary() == false");
1594 const FieldVector<int, dim> & size = _inside.
gridlevel().mg()->begin().cell_overlap().size();
1595 const FieldVector<int, dim> & origin = _inside.
gridlevel().mg()->begin().cell_overlap().origin();
1596 FieldVector<int, dim> sides;
1598 for (
int i=0; i<dim; i++)
1601 ((_inside.
gridlevel().mg()->begin().cell_overlap().origin(i)
1602 == _inside.
gridlevel().mg()->begin().cell_global().origin(i))+
1603 (_inside.
gridlevel().mg()->begin().cell_overlap().origin(i) +
1604 _inside.
gridlevel().mg()->begin().cell_overlap().size(i)
1605 == _inside.
gridlevel().mg()->begin().cell_global().origin(i) +
1606 _inside.
gridlevel().mg()->begin().cell_global().size(i)));
1611 pos /= (1<<_inside.
level());
1614 FieldVector<int, dim> fsize;
1617 for (
int k=0; k<dim; k++)
1619 for (
int k=0; k<dim; k++)
1620 fsize[k] = vol/size[k];
1625 int localoffset = 1;
1626 for (
int k=dim-1; k>=0; k--)
1628 if (k == _dir)
continue;
1629 index += (pos[k]) * localoffset;
1630 localoffset *= size[k];
1635 for (
int k=0; k<_dir; k++)
1636 index += sides[k] * fsize[k];
1638 index += _face * (sides[_dir]>1) * fsize[_dir];
1652 FieldVector<ctype, dimworld>
outerNormal (
const FieldVector<ctype, dim-1>& local)
const
1654 return _faceInfo[_count].normal;
1658 FieldVector<ctype, dimworld>
unitOuterNormal (
const FieldVector<ctype, dim-1>& local)
const
1660 return _faceInfo[_count].normal;
1666 return _faceInfo[_count].normal;
1674 FieldVector<ctype, dimworld> n = _faceInfo[_count].normal;
1727 _inside(myself.yaspgrid(), myself.gridlevel(),
1728 myself.transformingsubiterator()),
1729 _outside(myself.yaspgrid(), myself.gridlevel(),
1730 myself.transformingsubiterator()),
1735 _pos_world(myself.transformingsubiterator().position())
1754 _inside(it._inside),
1755 _outside(it._outside),
1759 _pos_world(it._pos_world)
1765 _inside = it._inside;
1766 _outside = it._outside;
1770 _pos_world = it._pos_world;
1782 mutable FieldVector<ctype, dimworld> _pos_world;
1787 FieldVector<ctype, dimworld> normal;
1788 LocalGeometryImpl geom_inside;
1789 LocalGeometryImpl geom_outside;
1793 static const array<faceInfo, 2*GridImp::dimension> _faceInfo;
1795 static array<faceInfo, 2*dim> initFaceInfo()
1797 const FieldVector<typename GridImp::ctype, GridImp::dimension> ext_local(1.0);
1798 array<faceInfo, 2*dim> I;
1802 FieldVector<ctype, dim> a(0.5); a[i] = 0.0;
1803 FieldVector<ctype, dim> b(0.5); b[i] = 1.0;
1805 I[2*i].normal = 0.0;
1806 I[2*i+1].normal = 0.0;
1807 I[2*i].normal[i] = -1.0;
1808 I[2*i+1].normal[i] = +1.0;
1810 I[2*i].geom_inside =
1811 LocalGeometryImpl(a, ext_local, i);
1812 I[2*i].geom_outside =
1813 LocalGeometryImpl(b, ext_local, i);
1814 I[2*i+1].geom_inside =
1815 LocalGeometryImpl(b, ext_local, i);
1816 I[2*i+1].geom_outside =
1817 LocalGeometryImpl(a, ext_local, i);
1823 template<
class Gr
idImp>
1824 const array<typename YaspIntersection<GridImp>::faceInfo, 2*GridImp::dimension>
1825 YaspIntersection<GridImp>::_faceInfo =
1826 YaspIntersection<GridImp>::initFaceInfo();
1835 template<
class Gr
idImp>
1838 enum { dim=GridImp::dimension };
1839 enum { dimworld=GridImp::dimensionworld };
1840 typedef typename GridImp::ctype ctype;
1850 GridImp::getRealImplementation(*this).increment();
1856 return GridImp::getRealImplementation(*this).equals(
1857 GridImp::getRealImplementation(other));
1879 GridImp::getRealImplementation(*this).assign(
1880 GridImp::getRealImplementation(it));
1891 template<
class Gr
idImp>
1895 enum { dim=GridImp::dimension };
1896 enum { dimworld=GridImp::dimensionworld };
1897 typedef typename GridImp::ctype ctype;
1900 typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator
YGLI;
1901 typedef typename SubYGrid<dim,ctype>::TransformingSubIterator
TSI;
1902 typedef typename GridImp::template Codim<0>::Entity
Entity;
1906 typedef typename YGrid<dim,ctype>::iTupel
iTupel;
1913 StackElem se(this->
_g);
1914 se.coord = this->
_it.coord();
1918 _maxlevel =
std::min(maxlevel,this->
_g.mg()->maxlevel());
1921 if (this->
_g.level()<_maxlevel)
1934 _maxlevel(it._maxlevel), stack(it.stack)
1941 if (stack.empty())
return;
1944 if (this->
_g.level()<_maxlevel)
1953 s <<
"HIER: " <<
"level=" << this->
_g.level()
1954 <<
" position=" << this->
_it.coord()
1955 <<
" superindex=" << this->
_it.superindex()
1956 <<
" maxlevel=" << this->_maxlevel
1957 <<
" stacksize=" << stack.size()
1967 StackElem(
YGLI gg) : g(gg) {}
1969 std::stack<StackElem> stack;
1975 StackElem se(this->
_g.finer());
1976 for (
int i=0; i<(1<<dim); i++)
1978 for (
int k=0; k<dim; k++)
1980 se.coord[k] = this->_it.coord(k)*2+1;
1982 se.coord[k] = this->
_it.coord(k)*2;
1990 StackElem se = stack.top();
1993 this->
_it.reinit(this->
_g.cell_overlap(),se.coord);
2002 template<
int codim,
class Gr
idImp>
2006 enum { dim=GridImp::dimension };
2014 :
_l(level),
_c(coord)
2023 const FieldVector<int, dim> &
coord()
const {
return _c; }
2027 FieldVector<int, dim>
_c;
2041 template<
int codim,
class Gr
idImp>
2045 enum { dim=GridImp::dimension };
2047 enum { dimworld=GridImp::dimensionworld };
2048 typedef typename GridImp::ctype ctype;
2050 typedef typename GridImp::template Codim<codim>::Entity
Entity;
2051 typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator
YGLI;
2052 typedef typename SubYGrid<dim,ctype>::TransformingSubIterator
TSI;
2066 if (codim>0 && codim<dim)
2068 DUNE_THROW(
GridError,
"YaspEntityPointer: codim not implemented");
2078 if (codim>0 && codim<dim)
2080 DUNE_THROW(
GridError,
"YaspEntityPointer: codim not implemented");
2088 if (codim>0 && codim<dim)
2090 DUNE_THROW(
GridError,
"YaspEntityPointer: codim not implemented");
2157 template<
int codim, PartitionIteratorType pitype,
class Gr
idImp>
2162 enum { dim=GridImp::dimension };
2164 enum { dimworld=GridImp::dimensionworld };
2165 typedef typename GridImp::ctype ctype;
2167 typedef typename GridImp::template Codim<codim>::Entity
Entity;
2168 typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator
YGLI;
2169 typedef typename SubYGrid<dim,ctype>::TransformingSubIterator
TSI;
2195 template<
class Gr
idImp>
2197 :
public IndexSet< GridImp, YaspLevelIndexSet< GridImp >, unsigned int >
2219 IndexType index (
const typename GridImp::Traits::template Codim<cc>::Entity& e)
const
2221 assert( cc == 0 || cc == GridImp::dimension );
2222 return grid.getRealImplementation(e).compressedIndex();
2227 IndexType subIndex (
const typename remove_const< GridImp >::type::Traits::template Codim< cc >::Entity &e,
2228 int i,
unsigned int codim )
const
2230 assert( cc == 0 || cc == GridImp::dimension );
2231 if( cc == GridImp::dimension )
2232 return grid.getRealImplementation(e).compressedIndex();
2234 return grid.getRealImplementation(e).subCompressedIndex(i,codim);
2240 return grid.size( level, type );
2246 return grid.size( level, codim );
2250 template<
class EntityType>
2253 return e.level() == level;
2257 const std::vector<GeometryType>&
geomTypes (
int codim)
const
2259 return mytypes[codim];
2263 const GridImp& grid;
2265 std::vector<GeometryType> mytypes[GridImp::dimension+1];
2271 template<
class Gr
idImp>
2273 :
public IndexSet< GridImp, YaspLeafIndexSet< GridImp >, unsigned int >
2294 IndexType index (
const typename GridImp::Traits::template Codim<cc>::Entity& e)
const
2296 assert( cc == 0 || cc == GridImp::dimension );
2297 return grid.getRealImplementation(e).compressedIndex();
2302 IndexType subIndex (
const typename remove_const< GridImp >::type::Traits::template Codim< cc >::Entity &e,
2303 int i,
unsigned int codim )
const
2305 assert( cc == 0 || cc == GridImp::dimension );
2306 if( cc == GridImp::dimension )
2307 return grid.getRealImplementation(e).compressedIndex();
2309 return grid.getRealImplementation(e).subCompressedIndex(i,codim);
2315 return grid.size( grid.maxLevel(), type );
2321 return grid.size( grid.maxLevel(), codim );
2325 template<
class EntityType>
2329 return (e.level() == grid.maxLevel());
2333 const std::vector<GeometryType>&
geomTypes (
int codim)
const
2335 return mytypes[codim];
2339 const GridImp& grid;
2340 enum { ncodim = remove_const<GridImp>::type::dimension+1 };
2341 std::vector<GeometryType> mytypes[ncodim];
2354 template<
class Gr
idImp>
2356 typename remove_const<GridImp>::type::PersistentIndexType >
2366 typedef typename remove_const<GridImp>::type::PersistentIndexType
IdType;
2381 IdType id (
const typename remove_const<GridImp>::type::Traits::template Codim<cd>::Entity& e)
const
2383 return grid.getRealImplementation(e).persistentIndex();
2391 IdType subId (
const typename remove_const<GridImp>::type::Traits::template Codim< 0 >::Entity &e,
2392 int i,
unsigned int codim )
const
2394 return grid.getRealImplementation(e).subPersistentIndex(i,codim);
2398 const GridImp& grid;
2402 template<
int dim,
int dimworld>
2406 typedef CollectiveCommunication<MPI_Comm>
CCType;
2408 typedef CollectiveCommunication<Dune::YaspGrid<dim> >
CCType;
2423 bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+yaspgrid_codim_bits>,
2425 bigunsignedint<dim*yaspgrid_dim_bits+yaspgrid_level_bits+yaspgrid_codim_bits>,
2432 template<
int dim,
int codim>
2434 template<
class G,
class DataHandle>
2437 if (data.contains(dim,codim))
2439 DUNE_THROW(
GridError,
"interface communication not implemented");
2447 template<
class G,
class DataHandle>
2450 if (data.contains(dim,dim))
2451 g.template communicateCodim<DataHandle,dim>(data,iftype,dir,level);
2458 template<
class G,
class DataHandle>
2461 if (data.contains(dim,0))
2462 g.template communicateCodim<DataHandle,0>(data,iftype,dir,level);
2486 public MultiYGrid<dim,yaspgrid_ctype>
2496 boundarysegmentssize();
2499 void boundarysegmentssize()
2502 const FieldVector<int, dim> & size = MultiYGrid<dim,ctype>::begin().cell_overlap().size();
2503 FieldVector<int, dim> sides;
2505 for (
int i=0; i<dim; i++)
2508 ((MultiYGrid<dim,ctype>::begin().cell_overlap().origin(i)
2509 == MultiYGrid<dim,ctype>::begin().cell_global().origin(i))+
2510 (MultiYGrid<dim,ctype>::begin().cell_overlap().origin(i) +
2511 MultiYGrid<dim,ctype>::begin().cell_overlap().size(i)
2512 == MultiYGrid<dim,ctype>::begin().cell_global().origin(i) +
2513 MultiYGrid<dim,ctype>::begin().cell_global().size(i)));
2517 for (
int k=0; k<dim; k++)
2520 for (
int l=0; l<dim; l++)
2525 nBSegments += sides[k]*offset;
2531 using MultiYGrid<dim,yaspgrid_ctype>::defaultLoadbalancer;
2553 typedef MultiYGrid<dim,ctype>
YMG;
2554 typedef typename MultiYGrid<dim,ctype>::YGridLevelIterator
YGLI;
2555 typedef typename SubYGrid<dim,ctype>::TransformingSubIterator
TSI;
2556 typedef typename MultiYGrid<dim,ctype>::Intersection
IS;
2557 typedef typename std::deque<IS>::const_iterator
ISIT;
2568 Dune::FieldVector<ctype, dim> L,
2569 Dune::FieldVector<int, dim> s,
2570 Dune::FieldVector<bool, dim> periodic,
int overlap,
2571 const YLoadBalance<dim>* lb = defaultLoadbalancer())
2573 :
YMG(comm,L,s,periodic,overlap,lb), ccobj(comm),
2574 keep_ovlp(
true), adaptRefCount(0), adaptActive(
false)
2576 :
YMG(L,s,periodic,overlap,lb),
2577 keep_ovlp(
true), adaptRefCount(0), adaptActive(
false)
2597 Dune::FieldVector<int, dim> s,
2598 Dune::FieldVector<bool, dim> periodic,
int overlap,
2599 const YLoadBalance<dim>* lb = YMG::defaultLoadbalancer())
2601 :
YMG(MPI_COMM_SELF,L,s,periodic,overlap,lb), ccobj(MPI_COMM_SELF),
2602 keep_ovlp(
true), adaptRefCount(0), adaptActive(
false)
2604 :
YMG(L,s,periodic,overlap,lb),
2605 keep_ovlp(
true), adaptRefCount(0), adaptActive(
false)
2613 deallocatePointers(indexsets);
2614 deallocatePointers(theleafindexset);
2615 deallocatePointers(theglobalidset);
2627 int maxLevel()
const {
return MultiYGrid<dim,ctype>::maxlevel();}
2634 "Coarsening " << -refCount <<
" levels requested!");
2635 for (
int k=refCount; k<0; k++)
2637 MultiYGrid<dim,ctype>::coarsen();
2639 indexsets.pop_back();
2641 for (
int k=0; k<refCount; k++)
2643 MultiYGrid<dim,ctype>::refine(keep_ovlp);
2655 keep_ovlp = keepPhysicalOverlap;
2671 assert(adaptActive ==
false);
2672 if (e.level() !=
maxLevel())
return false;
2673 adaptRefCount =
std::max(adaptRefCount, refCount);
2685 return ( e.level() ==
maxLevel() ) ? adaptRefCount : 0;
2692 return (adaptRefCount > 0);
2699 adaptRefCount =
comm().max(adaptRefCount);
2700 return (adaptRefCount < 0);
2706 adaptActive =
false;
2711 template<
int cd, PartitionIteratorType pitype>
2714 return levelbegin<cd,pitype>(level);
2718 template<
int cd, PartitionIteratorType pitype>
2721 return levelend<cd,pitype>(level);
2728 return levelbegin<cd,All_Partition>(level);
2735 return levelend<cd,All_Partition>(level);
2739 template<
int cd, PartitionIteratorType pitype>
2742 return levelbegin<cd,pitype>(
maxLevel());
2746 template<
int cd, PartitionIteratorType pitype>
2749 return levelend<cd,pitype>(
maxLevel());
2756 return levelbegin<cd,All_Partition>(
maxLevel());
2763 return levelend<cd,All_Partition>(
maxLevel());
2767 template <
typename Seed>
2768 typename Traits::template Codim<Seed::codimension>::EntityPointer
2771 const int codim = Seed::codimension;
2772 YGLI g = MultiYGrid<dim,ctype>::begin(seed.level());
2777 TSI(g.cell_overlap(), seed.coord()));
2780 TSI(g.vertex_overlap(), seed.coord()));
2782 DUNE_THROW(
GridError,
"YaspEntityPointer: codim not implemented");
2789 YGLI g = MultiYGrid<dim,ctype>::begin(level);
2813 int size (
int level,
int codim)
const
2815 return sizes[level][codim];
2827 return (type.isCube()) ? sizes[level][dim-type.dim()] : 0;
2846 template<
class DataHandleImp,
class DataType>
2856 template<
class DataHandleImp,
class DataType>
2866 template<
class DataHandle,
int codim>
2870 if (!data.contains(dim,codim))
return;
2873 typedef typename DataHandle::DataType DataType;
2876 YGLI g = MultiYGrid<dim,ctype>::begin(level);
2879 const std::deque<IS>* sendlist=0;
2880 const std::deque<IS>* recvlist=0;
2887 sendlist = &g.send_cell_interior_overlap();
2888 recvlist = &g.recv_cell_overlap_interior();
2892 sendlist = &g.send_cell_overlap_overlap();
2893 recvlist = &g.recv_cell_overlap_overlap();
2900 sendlist = &g.send_vertex_interiorborder_interiorborder();
2901 recvlist = &g.recv_vertex_interiorborder_interiorborder();
2906 sendlist = &g.send_vertex_interiorborder_overlapfront();
2907 recvlist = &g.recv_vertex_overlapfront_interiorborder();
2911 sendlist = &g.send_vertex_overlap_overlapfront();
2912 recvlist = &g.recv_vertex_overlapfront_overlap();
2916 sendlist = &g.send_vertex_overlapfront_overlapfront();
2917 recvlist = &g.recv_vertex_overlapfront_overlapfront();
2923 std::swap(sendlist,recvlist);
2928 std::vector<int> send_size(sendlist->size(),-1);
2929 std::vector<int> recv_size(recvlist->size(),-1);
2930 std::vector<size_t*> send_sizes(sendlist->size(),
static_cast<size_t*
>(0));
2931 std::vector<size_t*> recv_sizes(recvlist->size(),
static_cast<size_t*
>(0));
2932 if (data.fixedsize(dim,codim))
2936 for (
ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
2940 send_size[cnt] = is->grid.totalsize() * data.size(*it);
2944 for (
ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
2948 recv_size[cnt] = is->grid.totalsize() * data.size(*it);
2956 for (
ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
2959 size_t *buf =
new size_t[is->grid.totalsize()];
2960 send_sizes[cnt] = buf;
2963 int i=0;
size_t n=0;
2968 for ( ; it!=tsubend; ++it)
2970 buf[i] = data.size(*it);
2979 MultiYGrid<dim,ctype>::torus().send(is->rank,buf,is->grid.totalsize()*
sizeof(size_t));
2985 for (
ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
2988 size_t *buf =
new size_t[is->grid.totalsize()];
2989 recv_sizes[cnt] = buf;
2992 MultiYGrid<dim,ctype>::torus().recv(is->rank,buf,is->grid.totalsize()*
sizeof(size_t));
2997 MultiYGrid<dim,ctype>::torus().exchange();
3001 for (
ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
3003 delete[] send_sizes[cnt];
3004 send_sizes[cnt] = 0;
3010 for (
ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
3013 size_t *buf = recv_sizes[cnt];
3017 for (
int i=0; i<is->grid.totalsize(); ++i)
3028 std::vector<DataType*> sends(sendlist->size(),
static_cast<DataType*
>(0));
3030 for (
ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
3038 DataType *buf =
new DataType[send_size[cnt]];
3044 MessageBuffer<DataType> mb(buf);
3051 for ( ; it!=tsubend; ++it)
3052 data.gather(mb,*it);
3055 MultiYGrid<dim,ctype>::torus().send(is->rank,buf,send_size[cnt]*
sizeof(DataType));
3060 std::vector<DataType*> recvs(recvlist->size(),
static_cast<DataType*
>(0));
3062 for (
ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
3070 DataType *buf =
new DataType[recv_size[cnt]];
3076 MultiYGrid<dim,ctype>::torus().recv(is->rank,buf,recv_size[cnt]*
sizeof(DataType));
3081 MultiYGrid<dim,ctype>::torus().exchange();
3085 for (
ISIT is=sendlist->begin(); is!=sendlist->end(); ++is)
3087 delete[] sends[cnt];
3094 for (
ISIT is=recvlist->begin(); is!=recvlist->end(); ++is)
3097 DataType *buf = recvs[cnt];
3100 MessageBuffer<DataType> mb(buf);
3103 if (data.fixedsize(dim,codim))
3107 size_t n=data.size(*it);
3110 for ( ; it!=tsubend; ++it)
3111 data.scatter(mb,*it,n);
3116 size_t *sbuf = recv_sizes[cnt];
3121 for ( ; it!=tsubend; ++it)
3122 data.scatter(mb,*it,sbuf[i++]);
3135 return *(theglobalidset[0]);
3140 return *(theglobalidset[0]);
3145 if (level<0 || level>
maxLevel()) DUNE_THROW(RangeError,
"level out of range");
3146 return *(indexsets[level]);
3151 return *(theleafindexset[0]);
3157 const CollectiveCommunication<MPI_Comm>&
comm ()
const
3164 const CollectiveCommunication<YaspGrid>&
comm ()
const
3170 YaspIntersectionIterator<const YaspGrid<dim> >&
3185 CollectiveCommunication<MPI_Comm> ccobj;
3187 CollectiveCommunication<YaspGrid> ccobj;
3190 std::vector<YaspLevelIndexSet<const YaspGrid<dim> >*> indexsets;
3191 std::vector<YaspLeafIndexSet<const YaspGrid<dim> >*> theleafindexset;
3192 std::vector<YaspGlobalIdSet<const YaspGrid<dim> >*> theglobalidset;
3205 void deallocatePointers(T& container)
3207 typedef typename T::iterator Iterator;
3209 for(Iterator entry=container.begin(); entry != container.end(); ++entry)
3213 template<
int codim_,
int dim_,
class Gr
idImp_,
template<
int,
int,
class>
class EntityImp_>
3217 class MessageBuffer {
3220 MessageBuffer (DT *p)
3229 void write (
const Y& data)
3231 dune_static_assert(( is_same<DT,Y>::value ),
"DataType missmatch");
3237 void read (Y& data)
const
3239 dune_static_assert(( is_same<DT,Y>::value ),
"DataType missmatch");
3251 for (
YGLI g=MultiYGrid<dim,ctype>::begin(); g!=MultiYGrid<dim,ctype>::end(); ++g)
3254 sizes[g.level()][0] = 1;
3255 for (
int i=0; i<dim; ++i)
3256 sizes[g.level()][0] *= g.cell_overlap().size(i);
3261 sizes[g.level()][1] = 0;
3262 for (
int i=0; i<dim; ++i)
3264 int s=g.cell_overlap().size(i)+1;
3265 for (
int j=0; j<dim; ++j)
3267 s *= g.cell_overlap().size(j);
3268 sizes[g.level()][1] += s;
3275 sizes[g.level()][dim-1] = 0;
3276 for (
int i=0; i<dim; ++i)
3278 int s=g.cell_overlap().size(i);
3279 for (
int j=0; j<dim; ++j)
3281 s *= g.cell_overlap().size(j)+1;
3282 sizes[g.level()][dim-1] += s;
3287 sizes[g.level()][dim] = 1;
3288 for (
int i=0; i<dim; ++i)
3289 sizes[g.level()][dim] *= g.vertex_overlapfront().size(i);
3294 template<
int cd, PartitionIteratorType pitype>
3295 YaspLevelIterator<cd,pitype,GridImp> levelbegin (
int level)
const
3297 dune_static_assert( cd == dim || cd == 0 ,
3298 "YaspGrid only supports Entities with codim=dim and codim=0");
3299 YGLI g = MultiYGrid<dim,ctype>::begin(level);
3300 if (level<0 || level>
maxLevel()) DUNE_THROW(RangeError,
"level out of range");
3302 return levelend <cd, pitype> (level);
3306 return YaspLevelIterator<cd,pitype,GridImp>(
this,g,g.cell_interior().tsubbegin());
3308 return YaspLevelIterator<cd,pitype,GridImp>(
this,g,g.cell_overlap().tsubbegin());
3313 return YaspLevelIterator<cd,pitype,GridImp>(
this,g,g.vertex_interior().tsubbegin());
3315 return YaspLevelIterator<cd,pitype,GridImp>(
this,g,g.vertex_interiorborder().tsubbegin());
3317 return YaspLevelIterator<cd,pitype,GridImp>(
this,g,g.vertex_overlap().tsubbegin());
3319 return YaspLevelIterator<cd,pitype,GridImp>(
this,g,g.vertex_overlapfront().tsubbegin());
3321 DUNE_THROW(GridError,
"YaspLevelIterator with this codim or partition type not implemented");
3325 template<
int cd, PartitionIteratorType pitype>
3326 YaspLevelIterator<cd,pitype,GridImp> levelend (
int level)
const
3328 dune_static_assert( cd == dim || cd == 0 ,
3329 "YaspGrid only supports Entities with codim=dim and codim=0");
3330 YGLI g = MultiYGrid<dim,ctype>::begin(level);
3331 if (level<0 || level>
maxLevel()) DUNE_THROW(RangeError,
"level out of range");
3335 return YaspLevelIterator<cd,pitype,GridImp>(
this,g,g.cell_interior().tsubend());
3337 return YaspLevelIterator<cd,pitype,GridImp>(
this,g,g.cell_overlap().tsubend());
3342 return YaspLevelIterator<cd,pitype,GridImp>(
this,g,g.vertex_interior().tsubend());
3344 return YaspLevelIterator<cd,pitype,GridImp>(
this,g,g.vertex_interiorborder().tsubend());
3346 return YaspLevelIterator<cd,pitype,GridImp>(
this,g,g.vertex_overlap().tsubend());
3348 return YaspLevelIterator<cd,pitype,GridImp>(
this,g,g.vertex_overlapfront().tsubend());
3350 DUNE_THROW(GridError,
"YaspLevelIterator with this codim or partition type not implemented");
3353 int sizes[
MAXL][dim+1];
3359 namespace Capabilities
3376 static const bool v =
true;
3377 static const unsigned int topologyId = GenericGeometry :: CubeTopology< dim > :: type :: id ;
3386 static const bool v =
true;
3395 static const bool v =
true;
3404 static const bool v =
true;
3410 static const bool v =
true;
3416 static const bool v =
true;
3425 static const bool v =
true;
3434 static const bool v =
true;
3443 static const bool v =
true;