dune-localfunctions 2.10
Loading...
Searching...
No Matches
raviartthomas0pyramidlocalbasis.hh
Go to the documentation of this file.
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_RAVIARTTHOMAS0_PYRAMID_LOCALBASIS_HH
6#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_PYRAMID_LOCALBASIS_HH
7
8#include <numeric>
9#include <vector>
10
11#include <dune/common/fmatrix.hh>
12
13#include "../../common/localbasis.hh"
14
15namespace Dune
16{
26 template<class D, class R>
28 {
29
30 public:
31 typedef LocalBasisTraits<D,3,Dune::FieldVector<D,3>,R,3,Dune::FieldVector<R,3>,
32 Dune::FieldMatrix<R,3,3> > Traits;
33
39 RT0PyramidLocalBasis (std::bitset<5> s = 0)
40 {
41 for (size_t i=0; i<size(); i++)
42 sign[i] = s[i] ? -1.0 : 1.0;
43 }
44
46 unsigned int size () const
47 {
48 return 5;
49 }
50
57 inline void evaluateFunction (const typename Traits::DomainType& in,
58 std::vector<typename Traits::RangeType>& out) const
59 {
60 out.resize(5);
61 for (std::size_t i=0; i<out.size(); i++)
62 out[i] = {0.0,0.0,0.0};
63
64 out[0][0] = 1.5*in[0];
65 out[0][1] = 1.5*in[1];
66 out[0][2] = -1.0;
67
68 out[1][0] = -2.0 + 3.0*in[0];
69
70 out[2][0] = 3.0*in[0];
71
72 out[3][1] = -2.0 + 3.0*in[1];
73
74 out[4][1] = 3.0*in[1];
75
76 for (std::size_t i=0; i<out.size(); i++)
77 out[i] *= sign[i];
78
79 }
80
87 inline void evaluateJacobian (const typename Traits::DomainType& in,
88 std::vector<typename Traits::JacobianType>& out) const
89 {
90 out.resize(5);
91
92 for(int i=0; i<size(); i++)
93 for(int j=0; j<3; j++)
94 out[i][j] = {0.0, 0.0, 0.0};
95
96 out[0][0][0] = sign[0]*(1.5);
97 out[0][1][1] = sign[0]*(1.5);
98
99 out[1][0][0] = sign[1]*(3.0);
100
101 out[2][0][0] = sign[2]*(3.0);
102
103 out[3][1][1] = sign[3]*(3.0);
104
105 out[4][1][1] = sign[4]*(3.0);
106 }
107
109 void partial (const std::array<unsigned int, 3>& order,
110 const typename Traits::DomainType& in, // position
111 std::vector<typename Traits::RangeType>& out) const // return value
112 {
113 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
114 if (totalOrder == 0) {
115 evaluateFunction(in, out);
116 } else {
117 DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
118 }
119 }
120
122 unsigned int order () const
123 {
124 return 1;
125 }
126
127 private:
128 std::array<R,5> sign;
129 };
130}
131#endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_PYRAMID_LOCALBASIS_HH
Definition bdfmcube.hh:18
Type traits for LocalBasisVirtualInterface.
Definition common/localbasis.hh:35
D DomainType
domain type
Definition common/localbasis.hh:43
First order Raviart-Thomas shape functions on the reference pyramid.
Definition raviartthomas0pyramidlocalbasis.hh:28
RT0PyramidLocalBasis(std::bitset< 5 > s=0)
Make set number s, where 0 <= s < 32.
Definition raviartthomas0pyramidlocalbasis.hh:39
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition raviartthomas0pyramidlocalbasis.hh:57
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition raviartthomas0pyramidlocalbasis.hh:87
LocalBasisTraits< D, 3, Dune::FieldVector< D, 3 >, R, 3, Dune::FieldVector< R, 3 >, Dune::FieldMatrix< R, 3, 3 > > Traits
Definition raviartthomas0pyramidlocalbasis.hh:32
unsigned int order() const
Polynomial order of the shape functions.
Definition raviartthomas0pyramidlocalbasis.hh:122
unsigned int size() const
number of shape functions
Definition raviartthomas0pyramidlocalbasis.hh:46
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 raviartthomas0pyramidlocalbasis.hh:109