1#ifndef DUNE_FEM_INTEGRAL_HH 
    2#define DUNE_FEM_INTEGRAL_HH 
    7#include <dune/grid/common/rangegenerators.hh> 
    9#include <dune/fem/quadrature/cachingquadrature.hh> 
   10#include <dune/fem/quadrature/integrator.hh> 
   12#include <dune/fem/function/common/gridfunctionadapter.hh> 
   15#include <dune/fem/misc/bartonnackmaninterface.hh> 
   17#include <dune/fem/misc/mpimanager.hh> 
   18#include <dune/fem/misc/threads/threaditerator.hh> 
   19#include <dune/fem/common/bindguard.hh> 
   29    template< 
class Gr
idPart, 
class NormImplementation >
 
   31      : 
public BartonNackmanInterface< IntegralBase< GridPart, NormImplementation >,
 
   34      typedef BartonNackmanInterface< IntegralBase< GridPart, NormImplementation >,
 
   35                                      NormImplementation >  BaseType ;
 
   36      typedef IntegralBase< GridPart, NormImplementation > ThisType;
 
   39      typedef GridPart GridPartType;
 
   42      using BaseType :: asImp ;
 
   44      typedef typename GridPartType::template Codim< 0 >::EntityType EntityType;
 
   46      template <
bool uDiscrete, 
bool vDiscrete>
 
   50        template <
class IteratorRange, 
class UDiscreteFunctionType, 
class VDiscreteFunctionType, 
class ReturnType>
 
   51        static ReturnType forEachLocal ( 
const ThisType &norm, 
const IteratorRange& iterators,
 
   52                                         const UDiscreteFunctionType &u, 
const VDiscreteFunctionType &v,
 
   53                                         const ReturnType &initialValue,
 
   56          static_assert( uDiscrete && vDiscrete, 
"Distance can only be calculated between GridFunctions" );
 
   60            ConstLocalFunction< UDiscreteFunctionType > uLocal( u );
 
   61            ConstLocalFunction< VDiscreteFunctionType > vLocal( v );
 
   62            for( 
const EntityType &entity : iterators )
 
   64              auto uGuard = bindGuard( uLocal, entity );
 
   65              auto vGuard = bindGuard( vLocal, entity );
 
   67              const unsigned int uOrder = uLocal.order();
 
   68              const unsigned int vOrder = vLocal.order();
 
   69              const unsigned int orderLocal = (order == 0 ? 2*std::max( uOrder, vOrder ) : order);
 
   70              norm.distanceLocal( entity, orderLocal, uLocal, vLocal, sum );
 
   77        template <
class UDiscreteFunctionType, 
class VDiscreteFunctionType, 
class ReturnType, 
class PartitionSet>
 
   78        static ReturnType 
forEach ( 
const ThisType &norm, 
const UDiscreteFunctionType &u, 
const VDiscreteFunctionType &v,
 
   79                                    const ReturnType &initialValue,
 
   80                                    const PartitionSet& partitionSet,
 
   83          const int nThreads = MPIManager::numThreads();
 
   85            return forEachLocal( norm, elements( norm.gridPart_, partitionSet ), u, v, initialValue, order );
 
   87          std::vector< ReturnType > sums( nThreads, ReturnType(0) );
 
   90            iterators( norm.gridPart_ );
 
   92          auto doIntegrate = [ &norm, &iterators, &sums, &u, &v, &initialValue, &order ] ()
 
   94            sums[ MPIManager::thread() ] = forEachLocal( norm, iterators, u, v, initialValue, order );
 
   99            MPIManager::run( doIntegrate );
 
  101          catch ( 
const SingleThreadModeError& e )
 
  104            return forEachLocal( norm, elements( norm.gridPart_, partitionSet ), u, v, initialValue, order );
 
  112      template <
bool vDiscrete>
 
  113      struct ForEachCaller<false, vDiscrete>
 
  116                  class VDiscreteFunctionType,
 
  120        ReturnType 
forEach ( 
const ThisType& norm,
 
  121                             const F& f, 
const VDiscreteFunctionType& v,
 
  122                             const ReturnType& initialValue,
 
  123                             const PartitionSet& partitionSet,
 
  124                             const unsigned int order )
 
  126          typedef GridFunctionAdapter< F, GridPartType>  GridFunction ;
 
  127          GridFunction u( 
"Integral::adapter" , f , v.space().gridPart(), v.space().order() );
 
  130            forEach( norm, u, v, initialValue, partitionSet, order );
 
  135      template <
bool uDiscrete>
 
  136      struct ForEachCaller<uDiscrete, false>
 
  138        template <
class UDiscreteFunctionType,
 
  143        ReturnType 
forEach ( 
const ThisType& norm,
 
  144                             const UDiscreteFunctionType& u,
 
  146                             const ReturnType& initialValue,
 
  147                             const PartitionSet& partitionSet,
 
  148                             const unsigned int order )
 
  151            forEach( norm, f, u, initialValue, partitionSet, order );
 
  155      template <
class IteratorRange, 
class UDiscreteFunctionType, 
class ReturnType>
 
  156      ReturnType forEachLocal ( 
const IteratorRange& iterators,
 
  157                                const UDiscreteFunctionType &u,
 
  158                                const ReturnType &initialValue,
 
  159                                unsigned int order )
 const 
  161        ReturnType sum( initialValue );
 
  163          ConstLocalFunction< UDiscreteFunctionType > uLocal( u );
 
  164          for( 
const EntityType &entity : iterators )
 
  166            auto uGuard = bindGuard( uLocal, entity );
 
  168            const unsigned int uOrder = uLocal.order();
 
  169            const unsigned int orderLocal = (order == 0 ? 2*uOrder : order);
 
  170            normLocal( entity, orderLocal, uLocal, sum );
 
  176      template< 
class DiscreteFunctionType, 
class ReturnType, 
class PartitionSet >
 
  178                           const PartitionSet& partitionSet,
 
  179                           unsigned int order = 0 )
 const 
  181        static_assert( (std::is_base_of<Fem::HasLocalFunction, DiscreteFunctionType>::value),
 
  182                            "Norm only implemented for quantities implementing a local function!" );
 
  183        const int nThreads = MPIManager::numThreads();
 
  185          return forEachLocal( elements( gridPart_, partitionSet ), u, initialValue, order );
 
  187        std::vector< ReturnType > sums( nThreads, ReturnType(0) );
 
  190          iterators( gridPart_ );
 
  192        auto doIntegrate = [ 
this, &iterators, &sums, &u, &initialValue, &order ] ()
 
  194          sums[ MPIManager::thread() ] = forEachLocal( iterators, u, initialValue, order );
 
  199          MPIManager::run( doIntegrate );
 
  201        catch ( 
const SingleThreadModeError& e )
 
  204          return forEachLocal( elements( gridPart_, partitionSet ), u, initialValue, order );
 
  210      template< 
class UDiscreteFunctionType, 
class VDiscreteFunctionType, 
class ReturnType, 
class PartitionSet >
 
  211      ReturnType 
forEach ( 
const UDiscreteFunctionType &u, 
const VDiscreteFunctionType &v,
 
  212                           const ReturnType &initialValue, 
const PartitionSet& partitionSet,
 
  213                           unsigned int order = 0 )
 const 
  215        enum { uDiscrete = std::is_convertible<UDiscreteFunctionType, HasLocalFunction>::value };
 
  216        enum { vDiscrete = std::is_convertible<VDiscreteFunctionType, HasLocalFunction>::value };
 
  221                  forEach( *
this, u, v, initialValue, partitionSet, order );
 
  225      explicit IntegralBase ( 
const GridPartType &gridPart )
 
  226        : gridPart_( gridPart )
 
  230      template< 
class ULocalFunctionType, 
class VLocalFunctionType, 
class ReturnType >
 
  231      void distanceLocal ( 
const EntityType &entity, 
unsigned int order, 
const ULocalFunctionType &uLocal, 
const VLocalFunctionType &vLocal, ReturnType &sum )
 const 
  236      template< 
class LocalFunctionType, 
class ReturnType >
 
  237      void normLocal ( 
const EntityType &entity, 
unsigned int order, 
const LocalFunctionType &uLocal, ReturnType &sum )
 const 
  242      const GridPartType &gridPart ()
 const { 
return gridPart_; }
 
  246        return gridPart().comm();
 
  249      bool checkCommunicateFlag( 
bool communicate )
 const 
  252        bool commMax = communicate;
 
  253        assert( communicate == comm().
max( commMax ) );
 
  259      const GridPartType &gridPart_;
 
  266    template< 
class Gr
idPart >
 
  267    class Integral : 
public IntegralBase< GridPart, Integral< GridPart > >
 
  269      typedef IntegralBase< GridPart, Integral< GridPart > > BaseType ;
 
  270      typedef Integral< GridPart > ThisType;
 
  273      typedef GridPart GridPartType;
 
  275      using BaseType :: gridPart ;
 
  276      using BaseType :: comm ;
 
  280      template< 
class UFunction, 
class VFunction >
 
  281      struct FunctionDistance;
 
  283      typedef typename BaseType::EntityType EntityType;
 
  284      typedef CachingQuadrature< GridPartType, 0 > QuadratureType;
 
  286      const unsigned int order_;
 
  287      const bool communicate_;
 
  294      explicit Integral ( 
const GridPartType &gridPart,
 
  295                          const unsigned int order = 0,
 
  296                          const bool communicate = 
true );
 
  300      template< 
class DiscreteFunctionType, 
class PartitionSet >
 
  305      template< 
class DiscreteFunctionType >
 
  313      template< 
class UDiscreteFunctionType, 
class VDiscreteFunctionType, 
class PartitionSet >
 
  314      typename UDiscreteFunctionType::RangeType
 
  315      distance ( 
const UDiscreteFunctionType &u, 
const VDiscreteFunctionType &v, 
const PartitionSet& partitionSet ) 
const;
 
  318      template< 
class UDiscreteFunctionType, 
class VDiscreteFunctionType >
 
  319      typename UDiscreteFunctionType::RangeType
 
  320      distance ( 
const UDiscreteFunctionType &u, 
const VDiscreteFunctionType &v )
 const 
  325      template< 
class LocalFunctionType, 
class ReturnType >
 
  326      void normLocal ( 
const EntityType &entity, 
unsigned int order, 
const LocalFunctionType &uLocal, ReturnType &sum ) 
const;
 
  328      template< 
class ULocalFunctionType, 
class VLocalFunctionType, 
class ReturnType >
 
  329      void distanceLocal ( 
const EntityType &entity, 
unsigned int order, 
const ULocalFunctionType &uLocal, 
const VLocalFunctionType &vLocal, ReturnType &sum ) 
const;
 
  337    template< 
class Gr
idPart >
 
  338    inline Integral< GridPart >::Integral ( 
const GridPartType &gridPart, 
const unsigned int order, 
const bool communicate )
 
  339    : BaseType( gridPart ),
 
  341      communicate_( BaseType::checkCommunicateFlag( communicate ) )
 
  345    template< 
class Gr
idPart >
 
  346    template< 
class DiscreteFunctionType, 
class PartitionSet >
 
  351      typedef RangeType ReturnType ;
 
  359        sum = comm().sum( sum );
 
  366    template< 
class Gr
idPart >
 
  367    template< 
class UDiscreteFunctionType, 
class VDiscreteFunctionType, 
class PartitionSet >
 
  368    inline typename UDiscreteFunctionType::RangeType
 
  370      ::distance ( 
const UDiscreteFunctionType &u, 
const VDiscreteFunctionType &v, 
const PartitionSet& 
partitionSet )
 const 
  372      typedef typename UDiscreteFunctionType::RangeType RangeType;
 
  373      typedef RangeType ReturnType ;
 
  381        sum = comm().sum( sum );
 
  387    template< 
class Gr
idPart >
 
  388    template< 
class LocalFunctionType, 
class ReturnType >
 
  390    Integral< GridPart >::normLocal ( 
const EntityType &entity, 
unsigned int order, 
const LocalFunctionType &uLocal, ReturnType &sum )
 const 
  392      Integrator< QuadratureType > integrator( order );
 
  394      integrator.integrateAdd( entity, uLocal, sum );
 
  397    template< 
class Gr
idPart >
 
  398    template< 
class ULocalFunctionType, 
class VLocalFunctionType, 
class ReturnType >
 
  400    Integral< GridPart >::distanceLocal ( 
const EntityType &entity, 
unsigned int order, 
const ULocalFunctionType &uLocal, 
const VLocalFunctionType &vLocal, ReturnType &sum )
 const 
  402      Integrator< QuadratureType > integrator( order );
 
  404      typedef FunctionDistance< ULocalFunctionType, VLocalFunctionType > LocalDistanceType;
 
  406      LocalDistanceType dist( uLocal, vLocal );
 
  408      integrator.integrateAdd( entity, dist, sum );
 
  411    template< 
class Gr
idPart >
 
  412    template< 
class UFunction, 
class VFunction >
 
  413    struct Integral< GridPart >::FunctionDistance
 
  415      typedef UFunction UFunctionType;
 
  416      typedef VFunction VFunctionType;
 
  418      typedef typename UFunctionType::RangeFieldType RangeFieldType;
 
  419      typedef typename UFunctionType::RangeType RangeType;
 
  420      typedef typename UFunctionType::JacobianRangeType JacobianRangeType;
 
  422      FunctionDistance ( 
const UFunctionType &u, 
const VFunctionType &v )
 
  426      template< 
class Po
int >
 
  427      void evaluate ( 
const Point &x, RangeType &ret )
 const 
  430        u_.evaluate( x, ret );
 
  431        v_.evaluate( x, phi );
 
  435      template< 
class Po
int >
 
  436      void jacobian ( 
const Point &x, JacobianRangeType &ret )
 const 
  438        JacobianRangeType phi;
 
  439        u_.jacobian( x, ret );
 
  440        v_.jacobian( x, phi );
 
  445      const UFunctionType &u_;
 
  446      const VFunctionType &v_;
 
Provides check for implementation of interface methods when using static polymorphism,...
 
#define CHECK_AND_CALL_INTERFACE_IMPLEMENTATION(__interface_method_to_call__)
Definition: bartonnackmanifcheck.hh:61
 
DiscreteFunctionSpaceType::RangeType RangeType
type of range vector
Definition: discretefunction.hh:614
 
Traits::CommunicationType CommunicationType
Collective communication.
Definition: gridpart.hh:97
 
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 T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:280
 
constexpr Interior interior
PartitionSet for the interior partition.
Definition: partitionset.hh:271
 
Dune namespace.
Definition: alignedallocator.hh:13
 
static constexpr PartitionIteratorType partitionIterator()
Returns the PartitionIteratorType that can be used to iterate over the partitions in the set.
Definition: partitionset.hh:182