3#ifndef DUNE_HIERARCHICAL_SIMPLEX_P2_WITH_ELEMENT_BUBBLE_LOCALBASIS_HH
4#define DUNE_HIERARCHICAL_SIMPLEX_P2_WITH_ELEMENT_BUBBLE_LOCALBASIS_HH
16#include <dune/localfunctions/common/localbasis.hh>
17#include <dune/localfunctions/common/localkey.hh>
18#include <dune/localfunctions/common/localinterpolation.hh>
22 template<
class D,
class R,
int dim>
23 class HierarchicalSimplexP2WithElementBubbleLocalBasis
26 HierarchicalSimplexP2WithElementBubbleLocalBasis()
46 template<
class D,
class R>
47 class HierarchicalSimplexP2WithElementBubbleLocalBasis<D,R,1>
62 std::vector<typename Traits::RangeType>& out)
const
68 out[2] = 1-4*(in[0]-0.5)*(in[0]-0.5);
74 std::vector<typename Traits::JacobianType>& out)
const
80 out[2][0][0] = 4-8*in[0];
84 void partial (
const std::array<unsigned int, 1>& order,
86 std::vector<typename Traits::RangeType>& out)
const
88 auto totalOrder = order[0];
89 if (totalOrder == 0) {
90 evaluateFunction(in, out);
91 }
else if (totalOrder == 1) {
96 }
else if (totalOrder == 2) {
103 out[0] = out[1] = out[2] = 0;
136 template<
class D,
class R>
137 class HierarchicalSimplexP2WithElementBubbleLocalBasis<D,R,2>
152 std::vector<typename Traits::RangeType>& out)
const
156 out[0] = 1 - in[0] - in[1];
157 out[1] = 4*in[0]*(1-in[0]-in[1]);
159 out[3] = 4*in[1]*(1-in[0]-in[1]);
160 out[4] = 4*in[0]*in[1];
162 out[6] = 27*in[0]*in[1]*(1-in[0]-in[1]);
169 std::vector<typename Traits::JacobianType>& out)
const
173 out[0][0][0] = -1; out[0][0][1] = -1;
174 out[1][0][0] = 4-8*in[0]-4*in[1]; out[1][0][1] = -4*in[0];
175 out[2][0][0] = 1; out[2][0][1] = 0;
176 out[3][0][0] = -4*in[1]; out[3][0][1] = 4-4*in[0]-8*in[1];
177 out[4][0][0] = 4*in[1]; out[4][0][1] = 4*in[0];
178 out[5][0][0] = 0; out[5][0][1] = 1;
181 out[6][0][0] = 27 * in[1] * (1 - 2*in[0] - in[1]);
182 out[6][0][1] = 27 * in[0] * (1 - 2*in[1] - in[0]);
187 void partial (
const std::array<unsigned int, 2>& order,
189 std::vector<typename Traits::RangeType>& out)
const
192 if (totalOrder == 0) {
193 evaluateFunction(in, out);
194 }
else if (totalOrder == 1) {
196 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
201 out[1] = 4-8*in[0]-4*in[1];
206 out[6] = 27 * in[1] * (1 - 2*in[0] - in[1]);
212 out[3] = 4-4*in[0]-8*in[1];
215 out[6] = 27 * in[0] * (1 - 2*in[1] - in[0]);
258 template<
class D,
class R>
259 class HierarchicalSimplexP2WithElementBubbleLocalBasis<D,R,3>
274 std::vector<typename Traits::RangeType>& out)
const
278 out[0] = 1 - in[0] - in[1] - in[2];
279 out[1] = 4 * in[0] * (1 - in[0] - in[1] - in[2]);
281 out[3] = 4 * in[1] * (1 - in[0] - in[1] - in[2]);
282 out[4] = 4 * in[0] * in[1];
284 out[6] = 4 * in[2] * (1 - in[0] - in[1] - in[2]);
285 out[7] = 4 * in[0] * in[2];
286 out[8] = 4 * in[1] * in[2];
290 out[10] = 81*in[0]*in[1]*in[2]*(1-in[0]-in[1]-in[2]);
295 std::vector<typename Traits::JacobianType>& out)
const
299 out[0][0][0] = -1; out[0][0][1] = -1; out[0][0][2] = -1;
300 out[1][0][0] = 4-8*in[0]-4*in[1]-4*in[2]; out[1][0][1] = -4*in[0]; out[1][0][2] = -4*in[0];
301 out[2][0][0] = 1; out[2][0][1] = 0; out[2][0][2] = 0;
302 out[3][0][0] = -4*in[1]; out[3][0][1] = 4-4*in[0]-8*in[1]-4*in[2]; out[3][0][2] = -4*in[1];
303 out[4][0][0] = 4*in[1]; out[4][0][1] = 4*in[0]; out[4][0][2] = 0;
304 out[5][0][0] = 0; out[5][0][1] = 1; out[5][0][2] = 0;
305 out[6][0][0] = -4*in[2]; out[6][0][1] = -4*in[2]; out[6][0][2] = 4-4*in[0]-4*in[1]-8*in[2];
306 out[7][0][0] = 4*in[2]; out[7][0][1] = 0; out[7][0][2] = 4*in[0];
307 out[8][0][0] = 0; out[8][0][1] = 4*in[2]; out[8][0][2] = 4*in[1];
308 out[9][0][0] = 0; out[9][0][1] = 0; out[9][0][2] = 1;
310 out[10][0][0] = 81 * in[1] * in[2] * (1 - 2*in[0] - in[1] - in[2]);
311 out[10][0][1] = 81 * in[0] * in[2] * (1 - in[0] - 2*in[1] - in[2]);
312 out[10][0][2] = 81 * in[0] * in[1] * (1 - in[0] - in[1] - 2*in[2]);
316 void partial (
const std::array<unsigned int, 3>& order,
318 std::vector<typename Traits::RangeType>& out)
const
321 if (totalOrder == 0) {
322 evaluateFunction(in, out);
323 }
else if (totalOrder == 1) {
325 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
330 out[1] = 4-8*in[0]-4*in[1]-4*in[2];
339 out[10] = 81 * in[1] * in[2] * (1 - 2*in[0] - in[1] - in[2]);
345 out[3] = 4-4*in[0]-8*in[1]-4*in[2];
352 out[10] = 81 * in[0] * in[2] * (1 - in[0] - 2*in[1] - in[2]);
361 out[6] = 4-4*in[0]-4*in[1]-8*in[2];
365 out[10] = 81 * in[0] * in[1] * (1 - in[0] - in[1] - 2*in[2]);
414 static const int numVertices = dim+1;
417 static const int numEdges = (dim+1)*dim / 2;
422 : li(numVertices+numEdges + 1)
439 return numVertices+numEdges + 1;
449 std::vector<Dune::LocalKey> li;
453 class HierarchicalSimplexP2WithElementBubbleLocalInterpolation
458 template<
typename F,
typename C>
459 void interpolate (
const F& ff, std::vector<C>& out)
const
461 typename LB::Traits::DomainType x;
462 typename LB::Traits::RangeType y;
466 auto&& f = Impl::makeFunctionWithCallOperator<decltype(x)>(ff);
469 x[0] = 0.0; x[1] = 0.0; out[0] = f(x);
470 x[0] = 1.0; x[1] = 0.0; out[2] = f(x);
471 x[0] = 0.0; x[1] = 1.0; out[5] = f(x);
474 x[0] = 0.5; x[1] = 0.0; y = f(x);
475 out[1] = y - out[0]*(1-x[0]) - out[2]*x[0];
477 x[0] = 0.0; x[1] = 0.5; y = f(x);
478 out[3] = y - out[0]*(1-x[1]) - out[5]*x[1];
480 x[0] = 0.5; x[1] = 0.5; y = f(x);
481 out[4] = y - out[2]*x[0] - out[5]*x[1];
484 x[0] = 1.0/3; x[1] = 1.0/3; y = f(x);
487 HierarchicalSimplexP2WithElementBubbleLocalBasis<double,double,2> shapeFunctions;
488 std::vector<typename LB::Traits::RangeType> sfValues;
489 shapeFunctions.evaluateFunction(x, sfValues);
492 for (
int i=0; i<6; i++)
493 out[6] -= out[i]*sfValues[i];
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:95
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: hierarchicalsimplexp2withelementbubble.hh:73
LocalBasisTraits< D, 1, Dune::FieldVector< D, 1 >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, 1 > > Traits
export type traits for function signature
Definition: hierarchicalsimplexp2withelementbubble.hh:52
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: hierarchicalsimplexp2withelementbubble.hh:61
unsigned int order() const
Polynomial order of the shape functions (2, in this case)
Definition: hierarchicalsimplexp2withelementbubble.hh:109
unsigned int size() const
number of shape functions
Definition: hierarchicalsimplexp2withelementbubble.hh:55
void partial(const std::array< unsigned int, 1 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: hierarchicalsimplexp2withelementbubble.hh:84
unsigned int size() const
number of shape functions
Definition: hierarchicalsimplexp2withelementbubble.hh:145
unsigned int order() const
Polynomial order of the shape functions (3 in this case)
Definition: hierarchicalsimplexp2withelementbubble.hh:227
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: hierarchicalsimplexp2withelementbubble.hh:151
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: hierarchicalsimplexp2withelementbubble.hh:168
void partial(const std::array< unsigned int, 2 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: hierarchicalsimplexp2withelementbubble.hh:187
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, 2 > > Traits
export type traits for function signature
Definition: hierarchicalsimplexp2withelementbubble.hh:142
unsigned int size() const
number of shape functions
Definition: hierarchicalsimplexp2withelementbubble.hh:267
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: hierarchicalsimplexp2withelementbubble.hh:273
unsigned int order() const
Polynomial order of the shape functions (4 in this case)
Definition: hierarchicalsimplexp2withelementbubble.hh:377
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: hierarchicalsimplexp2withelementbubble.hh:294
void partial(const std::array< unsigned int, 3 > &order, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate partial derivatives of all shape functions.
Definition: hierarchicalsimplexp2withelementbubble.hh:316
LocalBasisTraits< D, 3, Dune::FieldVector< D, 3 >, R, 1, Dune::FieldVector< R, 1 >, Dune::FieldMatrix< R, 1, 3 > > Traits
export type traits for function signature
Definition: hierarchicalsimplexp2withelementbubble.hh:264
The local finite element needed for the Zou-Kornhuber estimator for Signorini problems.
Definition: hierarchicalsimplexp2withelementbubble.hh:412
size_t size() const
number of coefficients
Definition: hierarchicalsimplexp2withelementbubble.hh:437
const Dune::LocalKey & localKey(size_t i) const
get i'th index
Definition: hierarchicalsimplexp2withelementbubble.hh:443
HierarchicalSimplexP2WithElementBubbleLocalCoefficients()
Standard constructor.
Definition: hierarchicalsimplexp2withelementbubble.hh:421
Describe position of one degree of freedom.
Definition: localkey.hh:21
Default exception for dummy implementations.
Definition: exceptions.hh:261
Default exception class for range errors.
Definition: exceptions.hh:252
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
void interpolate(const F &f, const GFS &gfs, XG &xg)
interpolation from a given grid function
Definition: interpolate.hh:177
constexpr T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:289
Dune namespace.
Definition: alignedallocator.hh:11
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:32
D DomainType
domain type
Definition: localbasis.hh:43