1#ifndef DUNE_FEM_VTXPROJECTION_HH 
    2#define DUNE_FEM_VTXPROJECTION_HH 
   11#include <dune/grid/common/rangegenerators.hh> 
   13#include <dune/fem/common/coordinate.hh> 
   14#include <dune/fem/operator/common/operator.hh> 
   15#include <dune/fem/space/common/communicationmanager.hh> 
   16#include <dune/fem/space/common/interpolate.hh> 
   17#include <dune/fem/space/lagrange.hh> 
   24    template <
class Gr
idPartType>
 
   25    struct WeightDefault {
 
   26      typedef typename GridPartType::template Codim<0>::EntityType EntityType;
 
   28      typedef typename EntityType::Geometry::GlobalCoordinate WorldDomainType;
 
   29      typedef typename EntityType::Geometry::LocalCoordinate  DomainType;
 
   33      void setEntity(
const EntityType& en) {
 
   35        volume_ = en.geometry().volume();
 
   36        baryCenter_ = en.geometry().center();
 
   38      double operator()(
const DomainType& point) {
 
   40        auto tau = en_.geometry().global(point);
 
   42        return tau.two_norm() / volume_;
 
   45      WorldDomainType baryCenter_;
 
   49    struct VtxProjectionImpl
 
   51      template< 
class DiscreteFunction, 
class Intersection >
 
   52      struct OutsideLocalFunction
 
   54        typedef typename DiscreteFunction::LocalFunctionType LocalFunctionType;
 
   56        typedef typename LocalFunctionType::EntityType EntityType;
 
   68        typedef typename EntityType::Geometry::LocalCoordinate LocalCoordinateType;
 
   75        explicit OutsideLocalFunction ( 
const DiscreteFunction &df ) : order_(df.order()), localFunction_( df ) {}
 
   77        void init ( 
const EntityType &outside, 
const GeometryType &geoIn, 
const GeometryType &geoOut )
 
   79          localFunction_.init( outside );
 
   85        const EntityType &entity()
 const { 
return entity_; }
 
   86        bool valid ()
 const { 
return geoIn_ != 
nullptr; }
 
   88        int order()
 const { 
return order_; }
 
   90        template< 
class Po
int >
 
   91        void evaluate ( 
const Point &x, RangeType &value )
 const 
   93          localFunction_.evaluate( geoOut_->global( geoIn_->local( coordinate( x ) ) ), value );
 
   96        template< 
class Po
int >
 
   97        void jacobian ( 
const Point &x, JacobianRangeType &jacobian )
 const 
   99          DUNE_THROW( NotImplemented, 
"Vertex projection weights do not provide Jacobians." );
 
  102        template< 
class Po
int >
 
  103        void hessian ( 
const Point &x, HessianRangeType &hessian )
 const 
  105          DUNE_THROW( NotImplemented, 
"Vertex projection weights do not provide Hessians." );
 
  108        template< 
class Quadrature, 
class Values >
 
  109        void evaluateQuadrature ( 
const Quadrature &quadrature, Values &values )
 const 
  111          for( 
const auto &qp : quadrature )
 
  112            doEvaluate( qp, values[ qp.index() ] );
 
  116        template< 
class Po
int >
 
  117        void doEvaluate ( 
const Point &x, RangeType &value )
 const 
  119          evaluate( x, value );
 
  122        template< 
class Po
int >
 
  123        void doEvaluate ( 
const Point &x, JacobianRangeType &value )
 const 
  125          jacobian( x, value );
 
  128        template< 
class Po
int >
 
  129        void doEvaluate ( 
const Point &x, HessianRangeType &value )
 const 
  135        LocalFunctionType localFunction_;
 
  142      template< 
class DiscreteFunction >
 
  143      static void makeConforming ( DiscreteFunction &u )
 
  145        typedef typename DiscreteFunction::GridPartType GridPartType;
 
  148        typedef typename GridPartType::template Codim< 0 >::EntityType EntityType;
 
  149        typedef typename DiscreteFunction::DiscreteFunctionSpaceType
 
  150          DiscreteFunctionSpaceType;
 
  152        if( u.space().gridPart().isConforming() || !Fem::GridPartCapabilities::hasGrid< GridPartType >::v )
 
  155        const auto &blockMapper = u.space().blockMapper();
 
  157        std::vector< typename DiscreteFunction::DofType > ldu, ldw;
 
  158        const std::size_t localBlockSize = DiscreteFunctionSpaceType::localBlockSize;
 
  159        const std::size_t maxNumBlocks = blockMapper.maxNumDofs();
 
  160        ldu.reserve( maxNumBlocks * localBlockSize );
 
  161        ldw.reserve( maxNumBlocks * localBlockSize );
 
  163        std::vector< char > onSubEntity;
 
  164        onSubEntity.reserve( maxNumBlocks );
 
  166        OutsideLocalFunction< DiscreteFunction, IntersectionType > uOutside( u );
 
  168        LocalInterpolation< DiscreteFunctionSpaceType > interpolation( u.space() );
 
  170        for( 
const EntityType &inside : u.space() )
 
  172          for( 
const IntersectionType &intersection : intersections( u.gridPart(), inside ) )
 
  174            if( intersection.conforming() || !intersection.neighbor() )
 
  178            const EntityType outside = intersection.outside();
 
  179            if( inside.level() <= outside.level() )
 
  184            const auto geoIn = intersection.geometryInInside();
 
  185            const auto geoOut = intersection.geometryInOutside();
 
  186            uOutside.init( outside, geoIn, geoOut );
 
  189            const std::size_t numBlocks = blockMapper.numDofs( inside );
 
  190            ldu.resize( numBlocks * localBlockSize );
 
  191            ldw.resize( numBlocks * localBlockSize );
 
  194            auto guard = bindGuard( interpolation, inside );
 
  195            interpolation( uOutside, ldw );
 
  198            u.getLocalDofs( inside, ldu );
 
  201            blockMapper.onSubEntity( inside, intersection.indexInInside(), 1, onSubEntity );
 
  202            for( std::size_t i = 0; i < numBlocks; ++i )
 
  204              if( !onSubEntity[ i ] )
 
  206              for( std::size_t j = 0; j < localBlockSize; ++j )
 
  207                ldu[ i*localBlockSize + j ] = ldw[ i*localBlockSize + j ];
 
  211            u.setLocalDofs( inside, ldu );
 
  216      template< 
class Function, 
class DiscreteFunction, 
class Weight >
 
  217      static auto project ( 
const Function &f, DiscreteFunction &u, Weight &weight )
 
  218      -> std::enable_if_t<std::is_void< 
decltype( 
interpolate(f,u,weight,u )) >::value>
 
  221        DiscreteFunction w ( 
"weight", u.space() );
 
  225      template< 
class Function, 
class DiscreteFunction >
 
  226      static auto project ( 
const Function &f, DiscreteFunction &u )
 
  227      -> std::enable_if_t< std::is_void< 
decltype( project( f,u,std::declval<WeightDefault<typename DiscreteFunction::GridPartType>&>() ) ) >::value>
 
  229        WeightDefault<typename DiscreteFunction::GridPartType> weight;
 
  230        project(f, u, weight);
 
  242    template < 
typename DType, 
typename RType >
 
  246      typedef DType DomainType;
 
  247      typedef RType RangeType;
 
  248      typedef typename DomainType :: RangeFieldType DomainFieldType;
 
  249      typedef typename RType :: RangeFieldType  RangeFieldType;
 
  250      typedef typename Dune::FieldTraits< RangeFieldType >::real_type RealType;
 
  251      typedef typename RType :: DiscreteFunctionSpaceType :: GridPartType GridPartType;
 
  256      template <
class WeightType>
 
  258                       WeightType& weight)
 const 
  260        VtxProjectionImpl::project(f,discFunc,weight);
 
  263      void operator() (
const DomainType& f, RangeType& discFunc)
 const 
  265        VtxProjectionImpl::project(f,discFunc);
 
IntersectionIteratorType::Intersection IntersectionType
type of intersection
Definition: adaptiveleafgridpart.hh:95
 
Definition: explicitfieldvector.hh:75
 
FunctionSpaceTraits::DomainFieldType DomainFieldType
Intrinsic type used for values in the domain field (usually a double)
Definition: functionspaceinterface.hh:60
 
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::RangeFieldType RangeFieldType
Intrinsic type used for values in the range field (usually a double)
Definition: functionspaceinterface.hh:63
 
A vector valued function space.
Definition: functionspace.hh:60
 
The Projection class which average discontinuous function using the interpolation of the space (e....
Definition: vtxprojection.hh:244
 
VtxProjection()
Constructor.
Definition: vtxprojection.hh:253
 
void operator()(const DomainType &f, RangeType &discFunc, WeightType &weight) const
apply projection
Definition: vtxprojection.hh:257
 
GridImp::template Codim< 1 >::LocalGeometry LocalGeometry
Codim 1 geometry returned by geometryInInside() and geometryInOutside()
Definition: intersection.hh:204
 
actual interface class for quadratures
 
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
 
static void interpolate(const GridFunction &u, DiscreteFunction &v)
perform native interpolation of a discrete function space
Definition: interpolate.hh:55
 
#define DUNE_THROW(E,...)
Definition: exceptions.hh:314
 
Dune namespace.
Definition: alignedallocator.hh:13
 
abstract operator
Definition: operator.hh:34