Skip to content
Snippets Groups Projects
Commit b2b761ea authored by Dennis Gläser's avatar Dennis Gläser
Browse files

Merge branch 'fix/gridmanager-subgrid-spgrid' into 'master'

[gridmanager][bugfix] Move subgrid code to correct header

See merge request !1583
parents 08026235 a9b1f776
No related branches found
No related tags found
1 merge request!1583[gridmanager][bugfix] Move subgrid code to correct header
......@@ -34,8 +34,6 @@
#include <dumux/io/grid/gridmanager_base.hh>
#endif
#include <dumux/common/boundaryflag.hh>
namespace Dumux {
#if HAVE_DUNE_SPGRID
......@@ -99,25 +97,6 @@ public:
}
};
//! dune-subgrid doesn't have this implemented
template<int dim, class HostGrid>
class BoundaryFlag<Dune::SubGrid<dim, HostGrid>>
{
public:
BoundaryFlag() : flag_(-1) {}
template<class Intersection>
BoundaryFlag(const Intersection& i) : flag_(-1) {}
using value_type = int;
value_type get() const
{ DUNE_THROW(Dune::NotImplemented, "Sub-grid doesn't implement boundary segment indices!"); }
private:
int flag_;
};
#endif // HAVE_DUNE_SPGRID
} // end namespace Dumux
......
......@@ -33,6 +33,7 @@
#include <dune/grid/io/file/dgfparser/dgfwriter.hh>
#include <dumux/common/parameters.hh>
#include <dumux/common/boundaryflag.hh>
namespace Dumux {
......@@ -99,6 +100,25 @@ public:
}
};
//! dune-subgrid doesn't have this implemented
template<int dim, class HostGrid>
class BoundaryFlag<Dune::SubGrid<dim, HostGrid>>
{
public:
BoundaryFlag() : flag_(-1) {}
template<class Intersection>
BoundaryFlag(const Intersection& i) : flag_(-1) {}
using value_type = int;
value_type get() const
{ DUNE_THROW(Dune::NotImplemented, "Sub-grid doesn't implement boundary segment indices!"); }
private:
int flag_;
};
} // end namespace Dumux
#endif // HAVE_DUNE_SUBGRID
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment