1#ifndef DUNE_FEM_SPACE_LOCALFINITEELEMENT_INTERPOLATION_HH
2#define DUNE_FEM_SPACE_LOCALFINITEELEMENT_INTERPOLATION_HH
12#include <dune/fem/misc/threads/threadsafevalue.hh>
14#include <dune/fem/function/common/localcontribution.hh>
15#include <dune/fem/space/basisfunctionset/transformed.hh>
16#include <dune/fem/space/basisfunctionset/vectorial.hh>
17#include <dune/fem/space/combinedspace/interpolation.hh>
18#include <dune/fem/space/common/localinterpolation.hh>
30 template<
class LocalFunction,
class BasisFunctionSet >
31 struct LocalFunctionWrapper
41 LocalFunctionWrapper (
const LocalFunction &lf,
const BasisFunctionSet &bset ) : lf_( lf ) {}
44 void evaluate (
const Arg &x,
typename Traits::RangeType &y )
const
49 typename Traits::RangeType operator()(
const Arg &x)
const
51 typename Traits::RangeType y;
57 const LocalFunction &lf_;
61 template<
class LocalFunction,
class Entity,
class ShapeFunctionSet,
class Transformation >
62 struct LocalFunctionWrapper< LocalFunction, TransformedBasisFunctionSet< Entity, ShapeFunctionSet, Transformation > >
64 typedef TransformedBasisFunctionSet< Entity, ShapeFunctionSet, Transformation > BasisFunctionSetType;
74 LocalFunctionWrapper (
const LocalFunction &lf,
const BasisFunctionSetType &bset )
75 : lf_( lf ), geometry_( lf_.entity().geometry() ), bset_( bset ) {}
78 void evaluate (
const Arg &x,
typename Traits::RangeType &y )
const
80 typename Traits::RangeType help;
81 lf_.evaluate( x, help );
82 typename Transformation::InverseTransformationType transf( geometry_, x );
83 y = transf.apply( help );
86 typename Traits::RangeType operator() (
const Arg &x )
const
88 typename Traits::RangeType y;
94 const LocalFunction &lf_;
96 const BasisFunctionSetType &bset_;
105 template<
class Space,
class LocalInterpolation,
bool scalarSFS >
106 class LocalFiniteElementInterpolation;
111 template<
class Space,
class LocalInterpolation >
112 class LocalFiniteElementInterpolation<Space,LocalInterpolation,false>
114 typedef LocalFiniteElementInterpolation< Space, LocalInterpolation,false > ThisType;
117 typedef typename Space::BasisFunctionSetType BasisFunctionSetType;
118 typedef LocalInterpolation LocalInterpolationType;
127 typedef std::size_t size_type;
129 template<
class LocalFunction >
130 using LocalFunctionWrapper = Impl::LocalFunctionWrapper< LocalFunction, BasisFunctionSetType >;
133 LocalFiniteElementInterpolation()
134 : basisFunctionSet_()
135 , localInterpolation_()
139 explicit LocalFiniteElementInterpolation (
const BasisFunctionSetType &basisFunctionSet,
140 const LocalInterpolationType &localInterpolation = LocalInterpolationType() )
141 : basisFunctionSet_( basisFunctionSet ),
142 localInterpolation_( localInterpolation )
145 ThisType& operator=(
const ThisType& other )
147 basisFunctionSet_ = other.basisFunctionSet_;
148 localInterpolation_.emplace( other.localInterpolation() );
157 template<
class LocalFunction,
class Dof >
158 void operator() (
const LocalFunction &localFunction, std::vector< Dof > &localDofVector )
const
160 LocalFunctionWrapper< LocalFunction > wrapper( localFunction, basisFunctionSet() );
161 localInterpolation().interpolate( wrapper, localDofVector );
164 template<
class LocalFunction,
class Dof,
class Allocator >
167 (*this)(localFunction,localDofVector.container() );
170 template<
class LocalFunction,
class DiscreteFunction,
template<
class >
class Assembly >
171 void operator() (
const LocalFunction &localFunction, LocalContribution< DiscreteFunction, Assembly > &localContribution )
const
173 (*this)(localFunction,localContribution.localDofVector());
176 BasisFunctionSetType basisFunctionSet ()
const {
return basisFunctionSet_; }
177 const LocalInterpolationType &localInterpolation ()
const
179 assert( localInterpolation_.has_value() );
180 return *localInterpolation_;
184 BasisFunctionSetType basisFunctionSet_;
185 std::optional< LocalInterpolationType > localInterpolation_;
190 template <
int dimRange>
191 struct RangeConverter
193 RangeConverter ( std::size_t range ) : range_( range ) {}
196 FieldVector< T, 1 > operator() (
const FieldVector< T, dimRange > &in )
const
201 template<
class T,
int j >
202 FieldMatrix< T, 1, j > operator() (
const FieldMatrix< T, dimRange, j > &in )
const
205 FieldMatrix<T,1,j> ret;
213 template <
class DofVector,
class DofAlignment>
214 struct SubDofVectorWrapper
215 :
public SubDofVector< DofVector, DofAlignment >
217 typedef SubDofVector< DofVector, DofAlignment > BaseType;
219 SubDofVectorWrapper( DofVector& dofs,
int coordinate,
const DofAlignment &dofAlignment )
220 : BaseType( dofs, coordinate, dofAlignment )
225 void resize(
const unsigned int) {}
233 template<
class Space,
class LocalInterpolation >
234 class LocalFiniteElementInterpolation<Space,LocalInterpolation,true>
236 typedef LocalFiniteElementInterpolation< Space, LocalInterpolation,true > ThisType;
239 typedef typename Space::BasisFunctionSetType BasisFunctionSetType;
240 typedef LocalInterpolation LocalInterpolationType;
249 static const int dimR = Space::FunctionSpaceType::dimRange;
251 typedef std::size_t size_type;
253 typedef VerticalDofAlignment< BasisFunctionSetType, RangeType> DofAlignmentType;
256 LocalFiniteElementInterpolation ()
257 : basisFunctionSet_(),
258 localInterpolation_(),
262 explicit LocalFiniteElementInterpolation (
const BasisFunctionSetType &basisFunctionSet,
263 const LocalInterpolationType &localInterpolation = LocalInterpolationType() )
264 : basisFunctionSet_( basisFunctionSet ),
265 localInterpolation_( localInterpolation ),
266 dofAlignment_( basisFunctionSet_ )
269 LocalFiniteElementInterpolation (
const ThisType& other )
270 : basisFunctionSet_( other.basisFunctionSet_ ),
271 localInterpolation_( other.localInterpolation_ ),
272 dofAlignment_( basisFunctionSet_ )
275 ThisType& operator=(
const ThisType& other )
277 basisFunctionSet_ = other.basisFunctionSet_;
278 localInterpolation_.emplace( other.localInterpolation() );
287 template<
class LocalFunction,
class Vector>
288 void operator() (
const LocalFunction &localFunction, Vector &localDofVector )
const
292 std::fill(localDofVector.begin(),localDofVector.end(),0);
293 for( std::size_t i = 0; i < dimR; ++i )
295 Impl::SubDofVectorWrapper< Vector, DofAlignmentType > subLdv( localDofVector, i, dofAlignment_ );
297 localFunctionConverter( localFunction, Impl::RangeConverter<dimR>(i) ),
298 subLdv, PriorityTag<42>()
303 template<
class LocalFunction,
class Dof,
class Allocator >
306 (*this)(localFunction,localDofVector.container() );
309 template<
class LocalFunction,
class DiscreteFunction,
template<
class >
class Assembly >
310 void operator() (
const LocalFunction &localFunction, LocalContribution< DiscreteFunction, Assembly > &localContribution )
const
312 (*this)(localFunction,localContribution.localDofVector());
315 BasisFunctionSetType basisFunctionSet ()
const {
return basisFunctionSet_; }
316 const LocalInterpolationType &localInterpolation ()
const
318 assert( localInterpolation_.has_value() );
319 return *localInterpolation_;
323 template<
class LocalFunction,
class Vector>
324 auto operator() (
const LocalFunction &localFunction, Vector &localDofVector, PriorityTag<1> )
const
326 std::declval<LocalInterpolationType>().interpolate(
327 localFunction, localDofVector)) >
329 localInterpolation().interpolate( localFunction, localDofVector );
331 template<
class LocalFunction,
class Vector>
332 void operator() (
const LocalFunction &localFunction, Vector &localDofVector, PriorityTag<0> )
const
334 std::vector<typename Vector::value_type> tmp(basisFunctionSet_.size());
335 localInterpolation().interpolate( localFunction, tmp);
336 for (
unsigned int i=0;i<tmp.size();++i)
337 localDofVector[i] = tmp[i];
339 BasisFunctionSetType basisFunctionSet_;
340 std::optional< LocalInterpolationType > localInterpolation_;
341 DofAlignmentType dofAlignment_;
344 template <
class DiscreteFunctionSpace >
345 class LocalFEInterpolationWrapper
346 :
public LocalInterpolationWrapper< DiscreteFunctionSpace >
348 typedef LocalInterpolationWrapper< DiscreteFunctionSpace > BaseType;
350 using BaseType :: interpolation_;
351 typedef std::vector< typename DiscreteFunctionSpace::RangeFieldType >
352 TemporarayDofVectorType;
355 using BaseType::interpolation;
361 using BaseType :: bind;
362 using BaseType :: unbind;
364 template<
class LocalFunction,
class Dof >
365 void operator() (
const LocalFunction &localFunction, std::vector< Dof > &localDofVector )
const
368 interpolation()( localFunction, localDofVector );
371 template<
class LocalFunction,
class Dof,
class Allocator >
375 (*this)(localFunction, localDofVector.container() );
379 template<
class LocalFunction,
class DiscreteFunction,
template<
class >
class Assembly >
380 void operator() (
const LocalFunction &localFunction, LocalContribution< DiscreteFunction, Assembly > &localContribution )
const
383 (*this)(localFunction, localContribution.localDofVector());
387 template<
class LocalFunction,
class LocalDofVector >
388 void operator () (
const LocalFunction &localFunction, LocalDofVector &dofs )
const
390 const int size = dofs.size();
391 TemporarayDofVectorType& tmpDofs = *tmpDofs_;
392 tmpDofs.resize(
size );
393 (*this)(localFunction, tmpDofs );
396 for(
int i=0; i<
size; ++i )
398 dofs[ i ] = tmpDofs[ i ];
Construct a vector with a dynamic size.
Definition: dynvector.hh:59
GridImp::template Codim< cd >::Geometry Geometry
The corresponding geometry type.
Definition: entity.hh:100
FunctionSpaceTraits::RangeType RangeType
Type of range vector (using type of range field) has a Dune::FieldVector type interface.
Definition: functionspaceinterface.hh:71
@ dimRange
dimension of range vector space
Definition: functionspaceinterface.hh:48
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
FunctionSpaceType::DomainType DomainType
type of domain vectors, i.e., type of coordinates
Definition: localfunction.hh:105
FunctionSpaceType::RangeType RangeType
type of range vectors, i.e., type of function values
Definition: localfunction.hh:107
Implements a vector constructed from a given type representing a field and a compile-time given size.
typename Impl::voider< Types... >::type void_t
Is void for all valid input types. The workhorse for C++11 SFINAE-techniques.
Definition: typetraits.hh:40
Dune namespace.
Definition: alignedallocator.hh:13
constexpr std::integral_constant< std::size_t, sizeof...(II)> size(std::integer_sequence< T, II... >)
Return the size of the sequence.
Definition: integersequence.hh:75