1#ifndef DUNE_FEMPY_FUNCTION_GRIDFUNCTIONVIEW_HH 
    2#define DUNE_FEMPY_FUNCTION_GRIDFUNCTIONVIEW_HH 
    4#include <dune/fem/function/common/discretefunction.hh> 
    5#include <dune/fem/function/localfunction/bindable.hh> 
    6#include <dune/fem/function/localfunction/const.hh> 
    7#include <dune/fem/common/intersectionside.hh> 
   17    struct IsGridFunctionView {};
 
   19    template< class GF, bool isDiscreteFunction = std::is_base_of< Fem::IsDiscreteFunction, GF >::value >
 
   20    struct GridFunctionView;
 
   23    struct GridFunctionView< GF, false >
 
   24    : 
public BindableGridFunction<typename GF::GridPartType, typename GF::RangeType>,
 
   25             public IsGridFunctionView
 
   27      using Base = BindableGridFunction<typename GF::GridPartType, typename GF::RangeType>;
 
   28      typedef typename GF::EntityType Entity;
 
   29      typedef typename GF::RangeType Value;
 
   31      typedef typename Entity::Geometry::LocalCoordinate LocalCoordinate;
 
   33      GridFunctionView ( 
const GF &gf )
 
   35      , localFunction_( gf ) {}
 
   37      Value operator() ( 
const LocalCoordinate &x )
 const 
   40        localFunction_.evaluate( x, value );
 
   44      void bind ( 
const Entity &entity ) { localFunction_.bind( entity ); }
 
   46      template <
class IntersectionType>
 
   47      void bind(
const IntersectionType &intersection, IntersectionSide side)
 
   48      { localFunction_.bind(intersection, side); }
 
   51      Dune::Fem::ConstLocalFunction<GF> localFunction_;
 
   55    struct GridFunctionView< GF, true >
 
   56    : 
public BindableGridFunctionWithSpace<typename GF::GridPartType, typename GF::RangeType>,
 
   57             public IsGridFunctionView
 
   59      using Base = BindableGridFunctionWithSpace<typename GF::GridPartType, typename GF::RangeType>;
 
   60      typedef typename GF::EntityType Entity;
 
   61      typedef typename GF::RangeType Value;
 
   63      typedef typename Entity::Geometry::LocalCoordinate LocalCoordinate;
 
   67      typedef typename DiscreteFunctionSpace::BasisFunctionSetType BasisFunctionSet;
 
   70      GridFunctionView ( 
const GF &gf )
 
   71        : Base(gf.gridPart(), gf.name(), gf.order())
 
   74        localDofVector_.reserve( DiscreteFunctionSpace::localBlockSize * space().blockMapper().maxNumDofs() );
 
   77      Value operator() ( 
const LocalCoordinate &x )
 const 
   80        basisFunctionSet_.evaluateAll( x, localDofVector_, value );
 
   84      void bind ( 
const Entity &entity )
 
   86        basisFunctionSet_ = space().basisFunctionSet( entity );
 
   87        localDofVector_.resize( basisFunctionSet_.size() );
 
   88        gf_.getLocalDofs( entity, localDofVector_ );
 
   93        basisFunctionSet_ = BasisFunctionSet();
 
   94        localDofVector_.resize( 0u );
 
   97      template <
class IntersectionType>
 
   98      void bind(
const IntersectionType &intersection, IntersectionSide side)
 
   99      { defaultIntersectionBind(gf_, intersection, side); }
 
  105      BasisFunctionSet basisFunctionSet_;
 
  115      std::enable_if_t< !std::is_base_of< Fem::IsGridFunctionView, GF >::value &&
 
  116                         std::is_base_of< Fem::HasLocalFunction, GF >::value, 
int > = 0
 
  118    inline static GridFunctionView< GF > localFunction ( 
const GF &gf )
 
  120      return GridFunctionView< GF >( gf );
 
Dune namespace.
Definition: alignedallocator.hh:13