1#ifndef DUNE_FEM_LPNORM_HH 
    2#define DUNE_FEM_LPNORM_HH 
    4#include <dune/fem/quadrature/integrator.hh> 
    6#include "domainintegral.hh" 
   25      virtual int operator() (
const double p)=0;
 
   34      int operator() (
const double p)
 
   37        const double q = p / (p-1);
 
   38        double max = std::max(p,q);
 
   45    template< 
class Gr
idPart, 
class OrderCalculator = DefaultOrderCalculator >
 
   46    class LPNorm : 
public IntegralBase < GridPart, LPNorm< GridPart, OrderCalculator > >
 
   48      typedef LPNorm< GridPart, OrderCalculator > ThisType;
 
   49      typedef IntegralBase< GridPart, ThisType > BaseType ;
 
   52      typedef GridPart GridPartType;
 
   54      using BaseType::gridPart;
 
   58      template< 
class Function >
 
   59      struct FunctionMultiplicator;
 
   61      template< 
class UFunction, 
class VFunction >
 
   62      struct FunctionDistance;
 
   64      typedef typename BaseType::EntityType EntityType;
 
   74      explicit LPNorm ( 
const GridPartType &gridPart, 
const double p, 
const bool communicate = 
true  );
 
   77      template< 
class DiscreteFunctionType, 
class PartitionSet >
 
   78      typename Dune::FieldTraits< typename DiscreteFunctionType::RangeFieldType >::real_type
 
   82      template< 
class DiscreteFunctionType >
 
   83      typename Dune::FieldTraits< typename DiscreteFunctionType::RangeFieldType >::real_type
 
   90      template< 
class UDiscreteFunctionType, 
class VDiscreteFunctionType, 
class PartitionSet >
 
   91      typename Dune::FieldTraits< typename UDiscreteFunctionType::RangeFieldType >::real_type
 
   92      distance ( 
const UDiscreteFunctionType &u, 
const VDiscreteFunctionType &v, 
const PartitionSet& partitionSet ) 
const;
 
   95      template< 
class UDiscreteFunctionType, 
class VDiscreteFunctionType >
 
   96      typename Dune::FieldTraits< typename UDiscreteFunctionType::RangeFieldType >::real_type
 
   97      distance ( 
const UDiscreteFunctionType &u, 
const VDiscreteFunctionType &v )
 const 
  102      template< 
class ULocalFunctionType, 
class VLocalFunctionType, 
class ReturnType >
 
  103      void distanceLocal ( 
const EntityType &entity, 
unsigned int order, 
const ULocalFunctionType &uLocal, 
const VLocalFunctionType &vLocal, ReturnType &sum ) 
const;
 
  105      template< 
class LocalFunctionType, 
class ReturnType >
 
  106      void normLocal ( 
const EntityType &entity, 
unsigned int order, 
const LocalFunctionType &uLocal, ReturnType &sum ) 
const;
 
  108      int order ( 
const int spaceOrder ) 
const ;
 
  112      OrderCalculator *orderCalc_;
 
  113      const bool communicate_;
 
  122    template< 
class WeightFunction, 
class OrderCalculator = DefaultOrderCalculator >
 
  124    : 
public LPNorm< typename WeightFunction::DiscreteFunctionSpaceType::GridPartType,
 
  127      typedef WeightedLPNorm< WeightFunction, OrderCalculator > ThisType;
 
  128      typedef LPNorm< typename WeightFunction::DiscreteFunctionSpaceType::GridPartType, OrderCalculator> BaseType;
 
  131      typedef WeightFunction WeightFunctionType;
 
  133      typedef typename WeightFunctionType::DiscreteFunctionSpaceType WeightFunctionSpaceType;
 
  134      typedef typename WeightFunctionSpaceType::GridPartType GridPartType;
 
  137      template< 
class Function >
 
  138      struct WeightedFunctionMultiplicator;
 
  140      typedef ConstLocalFunction< WeightFunctionType > LocalWeightFunctionType;
 
  141      typedef typename WeightFunctionType::RangeType WeightType;
 
  143      typedef typename BaseType::EntityType EntityType;
 
  144      typedef typename BaseType::IntegratorType IntegratorType;
 
  146      using BaseType::gridPart;
 
  147      using BaseType::comm;
 
  150      using BaseType::norm;
 
  151      using BaseType::distance;
 
  153      explicit WeightedLPNorm ( 
const WeightFunctionType &weightFunction, 
const double p, 
const bool communicate = 
true );
 
  155      template< 
class LocalFunctionType, 
class ReturnType >
 
  156      void normLocal ( 
const EntityType &entity, 
unsigned int order, 
const LocalFunctionType &uLocal, ReturnType &sum ) 
const;
 
  158      template< 
class ULocalFunctionType, 
class VLocalFunctionType, 
class ReturnType >
 
  159      void distanceLocal ( 
const EntityType &entity, 
unsigned int order, 
const ULocalFunctionType &uLocal, 
const VLocalFunctionType &vLocal, ReturnType &sum ) 
const;
 
  162      LocalWeightFunctionType wfLocal_;
 
  170    template< 
class Gr
idPart, 
class OrderCalculator >
 
  171    inline LPNorm< GridPart, OrderCalculator >::LPNorm ( 
const GridPartType &gridPart, 
const double p, 
const bool communicate )
 
  172      :  BaseType( gridPart ),
 
  174         orderCalc_( new OrderCalculator() ),
 
  175         communicate_( BaseType::checkCommunicateFlag( communicate ) )
 
  179    template< 
class Gr
idPart, 
class OrderCalculator>
 
  180    inline int LPNorm< GridPart, OrderCalculator>::order(
const int spaceOrder)
 const 
  182      return spaceOrder * orderCalc_->operator() (p_);
 
  186    template< 
class Gr
idPart, 
class OrderCalculator>
 
  187    template< 
class DiscreteFunctionType, 
class PartitionSet >
 
  188    inline typename Dune::FieldTraits< typename DiscreteFunctionType::RangeFieldType >::real_type
 
  189    LPNorm< GridPart, OrderCalculator >::norm ( 
const DiscreteFunctionType &u, 
const PartitionSet& partitionSet )
 const 
  192      typedef typename Dune::FieldTraits< RangeFieldType >::real_type RealType;
 
  193      typedef FieldVector< RealType, 1 > ReturnType ;
 
  201        sum[ 0 ] = comm().sum( sum[ 0 ] );
 
  205      return std::pow ( sum[ 0 ], (1.0 / p_) );
 
  208    template< 
class Gr
idPart, 
class OrderCalculator >
 
  209    template< 
class UDiscreteFunctionType, 
class VDiscreteFunctionType, 
class PartitionSet >
 
  210    inline typename Dune::FieldTraits< typename UDiscreteFunctionType::RangeFieldType >::real_type
 
  211    LPNorm< GridPart, OrderCalculator >
 
  212      ::distance ( 
const UDiscreteFunctionType &u, 
const VDiscreteFunctionType &v, 
const PartitionSet& partitionSet )
 const 
  214      typedef typename UDiscreteFunctionType::RangeFieldType RangeFieldType;
 
  215      typedef typename Dune::FieldTraits< RangeFieldType >::real_type RealType;
 
  216      typedef FieldVector< RealType, 1 > ReturnType ;
 
  224        sum[ 0 ] = comm().sum( sum[ 0 ] );
 
  228      return std::pow( sum[ 0 ], (1.0/p_) );
 
  231    template< 
class Gr
idPart, 
class OrderCalculator >
 
  232    template< 
class LocalFunctionType, 
class ReturnType >
 
  234    LPNorm< GridPart, OrderCalculator >::normLocal ( 
const EntityType &entity, 
unsigned int order, 
const LocalFunctionType &uLocal, ReturnType &sum )
 const 
  236      IntegratorType integrator( order );
 
  238      FunctionMultiplicator< LocalFunctionType > uLocalp( uLocal, p_ );
 
  240      integrator.integrateAdd( entity, uLocalp, sum );
 
  243    template< 
class Gr
idPart, 
class OrderCalculator >
 
  244    template< 
class ULocalFunctionType, 
class VLocalFunctionType, 
class ReturnType >
 
  246    LPNorm< GridPart, OrderCalculator >::distanceLocal ( 
const EntityType &entity, 
unsigned int order, 
const ULocalFunctionType &uLocal, 
const VLocalFunctionType &vLocal, ReturnType &sum )
 const 
  248      typedef FunctionDistance< ULocalFunctionType, VLocalFunctionType > LocalDistanceType;
 
  250      IntegratorType integrator( order );
 
  252      LocalDistanceType dist( uLocal, vLocal );
 
  253      FunctionMultiplicator< LocalDistanceType > distp( dist, p_ );
 
  255      integrator.integrateAdd( entity, distp, sum );
 
  259    template< 
class Gr
idPart, 
class OrderCalculator >
 
  260    template< 
class Function >
 
  261    struct LPNorm< GridPart, OrderCalculator >::FunctionMultiplicator
 
  263      typedef Function FunctionType;
 
  265      typedef typename FunctionType::RangeFieldType RangeFieldType;
 
  266      typedef typename Dune::FieldTraits< RangeFieldType >::real_type RealType;
 
  267      typedef FieldVector< RealType, 1 > RangeType ;
 
  269      explicit FunctionMultiplicator ( 
const FunctionType &function, 
double p )
 
  270      : function_( function ),
 
  274      template< 
class Po
int >
 
  275      void evaluate ( 
const Point &x, RangeType &ret )
 const 
  277        typename FunctionType::RangeType phi;
 
  278        function_.evaluate( x, phi );
 
  279        ret = std :: pow ( phi.two_norm(), p_);
 
  283      const FunctionType &function_;
 
  288    template< 
class Gr
idPart, 
class OrderCalculator >
 
  289    template< 
class UFunction, 
class VFunction >
 
  290    struct LPNorm< GridPart, OrderCalculator >::FunctionDistance
 
  292      typedef UFunction UFunctionType;
 
  293      typedef VFunction VFunctionType;
 
  295      typedef typename UFunctionType::RangeFieldType RangeFieldType;
 
  296      typedef typename Dune::FieldTraits< RangeFieldType >::real_type RealType;
 
  297      typedef typename UFunctionType::RangeType RangeType;
 
  298      typedef typename UFunctionType::JacobianRangeType JacobianRangeType;
 
  300      FunctionDistance ( 
const UFunctionType &u, 
const VFunctionType &v )
 
  304      template< 
class Po
int >
 
  305      void evaluate ( 
const Point &x, RangeType &ret )
 const 
  308        u_.evaluate( x, ret );
 
  309        v_.evaluate( x, phi );
 
  313      template< 
class Po
int >
 
  314      void jacobian ( 
const Point &x, JacobianRangeType &ret )
 const 
  316        JacobianRangeType phi;
 
  317        u_.jacobian( x, ret );
 
  318        v_.jacobian( x, phi );
 
  323      const UFunctionType &u_;
 
  324      const VFunctionType &v_;
 
  331    template< 
class WeightFunction, 
class OrderCalculator >
 
  332    inline WeightedLPNorm< WeightFunction, OrderCalculator >
 
  333      ::WeightedLPNorm ( 
const WeightFunctionType &weightFunction, 
double p, 
const bool communicate )
 
  334    : BaseType( weightFunction.space().gridPart(), p, communicate ),
 
  335      wfLocal_( weightFunction ),
 
  338      static_assert( (WeightFunctionSpaceType::dimRange == 1),
 
  339                          "Weight function must be scalar." );
 
  343    template< 
class WeightFunction, 
class OrderCalculator >
 
  344    template< 
class LocalFunctionType, 
class ReturnType >
 
  346    WeightedLPNorm< WeightFunction, OrderCalculator >::normLocal ( 
const EntityType &entity, 
unsigned int order, 
const LocalFunctionType &uLocal, ReturnType &sum )
 const 
  349      IntegratorType integrator( order );
 
  351      wfLocal_.bind( entity );
 
  352      WeightedFunctionMultiplicator< LocalFunctionType > uLocal2( wfLocal_, uLocal, p_ );
 
  353      integrator.integrateAdd( entity, uLocal2, sum );
 
  358    template< 
class WeightFunction,
class OrderCalculator >
 
  359    template< 
class ULocalFunctionType, 
class VLocalFunctionType, 
class ReturnType >
 
  361    WeightedLPNorm< WeightFunction, OrderCalculator >::distanceLocal ( 
const EntityType &entity, 
unsigned int order, 
const ULocalFunctionType &uLocal, 
const VLocalFunctionType &vLocal, ReturnType &sum )
 const 
  363      typedef typename BaseType::template FunctionDistance< ULocalFunctionType, VLocalFunctionType > LocalDistanceType;
 
  366      IntegratorType integrator( order );
 
  368      wfLocal_.bind( entity );
 
  369      LocalDistanceType dist( uLocal, vLocal );
 
  370      WeightedFunctionMultiplicator< LocalDistanceType > dist2( wfLocal_, dist );
 
  372      integrator.integrateAdd( entity, dist2, sum );
 
  377    template< 
class WeightFunction, 
class OrderCalculator>
 
  378    template< 
class Function >
 
  379    struct WeightedLPNorm< WeightFunction, OrderCalculator>::WeightedFunctionMultiplicator
 
  381      typedef Function FunctionType;
 
  383      typedef typename FunctionType::RangeFieldType RangeFieldType;
 
  384      typedef typename Dune::FieldTraits< RangeFieldType >::real_type RealType;
 
  385      typedef FieldVector< RealType, 1 > RangeType;
 
  387      WeightedFunctionMultiplicator ( 
const LocalWeightFunctionType &weightFunction,
 
  388                                      const FunctionType &function,
 
  390      : weightFunction_( weightFunction ),
 
  391        function_( function ),
 
  395      template< 
class Po
int >
 
  396      void evaluate ( 
const Point &x, RangeType &ret )
 const 
  399        weightFunction_.evaluate( x, weight );
 
  401        typename FunctionType::RangeType phi;
 
  402        function_.evaluate( x, phi );
 
  403        ret = weight[ 0 ] * std::pow ( phi.two_norm(), p_);
 
  407      const LocalWeightFunctionType &weightFunction_;
 
  408      const FunctionType &function_;
 
quadrature class supporting base function caching
Definition: cachingquadrature.hh:106
 
DiscreteFunctionSpaceType::RangeFieldType RangeFieldType
type of range field (usually a float type)
Definition: discretefunction.hh:623
 
integrator for arbitrary functions providing evaluate
Definition: integrator.hh:28
 
forward declaration
Definition: discretefunction.hh:51
 
constexpr void forEach(Range &&range, F &&f)
Range based for loop.
Definition: hybridutilities.hh:257
 
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:485
 
constexpr Interior interior
PartitionSet for the interior partition.
Definition: partitionset.hh:271
 
Dune namespace.
Definition: alignedallocator.hh:13
 
default Quadrature Order Calculator
Definition: lpnorm.hh:33
 
Quadrature Order Interface.
Definition: lpnorm.hh:24
 
A set of PartitionType values.
Definition: partitionset.hh:137