8 #include <dune/common/fvector.hh>
9 #include <dune/common/fmatrix.hh>
10 #include <dune/common/bigunsignedint.hh>
11 #include <dune/common/collectivecommunication.hh>
12 #include <dune/common/reservedvector.hh>
13 #include <dune/geometry/genericgeometry/topologytypes.hh>
16 #include <dune/grid/sgrid/numbering.hh>
40 template<
int dim,
int dimworld,
class Gr
idImp>
class SGeometry;
41 template<
int codim,
int dim,
class Gr
idImp>
class SEntity;
43 template<
int codim,
class Gr
idImp>
class SEntitySeed;
44 template<
int codim, PartitionIteratorType,
class Gr
idImp>
class SLevelIterator;
45 template<
int dim,
int dimworld,
class ctype>
class SGrid;
50 namespace FacadeOptions
53 template<
int dim,
int dimworld,
class ctype,
int mydim,
int cdim>
57 static const bool v =
false;
60 template<
int dim,
int dimworld,
class ctype,
int mydim,
int cdim>
64 static const bool v =
false;
100 template<
int mydim,
int cdim,
class Gr
idImp>
106 typedef typename GridImp::ctype
ctype;
125 FieldVector< ctype, cdim >
corner (
const int i )
const
131 FieldVector<ctype, cdim >
center ( )
const
137 FieldVector<ctype, cdim>
global (
const FieldVector<ctype, mydim>&
local)
const;
140 FieldVector<ctype, mydim>
local (
const FieldVector<ctype, cdim>&
global)
const;
173 void print (std::ostream& ss,
int indent)
const;
179 void make (FieldMatrix<ctype,mydim+1,cdim>& __As);
185 FieldVector<ctype, cdim> s;
186 FieldVector<ctype, cdim> centroid;
187 FieldMatrix<ctype,mydim,cdim> A;
188 array<FieldVector<ctype, cdim>, 1<<mydim> c;
189 mutable FieldMatrix<ctype,cdim,mydim> Jinv;
190 mutable bool builtinverse;
194 template<
int cdim,
class Gr
idImp>
200 typedef typename GridImp::ctype
ctype;
219 FieldVector<ctype, cdim >
corner (
const int i )
const
225 FieldVector<ctype, cdim >
center ( )
const
231 void print (std::ostream& ss,
int indent)
const;
234 void make (FieldMatrix<ctype,1,cdim>& __As);
237 FieldVector<ctype, cdim>
global (
const FieldVector<ctype, 0>&
local)
const {
return corner(0); }
240 FieldVector<ctype, 0>
local (
const FieldVector<ctype, cdim>&
global)
const {
return FieldVector<ctype,0> (0.0); }
287 static const FieldMatrix<ctype, 0, cdim > dummy (
ctype( 0 ) );
293 static const FieldMatrix<ctype,cdim,0> dummy(
ctype(0) );
298 FieldVector<ctype, cdim>
s;
301 template <
int mydim,
int cdim,
class Gr
idImp>
302 inline std::ostream& operator<< (std::ostream& s, SGeometry<mydim,cdim,GridImp>& e)
313 template<
int codim,
int dim,
class Gr
idImp,
template<
int,
int,
class>
class EntityImp>
319 enum { dimworld = GridImp::dimensionworld };
321 typedef typename GridImp::Traits::template Codim< codim >::GeometryImpl GeometryImpl;
324 typedef typename GridImp::ctype
ctype;
325 typedef typename GridImp::template Codim<codim>::Geometry
Geometry;
386 void make (GridImp* _grid,
int _l,
int _id);
389 void make (
int _l,
int _id);
397 return grid->persistentIndex(
l, codim,
z);
411 if (codim<dim || l==grid->maxLevel())
416 array<int,dim> coord;
417 for (
int k=0; k<dim; k++)
418 coord[k] =
z[k]*(1<<(
grid->maxLevel()-
l));
421 return grid->n(
grid->maxLevel(),coord);
440 template<
int codim,
int dim,
class Gr
idImp>
478 template<
int dim,
class Gr
idImp>
481 enum { dimworld = GridImp::dimensionworld };
488 typedef typename GridImp::Traits::template Codim< 0 >::GeometryImpl GeometryImpl;
489 typedef typename GridImp::Traits::template Codim< 0 >::LocalGeometryImpl LocalGeometryImpl;
495 typedef typename GridImp::ctype
ctype;
496 typedef typename GridImp::template Codim<0>::Geometry
Geometry;
515 template<
int cc>
int count ()
const;
524 int subCompressedIndex (
int codim,
int i)
const
528 return (this->
grid)->n(this->
l, this->
grid->subz(this->z,i,codim));
534 int subCompressedLeafIndex (
int codim,
int i)
const
538 assert(this->
l == this->
grid->maxLevel());
540 return (this->
grid)->n(this->
l, this->
grid->subz(this->z,i,codim));
548 return this->
grid->persistentIndex(this->
l, codim, this->
grid->subz(this->z,i,codim));
574 bool hasFather ()
const
576 return (this->
level()>0);
582 return ( this->
grid->maxLevel() == this->
level() );
596 LocalGeometry geometryInFather ()
const;
611 SEntity (GridImp* _grid,
int _l,
int _index) :
622 void make (GridImp* _grid,
int _l,
int _id)
625 built_father =
false;
632 built_father =
false;
639 mutable bool built_father;
640 mutable int father_index;
641 mutable LocalGeometryImpl in_father_local;
642 void make_father()
const;
663 template<
class Gr
idImp>
668 enum { dim = GridImp::dimension };
669 enum { dimworld = GridImp::dimensionworld };
676 typedef typename GridImp::template Codim<0>::Entity
Entity;
677 typedef typename GridImp::ctype
ctype;
690 int _maxLevel,
bool makeend) :
698 orig_l = this->
entity().level();
699 orig_index = _grid->getRealImplementation(this->
entity()).compressedIndex();
703 stack.push(originalElement);
709 push_sons(orig_l,orig_index);
717 int orig_l, orig_index;
720 std::stack<SHierarchicStackElem, Dune::ReservedVector<SHierarchicStackElem,GridImp::MAXL> > stack;
722 void push_sons (
int level,
int fatherid);
732 template<
class Gr
idImp>
735 enum { dim=GridImp::dimension };
736 enum { dimworld=GridImp::dimensionworld };
738 typedef typename GridImp::Traits::template Codim< 1 >::GeometryImpl GeometryImpl;
739 typedef typename GridImp::Traits::template Codim< 1 >::LocalGeometryImpl LocalGeometryImpl;
742 typedef typename GridImp::template Codim<0>::Entity
Entity;
744 typedef typename GridImp::template Codim<1>::Geometry
Geometry;
753 typedef typename GridImp::ctype
ctype;
787 return grid->boundarySegmentIndex(
self.level(), count, zred);
795 FieldVector<ctype, GridImp::dimensionworld>
outerNormal (
const FieldVector<ctype, GridImp::dimension-1>& local)
const
800 FieldVector<ctype, GridImp::dimensionworld>
unitOuterNormal (
const FieldVector<ctype, GridImp::dimension-1>& local)
const
808 FieldVector<ctype, dimworld> normal(0.0);
810 normal[count/2] = 1.0;
812 normal[count/2] = -1.0;
817 FieldVector<ctype, GridImp::dimensionworld>
integrationOuterNormal (
const FieldVector<ctype, GridImp::dimension-1>& local)
const
820 n *=
geometry().integrationElement(local);
851 self(*_self), ne(self), grid(_grid),
852 partition(_grid->partition(grid->getRealImplementation(ne).l,_self->z)),
853 zred(_grid->compress(grid->getRealImplementation(ne).l,_self->z)),
861 self(other.self), ne(other.ne), grid(other.grid),
862 partition(other.partition), zred(other.zred),
863 count(other.count), valid_count(other.valid_count),
864 valid_nb(other.valid_nb), is_on_boundary(other.is_on_boundary),
865 built_intersections(false),
874 assert(grid == other.grid);
879 partition = other.partition;
882 valid_count = other.valid_count;
883 valid_nb = other.valid_nb;
884 is_on_boundary = other.is_on_boundary;
887 built_intersections =
false;
893 void make (
int _count)
const;
894 void makeintersections ()
const;
897 const GridImp * grid;
901 mutable bool valid_count;
902 mutable bool valid_nb;
903 mutable bool is_on_boundary;
904 mutable bool built_intersections;
905 mutable LocalGeometryImpl is_self_local;
906 mutable GeometryImpl is_global;
907 mutable LocalGeometryImpl is_nb_local;
911 template<
class Gr
idImp>
914 enum { dim=GridImp::dimension };
915 enum { dimworld=GridImp::dimensionworld };
917 typedef typename GridImp::template Codim<0>::Entity
Entity;
919 typedef typename GridImp::template Codim<1>::Geometry
Geometry;
929 typedef typename GridImp::ctype
ctype;
933 return is.boundary();
939 return is.boundaryId();
945 return is.boundarySegmentIndex();
951 return is.neighbor();
969 return is.conforming();
975 return is.geometryInInside();
981 return is.geometryInOutside();
987 return is.geometry();
999 return is.indexInInside();
1005 return is.indexInOutside();
1011 return is.outerNormal(local);
1017 return is.integrationOuterNormal(local);
1023 return is.unitOuterNormal(local);
1029 return is.centerUnitOuterNormal();
1036 #ifndef DOXYGEN // doxygen can't handle this recursive usage
1052 while(! this->empty())
1063 template<
int codim,
class Gr
idImp>
1066 enum { dim = GridImp::dimension };
1070 typedef typename GridImp::template Codim<codim>::Entity
Entity;
1130 return grid->getRealImplementation(
entity());
1162 grid->getRealImplementation(*e).make(_grid, _l,_id);
1175 template<
int codim,
class Gr
idImp>
1178 enum { dim = GridImp::dimension };
1184 _l(l), _index(index)
1188 int index ()
const {
return this->_index; }
1200 template<
int codim, PartitionIteratorType pitype,
class Gr
idImp>
1205 enum { dim = GridImp::dimension };
1211 typedef typename GridImp::template Codim<codim>::Entity
Entity;
1229 template<
class Gr
idImp>
1235 enum { dim = GridImp::dimension };
1252 int index (
const typename GridImp::Traits::template Codim<cd>::Entity& e)
const
1254 return grid.getRealImplementation(e).compressedIndex();
1258 int subIndex (
const typename GridImp::Traits::template Codim< cc >::Entity &e,
1259 int i,
unsigned int codim )
const
1262 return grid.getRealImplementation(e).subCompressedIndex(codim, i);
1264 DUNE_THROW( NotImplemented,
"subIndex for higher codimension entity not implemented for SGrid." );
1268 template<
class EntityType >
1271 return (e.level() == level);
1277 return grid.size( level, type );
1283 return grid.size( level, codim );
1287 const std::vector<GeometryType>&
geomTypes (
int codim)
const
1289 return mytypes[codim];
1293 const GridImp& grid;
1295 std::vector<GeometryType> mytypes[GridImp::dimension+1];
1300 template<
class Gr
idImp>
1306 enum { dim = GridImp::dimension };
1316 for (
int codim=0; codim<=dim; codim++)
1326 int index (
const typename remove_const<GridImp>::type::Traits::template Codim<cd>::Entity& e)
const
1328 return grid.getRealImplementation(e).compressedLeafIndex();
1332 int subIndex (
const typename GridImp::Traits::template Codim< cc >::Entity &e,
1333 int i,
unsigned int codim )
const
1336 return grid.getRealImplementation(e).subCompressedIndex(codim, i);
1338 DUNE_THROW( NotImplemented,
"subIndex for higher codimension entity not implemented for SGrid." );
1344 return grid.size( grid.maxLevel(), type );
1350 return grid.size( grid.maxLevel(), codim );
1354 template<
class EntityType >
1357 return (e.level() == grid.maxLevel());
1361 const std::vector<GeometryType>&
geomTypes (
int codim)
const
1363 return mytypes[codim];
1367 const GridImp& grid;
1368 std::vector<GeometryType> mytypes[dim+1];
1379 template<
class Gr
idImp>
1381 public IdSet<GridImp,SGridGlobalIdSet<GridImp>, typename remove_const<GridImp>::type::PersistentIndexType>
1400 typedef typename remove_const<GridImp>::type::PersistentIndexType
IdType;
1413 IdType id (
const typename remove_const<GridImp>::type::Traits::template Codim<cd>::Entity& e)
const
1415 return grid.getRealImplementation(e).persistentIndex();
1423 IdType subId (
const typename remove_const< GridImp >::type::Traits::template Codim< 0 >::Entity &e,
1424 int i,
unsigned int codim )
const
1426 return grid.getRealImplementation(e).subPersistentIndex(codim, i);
1430 const GridImp& grid;
1434 template<
int dim,
int dimworld,
class ctype>
1449 bigunsignedint<dim*sgrid_dim_bits+sgrid_level_bits+sgrid_codim_bits>,
1451 bigunsignedint<dim*sgrid_dim_bits+sgrid_level_bits+sgrid_codim_bits>,
1452 CollectiveCommunication<Dune::SGrid<dim,dimworld,ctype> >,
1513 template<
int dim,
int dimworld,
typename _ctype = sgr
id_ctype>
1542 SGrid (
const int *
const N_,
const ctype *
const H_);
1551 SGrid (
const int *
const N_,
const ctype *
const L_,
const ctype *
const H_);
1562 SGrid (FieldVector<int,dim> N_, FieldVector<ctype,dim> L_, FieldVector<ctype,dim> H_);
1575 template<
int cd, PartitionIteratorType pitype>
1579 template<
int cd, PartitionIteratorType pitype>
1586 return lbegin<cd,All_Partition>(level);
1593 return lend<cd,All_Partition>(level);
1597 template<
int cd, PartitionIteratorType pitype>
1598 typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator
leafbegin ()
const;
1601 template<
int cd, PartitionIteratorType pitype>
1602 typename Traits::template Codim<cd>::template Partition<pitype>::LeafIterator
leafend ()
const;
1608 return leafbegin<cd,All_Partition>();
1615 return leafend<cd,All_Partition>();
1619 template <
typename Seed>
1620 typename Traits::template Codim<Seed::codimension>::EntityPointer
1623 enum { codim = Seed::codimension };
1640 template<
class T,
template<
class>
class P,
int codim>
1648 int size (
int level,
int codim)
const;
1659 return (type.isCube()) ?
size(level,dim-type.dim()) : 0;
1671 return boundarysize;
1695 const array<int, dim>&
dims(
int level)
const {
1719 return *theglobalidset;
1724 return *theglobalidset;
1729 assert(level>=0 && level<=
maxLevel());
1730 return *(indexsets[level]);
1735 return *theleafindexset;
1742 template<
class DataHandle>
1747 template<
class DataHandle>
1803 template<int codim_, int dim_, class GridImp_, template<int,int,class> class EntityImp_>
1806 template<int codim_, class GridImp_>
1809 template<int codim_, int dim_, class GridImp_, template<int,int,class> class EntityImp_>
1813 FieldVector<ctype, dimworld> pos (int level, array<int,dim>& z) const;
1816 int calc_codim (int level, const array<int,dim>& z) const;
1819 int n (int level, const array<int,dim>& z) const;
1822 array<int,dim> z (int level, int i, int codim) const;
1825 array<int,dim> subz (const array<int,dim> & z, int i, int codim) const;
1828 array<int,dim> compress (int level, const array<int,dim>& z) const;
1831 array<int,dim> expand (int level, const array<int,dim>& r, int b) const;
1836 int partition (int level, const array<int,dim>& z) const;
1839 bool exists (int level, const array<int,dim>& zred) const;
1842 int boundarySegmentIndex (int l, int face, const array<int,dim> & zentity) const
1844 array<int,dim-1> zface;
1848 for (
int i=0; i<dir; i++) zface[i] = zentity[i]/(1<<l);
1849 for (
int i=dir+1; i<dim; i++) zface[i-1] = zentity[i]/(1<<l);
1850 zface = boundarymapper[dir].expand(zface, 0);
1852 int index = boundarymapper[dir].n(zface);
1854 for (
int i=0; i<dir; i++)
1855 index += 2*boundarymapper[i].elements(0);
1856 index += side*boundarymapper[dir].elements(0);
1861 PersistentIndexType persistentIndex (
int l,
int codim,
const array<int,dim> & zentity)
const
1874 for (
int i=dim-1; i>=0; i--)
1886 int trailing = 1000;
1887 for (
int i=0; i<dim; i++)
1891 for (
int j=0; j<l; j++)
1892 if (zentity[i]&(1<<(j+1)))
1896 trailing =
std::min(trailing,zeros);
1900 int level = l-trailing;
1910 for (
int i=dim-1; i>=0; i--)
1922 SGrid & operator = (
const SGrid &) {
return *
this; };
1924 void makeSGrid (
const array<int,dim>& N_,
const FieldVector<ctype, dim>& L_,
const FieldVector<ctype, dim>& H_);
1929 CollectiveCommunication<SGrid> ccobj;
1931 ReservedVector<SGridLevelIndexSet<const SGrid<dim,dimworld> >*,
MAXL> indexsets;
1932 SGridLeafIndexSet<const SGrid<dim,dimworld> > *theleafindexset;
1933 SGridGlobalIdSet<const SGrid<dim,dimworld> > *theglobalidset;
1936 FieldVector<ctype, dim> low;
1937 FieldVector<ctype, dim> H;
1939 FieldVector<ctype, dim> *h;
1940 mutable CubeMapper<dim> *mapper;
1943 array<CubeMapper<dim-1>, dim> boundarymapper;
1947 mutable array <int,dim> zrefStatic;
1948 mutable array <int,dim> zentityStatic;
1951 namespace Capabilities
1965 template<
int dim,
int dimw>
1968 static const bool v =
true;
1969 static const unsigned int topologyId = GenericGeometry :: CubeTopology< dim > :: type :: id ;
1975 template<
int dim,
int dimw>
1978 static const bool v =
true;
1984 template<
int dim,
int dimw,
int cdim>
1987 static const bool v =
true;
1993 template<
int dim,
int dimw>
1996 static const bool v =
true;
2002 template<
int dim,
int dimw>
2005 static const bool v =
true;
2012 #include"sgrid/sgrid.cc"