Dune Core Modules (2.9.1)

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 (C) 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/localinterpolation.hh>
17#include <dune/localfunctions/common/localkey.hh>
18#include <dune/localfunctions/common/localfiniteelementtraits.hh>
19
20namespace Dune
21{
22
23 // forward declaration needed by the helper traits
24 template<class DomainType, class RangeType>
25 class LocalInterpolationVirtualInterface;
26
27 // -----------------------------------------------------------------
28 // Helper traits classes
29 // -----------------------------------------------------------------
30
43 template<class FE>
44 class
45 [[deprecated("Dune::LocalFiniteElementFunctionBase is deprecated after Dune 2.7. You can now pass functions providing operator() to interpolate.")]]
47 {
48 typedef typename FE::Traits::LocalBasisType::Traits::DomainType Domain;
49 typedef typename FE::Traits::LocalBasisType::Traits::RangeType Range;
50
51 // Hack: Keep a copy of Dune::Function here. This allows to avoid depending
52 // on the deprecated dune-common header while still keeping the LocalFiniteElementFunctionBase
53 // mechanism working during its deprecation period.
54 class FunctionBaseDummy
55 {
56 public:
57
58 using RangeType = Range;
59 using DomainType = Domain;
60
61 struct Traits
62 {
63 using RangeType = Range;
64 using DomainType = Domain;
65 };
66
67 void evaluate(const DomainType& x, RangeType& y) const;
68 };
69
70 public:
71
72 using VirtualFunctionBase = FunctionBaseDummy;
73 using FunctionBase = FunctionBaseDummy;
74
80 using type = FunctionBaseDummy;
81 };
82
83
84
85 // -----------------------------------------------------------------
86 // Basis
87 // -----------------------------------------------------------------
88
95 template<class T>
97 {
98 public:
99 using Traits = T;
100
101
102 virtual ~LocalBasisVirtualInterface() {}
103
105 virtual unsigned int size () const = 0;
106
108 virtual unsigned int order () const = 0;
109
115 virtual void evaluateFunction (const typename Traits::DomainType& in,
116 std::vector<typename Traits::RangeType>& out) const = 0;
117
126 virtual void evaluateJacobian(const typename Traits::DomainType& in, // position
127 std::vector<typename Traits::JacobianType>& out) const = 0;
128
134 virtual void partial(const std::array<unsigned int,Traits::dimDomain>& order,
135 const typename Traits::DomainType& in,
136 std::vector<typename Traits::RangeType>& out) const = 0;
137 };
138
139
140
141 // -----------------------------------------------------------------
142 // Interpolation
143 // -----------------------------------------------------------------
144
157 template<class DomainType, class RangeType>
159 {
160 public:
161
163 using FunctionType = std::function<RangeType(DomainType)>;
164
166 typedef typename RangeType::field_type CoefficientType;
167
169
177 virtual void interpolate (const FunctionType& f, std::vector<CoefficientType>& out) const = 0;
178 };
179
187 template<class DomainType, class RangeType>
189 : public LocalInterpolationVirtualInterfaceBase<DomainType, RangeType>
190 {
191 public:
192
194 using FunctionType = std::function<RangeType(DomainType)>;
195
197 typedef typename RangeType::field_type CoefficientType;
198
199
201
202 // This method is only noted again for to make the documentation complete.
203
211 virtual void interpolate (const FunctionType& f, std::vector<CoefficientType>& out) const = 0;
212
218 template<class F,
219 std::enable_if_t<not std::is_base_of<FunctionType, F>::value, int> = 0>
220 void interpolate (const F& ff, std::vector<CoefficientType>& out) const
221 {
222 const auto& f = Impl::makeFunctionWithCallOperator<DomainType>(ff);
223
225 asBase.interpolate(FunctionType(std::cref(f)),out);
226 }
227
233 template<class F, class C>
234 void interpolate (const F& ff, std::vector<C>& out) const
235 {
236 const auto& f = Impl::makeFunctionWithCallOperator<DomainType>(ff);
237
238 std::vector<CoefficientType> outDummy;
240 asBase.interpolate(FunctionType(std::cref(f)),outDummy);
241 out.resize(outDummy.size());
242 for(typename std::vector<CoefficientType>::size_type i=0; i<outDummy.size(); ++i)
243 out[i] = outDummy[i];
244 }
245 };
246
247
248
249 // -----------------------------------------------------------------
250 // Coefficients
251 // -----------------------------------------------------------------
252
259 {
260 public:
261
263
265 virtual std::size_t size () const = 0;
266
268 const virtual LocalKey& localKey (std::size_t i) const = 0;
269
270 };
271
272
273
274 // -----------------------------------------------------------------
275 // Finite Element
276 // -----------------------------------------------------------------
277
278
284 template<class T>
286 {
287 using LocalBasisTraits = T;
288 public:
295
297
299 virtual const typename Traits::LocalBasisType& localBasis () const = 0;
300
302 virtual const typename Traits::LocalCoefficientsType& localCoefficients () const = 0;
303
305 virtual const typename Traits::LocalInterpolationType& localInterpolation () const = 0;
306
308 virtual unsigned int size () const = 0;
309
311 virtual const GeometryType type () const = 0;
312
313 virtual LocalFiniteElementVirtualInterface<T>* clone() const = 0;
314 };
315}
316#endif
Unique label for each type of entities that can occur in DUNE grids.
Definition: type.hh:126
virtual base class for a local basis
Definition: virtualinterface.hh:97
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:259
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:47
FunctionBaseDummy type
Base class type for functions to use with LocalInterpolation.
Definition: virtualinterface.hh:80
virtual base class for local finite elements with functions
Definition: virtualinterface.hh:286
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:159
std::function< RangeType(DomainType)> FunctionType
type of function to interpolate
Definition: virtualinterface.hh:163
RangeType::field_type CoefficientType
type of the coefficient vector in the interpolate method
Definition: virtualinterface.hh:166
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:190
std::function< RangeType(DomainType)> FunctionType
type of function to interpolate
Definition: virtualinterface.hh:194
RangeType::field_type CoefficientType
type of the coefficient vector in the interpolate method
Definition: virtualinterface.hh:197
void interpolate(const F &ff, std::vector< CoefficientType > &out) const
determine coefficients interpolating a given function
Definition: virtualinterface.hh:220
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:234
Describe position of one degree of freedom.
Definition: localkey.hh:23
Dune namespace.
Definition: alignedallocator.hh:13
D DomainType
domain type
Definition: localbasis.hh:42
R RangeType
range type
Definition: localbasis.hh:51
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 (Jul 15, 22:36, 2024)