1#ifndef DUNE_FEM_MISC_BOUNDARYIDPROVIDER_HH 
    2#define DUNE_FEM_MISC_BOUNDARYIDPROVIDER_HH 
    6#include <dune/fem/misc/griddeclaration.hh> 
    7#include <dune/fem/function/common/localcontribution.hh> 
   15  template< 
class Gr
id >
 
   16  struct HostGridAccess;
 
   26    template< 
class Gr
id >
 
   27    struct BoundaryIdProvider;
 
   35    template< 
int dim, 
int dimW >
 
   36    struct BoundaryIdProvider< AlbertaGrid< dim, dimW > >
 
   38      typedef AlbertaGrid< dim, dimW > GridType;
 
   40      template< 
class Intersection >
 
   41      static int boundaryId ( 
const Intersection &intersection )
 
   43        return intersection.impl().boundaryId();
 
   54    template< 
int dim, 
int dimw, ALUGr
idElementType elType, ALUGr
idRefinementType refineType, 
class Comm >
 
   55    struct BoundaryIdProvider< ALUGrid< dim, dimw, elType, refineType, Comm > >
 
   57      typedef ALUGrid< dim, dimw, elType, refineType, Comm > GridType;
 
   59      template< 
class Intersection >
 
   60      static int boundaryId ( 
const Intersection &intersection )
 
   62        return intersection.impl().boundaryId();
 
   73    template< 
class HostGr
id >
 
   74    struct BoundaryIdProvider< CacheItGrid< HostGrid > >
 
   76      typedef CacheItGrid< HostGrid > GridType;
 
   78      template< 
class Intersection >
 
   79      static int boundaryId ( 
const Intersection &intersection )
 
   81        return BoundaryIdProvider< HostGrid >
 
   82          ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
 
   93    template< 
class HostGr
id >
 
   94    struct BoundaryIdProvider< CartesianGrid< HostGrid > >
 
   96      typedef CartesianGrid< HostGrid > GridType;
 
   98      template< 
class Intersection >
 
   99      static int boundaryId ( 
const Intersection &intersection )
 
  101        return BoundaryIdProvider< HostGrid >
 
  102          ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
 
  112#if HAVE_DUNE_METAGRID 
  113    template< 
class HostGr
id >
 
  114    struct BoundaryIdProvider< FilteredGrid< HostGrid > >
 
  116      typedef FilteredGrid< HostGrid > GridType;
 
  120      template< 
class Intersection >
 
  121      static int boundaryId ( 
const Intersection &intersection )
 
  123        if( !HostGridAccess< GridType >::hostIntersection( intersection ).boundary() )
 
  124          DUNE_THROW( NotImplemented, 
"BoundaryIdProvider for artificial boundaries of FilteredGrid not implemented." );
 
  125        return BoundaryIdProvider< HostGrid >
 
  126          ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
 
  136    template< 
class HostGr
id, 
class CoordFunction, 
class Allocator >
 
  137    struct BoundaryIdProvider< GeometryGrid< HostGrid, CoordFunction, Allocator > >
 
  139      typedef GeometryGrid< HostGrid, CoordFunction, Allocator > GridType;
 
  141      template< 
class Intersection >
 
  142      static int boundaryId ( 
const Intersection &intersection )
 
  144        return BoundaryIdProvider< HostGrid >
 
  145          ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
 
  154#if HAVE_DUNE_METAGRID 
  155    template< 
class HostGr
id >
 
  156    struct BoundaryIdProvider< IdGrid< HostGrid > >
 
  158      typedef IdGrid< HostGrid > GridType;
 
  160      template< 
class Intersection >
 
  161      static int boundaryId ( 
const Intersection &intersection )
 
  163        return BoundaryIdProvider< HostGrid >
 
  164          ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
 
  175    struct BoundaryIdProvider< OneDGrid >
 
  177      typedef OneDGrid GridType;
 
  179      template< 
class Intersection >
 
  180      static int boundaryId ( 
const Intersection &intersection )
 
  182        return intersection.boundarySegmentIndex();
 
  191#if HAVE_DUNE_METAGRID 
  192    template< 
class HostGr
id >
 
  193    struct BoundaryIdProvider< ParallelGrid< HostGrid > >
 
  195      typedef ParallelGrid< HostGrid > GridType;
 
  197      template< 
class Intersection >
 
  198      static int boundaryId ( 
const Intersection &intersection )
 
  200        return BoundaryIdProvider< HostGrid >
 
  201          ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
 
  211#if HAVE_DUNE_METAGRID 
  212    template< 
class HostGr
id, 
class MapToSphere >
 
  213    struct BoundaryIdProvider< SphereGrid< HostGrid, MapToSphere > >
 
  215      typedef SphereGrid< HostGrid, MapToSphere > GridType;
 
  217      template< 
class Intersection >
 
  218      static int boundaryId ( 
const Intersection &intersection )
 
  220        return BoundaryIdProvider< HostGrid >
 
  221          ::boundaryId ( HostGridAccess< GridType >::hostIntersection( intersection ) );
 
  232    template< 
class ct, 
int dim, 
template< 
int > 
class Strategy, 
class Comm >
 
  233    struct BoundaryIdProvider< SPGrid< ct, dim, Strategy, Comm > >
 
  235      typedef SPGrid< ct, dim, Strategy, Comm > GridType;
 
  237      template< 
class Intersection >
 
  238      static int boundaryId ( 
const Intersection &intersection )
 
  240        return (intersection.boundary() ? (intersection.indexInInside()+1) : 0);
 
  249#if HAVE_DUNE_POLYGONGRID 
  251    struct BoundaryIdProvider< PolygonGrid< ct > >
 
  253      typedef PolygonGrid< ct > GridType;
 
  255      template< 
class Intersection >
 
  256      static int boundaryId ( 
const Intersection &intersection )
 
  258        return (intersection.boundary() ? (intersection.impl().boundaryId()) : 0);
 
  266#if HAVE_DUNE_P4ESTGRID 
  267    template< 
int dim, 
int dimworld, P4estType elType, 
class ct >
 
  268    struct BoundaryIdProvider< P4estGrid<dim, dimworld, elType, ct > >
 
  270      typedef P4estGrid<dim, dimworld, elType, ct > GridType;
 
  272      template< 
class Intersection >
 
  273      static int boundaryId ( 
const Intersection &intersection )
 
  275        return (intersection.boundary() ? (intersection.impl().boundaryId()) : 0);
 
  284    template< 
int dim, 
int dimworld, 
class ct >
 
  285    struct BoundaryIdProvider< PolyhedralGrid< dim, dimworld, ct > >
 
  287      typedef PolyhedralGrid< dim, dimworld, ct > GridType;
 
  289      template< 
class Intersection >
 
  290      static int boundaryId ( 
const Intersection &intersection )
 
  292        return (intersection.boundary() ? (intersection.impl().boundaryId()) : 0);
 
  302    struct BoundaryIdProvider< UGGrid< dim > >
 
  304      typedef UGGrid< dim > GridType;
 
  306      template< 
class Intersection >
 
  307      static int boundaryId ( 
const Intersection &intersection )
 
  309        return intersection.boundarySegmentIndex();
 
  318    template< 
int dim, 
class CoordCont >
 
  319    struct BoundaryIdProvider< YaspGrid< dim, CoordCont > >
 
  323      template< 
class Intersection >
 
  324      static int boundaryId ( 
const Intersection &intersection )
 
  326        return (intersection.boundary() ? (intersection.indexInInside()+1) : 0);
 
  336    template< 
class Gr
idView, 
class Intersection>
 
  337    inline static int boundaryId ( 
const Intersection &intersection )
 
  339      return Dune::Fem::BoundaryIdProvider< typename GridView::Grid > ::
 
  340             boundaryId( intersection );
 
  344    template< 
class Gr
idView, 
class Intersection>
 
  345    inline static int boundaryId ( 
const GridView&, 
const Intersection &intersection )
 
  347      return Dune::Fem::BoundaryIdProvider< typename GridView::Grid > ::
 
  348             boundaryId( intersection );
 
  356    template <
class DiscreteFunction>
 
  357    void projectBoundaryIds( DiscreteFunction& df )
 
  359      if( df.space().order() != 0 ) 
 
  361        DUNE_THROW(InvalidStateException,
"projectBoundaryIds: expect piecewise constant discrete function");
 
  366      const auto& gridPart = df.space().gridPart();
 
  368      typedef AddLocalContribution< DiscreteFunction > AddLocalContributionType;
 
  369      AddLocalContributionType localDf( df );
 
  370      for( 
const auto& element : df.space() )
 
  373        auto guard = bindGuard( localDf, element );
 
  374        for( 
const auto& intersection : intersections( gridPart, element ) )
 
  376          if( intersection.boundary() )
 
  378            localDf[ 0 ] += boundaryId( gridPart, intersection );
 
  383          localDf[ 0 ] /= double(count);
 
A few common exception classes.
 
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
 
Dune namespace.
Definition: alignedallocator.hh:13