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>
CellMapper;
103 const typename GridView::IndexSet& indexSet =
grid.leafGridView().indexSet();
105 node_info = std::vector<NodeInfo>(indexSet.size(
dim));
110 for(
const auto& cell : elements(gv)) {
112 auto reference_element = referenceElement(cell.geometry());
115 const unsigned short level = cell.level();
118 const IndexType v_size = reference_element.size(
dim);
123 for(IndexType i=0; i<v_size; ++i){
124 const IndexType v_globalindex = indexSet.subIndex(cell,i,
dim);
125 NodeInfo& ni = node_info[v_globalindex];
130 std::cout <<
" cell-id=" <<
cell_mapper.index(cell);
131 std::cout <<
" level=" << level;
132 std::cout <<
" v_size=" << v_size;
133 std::cout <<
" v_globalindex = " << v_globalindex;
134 std::cout <<
" maximum_level = " << ni.maximum_level;
135 std::cout <<
" minimum_level = " << ni.minimum_level;
136 std::cout << std::endl;
144 typedef typename IntersectionIterator::Intersection Intersection;
147 unsigned int intersection_index = 0;
148 for(
const auto& intersection : intersections(gv,cell)) {
149 ++intersection_index;
151 auto reference_face_element = referenceElement(intersection.geometry());
153 const int eLocalIndex = intersection.indexInInside();
154 const int e_level = intersection.inside().level();
157 const int e_v_size = reference_element.size(eLocalIndex,1,
dim);
159 if(intersection.boundary()) {
162 for(
int i=0; i<e_v_size;++i){
163 const int e_v_index = reference_element.subEntity(eLocalIndex,1,i,
dim);
164 const IndexType v_globalindex = indexSet.subIndex(cell,e_v_index,
dim);
166 const FacePoint facelocal_position = reference_face_element.position(i,
dim-1);
178 NodeInfo& ni = node_info[v_globalindex];
180 ni.addTouchingLevel(e_level);
187 const int f_level = intersection.outside().level();
190 if(intersection.conforming())
194 assert(e_level != f_level);
198 if(e_level < f_level)
202 for(
int i=0; i<e_v_size;++i){
203 const int e_v_index = reference_element.subEntity(eLocalIndex,1,i,
dim);
204 const IndexType v_globalindex = indexSet.subIndex(cell,e_v_index,
dim);
206 node_info[v_globalindex].addTouchingLevel(f_level);
222 const typename GridView::IndexSet& indexSet =
grid.leafGridView().indexSet();
223 std::vector<NodeState> is_hanging;
225 auto reference_element = referenceElement(e.geometry());
228 const IndexType v_size = reference_element.size(
dim);
231 is_hanging.resize(v_size);
236 for(IndexType i=0; i<v_size; ++i){
237 const IndexType v_globalindex = indexSet.subIndex(e,i,
dim);
242 const NodeInfo & v_info = node_info[v_globalindex];
243 if(v_info.minimum_touching_level < v_info.minimum_level){
244 is_hanging[i].is_hanging =
true;
245 is_hanging[i].is_boundary = v_info.is_boundary;
248 const Point & local = reference_element.position(i,
dim);
249 const Point global = e.geometry().global(local);
251 std::cout <<
"Found hanging node with id " << v_globalindex <<
" at " << global << std::endl;
256 is_hanging[i].is_hanging =
false;
257 is_hanging[i].is_boundary = v_info.is_boundary;
267 std::cout <<
"Begin isolation of hanging nodes" << std::endl;
269 const typename GridView::IndexSet& indexSet =
grid.leafGridView().indexSet();
271 size_t iterations(0);
277 size_t refinements(0);
283 for(
const auto& cell : elements(gv)) {
285 auto reference_element = referenceElement(cell.geometry());
290 const unsigned short level = cell.level();
295 const IndexType v_size = reference_element.size(
dim);
300 for(IndexType i=0; i<v_size; ++i){
302 const IndexType v_globalindex = indexSet.subIndex(cell,i,
dim);
304 const NodeInfo & v_info = node_info[v_globalindex];
308 const unsigned short level_diff = v_info.maximum_level - level;
316 std::cout <<
" cell-id=" <<
cell_mapper.index(cell);
317 std::cout <<
" level=" << level;
318 std::cout <<
" v_size=" << v_size;
319 std::cout <<
" v_globalindex = " << v_globalindex;
320 std::cout << std::endl;
321 std::cout <<
"Refining element nr " <<
cell_mapper.index(cell)
322 <<
" to isolate hanging nodes. Level diff = "
323 << v_info.maximum_level <<
" - " << level<< std::endl;
330 if( cell.geometry().type().isSimplex() ){
341 unsigned int intersection_index = 0;
343 bool bJumpOut =
false;
346 for(
const auto& intersection : intersections(gv,cell)) {
347 ++intersection_index;
350 if(!intersection.boundary()) {
352 const int e_level = intersection.inside().level();
353 const int f_level = intersection.outside().level();
355 if( f_level > e_level ){
361 const int eLocalIndex = intersection.indexInInside();
368 int nEdgeCenters = 0;
385 std::vector<Dune::FieldVector<ctype,dim> >
386 edgecenter( nEdgeCenters, Dune::FieldVector<ctype,dim>(0) );
390 for(
int counter=0; counter<nEdgeCenters; ++counter){
392 int cornerIndex1 = counter %
dim;
393 int cornerIndex2 = (counter+1) %
dim;
395 const int e_v_index_1 =
396 reference_element.subEntity(eLocalIndex,1,cornerIndex1,
dim);
398 const int e_v_index_2 =
399 reference_element.subEntity(eLocalIndex,1,cornerIndex2,
dim);
401 auto vertex1 = cell.template subEntity<dim>(e_v_index_1);
403 auto vertex2 = cell.template subEntity<dim>(e_v_index_2);
405 edgecenter[counter] += vertex1.geometry().center();
406 edgecenter[counter] += vertex2.geometry().center();
407 edgecenter[counter] /=
ctype( 2 );
414 auto nb_reference_element = referenceElement(intersection.outside().geometry());
417 const IndexType nb_v_size = nb_reference_element.size(
dim);
420 for(IndexType i=0; i<nb_v_size; ++i){
422 auto nb_vertex = intersection.outside().template subEntity<dim>(i);
424 bool doExtraCheck =
false;
426 Dune::FieldVector<ctype,dim> center_diff ( edgecenter[counter] );
428 center_diff -= nb_vertex.geometry().center();
432 if( center_diff.two_norm2() < 5e-12 ){
442 const IndexType nb_v_globalindex = indexSet.index(nb_vertex);
444 const NodeInfo & nb_v_info = node_info[nb_v_globalindex];
446 const unsigned short level_diff = nb_v_info.maximum_level - level;
456 std::cout <<
" cell-id=" <<
cell_mapper.index(cell);
457 std::cout <<
" level=" << level;
458 std::cout <<
" v_size=" << v_size;
459 std::cout << std::endl;
460 std::cout <<
" Extra refining for element nr " <<
cell_mapper.index(cell)
461 <<
" to isolate hanging nodes. Level diff = "
462 << nb_v_info.maximum_level <<
" - " << level<< std::endl;
469 if( bJumpOut )
break;
471 if( bJumpOut )
break;
477 if( bJumpOut )
break;
489 std::cout <<
"Re-adapt for isolation of hanging nodes..." << std::endl;
498 std::cout <<
"In iteration " << iterations <<
" there were "
499 << refinements <<
" grid cells to be refined additionally to isolate hanging nodes." << std::endl;
For backward compatibility – Do not use this!
Definition: adaptivity.hh:28
Wrap intersection.
Definition: geometrywrapper.hh:57
Definition: hangingnodemanager.hh:18
const BoundaryFunction & boundaryFunction
Definition: hangingnodemanager.hh:95
const std::vector< NodeState > hangingNodes(const Cell &e) const
Definition: hangingnodemanager.hh:220
Dune::MultipleCodimMultipleGeomTypeMapper< GridView > CellMapper
Definition: hangingnodemanager.hh:92
Dune::FieldVector< ctype, dim-1 > FacePoint
Definition: hangingnodemanager.hh:90
Grid::LeafGridView GridView
Definition: hangingnodemanager.hh:82
void adaptToIsolatedHangingNodes()
Definition: hangingnodemanager.hh:264
GridView::IntersectionIterator IntersectionIterator
Definition: hangingnodemanager.hh:87
GridView::template Codim< 0 >::Entity Cell
Definition: hangingnodemanager.hh:84
CellMapper cell_mapper
Definition: hangingnodemanager.hh:96
GridView::template Codim< 0 >::Iterator Iterator
Definition: hangingnodemanager.hh:86
Dune::FieldVector< ctype, dim > Point
Definition: hangingnodemanager.hh:89
HangingNodeManager(Grid &_grid, const BoundaryFunction &_boundaryFunction)
Definition: hangingnodemanager.hh:214
Grid & grid
Definition: hangingnodemanager.hh:94
@ dim
Definition: hangingnodemanager.hh:83
GridView::Grid::ctype ctype
Definition: hangingnodemanager.hh:88
void analyzeView()
Definition: hangingnodemanager.hh:100
Definition: hangingnodemanager.hh:65
NodeState()
Definition: hangingnodemanager.hh:76
bool isHanging() const
Definition: hangingnodemanager.hh:74
bool isBoundary() const
Definition: hangingnodemanager.hh:73