1#ifndef DUNE_FEM_L2NORM_HH 
    2#define DUNE_FEM_L2NORM_HH 
    4#include <dune/fem/quadrature/cachingquadrature.hh> 
    5#include <dune/fem/quadrature/integrator.hh> 
    7#include <dune/fem/misc/domainintegral.hh> 
   18    template< 
class Gr
idPart >
 
   19    class L2Norm : 
public IntegralBase< GridPart, L2Norm< GridPart > >
 
   21      typedef IntegralBase< GridPart, L2Norm< GridPart > > BaseType ;
 
   22      typedef L2Norm< GridPart > ThisType;
 
   25      typedef GridPart GridPartType;
 
   27      using BaseType :: gridPart ;
 
   28      using BaseType :: comm ;
 
   30      template< 
class Function >
 
   31      struct FunctionSquare;
 
   33      template< 
class UFunction, 
class VFunction >
 
   34      struct FunctionDistance;
 
   37      typedef typename BaseType::EntityType EntityType;
 
   38      typedef CachingQuadrature< GridPartType, 0 > QuadratureType;
 
   39      typedef Integrator< QuadratureType > IntegratorType;
 
   41      const unsigned int order_;
 
   42      const bool communicate_;
 
   49      explicit L2Norm ( 
const GridPartType &gridPart,
 
   50                        const unsigned int order = 0,
 
   51                        const bool communicate = 
true );
 
   54      template< 
class DiscreteFunctionType, 
class PartitionSet >
 
   55      typename Dune::FieldTraits< typename DiscreteFunctionType::RangeFieldType >::real_type
 
   59      template< 
class DiscreteFunctionType >
 
   60      typename Dune::FieldTraits< typename DiscreteFunctionType::RangeFieldType >::real_type
 
   67      template< 
class UDiscreteFunctionType, 
class VDiscreteFunctionType, 
class PartitionSet >
 
   68      typename Dune::FieldTraits< typename UDiscreteFunctionType::RangeFieldType >::real_type
 
   69      distance ( 
const UDiscreteFunctionType &u,
 
   70                 const VDiscreteFunctionType &v,
 
   71                 const PartitionSet& partitionSet ) 
const;
 
   74      template< 
class UDiscreteFunctionType, 
class VDiscreteFunctionType >
 
   75      typename Dune::FieldTraits< typename UDiscreteFunctionType::RangeFieldType >::real_type
 
   76      distance ( 
const UDiscreteFunctionType &u, 
const VDiscreteFunctionType &v )
 const 
   81      template< 
class ULocalFunctionType, 
class VLocalFunctionType, 
class ReturnType >
 
   82      void distanceLocal ( 
const EntityType &entity, 
unsigned int order, 
const ULocalFunctionType &uLocal, 
const VLocalFunctionType &vLocal, ReturnType &sum ) 
const;
 
   84      template< 
class LocalFunctionType, 
class ReturnType >
 
   85      void normLocal ( 
const EntityType &entity, 
unsigned int order, 
const LocalFunctionType &uLocal, ReturnType &sum ) 
const;
 
   92    template< 
class Gr
idPart >
 
   93    inline L2Norm< GridPart >::L2Norm ( 
const GridPartType &gridPart, 
const unsigned int order, 
const bool communicate )
 
   94    : BaseType( gridPart ),
 
   96      communicate_( BaseType::checkCommunicateFlag( communicate ) )
 
  101    template< 
class Gr
idPart >
 
  102    template< 
class DiscreteFunctionType, 
class PartitionSet >
 
  103    typename Dune::FieldTraits< typename DiscreteFunctionType::RangeFieldType >::real_type
 
  106      typedef typename Dune::FieldTraits< typename DiscreteFunctionType::RangeFieldType >::real_type RealType;
 
  107      typedef FieldVector< RealType, 1 > ReturnType ;
 
  115        sum[ 0 ] = comm().sum( sum[ 0 ] );
 
  119      return std::sqrt( sum[ 0 ] );
 
  123    template< 
class Gr
idPart >
 
  124    template< 
class UDiscreteFunctionType, 
class VDiscreteFunctionType, 
class PartitionSet >
 
  125    typename Dune::FieldTraits< typename UDiscreteFunctionType::RangeFieldType >::real_type
 
  127      ::distance ( 
const UDiscreteFunctionType &u,
 
  128                   const VDiscreteFunctionType &v,
 
  131      typedef typename Dune::FieldTraits< typename UDiscreteFunctionType::RangeFieldType >::real_type RealType;
 
  132      typedef FieldVector< RealType, 1 > ReturnType ;
 
  140        sum[ 0 ] = comm().sum( sum[ 0 ] );
 
  144      return std::sqrt( sum[ 0 ] );
 
  147    template< 
class Gr
idPart >
 
  148    template< 
class LocalFunctionType, 
class ReturnType >
 
  149    inline void L2Norm< GridPart >::normLocal ( 
const EntityType &entity, 
unsigned int order, 
const LocalFunctionType &uLocal, ReturnType &sum )
 const 
  152      IntegratorType integrator( order );
 
  154      FunctionSquare< LocalFunctionType > uLocal2( uLocal );
 
  156      integrator.integrateAdd( entity, uLocal2, sum );
 
  159    template< 
class Gr
idPart >
 
  160    template< 
class ULocalFunctionType, 
class VLocalFunctionType, 
class ReturnType >
 
  162    L2Norm< GridPart >::distanceLocal ( 
const EntityType &entity, 
unsigned int order, 
const ULocalFunctionType &uLocal, 
const VLocalFunctionType &vLocal, ReturnType &sum )
 const 
  165      IntegratorType integrator( order );
 
  167      typedef FunctionDistance< ULocalFunctionType, VLocalFunctionType > LocalDistanceType;
 
  169      LocalDistanceType dist( uLocal, vLocal );
 
  170      FunctionSquare< LocalDistanceType > dist2( dist );
 
  172      integrator.integrateAdd( entity, dist2, sum );
 
  176    template< 
class Gr
idPart >
 
  177    template< 
class Function >
 
  178    struct L2Norm< GridPart >::FunctionSquare
 
  180      typedef Function FunctionType;
 
  182      typedef typename FunctionType::RangeFieldType RangeFieldType;
 
  183      typedef typename Dune::FieldTraits< RangeFieldType >::real_type RealType;
 
  184      typedef FieldVector< RealType, 1 > RangeType;
 
  186      explicit FunctionSquare ( 
const FunctionType &function )
 
  187      : function_( function )
 
  190      template< 
class Po
int >
 
  191      void evaluate ( 
const Point &x, RangeType &ret )
 const 
  193        typename FunctionType::RangeType phi;
 
  194        function_.evaluate( x, phi );
 
  195        ret = phi.two_norm2();
 
  199      const FunctionType &function_;
 
  203    template< 
class Gr
idPart >
 
  204    template< 
class UFunction, 
class VFunction >
 
  205    struct L2Norm< GridPart >::FunctionDistance
 
  207      typedef UFunction UFunctionType;
 
  208      typedef VFunction VFunctionType;
 
  210      typedef typename UFunctionType::RangeFieldType RangeFieldType;
 
  211      typedef typename UFunctionType::RangeType RangeType;
 
  212      typedef typename UFunctionType::JacobianRangeType JacobianRangeType;
 
  214      FunctionDistance ( 
const UFunctionType &u, 
const VFunctionType &v )
 
  218      template< 
class Po
int >
 
  219      void evaluate ( 
const Point &x, RangeType &ret )
 const 
  222        u_.evaluate( x, ret );
 
  223        v_.evaluate( x, phi );
 
  227      template< 
class Po
int >
 
  228      void jacobian ( 
const Point &x, JacobianRangeType &ret )
 const 
  230        JacobianRangeType phi;
 
  231        u_.jacobian( x, ret );
 
  232        v_.jacobian( x, phi );
 
  237      const UFunctionType &u_;
 
  238      const VFunctionType &v_;
 
  246    template< 
class WeightFunction >
 
  248    : 
public L2Norm< typename WeightFunction::DiscreteFunctionSpaceType::GridPartType >
 
  250      typedef WeightedL2Norm< WeightFunction > ThisType;
 
  251      typedef L2Norm< typename WeightFunction::DiscreteFunctionSpaceType::GridPartType > BaseType;
 
  254      typedef WeightFunction WeightFunctionType;
 
  256      typedef typename WeightFunctionType::DiscreteFunctionSpaceType WeightFunctionSpaceType;
 
  257      typedef typename WeightFunctionSpaceType::GridPartType GridPartType;
 
  260      template< 
class Function >
 
  261      struct WeightedFunctionSquare;
 
  263      typedef ConstLocalFunction< WeightFunctionType > LocalWeightFunctionType;
 
  264      typedef typename WeightFunctionType::RangeType WeightType;
 
  266      typedef typename BaseType::EntityType EntityType;
 
  267      typedef typename BaseType::IntegratorType IntegratorType;
 
  269      using BaseType::gridPart;
 
  270      using BaseType::comm;
 
  274      using BaseType::norm;
 
  275      using BaseType::distance;
 
  277      explicit WeightedL2Norm ( 
const WeightFunctionType &weightFunction, 
const unsigned int order = 0, 
const bool communicate = 
true );
 
  279      template< 
class LocalFunctionType, 
class ReturnType >
 
  280      void normLocal ( 
const EntityType &entity, 
unsigned int order, 
const LocalFunctionType &uLocal, ReturnType &sum ) 
const;
 
  282      template< 
class ULocalFunctionType, 
class VLocalFunctionType, 
class ReturnType >
 
  283      void distanceLocal ( 
const EntityType &entity, 
unsigned int order, 
const ULocalFunctionType &uLocal, 
const VLocalFunctionType &vLocal, ReturnType &sum ) 
const;
 
  286      LocalWeightFunctionType wfLocal_;
 
  295    template< 
class WeightFunction >
 
  296    inline WeightedL2Norm< WeightFunction >
 
  297      ::WeightedL2Norm ( 
const WeightFunctionType &weightFunction, 
const unsigned int order, 
const bool communicate )
 
  298    : BaseType( weightFunction.space().gridPart(), order, communicate ),
 
  299      wfLocal_( weightFunction )
 
  301      static_assert( (WeightFunctionSpaceType::dimRange == 1),
 
  302                     "Wight function must be scalar." );
 
  306    template< 
class WeightFunction >
 
  307    template< 
class LocalFunctionType, 
class ReturnType >
 
  309    WeightedL2Norm< WeightFunction >::normLocal ( 
const EntityType &entity, 
unsigned int order, 
const LocalFunctionType &uLocal, ReturnType &sum )
 const 
  312      IntegratorType integrator( order );
 
  314      wfLocal_.bind( entity );
 
  315      WeightedFunctionSquare< LocalFunctionType > uLocal2( wfLocal_, uLocal );
 
  316      integrator.integrateAdd( entity, uLocal2, sum );
 
  321    template< 
class WeightFunction >
 
  322    template< 
class ULocalFunctionType, 
class VLocalFunctionType, 
class ReturnType >
 
  324    WeightedL2Norm< WeightFunction >::distanceLocal ( 
const EntityType& entity, 
unsigned int order, 
const ULocalFunctionType &uLocal, 
const VLocalFunctionType &vLocal, ReturnType& sum )
 const 
  326      typedef typename BaseType::template FunctionDistance< ULocalFunctionType, VLocalFunctionType > LocalDistanceType;
 
  329      IntegratorType integrator( order );
 
  331      wfLocal_.bind( entity );
 
  332      LocalDistanceType dist( uLocal, vLocal );
 
  333      WeightedFunctionSquare< LocalDistanceType > dist2( wfLocal_, dist );
 
  335      integrator.integrateAdd( entity, dist2, sum );
 
  340    template< 
class WeightFunction >
 
  341    template< 
class Function >
 
  342    struct WeightedL2Norm< WeightFunction >::WeightedFunctionSquare
 
  344      typedef Function FunctionType;
 
  346      typedef typename FunctionType::RangeFieldType RangeFieldType;
 
  347      typedef typename Dune::FieldTraits< RangeFieldType >::real_type RealType;
 
  348      typedef FieldVector< RealType, 1 > RangeType;
 
  350      WeightedFunctionSquare ( 
const LocalWeightFunctionType &weightFunction,
 
  351                               const FunctionType &function )
 
  352      : weightFunction_( weightFunction ),
 
  353        function_( function )
 
  356      template< 
class Po
int >
 
  357      void evaluate ( 
const Point &x, RangeType &ret )
 const 
  360        weightFunction_.evaluate( x, weight );
 
  362        typename FunctionType::RangeType phi;
 
  363        function_.evaluate( x, phi );
 
  364        ret = weight[ 0 ] * (phi * phi);
 
  368      const LocalWeightFunctionType &weightFunction_;
 
  369      const FunctionType &function_;
 
forward declaration
Definition: discretefunction.hh:51
 
constexpr void forEach(Range &&range, F &&f)
Range based for loop.
Definition: hybridutilities.hh:257
 
constexpr Interior interior
PartitionSet for the interior partition.
Definition: partitionset.hh:271
 
Dune namespace.
Definition: alignedallocator.hh:13