4 #ifndef DUNE_PDELAB_CONSTRAINTS_HANGINGNODEMANAGER_HH 5 #define DUNE_PDELAB_CONSTRAINTS_HANGINGNODEMANAGER_HH 7 #include<dune/grid/common/grid.hh> 8 #include<dune/grid/common/mcmgmapper.hh> 9 #include<dune/common/float_cmp.hh> 11 #include"../common/geometrywrapper.hh" 16 template<
class Gr
id,
class BoundaryFunction>
21 enum{ verbosity = 2 };
23 enum{ verbosity = 0 };
25 typedef typename Grid::LeafIndexSet::IndexType IndexType;
32 unsigned short minimum_level;
35 unsigned short maximum_level;
38 unsigned short minimum_touching_level;
42 void addLevel(
unsigned short level)
44 minimum_level = std::min(minimum_level,level);
45 maximum_level = std::max(maximum_level,level);
48 void addTouchingLevel(
unsigned short level)
50 minimum_touching_level = std::min(minimum_touching_level,level);
53 NodeInfo() : minimum_level(std::numeric_limits<unsigned short>::max()),
55 minimum_touching_level(std::numeric_limits<unsigned short>::max()),
60 std::vector<NodeInfo> node_info;
74 inline bool isHanging()
const {
return is_hanging; }
83 enum {
dim = GridView::dimension};
84 typedef typename GridView::template Codim<0>::Entity
Cell;
86 typedef typename GridView::template Codim<0>::Iterator
Iterator;
88 typedef typename GridView::Grid::ctype
ctype;
89 typedef typename Dune::FieldVector<ctype,dim>
Point;
92 typedef Dune::MultipleCodimMultipleGeomTypeMapper<
GridView,
103 cell_mapper.update();
104 const typename GridView::IndexSet& indexSet = grid.leafGridView().indexSet();
106 node_info = std::vector<NodeInfo>(indexSet.size(
dim));
108 GridView gv = grid.leafGridView();
111 for(
const auto& cell : elements(gv)) {
113 const Dune::ReferenceElement<double,dim> &
115 Dune::ReferenceElements<double,dim>::general(cell.geometry().type());
119 const unsigned short level = cell.level();
122 const IndexType v_size = reference_element.size(
dim);
127 for(IndexType i=0; i<v_size; ++i){
128 const IndexType v_globalindex = indexSet.subIndex(cell,i,
dim);
129 NodeInfo& ni = node_info[v_globalindex];
134 std::cout <<
" cell-id=" << cell_mapper.index(cell);
135 std::cout <<
" level=" << level;
136 std::cout <<
" v_size=" << v_size;
137 std::cout <<
" v_globalindex = " << v_globalindex;
138 std::cout <<
" maximum_level = " << ni.maximum_level;
139 std::cout <<
" minimum_level = " << ni.minimum_level;
140 std::cout << std::endl;
148 typedef typename IntersectionIterator::Intersection Intersection;
151 unsigned int intersection_index = 0;
152 for(
const auto& intersection : intersections(gv,cell)) {
153 ++intersection_index;
155 const Dune::ReferenceElement<double,
dim-1> &
156 reference_face_element =
157 Dune::ReferenceElements<double,dim-1>::general(intersection.geometry().type());
159 const int eLocalIndex = intersection.indexInInside();
160 const int e_level = intersection.inside().level();
163 const int e_v_size = reference_element.size(eLocalIndex,1,
dim);
165 if(intersection.boundary()) {
168 for(
int i=0; i<e_v_size;++i){
169 const int e_v_index = reference_element.subEntity(eLocalIndex,1,i,
dim);
170 const IndexType v_globalindex = indexSet.subIndex(cell,e_v_index,
dim);
172 const FacePoint facelocal_position = reference_face_element.position(i,
dim-1);
184 NodeInfo& ni = node_info[v_globalindex];
186 ni.addTouchingLevel(e_level);
193 const int f_level = intersection.outside().level();
196 if(intersection.conforming())
200 assert(e_level != f_level);
204 if(e_level < f_level)
208 for(
int i=0; i<e_v_size;++i){
209 const int e_v_index = reference_element.subEntity(eLocalIndex,1,i,
dim);
210 const IndexType v_globalindex = indexSet.subIndex(cell,e_v_index,
dim);
212 node_info[v_globalindex].addTouchingLevel(f_level);
222 boundaryFunction(_boundaryFunction),
223 cell_mapper(grid.leafGridView())
228 const typename GridView::IndexSet& indexSet = grid.leafGridView().indexSet();
229 std::vector<NodeState> is_hanging;
231 const Dune::ReferenceElement<double,dim> &
233 Dune::ReferenceElements<double,dim>::general(e.geometry().type());
236 const IndexType v_size = reference_element.size(
dim);
239 is_hanging.resize(v_size);
244 for(IndexType i=0; i<v_size; ++i){
245 const IndexType v_globalindex = indexSet.subIndex(e,i,
dim);
250 const NodeInfo & v_info = node_info[v_globalindex];
251 if(v_info.minimum_touching_level < v_info.minimum_level){
252 is_hanging[i].is_hanging =
true;
253 is_hanging[i].is_boundary = v_info.is_boundary;
256 const Point & local = reference_element.position(i,
dim);
257 const Point global = e.geometry().global(local);
259 std::cout <<
"Found hanging node with id " << v_globalindex <<
" at " << global << std::endl;
264 is_hanging[i].is_hanging =
false;
265 is_hanging[i].is_boundary = v_info.is_boundary;
275 std::cout <<
"Begin isolation of hanging nodes" << std::endl;
277 const typename GridView::IndexSet& indexSet = grid.leafGridView().indexSet();
279 size_t iterations(0);
285 size_t refinements(0);
288 GridView gv = grid.leafGridView();
291 for(
const auto& cell : elements(gv)) {
293 const Dune::ReferenceElement<double,dim> &
295 Dune::ReferenceElements<double,dim>::general(cell.geometry().type());
300 const unsigned short level = cell.level();
305 const IndexType v_size = reference_element.size(
dim);
310 for(IndexType i=0; i<v_size; ++i){
312 const IndexType v_globalindex = indexSet.subIndex(cell,i,
dim);
314 const NodeInfo & v_info = node_info[v_globalindex];
318 const unsigned short level_diff = v_info.maximum_level - level;
326 std::cout <<
" cell-id=" << cell_mapper.index(cell);
327 std::cout <<
" level=" << level;
328 std::cout <<
" v_size=" << v_size;
329 std::cout <<
" v_globalindex = " << v_globalindex;
330 std::cout << std::endl;
331 std::cout <<
"Refining element nr " << cell_mapper.index(cell)
332 <<
" to isolate hanging nodes. Level diff = " 333 << v_info.maximum_level <<
" - " << level<< std::endl;
340 if( cell.geometry().type().isSimplex() ){
351 unsigned int intersection_index = 0;
353 bool bJumpOut =
false;
356 for(
const auto& intersection : intersections(gv,cell)) {
357 ++intersection_index;
360 if(!intersection.boundary()) {
362 const int e_level = intersection.inside().level();
363 const int f_level = intersection.outside().level();
365 if( f_level > e_level ){
371 const int eLocalIndex = intersection.indexInInside();
378 int nEdgeCenters = 0;
395 std::vector<Dune::FieldVector<ctype,dim> >
396 edgecenter( nEdgeCenters, Dune::FieldVector<ctype,dim>(0) );
400 for(
int counter=0; counter<nEdgeCenters; ++counter){
402 int cornerIndex1 = counter %
dim;
403 int cornerIndex2 = (counter+1) %
dim;
405 const int e_v_index_1 =
406 reference_element.subEntity(eLocalIndex,1,cornerIndex1,
dim);
408 const int e_v_index_2 =
409 reference_element.subEntity(eLocalIndex,1,cornerIndex2,
dim);
411 auto vertex1 = cell.template subEntity<dim>(e_v_index_1);
413 auto vertex2 = cell.template subEntity<dim>(e_v_index_2);
415 edgecenter[counter] += vertex1.geometry().center();
416 edgecenter[counter] += vertex2.geometry().center();
417 edgecenter[counter] /=
ctype( 2 );
424 const Dune::ReferenceElement<double,dim> &
425 nb_reference_element =
426 Dune::ReferenceElements<double,dim>::general( intersection.outside().geometry().type() );
429 const IndexType nb_v_size = nb_reference_element.size(
dim);
432 for(IndexType i=0; i<nb_v_size; ++i){
434 auto nb_vertex = intersection.outside().template subEntity<dim>(i);
436 bool doExtraCheck =
false;
438 Dune::FieldVector<ctype,dim> center_diff ( edgecenter[counter] );
440 center_diff -= nb_vertex.geometry().center();
444 if( center_diff.two_norm2() < 5
e-12 ){
454 const IndexType nb_v_globalindex = indexSet.index(nb_vertex);
456 const NodeInfo & nb_v_info = node_info[nb_v_globalindex];
458 const unsigned short level_diff = nb_v_info.maximum_level - level;
468 std::cout <<
" cell-id=" << cell_mapper.index(cell);
469 std::cout <<
" level=" << level;
470 std::cout <<
" v_size=" << v_size;
471 std::cout << std::endl;
472 std::cout <<
" Extra refining for element nr " << cell_mapper.index(cell)
473 <<
" to isolate hanging nodes. Level diff = " 474 << nb_v_info.maximum_level <<
" - " << level<< std::endl;
481 if( bJumpOut )
break;
483 if( bJumpOut )
break;
489 if( bJumpOut )
break;
501 std::cout <<
"Re-adapt for isolation of hanging nodes..." << std::endl;
510 std::cout <<
"In iteration " << iterations <<
" there were " 511 << refinements <<
" grid cells to be refined additionally to isolate hanging nodes." << std::endl;
519 #endif // DUNE_PDELAB_CONSTRAINTS_HANGINGNODEMANAGER_HH Grid & grid
Definition: hangingnodemanager.hh:95
CellMapper cell_mapper
Definition: hangingnodemanager.hh:97
Definition: hangingnodemanager.hh:17
Dune::FieldVector< ctype, dim-1 > FacePoint
Definition: hangingnodemanager.hh:90
const Entity & e
Definition: localfunctionspace.hh:111
Wrap intersection.
Definition: geometrywrapper.hh:56
GridView::Grid::ctype ctype
Definition: hangingnodemanager.hh:88
void analyzeView()
Definition: hangingnodemanager.hh:101
Definition: hangingnodemanager.hh:83
Dune::MultipleCodimMultipleGeomTypeMapper< GridView, MCMGElementLayout > CellMapper
Definition: hangingnodemanager.hh:93
const BoundaryFunction & boundaryFunction
Definition: hangingnodemanager.hh:96
HangingNodeManager(Grid &_grid, const BoundaryFunction &_boundaryFunction)
Definition: hangingnodemanager.hh:220
const std::vector< NodeState > hangingNodes(const Cell &e) const
Definition: hangingnodemanager.hh:226
GridView::IntersectionIterator IntersectionIterator
Definition: hangingnodemanager.hh:87
For backward compatibility – Do not use this!
Definition: adaptivity.hh:27
GridView::template Codim< 0 >::Entity Cell
Definition: hangingnodemanager.hh:84
Grid::LeafGridView GridView
Definition: hangingnodemanager.hh:82
void adaptToIsolatedHangingNodes()
Definition: hangingnodemanager.hh:272
Definition: hangingnodemanager.hh:64
Dune::FieldVector< ctype, dim > Point
Definition: hangingnodemanager.hh:89
bool isHanging() const
Definition: hangingnodemanager.hh:74
NodeState()
Definition: hangingnodemanager.hh:76
GridView::template Codim< 0 >::Iterator Iterator
Definition: hangingnodemanager.hh:86
bool isBoundary() const
Definition: hangingnodemanager.hh:73