Dune::Copasi
Loading...
Searching...
No Matches
mark_stripes.hh
Go to the documentation of this file.
1#ifndef DUNE_COPASI_GRID_MARK_STRIPES_HH
2#define DUNE_COPASI_GRID_MARK_STRIPES_HH
3
4#include <dune-copasi-config.hh>
5
6#include <dune/grid/uggrid.hh>
7
8#include <list>
9
10namespace Dune::Copasi {
11
25template<int dim>
26void mark_stripes(UGGrid<dim>& grid, bool mark_others = true)
27{
28 using RuleType = typename UG_NS<dim>::RefinementRule;
29
30 auto grid_view = grid.leafGridView();
31 std::list<int> non_cube_side;
32
33 // Loop over the grid
34 for (auto&& entity : Dune::elements(grid_view))
35 {
36 if (entity.type().isCube())
37 {
38 // register side index with simplices (see DUNE cheatsheet for entity ids)
39 non_cube_side.clear();
40 for (auto&& ig : Dune::intersections(grid_view,entity))
41 if(ig.neighbor())
42 if (not ig.outside().type().isCube())
43 non_cube_side.push_back(ig.indexInInside());
44
45 bool is_stripe = false;
46
47 // opposite facets have consecutive indexing (e.g. [2,3] are opposite)
48 if (non_cube_side.size() == 2)
49 is_stripe = !(non_cube_side.front()/2 - non_cube_side.back()/2);
50
51 if (is_stripe)
52 {
53 // side orientation of the simplices
54 [[maybe_unused]] int orientation = *(non_cube_side.begin())/2;
55
56 if constexpr (dim == 2)
57 {
58 // mark entity with a blue type refinement
59 grid.mark(entity,RuleType::BLUE,!(bool)orientation);
60 }
61 else if constexpr (dim == 3)
62 {
63 DUNE_THROW(NotImplemented,"\tStripes on 3D is not available yet!");
64 // Need a mesh to correctly check which orientation needs which rule!
65 // if (orientation == 0)
66 // grid.mark(entity,RuleType::HEX_BISECT_0_1);
67 // if (orientation == 1)
68 // grid.mark(entity,RuleType::HEX_BISECT_0_2);
69 // if (orientation == 2)
70 // grid.mark(entity,RuleType::HEX_BISECT_0_3);
71 }
72 else
73 {
75 "\tStripe refinement not known for grids of dimension '"
76 << dim << "'");
77 }
78 }
79 else if (mark_others)
80 {
81 grid.mark(1,entity);
82 }
83 }
84 else if (mark_others)
85 {
86 grid.mark(1,entity);
87 }
88 }
89}
90
91} // namespace Dune::Copasi
92
93#endif // DUNE_COPASI_GRID_MARK_STRIPES_HH
Definition axis_names.hh:7
void mark_stripes(UGGrid< dim > &grid, bool mark_others=true)
Mark stripes for refinement.
Definition mark_stripes.hh:26
constexpr bool is_bitflags_v
Alias for Bitflag indicator.
Definition bit_flags.hh:24