5#ifndef DUNE_LOCALFUNCTIONS_NEDELEC_NEDELECSIMPLEX_NEDELECSIMPLEXINTERPOLATION_HH
6#define DUNE_LOCALFUNCTIONS_NEDELEC_NEDELECSIMPLEX_NEDELECSIMPLEXINTERPOLATION_HH
15#include <dune/geometry/referenceelements.hh>
18#include <dune/localfunctions/common/localkey.hh>
19#include <dune/localfunctions/utility/interpolationhelper.hh>
20#include <dune/localfunctions/utility/polynomialbasis.hh>
21#include <dune/localfunctions/orthonormal/orthonormalbasis.hh>
29 template <
unsigned int dim,
class Field >
30 struct NedelecL2InterpolationFactory;
37 class LocalCoefficientsContainer
39 typedef LocalCoefficientsContainer This;
42 template <
class Setter>
43 LocalCoefficientsContainer (
const Setter &setter )
45 setter.setLocalKeys(localKey_);
48 const LocalKey &localKey (
const unsigned int i )
const
51 return localKey_[ i ];
54 std::size_t size ()
const
56 return localKey_.size();
60 std::vector< LocalKey > localKey_;
68 template <
unsigned int dim >
69 struct NedelecCoefficientsFactory
71 typedef std::size_t Key;
72 typedef const LocalCoefficientsContainer Object;
74 template< GeometryType::Id geometryId >
75 static Object *create(
const Key &key )
77 typedef NedelecL2InterpolationFactory< dim, double > InterpolationFactory;
78 if( !supports< geometryId >( key ) )
80 typename InterpolationFactory::Object *interpolation = InterpolationFactory::template create< geometryId >( key );
81 Object *localKeys =
new Object( *interpolation );
82 InterpolationFactory::release( interpolation );
86 template< GeometryType::Id geometryId >
87 static bool supports (
const Key &key )
90 return gt.isTriangle() ||
gt.isTetrahedron() ;
92 static void release( Object *
object ) {
delete object; }
109 template<
unsigned int dim,
class Field >
110 struct NedelecL2InterpolationBuilder
112 static const unsigned int dimension = dim;
115 typedef OrthonormalBasisFactory< dimension, Field > TestBasisFactory;
116 typedef typename TestBasisFactory::Object TestBasis;
119 typedef OrthonormalBasisFactory< dimension-1, Field > TestFaceBasisFactory;
120 typedef typename TestFaceBasisFactory::Object TestFaceBasis;
123 typedef OrthonormalBasisFactory< 1, Field > TestEdgeBasisFactory;
124 typedef typename TestEdgeBasisFactory::Object TestEdgeBasis;
131 typedef std::array<FieldVector< Field, dimension >,dim-1> FaceTangents;
133 NedelecL2InterpolationBuilder () =
default;
135 NedelecL2InterpolationBuilder (
const NedelecL2InterpolationBuilder & ) =
delete;
136 NedelecL2InterpolationBuilder ( NedelecL2InterpolationBuilder && ) =
delete;
138 ~NedelecL2InterpolationBuilder ()
140 TestBasisFactory::release( testBasis_ );
141 for( FaceStructure &f : faceStructure_ )
142 TestFaceBasisFactory::release( f.basis_ );
143 for( EdgeStructure& e : edgeStructure_ )
144 TestEdgeBasisFactory::release( e.basis_ );
147 unsigned int topologyId ()
const
149 return geometry_.id();
157 std::size_t order ()
const
163 unsigned int faceSize ()
const
165 return numberOfFaces_;
169 unsigned int edgeSize ()
const
171 return numberOfEdges_;
175 TestBasis *testBasis ()
const
181 TestFaceBasis *testFaceBasis (
unsigned int f )
const
183 assert( f < faceSize() );
184 return faceStructure_[ f ].basis_;
188 TestEdgeBasis *testEdgeBasis (
unsigned int e )
const
190 assert( e < edgeSize() );
191 return edgeStructure_[ e ].basis_;
194 const Tangent& edgeTangent (
unsigned int e )
const
196 assert( e < edgeSize() );
197 return edgeStructure_[ e ].tangent_;
200 const FaceTangents& faceTangents (
unsigned int f )
const
202 assert( f < faceSize() );
203 return faceStructure_[ f ].faceTangents_;
206 const Normal &normal (
unsigned int f )
const
208 assert( f < faceSize() );
209 return faceStructure_[ f ].normal_;
212 template< GeometryType::Id geometryId >
213 void build ( std::size_t order )
217 geometry_ = geometry;
232 int requiredOrder =
static_cast<int>(dimension==3);
233 testBasis_ = (order > requiredOrder ? TestBasisFactory::template create< geometry >( order-1-requiredOrder ) : nullptr);
237 numberOfFaces_ = refElement.size( 1 );
238 faceStructure_.reserve( numberOfFaces_ );
241 for (std::size_t i=0; i<numberOfFaces_; i++)
243 FieldVector<Field,dimension> zero(0);
244 FaceTangents faceTangents;
245 faceTangents.fill(zero);
248 auto vertices = refElement.subEntities(i,1,dim).begin();
249 auto vertex1 = *vertices;
250 for(
int j=1; j<dim;j++)
252 auto vertex2 = vertices[j];
254 faceTangents[j-1] = refElement.position(vertex2,dim) - refElement.position(vertex1,dim);
259 faceTangents[j-1] *=-1;
275 TestFaceBasis *faceBasis = ( dim == 3 && order > 0 ? Impl::IfGeometryType< CreateFaceBasis, dimension-1 >::apply( refElement.type( i, 1 ), order-1 ) : nullptr);
276 faceStructure_.emplace_back( faceBasis, refElement.integrationOuterNormal(i), faceTangents );
278 assert( faceStructure_.size() == numberOfFaces_ );
280 numberOfEdges_ = refElement.size( dim-1 );
281 edgeStructure_.reserve( numberOfEdges_ );
284 for (std::size_t i=0; i<numberOfEdges_; i++)
286 auto vertexIterator = refElement.subEntities(i,dim-1,dim).begin();
287 auto v0 = *vertexIterator;
288 auto v1 = *(++vertexIterator);
294 auto tangent = std::move(refElement.position(v1,dim) - refElement.position(v0,dim));
296 TestEdgeBasis *edgeBasis = Impl::IfGeometryType< CreateEdgeBasis, 1 >::apply( refElement.type( i, dim-1 ), order );
297 edgeStructure_.emplace_back( edgeBasis, tangent );
299 assert( edgeStructure_.size() == numberOfEdges_ );
307 EdgeStructure( TestEdgeBasis *teb,
const Tangent &t )
308 : basis_( teb ), tangent_( t )
311 TestEdgeBasis *basis_;
315 template< GeometryType::Id edgeGeometryId >
316 struct CreateEdgeBasis
318 static TestEdgeBasis *apply ( std::size_t order ) {
return TestEdgeBasisFactory::template create< edgeGeometryId >( order ); }
324 FaceStructure( TestFaceBasis *tfb,
const Normal& normal,
const FaceTangents& faceTangents )
325 : basis_( tfb ), normal_(normal), faceTangents_( faceTangents )
328 TestFaceBasis *basis_;
330 const FaceTangents faceTangents_;
333 template< GeometryType::Id faceGeometryId >
334 struct CreateFaceBasis
336 static TestFaceBasis *apply ( std::size_t order ) {
return TestFaceBasisFactory::template create< faceGeometryId >( order ); }
339 TestBasis *testBasis_ =
nullptr;
340 std::vector< FaceStructure > faceStructure_;
341 unsigned int numberOfFaces_;
342 std::vector< EdgeStructure > edgeStructure_;
343 unsigned int numberOfEdges_;
358 template<
unsigned int dimension,
class F>
360 :
public InterpolationHelper< F ,dimension >
363 typedef InterpolationHelper<F,dimension> Base;
367 typedef NedelecL2InterpolationBuilder<dimension,Field> Builder;
368 typedef typename Builder::FaceTangents FaceTangents;
375 template<
class Function,
class Vector >
376 auto interpolate (
const Function &function, Vector &coefficients )
const
377 -> std::enable_if_t< std::is_same< decltype(std::declval<Vector>().resize(1) ),
void >::value,
void>
379 coefficients.resize(size());
380 typename Base::template Helper<Function,Vector,true> func( function,coefficients );
384 template<
class Basis,
class Matrix >
385 auto interpolate (
const Basis &basis,
Matrix &matrix )
const
386 -> std::enable_if_t< std::is_same<
387 decltype(std::declval<Matrix>().rowPtr(0)),
typename Matrix::Field* >::value,
void>
389 matrix.resize( size(), basis.size() );
390 typename Base::template Helper<Basis,Matrix,false> func( basis,matrix );
394 std::size_t order()
const
398 std::size_t size()
const
403 template <GeometryType::Id geometryId>
404 void build( std::size_t order )
408 builder_.template build<geometryId>(order_);
409 if (builder_.testBasis())
410 size_ += dimension*builder_.testBasis()->size();
412 for (
unsigned int f=0; f<builder_.faceSize(); ++f )
413 if (builder_.testFaceBasis(f))
414 size_ += (dimension-1)*builder_.testFaceBasis(f)->size();
416 for (
unsigned int e=0; e<builder_.edgeSize(); ++e )
417 if (builder_.testEdgeBasis(e))
418 size_ += builder_.testEdgeBasis(e)->size();
421 void setLocalKeys(std::vector< LocalKey > &keys)
const
424 unsigned int row = 0;
425 for (
unsigned int e=0; e<builder_.edgeSize(); ++e)
427 if (builder_.edgeSize())
428 for (
unsigned int i=0; i<builder_.testEdgeBasis(e)->size(); ++i,++row)
429 keys[row] =
LocalKey(e,dimension-1,i);
431 for (
unsigned int f=0; f<builder_.faceSize(); ++f)
433 if (builder_.faceSize())
434 for (
unsigned int i=0; i<builder_.testFaceBasis(f)->size()*(dimension-1); ++i,++row)
438 if (builder_.testBasis())
439 for (
unsigned int i=0; i<builder_.testBasis()->size()*dimension; ++i,++row)
441 assert( row == size() );
445 template<
class Func,
class Container,
bool type >
446 void interpolate (
typename Base::template Helper<Func,Container,type> &func )
const
450 std::vector<Field> testBasisVal;
452 for (
unsigned int i=0; i<size(); ++i)
453 for (
unsigned int j=0; j<func.size(); ++j)
456 unsigned int row = 0;
464 for (
unsigned int e=0; e<builder_.edgeSize(); ++e)
466 if (!builder_.testEdgeBasis(e))
468 testBasisVal.resize(builder_.testEdgeBasis(e)->size());
470 const auto &geometry = refElement.template geometry< dimension-1 >( e );
472 const EdgeQuadrature &edgeQuad = EdgeQuadratureRules::rule( subGeoType, 2*order_+2 );
474 const unsigned int quadratureSize = edgeQuad.size();
475 for(
unsigned int qi = 0; qi < quadratureSize; ++qi )
478 builder_.testEdgeBasis(e)->template evaluate<0>(edgeQuad[qi].position(),testBasisVal);
480 testBasisVal[0] = 1.;
483 func.evaluate( geometry.global( edgeQuad[qi].position() ) ),
484 builder_.edgeTangent(e),
485 edgeQuad[qi].weight(),
489 row += builder_.testEdgeBasis(e)->size();
496 for (
unsigned int f=0; f<builder_.faceSize(); ++f)
498 if (builder_.testFaceBasis(f))
500 testBasisVal.resize(builder_.testFaceBasis(f)->size());
502 const auto &geometry = refElement.template geometry< 1 >( f );
504 const FaceQuadrature &faceQuad = FaceQuadratureRules::rule( subGeoType, 2*order_+2 );
506 const unsigned int quadratureSize = faceQuad.size();
507 for(
unsigned int qi = 0; qi < quadratureSize; ++qi )
510 builder_.testFaceBasis(f)->template evaluate<0>(faceQuad[qi].position(),testBasisVal);
512 testBasisVal[0] = 1.;
514 computeFaceDofs( row,
516 func.evaluate( geometry.global( faceQuad[qi].position() ) ),
517 builder_.faceTangents(f),
519 faceQuad[qi].weight(),
523 row += builder_.testFaceBasis(f)->size()*(dimension-1);
528 if (builder_.testBasis())
530 testBasisVal.resize(builder_.testBasis()->size());
536 const unsigned int quadratureSize = elemQuad.size();
537 for(
unsigned int qi = 0; qi < quadratureSize; ++qi )
539 builder_.testBasis()->template evaluate<0>(elemQuad[qi].position(),testBasisVal);
540 computeInteriorDofs(row,
542 func.evaluate(elemQuad[qi].position()),
543 elemQuad[qi].weight(),
547 row += builder_.testBasis()->size()*dimension;
562 template <
class MVal,
class NedVal,
class Matrix>
563 void computeEdgeDofs (
unsigned int startRow,
565 const NedVal &nedVal,
570 const unsigned int endRow = startRow+mVal.size();
571 typename NedVal::const_iterator nedIter = nedVal.begin();
572 for (
unsigned int col = 0; col < nedVal.size() ; ++nedIter,++col)
574 Field cFactor = (*nedIter)*tangent;
575 typename MVal::const_iterator mIter = mVal.
begin();
576 for (
unsigned int row = startRow; row!=endRow; ++mIter, ++row )
577 matrix.add(row,col, (weight*cFactor)*(*mIter) );
579 assert( mIter == mVal.end() );
593 template <
class MVal,
class NedVal,
class Matrix>
594 void computeFaceDofs (
unsigned int startRow,
596 const NedVal &nedVal,
597 const FaceTangents& faceTangents,
602 const unsigned int endRow = startRow+mVal.size()*(dimension-1);
603 typename NedVal::const_iterator nedIter = nedVal.begin();
604 for (
unsigned int col = 0; col < nedVal.size() ; ++nedIter,++col)
606 auto const& u=*nedIter;
607 auto const& n=normal;
610 u[0]*n[1]-u[1]*n[0]};
611 typename MVal::const_iterator mIter = mVal.
begin();
612 for (
unsigned int row = startRow; row!=endRow; ++mIter)
614 for(
int i=0; i<dimension-1;i++)
616 auto test = *mIter*faceTangents[i];
617 matrix.add(row+i,col, weight*(nedTimesNormal*test) );
622 assert( mIter == mVal.end() );
634 template <
class MVal,
class NedVal,
class Matrix>
635 void computeInteriorDofs (
unsigned int startRow,
637 const NedVal &nedVal,
641 const unsigned int endRow = startRow+mVal.size()*dimension;
642 typename NedVal::const_iterator nedIter = nedVal.begin();
643 for (
unsigned int col = 0; col < nedVal.size() ; ++nedIter,++col)
645 typename MVal::const_iterator mIter = mVal.begin();
646 for (
unsigned int row = startRow; row!=endRow; ++mIter,row+=dimension )
647 for (
unsigned int i=0; i<dimension; ++i)
648 matrix.add(row+i,col, (weight*(*mIter))*(*nedIter)[i] );
650 assert( mIter == mVal.end() );
660 template <
unsigned int dim,
class Field >
661 struct NedelecL2InterpolationFactory
663 typedef NedelecL2InterpolationBuilder<dim,Field> Builder;
665 typedef std::size_t Key;
666 typedef typename std::remove_const<Object>::type NonConstObject;
668 template <GeometryType::Id geometryId>
669 static Object *create(
const Key &key )
671 if ( !supports<geometryId>(key) )
673 NonConstObject *interpol =
new NonConstObject();
674 interpol->template build<geometryId>(key);
678 template <GeometryType::Id geometryId>
679 static bool supports(
const Key &key )
681 GeometryType
gt = geometryId;
682 return gt.isTriangle() ||
gt.isTetrahedron() ;
684 static void release( Object *
object ) {
delete object; }
Iterator begin()
begin iterator
Definition: densevector.hh:347
Base class template for function classes.
Definition: function.hh:41
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:125
Describe position of one degree of freedom.
Definition: localkey.hh:23
A generic dynamic dense matrix.
Definition: matrix.hh:561
An L2-based interpolation for Nedelec.
Definition: nedelecsimplexinterpolation.hh:361
Abstract base class for quadrature rules.
Definition: quadraturerules.hh:154
A container for all quadrature rules of dimension dim
Definition: quadraturerules.hh:200
static const QuadratureRule & rule(const GeometryType &t, int p, QuadratureType::Enum qt=QuadratureType::GaussLegendre)
select the appropriate QuadratureRule for GeometryType t and order p
Definition: quadraturerules.hh:266
A few common exception classes.
GeometryType
Type representing VTK's entity geometry types.
Definition: common.hh:132
bool gt(const T &first, const T &second, typename EpsilonType< T >::Type epsilon)
test if first greater than second
Definition: float_cmp.cc:158
Dune namespace.
Definition: alignedallocator.hh:13
static const ReferenceElement & general(const GeometryType &type)
get general reference elements
Definition: referenceelements.hh:198
A unique label for each type of element that can occur in a grid.