13#ifndef DUMUX_IO_GRID_PERIODIC_GRID_TRAITS_HH
14#define DUMUX_IO_GRID_PERIODIC_GRID_TRAITS_HH
17#include <dune/common/exceptions.hh>
18#include <dune/common/std/type_traits.hh>
19#include <dune/grid/common/exceptions.hh>
22#include <dune/grid/spgrid.hh>
26#include <dune/subgrid/subgrid.hh>
31template<
typename Gr
id>
38 bool isPeriodic (
const typename Grid::LeafIntersection& intersection)
const
45template<
class ct,
int dim,
template<
int >
class Ref,
class Comm>
46struct PeriodicGridTraits<
Dune::SPGrid<ct, dim, Ref, Comm>>
49 using Grid = Dune::SPGrid<ct, dim, Ref, Comm>;
51 struct SupportsPeriodicity :
public std::true_type {};
55 bool isPeriodic (
const typename Grid::LeafIntersection& intersection)
const
57 return intersection.neighbor() && intersection.boundary();
64template<
int dim,
typename HostGr
id,
bool MapIndexStorage>
65struct PeriodicGridTraits<
Dune::SubGrid<dim, HostGrid, MapIndexStorage>>
71 const PeriodicGridTraits<HostGrid> hostTraits_;
74 struct SupportsPeriodicity :
public PeriodicGridTraits<HostGrid>::SupportsPeriodicity {};
77 : subGrid_(subGrid), hostTraits_(subGrid_.getHostGrid()) {};
79 bool isPeriodic (
const typename Grid::LeafIntersection& intersection)
const
81 const auto& hostElement = subGrid_.template getHostEntity<0>(intersection.inside());
82 for (
const auto& hostIntersection : intersections(subGrid_.getHostGrid().leafGridView(), hostElement))
84 if (hostIntersection.indexInInside() == intersection.indexInInside())
86 const bool periodicInHostGrid = hostTraits_.isPeriodic(hostIntersection);
87 return periodicInHostGrid && subGrid_.template contains<0>(hostIntersection.outside());
93 void verifyConformingPeriodicBoundary()
const
95 for (
const auto& element : elements(subGrid_.leafGridView()))
97 for (
const auto& intersection : intersections(subGrid_.leafGridView(), element))
99 const auto& hostElement = subGrid_.template getHostEntity<0>(intersection.inside());
100 for (
const auto& hostIntersection : intersections(subGrid_.getHostGrid().leafGridView(), hostElement))
102 if (hostIntersection.indexInInside() == intersection.indexInInside())
104 const bool periodicInHostGrid = hostTraits_.isPeriodic(hostIntersection);
105 if (periodicInHostGrid && !subGrid_.template contains<0>(hostIntersection.outside()))
106 DUNE_THROW(Dune::GridError,
"Periodic boundary in host grid but outside"
107 <<
" element not included in subgrid. If this is intentional,"
108 <<
" take additional care with boundary conditions and remove"
109 <<
" verification call.");
123 using SP =
typename G::SupportsPeriodicity;
125 using type =
typename Dune::Std::detected_or<std::false_type, SP, T>::type;
Definition: periodicgridtraits.hh:121
typename Dune::Std::detected_or< std::false_type, SP, T >::type type
Definition: periodicgridtraits.hh:125
Definition: consistentlyorientedgrid.hh:23
static constexpr bool supportsPeriodicity()
Definition: periodicgridtraits.hh:129
Definition: common/pdesolver.hh:24
Definition: periodicgridtraits.hh:34
Definition: periodicgridtraits.hh:33
bool isPeriodic(const typename Grid::LeafIntersection &intersection) const
Definition: periodicgridtraits.hh:38
PeriodicGridTraits(const Grid &grid)
Definition: periodicgridtraits.hh:36