1 #ifndef DUNE_FEM_CHECKGEOMETRYAFFINITY_HH 2 #define DUNE_FEM_CHECKGEOMETRYAFFINITY_HH 4 #include <dune/common/fvector.hh> 5 #include <dune/grid/common/gridenums.hh> 20 template <
class QuadratureType>
25 template<
class EntityType,
class ElementGeometryType >
27 const ElementGeometryType& geo,
35 QuadratureType volQuad( entity, quadOrd );
36 const int nop = volQuad.nop();
39 const double oldIntel = geo.integrationElement( volQuad.point(0) );
40 for(
int l=1; l<nop; ++l)
42 const double intel = geo.integrationElement( volQuad.point(l) );
43 if(
std::abs( oldIntel - intel ) > 1e-12 )
51 template<
class EntityType >
52 static bool checkGeometry(
const EntityType& entity,
const int quadOrd )
58 template <
class IteratorType>
60 const IteratorType& endit,
63 for(IteratorType it = begin; it != endit; ++it)
71 template <
class Gr
idPartType,
class Vector >
74 Vector& affineGeomtryVec )
76 typedef typename GridPartType :: template Codim< 0 > :: IteratorType IteratorType;
77 typedef typename GridPartType :: template Codim< 0 > :: EntityType EntityType;
78 const IteratorType endit = gridPart.template end<0> ();
79 affineGeomtryVec.resize( gridPart.indexSet().size( 0 ) );
80 for(IteratorType it = gridPart.template begin<0>(); it != endit; ++it)
82 const EntityType& entity = *it ;
83 const int index = gridPart.indexSet().index( entity );
93 template <
class Gr
idPartType>
96 typedef typename GridPartType :: GridType
GridType ;
99 template <
class IndexSetType>
100 static inline bool doCheck(
const GridType& grid,
const IndexSetType& indexSet)
102 typedef typename GridType :: template Codim<0> :: LevelIterator MacroIteratorType;
103 typedef typename GridType :: template Codim<0> :: Entity EntityType;
104 typedef typename GridType :: template Codim<0> :: Geometry Geometry;
106 typedef typename GridType :: template Partition< All_Partition > :: LevelGridView MacroGridView ;
109 MacroGridView macroView = grid.levelGridView( 0 );
111 const MacroIteratorType endit = macroView.template end<0> ();
112 MacroIteratorType it = macroView.template begin<0> ();
115 if( it == endit )
return true;
118 GeometryInformationType geoInfo( indexSet );
121 if( geoInfo.geomTypes(0).size() > 1 )
return false;
124 if( ! geoInfo.geomTypes(0)[0].isCube() )
return false;
126 typedef typename GridType :: ctype ctype;
127 enum { dimension = GridType :: dimension };
128 enum { dimworld = GridType :: dimensionworld };
131 FieldVector<ctype, dimension> h(0);
132 FieldVector<ctype, dimworld> enBary;
133 FieldVector<ctype, dimension-1> mid(0.5);
135 const int map[3] = {1, 2, 4};
137 const Geometry geo = it->geometry();
138 if ( ! geo.type().isCube() )
return false;
141 for(
int i=0; i<dimension; ++i)
143 h[i] = (geo.corner( 0 ) - geo.corner( map[i] )).two_norm();
148 for(MacroIteratorType it = macroView.template begin<0> ();
151 const EntityType& en = *it;
152 const Geometry geo = en.geometry();
154 const FieldVector<ctype, dimworld> enBary =
155 geo.global( geoInfo.localCenter( geo.type() ));
157 typedef typename MacroGridView :: IntersectionIterator IntersectionIteratorType;
160 if ( ! geo.type().isCube() )
return false;
162 for(
int i=0; i<dimension; ++i)
164 const ctype w = (geo.corner(0) - geo.corner( map[i] )).two_norm();
165 if(
std::abs( h[i] - w ) > 1e-15 )
return false;
168 IntersectionIteratorType endnit = macroView.iend( en );
169 for(IntersectionIteratorType nit = macroView.ibegin( en );
170 nit != endnit; ++nit)
172 const typename IntersectionIteratorType::Intersection& inter=*nit;
173 if( ! checkIntersection(inter) )
return false;
175 if( inter.neighbor() )
178 Geometry nbGeo = nb.geometry();
180 FieldVector<ctype, dimworld> diff = nbGeo.global( geoInfo.localCenter( nbGeo.type() ));
182 assert( diff.two_norm() > 1e-15 );
183 diff /= diff.two_norm();
186 if(
std::abs(
std::abs((diff * inter.unitOuterNormal(mid))) - 1.0) > 1e-12 )
return false;
193 template <
class ctype,
int dim>
196 const FieldVector<ctype,dim-1> mid_;
197 enum { numberOfNormals = 2 * dim };
198 FieldVector<ctype,dim> refNormal_[numberOfNormals];
202 for(
int i=0; i<numberOfNormals; ++i)
205 FieldVector<ctype,dim>& refNormal = refNormal_[i];
209 int comp = ((int) i/2);
210 refNormal[comp] = ((i % 2) == 0) ? -1 : 1;
216 assert( i >= 0 && i< numberOfNormals );
217 return refNormal_[i];
228 template <
class IntersectionType>
231 if ( ! nit.type().isCube() )
return false;
233 typedef typename IntersectionType :: Entity EntityType;
234 typedef typename EntityType :: Geometry :: ctype ctype;
235 enum { dimworld = EntityType :: Geometry :: coorddimension };
241 FieldVector<ctype,dimworld> unitNormal = nit.unitOuterNormal(normals.
faceMidPoint());
245 if( unitNormal.infinity_norm() > 1e-10 )
return false;
251 static inline bool check(
const GridPartType& gridPart)
253 bool cartesian = doCheck( gridPart.grid() , gridPart.indexSet() );
254 int val = (cartesian) ? 1 : 0;
256 val = gridPart.comm().min( val );
257 return (val == 1) ?
true :
false;
267 #endif // #ifndef DUNE_FEM_CHECKGEOMETRYAFFINITY_HH const FieldVector< ctype, dim-1 > & faceMidPoint() const
Definition: checkgeomaffinity.hh:220
static bool checkIntersection(const IntersectionType &nit)
Definition: checkgeomaffinity.hh:229
const FieldVector< ctype, dim > & referenceNormal(const int i) const
Definition: checkgeomaffinity.hh:214
ReferenceNormals()
Definition: checkgeomaffinity.hh:200
static bool checkGeometry(const EntityType &entity, const int quadOrd)
Definition: checkgeomaffinity.hh:52
Helper class to check whether grid is structured or not.
Definition: checkgeomaffinity.hh:94
static bool checkGeometry(const EntityType &entity, const ElementGeometryType &geo, const int quadOrd)
Definition: checkgeomaffinity.hh:26
GridPartType::GridType GridType
Definition: checkgeomaffinity.hh:96
Double abs(const Double &a)
Definition: double.hh:860
Dune::EntityPointer< Grid, Implementation >::Entity make_entity(const Dune::EntityPointer< Grid, Implementation > &entityPointer)
Definition: compatibility.hh:23
static bool check(const GridPartType &gridPart)
check whether all the is grid is a cartesian grid
Definition: checkgeomaffinity.hh:251
static bool doCheck(const GridType &grid, const IndexSetType &indexSet)
check whether this is a cartesian grid or not
Definition: checkgeomaffinity.hh:100
Definition: coordinate.hh:4
static void checkElementAffinity(const GridPartType &gridPart, const int quadOrd, Vector &affineGeomtryVec)
check whether all geometry mappings are affine
Definition: checkgeomaffinity.hh:72
Helper class to check affinity of the grid's geometries.
Definition: checkgeomaffinity.hh:22
default implementation uses method geomTypes of given index set. Used in DiscreteFunctionSpaces.
Definition: allgeomtypes.hh:89
static bool checkAffinity(const IteratorType &begin, const IteratorType &endit, const int quadOrd)
check whether all geometry mappings are affine
Definition: checkgeomaffinity.hh:59
Definition: checkgeomaffinity.hh:194