Dune Core Modules (unstable)

interpolation.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3// SPDX-FileCopyrightInfo: Copyright © DUNE Project contributors, see file LICENSE.md in module root
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
5#ifndef DUNE_LAGRANGEBASIS_INTERPOLATION_HH
6#define DUNE_LAGRANGEBASIS_INTERPOLATION_HH
7
8#include <type_traits>
9#include <utility>
10#include <vector>
11
13
14#include <dune/localfunctions/lagrange/lagrangecoefficients.hh>
15
16namespace Dune
17{
18
19 template< template <class,unsigned int> class LP,
20 unsigned int dim, class F >
21 struct LagrangeInterpolationFactory;
22
23 // LocalLagrangeInterpolation
24 // --------------------------
25
26 template< template <class,unsigned int> class LP, unsigned int dim, class F >
27 class LocalLagrangeInterpolation
28 {
29 typedef LocalLagrangeInterpolation< LP,dim,F > This;
30
31 public:
32 typedef LP<F,dim> LagrangePointSet;
33 typedef typename LagrangePointSet::Field Field;
34
35 static const unsigned int dimension = LagrangePointSet::dimension;
36
37 private:
38 friend struct LagrangeInterpolationFactory<LP,dim,F>;
39 const LagrangePointSet &lagrangePoints_;
40
41 explicit LocalLagrangeInterpolation ( const LagrangePointSet &lagrangePoints )
42 : lagrangePoints_( lagrangePoints )
43 {}
44
45 const LagrangePointSet *points () const { return &lagrangePoints_; }
46
47 template< class Fn, class Vector >
48 auto interpolate ( const Fn &fn, Vector &coefficients, PriorityTag< 1 > ) const
49 -> std::enable_if_t< std::is_invocable_v< const Fn &, decltype( this->lagrangePoints_.begin()->point() ) > >
50 {
51 unsigned int index = 0;
52 for( const auto &lp : lagrangePoints_ )
53 field_cast( fn( lp.point() ), coefficients[ index++ ] );
54 }
55
56 public:
57 template< class Fn, class Vector,
58 decltype(std::declval<Vector>().size(),bool{}) = true,
59 decltype(std::declval<Vector>().resize(0u),bool{}) = true>
60 void interpolate ( const Fn &fn, Vector &coefficients ) const
61 {
62 coefficients.resize( lagrangePoints_.size() );
63 interpolate( fn, coefficients, PriorityTag< 42 >() );
64 }
65
66 template< class Basis, class Matrix,
67 decltype(std::declval<Matrix>().rows(),bool{}) = true,
68 decltype(std::declval<Matrix>().cols(),bool{}) = true,
69 decltype(std::declval<Matrix>().resize(0u,0u),bool{}) = true>
70 void interpolate ( const Basis &basis, Matrix &coefficients ) const
71 {
72 coefficients.resize( lagrangePoints_.size(), basis.size( ) );
73
74 unsigned int index = 0;
75 for( const auto &lp : lagrangePoints_ )
76 basis.template evaluate< 0 >( lp.point(), coefficients[index++] );
77 }
78
79 const LagrangePointSet &lagrangePoints () const { return lagrangePoints_; }
80 };
81
82
83
84 // LocalLagrangeInterpolationFactory
85 // ---------------------------------
86 template< template <class,unsigned int> class LP,
87 unsigned int dim, class F >
88 struct LagrangeInterpolationFactory
89 {
90 typedef LagrangeCoefficientsFactory<LP,dim,F> LagrangePointSetFactory;
91 typedef typename LagrangePointSetFactory::Object LagrangePointSet;
92
93 typedef typename LagrangePointSetFactory::Key Key;
94 typedef const LocalLagrangeInterpolation< LP,dim,F > Object;
95
96 template< GeometryType::Id geometryId >
97 static Object *create ( const Key &key )
98 {
99 const LagrangePointSet *lagrangeCoeff
100 = LagrangePointSetFactory::template create< geometryId >( key );
101 if ( lagrangeCoeff == 0 )
102 return 0;
103 else
104 return new Object( *lagrangeCoeff );
105 }
106 template< GeometryType::Id geometryId >
107 static bool supports ( const Key &key )
108 {
109 return true;
110 }
111 static void release( Object *object)
112 {
113 LagrangePointSetFactory::release( object->points() );
114 delete object;
115 }
116 };
117
118}
119
120#endif // #ifndef DUNE_LAGRANGEBASIS_INTERPOLATION_HH
Dune namespace.
Definition: alignedallocator.hh:13
void field_cast(const F1 &f1, F2 &f2)
a helper class to cast from one field to another
Definition: field.hh:159
Utilities for type computations, constraining overloads, ...
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 22, 23:30, 2024)