Dune Core Modules (2.7.0)

virtualinterface.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_LOCALFUNCTIONS_COMMON_VIRTUALINTERFACE_HH
4#define DUNE_LOCALFUNCTIONS_COMMON_VIRTUALINTERFACE_HH
5
6#include <type_traits>
7#include <array>
8#include <vector>
9
11
12#include <dune/geometry/type.hh>
13
14#include <dune/localfunctions/common/localbasis.hh>
15#include <dune/localfunctions/common/localinterpolation.hh>
16#include <dune/localfunctions/common/localkey.hh>
17#include <dune/localfunctions/common/localfiniteelementtraits.hh>
18
19namespace Dune
20{
21
22 // forward declaration needed by the helper traits
23 template<class DomainType, class RangeType>
24 class LocalInterpolationVirtualInterface;
25
26 // -----------------------------------------------------------------
27 // Helper traits classes
28 // -----------------------------------------------------------------
29
35 template<class FE>
37 {
38 typedef typename FE::Traits::LocalBasisType::Traits::DomainType DomainType;
39 typedef typename FE::Traits::LocalBasisType::Traits::RangeType RangeType;
40
42 typedef typename FE::Traits::LocalInterpolationType Implementation;
43
44 public:
45
48
54 typedef typename std::conditional<std::is_base_of<Interface, Implementation>::value, VirtualFunctionBase, FunctionBase>::type type;
55 };
56
57
58
59 // -----------------------------------------------------------------
60 // Basis
61 // -----------------------------------------------------------------
62
69 template<class T>
71 {
72 public:
73 using Traits = T;
74
75
77
79 virtual unsigned int size () const = 0;
80
82 virtual unsigned int order () const = 0;
83
89 virtual void evaluateFunction (const typename Traits::DomainType& in,
90 std::vector<typename Traits::RangeType>& out) const = 0;
91
100 virtual void evaluateJacobian(const typename Traits::DomainType& in, // position
101 std::vector<typename Traits::JacobianType>& out) const = 0;
102
108 virtual void partial(const std::array<unsigned int,Traits::dimDomain>& order,
109 const typename Traits::DomainType& in,
110 std::vector<typename Traits::RangeType>& out) const = 0;
111 };
112
113
114
115 // -----------------------------------------------------------------
116 // Interpolation
117 // -----------------------------------------------------------------
118
131 template<class DomainType, class RangeType>
133 {
134 public:
135
138
140 typedef typename RangeType::field_type CoefficientType;
141
143
151 virtual void interpolate (const FunctionType& f, std::vector<CoefficientType>& out) const = 0;
152 };
153
161 template<class DomainType, class RangeType>
163 : public LocalInterpolationVirtualInterfaceBase<DomainType, RangeType>
164 {
165 public:
166
169
171 typedef typename RangeType::field_type CoefficientType;
172
173
175
176 // This method is only noted again for to make the documentation complete.
177
185 virtual void interpolate (const FunctionType& f, std::vector<CoefficientType>& out) const = 0;
186
192 template<class F,
193 std::enable_if_t<not std::is_base_of<FunctionType, F>::value, int> = 0>
194 void interpolate (const F& ff, std::vector<CoefficientType>& out) const
195 {
196 const auto& f = Impl::makeFunctionWithCallOperator<DomainType>(ff);
197
199 asBase.interpolate(makeVirtualFunction<DomainType, RangeType>(std::cref(f)),out);
200 }
201
207 template<class F, class C>
208 void interpolate (const F& ff, std::vector<C>& out) const
209 {
210 const auto& f = Impl::makeFunctionWithCallOperator<DomainType>(ff);
211
212 std::vector<CoefficientType> outDummy;
214 asBase.interpolate(makeVirtualFunction<DomainType, RangeType>(std::cref(f)),outDummy);
215 out.resize(outDummy.size());
216 for(typename std::vector<CoefficientType>::size_type i=0; i<outDummy.size(); ++i)
217 out[i] = outDummy[i];
218 }
219 };
220
221
222
223 // -----------------------------------------------------------------
224 // Coefficients
225 // -----------------------------------------------------------------
226
233 {
234 public:
235
237
239 virtual std::size_t size () const = 0;
240
242 const virtual LocalKey& localKey (std::size_t i) const = 0;
243
244 };
245
246
247
248 // -----------------------------------------------------------------
249 // Finite Element
250 // -----------------------------------------------------------------
251
252
258 template<class T>
260 {
261 using LocalBasisTraits = T;
262 public:
269
271
273 virtual const typename Traits::LocalBasisType& localBasis () const = 0;
274
276 virtual const typename Traits::LocalCoefficientsType& localCoefficients () const = 0;
277
279 virtual const typename Traits::LocalInterpolationType& localInterpolation () const = 0;
280
282 virtual unsigned int size () const = 0;
283
285 virtual const GeometryType type () const = 0;
286
287 virtual LocalFiniteElementVirtualInterface<T>* clone() const = 0;
288 };
289}
290#endif
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:279
virtual base class for a local basis
Definition: virtualinterface.hh:71
virtual unsigned int order() const =0
Polynomial order of the shape functions.
virtual void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const =0
Evaluate jacobian of all shape functions at given position.
virtual unsigned int size() const =0
Number of shape functions.
virtual void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const =0
Evaluate all basis function at given position.
virtual void partial(const std::array< unsigned int, Traits::dimDomain > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const =0
Evaluate partial derivatives of any order of all shape functions.
virtual base class for local coefficients
Definition: virtualinterface.hh:233
virtual std::size_t size() const =0
number of coefficients
virtual const LocalKey & localKey(std::size_t i) const =0
get i'th index
Return a proper base class for functions to use with LocalInterpolation.
Definition: virtualinterface.hh:37
std::conditional< std::is_base_of< Interface, Implementation >::value, VirtualFunctionBase, FunctionBase >::type type
Base class type for functions to use with LocalInterpolation.
Definition: virtualinterface.hh:54
virtual base class for local finite elements with functions
Definition: virtualinterface.hh:260
virtual const Traits::LocalInterpolationType & localInterpolation() const =0
virtual unsigned int size() const =0
virtual const Traits::LocalBasisType & localBasis() const =0
virtual const GeometryType type() const =0
virtual const Traits::LocalCoefficientsType & localCoefficients() const =0
virtual base class for a local interpolation
Definition: virtualinterface.hh:133
Dune::VirtualFunction< DomainType, RangeType > FunctionType
type of virtual function to interpolate
Definition: virtualinterface.hh:137
RangeType::field_type CoefficientType
type of the coefficient vector in the interpolate method
Definition: virtualinterface.hh:140
virtual void interpolate(const FunctionType &f, std::vector< CoefficientType > &out) const =0
determine coefficients interpolating a given function
virtual base class for a local interpolation
Definition: virtualinterface.hh:164
RangeType::field_type CoefficientType
type of the coefficient vector in the interpolate method
Definition: virtualinterface.hh:171
void interpolate(const F &ff, std::vector< CoefficientType > &out) const
determine coefficients interpolating a given function
Definition: virtualinterface.hh:194
Dune::VirtualFunction< DomainType, RangeType > FunctionType
type of virtual function to interpolate
Definition: virtualinterface.hh:168
virtual void interpolate(const FunctionType &f, std::vector< CoefficientType > &out) const =0
determine coefficients interpolating a given function
void interpolate(const F &ff, std::vector< C > &out) const
determine coefficients interpolating a given function
Definition: virtualinterface.hh:208
Describe position of one degree of freedom.
Definition: localkey.hh:21
Virtual base class template for function classes.
Definition: function.hh:72
Simple base class templates for functions.
Dune namespace.
Definition: alignedallocator.hh:14
D DomainType
domain type
Definition: localbasis.hh:43
R RangeType
range type
Definition: localbasis.hh:55
traits helper struct
Definition: localfiniteelementtraits.hh:11
LB LocalBasisType
Definition: localfiniteelementtraits.hh:14
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:18
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:22
A unique label for each type of element that can occur in a grid.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)