1#ifndef SPACE_P1BUBBLE_HH 
    2#define SPACE_P1BUBBLE_HH 
    6#include <dune/geometry/referenceelements.hh> 
    7#include <dune/grid/common/gridenums.hh> 
    9#include <dune/fem/common/hybrid.hh> 
   10#include <dune/fem/space/common/functionspace.hh> 
   11#include <dune/fem/space/common/defaultcommhandler.hh> 
   12#include <dune/fem/space/common/discretefunctionspace.hh> 
   13#include <dune/fem/space/common/localrestrictprolong.hh> 
   14#include <dune/fem/space/mapper/code.hh> 
   15#include <dune/fem/space/mapper/localkey.hh> 
   16#include <dune/fem/space/mapper/indexsetdofmapper.hh> 
   17#include <dune/fem/space/shapefunctionset/simple.hh> 
   18#include <dune/fem/space/shapefunctionset/selectcaching.hh> 
   19#include <dune/fem/space/shapefunctionset/vectorial.hh> 
   20#include <dune/fem/space/basisfunctionset/default.hh> 
   22#include <dune/fem/operator/matrix/istlmatrixadapter.hh> 
   30    template< 
class FunctionSpace >
 
   31    struct LocalBubbleElementInterpolation
 
   38      LocalBubbleElementInterpolation ()
 
   39        : points_( dimDomain + 2, DomainType( 0.0 ) )
 
   41        for( 
int i = 0; i < dimDomain; ++i )
 
   42          points_[ i + 1 ][ i ] = 1.0;
 
   44        points_[ dimDomain +1 ] = DomainType( 1.0 / ( dimDomain + 1.0 ) );
 
   47      LocalBubbleElementInterpolation ( 
const LocalBubbleElementInterpolation & ) = 
default;
 
   48      LocalBubbleElementInterpolation ( LocalBubbleElementInterpolation && ) = 
default;
 
   51      template< 
class LocalFunction, 
class LocalDofVector >
 
   52      void operator() ( 
const LocalFunction &lf, LocalDofVector &ldv ) 
const 
   56        for( 
const DomainType &x : points_ )
 
   59          lf.evaluate( x, phi );
 
   60          for( 
int i = 0; i < dimRange; ++i )
 
   61            ldv[ k++ ] = phi[ i ];
 
   65      template <
class Entity>
 
   66      void bind( 
const Entity & ) {}
 
   70      std::vector< DomainType > points_;
 
   75    struct BubbleElementLocalKeyMap
 
   78      BubbleElementLocalKeyMap ( 
int vertices )
 
   80        for( 
int i = 0; i < vertices; ++i )
 
   81          map_.emplace_back( i, dim, 0 );
 
   82        map_.emplace_back( 0, 0, 0 );
 
   86      std::size_t 
size()
 const { 
return map_.size(); }
 
   88      LocalKey& localKey ( std::size_t i ) { 
return map_[ i ]; }
 
   89      const LocalKey& localKey ( std::size_t i )
 const { 
return map_[ i ]; }
 
   92      std::vector< LocalKey > map_;
 
   95    struct BubbleElementDofMapperCodeFactory
 
   99      template< 
class RefElement,
 
  100                std::enable_if_t< std::is_same< std::decay_t< decltype( std::declval< const RefElement & >().size( 0 ) ) >, 
int >::value, 
int > = 0,
 
  101                std::enable_if_t< std::is_same< std::decay_t< 
decltype( std::declval< const RefElement & >().type( 0, 0 ) ) >, GeometryType >::value, 
int > = 0 >
 
  102      DofMapperCode 
operator() ( 
const RefElement &refElement ) 
const 
  104        static const int dim = RefElement::dimension;
 
  105        if( refElement.type().isSimplex() )
 
  106          return compile( refElement, BubbleElementLocalKeyMap< dim >(dim+1) );
 
  107        if( refElement.type().isCube() && refElement.type().dim() == 2)
 
  108          return compile( refElement, BubbleElementLocalKeyMap< dim >(pow(2,dim)) );
 
  110          return DofMapperCode();
 
  117    template< 
class FunctionSpace >
 
  118    class SimplexBubbleElementShapeFunctionSet
 
  120      typedef SimplexBubbleElementShapeFunctionSet< FunctionSpace > ThisType;
 
  123      static const int polynomialOrder = dimDomain + 1;
 
  124      static const int numShapeFunctions = dimDomain + 2;
 
  129      static_assert( RangeType::dimension == 1, 
"This is a scalar shapefunction set" );
 
  133      SimplexBubbleElementShapeFunctionSet () {}
 
  135      template<
class GeometryType >
 
  136      SimplexBubbleElementShapeFunctionSet ( 
const GeometryType& 
gt )
 
  138        if( !
gt.isSimplex() )
 
  139          DUNE_THROW( NotImplemented, 
"Wrong geometry type for Cube2DBubblelementShapeFunctionSet." );
 
  142      template< 
class Po
int, 
class Functor >
 
  143      void evaluateEach ( 
const Point &x, Functor functor )
 const 
  145        DomainType xRef = coordinate( x );
 
  146        RangeType phi(1), phi0(1);
 
  147        for( 
int i=0; i< dimDomain; ++i )
 
  149          functor( i+1, RangeType( xRef[ i ] ) );
 
  150          phi0[ 0 ] -= xRef[ i ];
 
  151          phi[ 0 ] *= xRef[ i ] ;
 
  154        phi[ 0 ] *= phi0[ 0 ] / std::pow( ( dimDomain + 1.0 ), dimDomain + 1.0 );
 
  156        functor( dimDomain +1, phi );
 
  159      template< 
class Po
int, 
class Functor >
 
  160      void jacobianEach ( 
const Point &x, Functor functor )
 const 
  162        DomainType xRef = coordinate( x );
 
  164        JacobianRangeType jac(0), jac0( -1 );
 
  169        for( 
int i=0; i< dimDomain; ++i )
 
  171          phi0[ 0 ] -= xRef[ i ];
 
  173          for( 
int j=1; j < dimDomain; ++j )
 
  174            jac0[ 0 ][ (i+j)%dimDomain ] *= xRef[ i ];
 
  181        for( 
int i=0; i< dimDomain; ++i )
 
  182          jac0[ 0 ][ i ] *= -(phi0[ 0 ] - xRef[ i ]);
 
  183        jac0[ 0 ] *= 1.0 / std::pow( dimDomain + 1.0, dimDomain + 1.0 );
 
  184        functor( dimDomain +1, jac0 );
 
  187      template< 
class Po
int, 
class Functor >
 
  188      void hessianEach ( 
const Point &x, Functor functor )
 const 
  190        DUNE_THROW( NotImplemented, 
"NotImplemented" );
 
  191        DomainType xRef = coordinate( x );
 
  192        HessianRangeType hes;
 
  199      int order ()
 const { 
return dimDomain + 1; }
 
  201      std::size_t 
size ()
 const { 
return dimDomain +2; }
 
  206    template< 
class FunctionSpace >
 
  207    class Cube2DBubbleElementShapeFunctionSet
 
  209      typedef Cube2DBubbleElementShapeFunctionSet< FunctionSpace > ThisType;
 
  212      static const int polynomialOrder = dimDomain + 1;
 
  213      static const int numShapeFunctions = dimDomain + 2;
 
  218      static_assert( RangeType::dimension == 1, 
"This is a scalar shapefunction set" );
 
  223      Cube2DBubbleElementShapeFunctionSet () {}
 
  225      template<
class GeometryType >
 
  226      Cube2DBubbleElementShapeFunctionSet ( 
const GeometryType& 
gt )
 
  229          DUNE_THROW( NotImplemented, 
"Wrong geometry type for Cube2DBubblelementShapeFunctionSet." );
 
  231          DUNE_THROW( NotImplemented, 
"2DCubeBubbleElementShapeFunctionSet only implemented for dimension = 2." );
 
  234      template< 
class Po
int, 
class Functor >
 
  235      void evaluateEach ( 
const Point &x, Functor functor ) 
const 
  238        DomainType xRef = coordinate( x );
 
  239        RangeType phi(1), phi0(1);
 
  241        phi0[0] = (1.-xRef[0])*(1.-xRef[1]);
 
  244        phi[0] = xRef[0]*(1.-xRef[1]);
 
  247        phi[0] = (1.-xRef[0])*xRef[1];
 
  250        phi[0] = xRef[0]*xRef[1];
 
  254        phi[0] = 64.*phi0[0]*phi0[0]*xRef[0]*xRef[1];
 
  258      template< 
class Po
int, 
class Functor >
 
  259      void jacobianEach ( 
const Point &x, Functor functor )
 const 
  261        DomainType xRef = coordinate( x );
 
  263        JacobianRangeType jac, jac0;
 
  265        phi0[0] = (1.-xRef[0])*(1.-xRef[1]);
 
  266        jac0[0] = {xRef[1]-1,xRef[0]-1};
 
  269        jac[0] = {1-xRef[1],xRef[0]};
 
  271        jac[0] = {xRef[1],1-xRef[0]};
 
  273        jac[0] = {xRef[1],xRef[0]};
 
  277        jac0[0] *= 128.*phi0*xRef[0]*xRef[1];
 
  278        jac[0] = {xRef[1],xRef[0]};
 
  279        jac[0] *= 64.*phi0*phi0;
 
  284      template< 
class Po
int, 
class Functor >
 
  285      void hessianEach ( 
const Point &x, Functor functor )
 const 
  287        DUNE_THROW( NotImplemented, 
"NotImplemented" );
 
  288        DomainType xRef = coordinate( x );
 
  289        HessianRangeType hes;
 
  296      int order ()
 const { 
return 4; }
 
  298      std::size_t 
size ()
 const { 
return 5; }
 
  301    template< 
class FunctionSpace, 
class Gr
idPart, 
class Storage >
 
  302    class BubbleElementSpace;
 
  307    template< 
class FunctionSpace, 
class Gr
idPart, 
class Storage >
 
  308    struct BubbleElementSpaceTraits
 
  310      typedef BubbleElementSpace< FunctionSpace, GridPart, Storage > DiscreteFunctionSpaceType;
 
  313      typedef GridPart GridPartType;
 
  315      static const int codimension = 0;
 
  318      typedef typename GridPartType::template Codim< codimension >::EntityType EntityType;
 
  328                     "bubble elements only implemented for grids with a single geometry type" );
 
  329      static const unsigned int topologyId =
 
  333      typedef SimplexBubbleElementShapeFunctionSet< ScalarFunctionSpaceType > ScalarSimplexShapeFunctionSetType;
 
  334      typedef Cube2DBubbleElementShapeFunctionSet< ScalarFunctionSpaceType > ScalarCubeShapeFunctionSetType;
 
  335      typedef typename std::conditional< topologyId == 0,
 
  336              ScalarSimplexShapeFunctionSetType, ScalarCubeShapeFunctionSetType >::type
 
  337                ScalarShapeFunctionSetType;
 
  339      typedef VectorialShapeFunctionSet< ScalarShapeFunctionSetType, typename FunctionSpaceType::RangeType > ShapeFunctionSetType;
 
  341      typedef DefaultBasisFunctionSet< EntityType, ShapeFunctionSetType > BasisFunctionSetType;
 
  345      static const int polynomialOrder = ScalarShapeFunctionSetType::polynomialOrder;
 
  347      typedef IndexSetDofMapper< GridPartType > BlockMapperType;
 
  348      typedef Hybrid::IndexRange< int, FunctionSpaceType::dimRange > LocalBlockIndices;
 
  350      template <
class DiscreteFunction, 
class Operation = DFCommunicationOperation::Add >
 
  351      struct CommDataHandle
 
  353        typedef Operation OperationType;
 
  354        typedef DefaultCommunicationHandler< DiscreteFunction, Operation > Type;
 
  362    template< 
class FunctionSpace, 
class Gr
idPart, 
class Storage = CachingStorage >
 
  371      typedef BubbleElementSpaceTraits< FunctionSpace, GridPart, Storage > Traits;
 
  372      typedef typename Traits::ShapeFunctionSetType ShapeFunctionSetType;
 
  373      static const int polynomialOrder = Traits::polynomialOrder;
 
  375      typedef typename BaseType::GridPartType GridPartType;
 
  379      typedef typename BaseType::EntityType EntityType;
 
  381      typedef typename BaseType::BlockMapperType BlockMapperType;
 
  384      typedef LocalBubbleElementInterpolation< FunctionSpace > InterpolationType;
 
  387      static const InterfaceType defaultInterface = GridPartType::indexSetInterfaceType;
 
  394        blockMapper_(  
gridPart, BubbleElementDofMapperCodeFactory() )
 
  402      bool contains ( 
const int codim )
 const 
  408      bool continuous ()
 const { 
return true; }
 
  411      bool continuous ( 
const IntersectionType & intersection )
 const 
  420        return polynomialOrder;
 
  424      template<
class Entity>
 
  427        return polynomialOrder;
 
  431      template< 
class EntityType >
 
  434        return BasisFunctionSetType( entity, ShapeFunctionSetType( entity.geometry().type() ) );
 
  448        return InterpolationType();
 
  463      mutable BlockMapperType blockMapper_;
 
  466    template< 
class FunctionSpace,
 
  469              class NewFunctionSpace >
 
  470    struct DifferentDiscreteFunctionSpace< BubbleElementSpace < FunctionSpace, GridPart, Storage >, NewFunctionSpace >
 
  472      typedef BubbleElementSpace< NewFunctionSpace, GridPart, Storage > Type;
 
  475    template< 
class FunctionSpace, 
class Gr
idPart, 
class Storage >
 
  476    class DefaultLocalRestrictProlong < BubbleElementSpace< FunctionSpace, GridPart, Storage > >
 
  477    : 
public EmptyLocalRestrictProlong< BubbleElementSpace< FunctionSpace, GridPart, Storage > >
 
  479      typedef EmptyLocalRestrictProlong< BubbleElementSpace< FunctionSpace, GridPart, Storage > > BaseType;
 
  481      DefaultLocalRestrictProlong( 
const BubbleElementSpace< FunctionSpace, GridPart, Storage > &space )
 
Wrapper class for entities.
Definition: entity.hh:66
 
[Class definition for new space]
Definition: p1bubble.hh:366
 
InterpolationType interpolation() const
return local interpolation object for LocalInterpolation
Definition: p1bubble.hh:446
 
BlockMapperType & blockMapper() const
obtain the DoF block mapper of this space
Definition: p1bubble.hh:440
 
bool continuous(const IntersectionType &intersection) const
returns true if the space contains only globally continuous functions
Definition: p1bubble.hh:411
 
int order() const
get global order of space
Definition: p1bubble.hh:418
 
BasisFunctionSetType basisFunctionSet(const EntityType &entity) const
get basis function set for given entity
Definition: p1bubble.hh:432
 
int order(const Entity &entity) const
get global order of space
Definition: p1bubble.hh:425
 
InterpolationType interpolation(const EntityType &entity) const
return local interpolation for given entity
Definition: p1bubble.hh:457
 
This is the class with default implementations for discrete function. The methods not marked with hav...
Definition: discretefunctionspace.hh:649
 
GridPartType & gridPart() const
Definition: discretefunctionspace.hh:766
 
Traits::BasisFunctionSetType BasisFunctionSetType
type of basis function set of this space
Definition: discretefunctionspace.hh:201
 
GridPartType::IntersectionType IntersectionType
type of the intersections
Definition: discretefunctionspace.hh:226
 
ExplicitFieldVector< FieldMatrix< RangeFieldType, dimDomain, dimDomain >, dimRange > HessianRangeType
Intrinsic type used for the hessian values has a Dune::FieldMatrix type interface.
Definition: functionspaceinterface.hh:79
 
FunctionSpaceTraits::RangeType RangeType
Type of range vector (using type of range field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:71
 
@ dimDomain
dimension of domain vector space
Definition: functionspaceinterface.hh:46
 
@ dimRange
dimension of range vector space
Definition: functionspaceinterface.hh:48
 
FunctionSpaceTraits::LinearMappingType JacobianRangeType
Intrinsic type used for the jacobian values has a Dune::FieldMatrix type interface.
Definition: functionspaceinterface.hh:75
 
FunctionSpaceTraits::DomainType DomainType
Type of domain vector (using type of domain field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:67
 
FunctionSpaceTraits::ScalarFunctionSpaceType ScalarFunctionSpaceType
corresponding scalar function space
Definition: functionspaceinterface.hh:83
 
A vector valued function space.
Definition: functionspace.hh:60
 
A few common exception classes.
 
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
 
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:158
 
Dune namespace.
Definition: alignedallocator.hh:13
 
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75
 
Specialize with 'true' for if the codimension 0 entity of the grid has only one possible geometry typ...
Definition: capabilities.hh:27