DUNE-FEM (unstable)

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// 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_LOCALFUNCTIONS_COMMON_VIRTUALINTERFACE_HH
6#define DUNE_LOCALFUNCTIONS_COMMON_VIRTUALINTERFACE_HH
7
8#include <type_traits>
9#include <array>
10#include <vector>
11#include <functional>
12
13#include <dune/geometry/type.hh>
14
15#include <dune/localfunctions/common/localbasis.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
28 // -----------------------------------------------------------------
29 // Basis
30 // -----------------------------------------------------------------
31
38 template<class T>
40 {
41 public:
42 using Traits = T;
43
44
46
48 virtual unsigned int size () const = 0;
49
51 virtual unsigned int order () const = 0;
52
58 virtual void evaluateFunction (const typename Traits::DomainType& in,
59 std::vector<typename Traits::RangeType>& out) const = 0;
60
69 virtual void evaluateJacobian(const typename Traits::DomainType& in, // position
70 std::vector<typename Traits::JacobianType>& out) const = 0;
71
77 virtual void partial(const std::array<unsigned int,Traits::dimDomain>& order,
78 const typename Traits::DomainType& in,
79 std::vector<typename Traits::RangeType>& out) const = 0;
80 };
81
82
83
84 // -----------------------------------------------------------------
85 // Interpolation
86 // -----------------------------------------------------------------
87
100 template<class DomainType, class RangeType>
102 {
103 public:
104
106 using FunctionType = std::function<RangeType(DomainType)>;
107
109 typedef typename RangeType::field_type CoefficientType;
110
112
120 virtual void interpolate (const FunctionType& f, std::vector<CoefficientType>& out) const = 0;
121 };
122
130 template<class DomainType, class RangeType>
132 : public LocalInterpolationVirtualInterfaceBase<DomainType, RangeType>
133 {
134 public:
135
137 using FunctionType = std::function<RangeType(DomainType)>;
138
140 typedef typename RangeType::field_type CoefficientType;
141
142
144
145 // This method is only noted again for to make the documentation complete.
146
154 virtual void interpolate (const FunctionType& f, std::vector<CoefficientType>& out) const = 0;
155
161 template<class F,
162 std::enable_if_t<not std::is_base_of<FunctionType, F>::value, int> = 0>
163 void interpolate (const F& f, std::vector<CoefficientType>& out) const
164 {
166 asBase.interpolate(FunctionType(std::cref(f)),out);
167 }
168
174 template<class F, class C>
175 void interpolate (const F& f, std::vector<C>& out) const
176 {
177 std::vector<CoefficientType> outDummy;
179 asBase.interpolate(FunctionType(std::cref(f)),outDummy);
180 out.resize(outDummy.size());
181 for(typename std::vector<CoefficientType>::size_type i=0; i<outDummy.size(); ++i)
182 out[i] = outDummy[i];
183 }
184 };
185
186
187
188 // -----------------------------------------------------------------
189 // Coefficients
190 // -----------------------------------------------------------------
191
198 {
199 public:
200
202
204 virtual std::size_t size () const = 0;
205
207 const virtual LocalKey& localKey (std::size_t i) const = 0;
208
209 };
210
211
212
213 // -----------------------------------------------------------------
214 // Finite Element
215 // -----------------------------------------------------------------
216
217
223 template<class T>
225 {
226 using LocalBasisTraits = T;
227 public:
234
236
238 virtual const typename Traits::LocalBasisType& localBasis () const = 0;
239
241 virtual const typename Traits::LocalCoefficientsType& localCoefficients () const = 0;
242
244 virtual const typename Traits::LocalInterpolationType& localInterpolation () const = 0;
245
247 virtual unsigned int size () const = 0;
248
250 virtual const GeometryType type () const = 0;
251
252 virtual LocalFiniteElementVirtualInterface<T>* clone() const = 0;
253 };
254}
255#endif
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:114
virtual base class for a local basis
Definition: virtualinterface.hh:40
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:198
virtual std::size_t size() const =0
number of coefficients
virtual const LocalKey & localKey(std::size_t i) const =0
get i'th index
virtual base class for local finite elements with functions
Definition: virtualinterface.hh:225
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:102
std::function< RangeType(DomainType)> FunctionType
type of function to interpolate
Definition: virtualinterface.hh:106
RangeType::field_type CoefficientType
type of the coefficient vector in the interpolate method
Definition: virtualinterface.hh:109
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:133
std::function< RangeType(DomainType)> FunctionType
type of function to interpolate
Definition: virtualinterface.hh:137
void interpolate(const F &f, std::vector< C > &out) const
determine coefficients interpolating a given function
Definition: virtualinterface.hh:175
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
void interpolate(const F &f, std::vector< CoefficientType > &out) const
determine coefficients interpolating a given function
Definition: virtualinterface.hh:163
Describe position of one degree of freedom.
Definition: localkey.hh:24
Dune namespace.
Definition: alignedallocator.hh:13
D DomainType
domain type
Definition: localbasis.hh:43
R RangeType
range type
Definition: localbasis.hh:52
traits helper struct
Definition: localfiniteelementtraits.hh:13
LB LocalBasisType
Definition: localfiniteelementtraits.hh:16
LC LocalCoefficientsType
Definition: localfiniteelementtraits.hh:20
LI LocalInterpolationType
Definition: localfiniteelementtraits.hh:24
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 (Nov 12, 23:30, 2024)