dune-localfunctions 3.0-git
raviartthomas0cube2dall.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#ifndef DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_CUBE2D_ALL_HH
4#define DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_CUBE2D_ALL_HH
5
6#include <cstddef>
7#include <numeric>
8#include <vector>
9
10#include <dune/common/fmatrix.hh>
11
14
15namespace Dune
16{
25 template<class D, class R>
27 {
28 public:
29 typedef LocalBasisTraits<D,2,Dune::FieldVector<D,2>,R,2,Dune::FieldVector<R,2>,
30 Dune::FieldMatrix<R,2,2> > Traits;
31
34 {
35 sign0 = sign1 = sign2 = sign3 = 1.0;
36 }
37
39 RT0Cube2DLocalBasis (unsigned int s)
40 {
41 sign0 = sign1 = sign2 = sign3 = 1.0;
42 if (s&1) sign0 = -1.0;
43 if (s&2) sign1 = -1.0;
44 if (s&4) sign2 = -1.0;
45 if (s&8) sign3 = -1.0;
46 }
47
49 unsigned int size () const
50 {
51 return 4;
52 }
53
55 inline void evaluateFunction (const typename Traits::DomainType& in,
56 std::vector<typename Traits::RangeType>& out) const
57 {
58 out.resize(4);
59 out[0][0] = sign0*(in[0]-1.0); out[0][1]=0.0;
60 out[1][0] = sign1*(in[0]); out[1][1]=0.0;
61 out[2][0] = 0.0; out[2][1]=sign2*(in[1]-1.0);
62 out[3][0] = 0.0; out[3][1]=sign3*(in[1]);
63 }
64
66 inline void
67 evaluateJacobian (const typename Traits::DomainType& in, // position
68 std::vector<typename Traits::JacobianType>& out) const // return value
69 {
70 out.resize(4);
71 out[0][0][0] = sign0; out[0][0][1] = 0;
72 out[0][1][0] = 0; out[0][1][1] = 0;
73
74 out[1][0][0] = sign1; out[1][0][1] = 0;
75 out[1][1][0] = 0; out[1][1][1] = 0;
76
77 out[2][0][0] = 0; out[2][0][1] = 0;
78 out[2][1][0] = 0; out[2][1][1] = sign2;
79
80 out[3][0][0] = 0; out[3][0][1] = 0;
81 out[3][1][0] = 0; out[3][1][1] = sign3;
82 }
83
85 void partial (const std::array<unsigned int, 2>& order,
86 const typename Traits::DomainType& in, // position
87 std::vector<typename Traits::RangeType>& out) const // return value
88 {
89 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
90 if (totalOrder == 0) {
91 evaluateFunction(in, out);
92 } else if (totalOrder == 1) {
93 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
94 out.resize(size());
95
96 for (std::size_t i = 0; i < size(); ++i)
97 out[i][0] = out[i][1] = 0;
98
99 switch (direction) {
100 case 0:
101 out[0][0] = sign0;
102 out[1][0] = sign1;
103 break;
104 case 1:
105 out[2][1] = sign2;
106 out[3][1] = sign3;
107 break;
108 default:
109 DUNE_THROW(RangeError, "Component out of range.");
110 }
111 } else {
112 out.resize(size());
113 for (std::size_t i = 0; i < size(); ++i)
114 for (std::size_t j = 0; j < 2; ++j)
115 out[i][j] = 0;
116 }
117
118 }
119
121 unsigned int order () const
122 {
123 return 1;
124 }
125
126 private:
127 R sign0, sign1, sign2, sign3;
128 };
129
130
138 template<class LB>
140 {
141 public:
142
145 {
146 sign0 = sign1 = sign2 = sign3 = 1.0;
147 }
148
151 {
152 sign0 = sign1 = sign2 = sign3 = 1.0;
153 if (s&1) sign0 *= -1.0;
154 if (s&2) sign1 *= -1.0;
155 if (s&4) sign2 *= -1.0;
156 if (s&8) sign3 *= -1.0;
157
158 m0[0] = 0.0; m0[1] = 0.5;
159 m1[0] = 1.0; m1[1] = 0.5;
160 m2[0] = 0.5; m2[1] = 0.0;
161 m3[0] = 0.5; m3[1] = 1.0;
162
163 n0[0] = -1.0; n0[1] = 0.0;
164 n1[0] = 1.0; n1[1] = 0.0;
165 n2[0] = 0.0; n2[1] = -1.0;
166 n3[0] = 0.0; n3[1] = 1.0;
167 }
168
169 template<typename F, typename C>
170 void interpolate (const F& f, std::vector<C>& out) const
171 {
172 // f gives v*outer normal at a point on the edge!
173 typename F::Traits::RangeType y;
174
175 out.resize(4);
176
177 f.evaluate(m0,y); out[0] = (y[0]*n0[0]+y[1]*n0[1])*sign0;
178 f.evaluate(m1,y); out[1] = (y[0]*n1[0]+y[1]*n1[1])*sign1;
179 f.evaluate(m2,y); out[2] = (y[0]*n2[0]+y[1]*n2[1])*sign2;
180 f.evaluate(m3,y); out[3] = (y[0]*n3[0]+y[1]*n3[1])*sign3;
181 }
182
183 private:
184 typename LB::Traits::RangeFieldType sign0,sign1,sign2,sign3;
185 typename LB::Traits::DomainType m0,m1,m2,m3;
186 typename LB::Traits::DomainType n0,n1,n2,n3;
187 };
188
196 {
197 public:
200 {
201 for (std::size_t i=0; i<4; i++)
202 li[i] = LocalKey(i,1,0);
203 }
204
206 std::size_t size () const
207 {
208 return 4;
209 }
210
212 const LocalKey& localKey (std::size_t i) const
213 {
214 return li[i];
215 }
216
217 private:
218 std::vector<LocalKey> li;
219 };
220
221}
222#endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_CUBE2D_ALL_HH
Definition brezzidouglasmarini1cube2d.hh:14
Type traits for LocalBasisVirtualInterface.
Definition localbasis.hh:38
D DomainType
domain type
Definition localbasis.hh:49
Describe position of one degree of freedom.
Definition localkey.hh:21
Lowest order Raviart-Thomas shape functions on the reference quadrilateral.
Definition raviartthomas0cube2dall.hh:27
LocalBasisTraits< D, 2, Dune::FieldVector< D, 2 >, R, 2, Dune::FieldVector< R, 2 >, Dune::FieldMatrix< R, 2, 2 > > Traits
Definition raviartthomas0cube2dall.hh:30
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition raviartthomas0cube2dall.hh:55
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 raviartthomas0cube2dall.hh:85
RT0Cube2DLocalBasis(unsigned int s)
Make set numer s, where 0<=s<16.
Definition raviartthomas0cube2dall.hh:39
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition raviartthomas0cube2dall.hh:67
unsigned int order() const
Polynomial order of the shape functions.
Definition raviartthomas0cube2dall.hh:121
RT0Cube2DLocalBasis()
Standard constructor.
Definition raviartthomas0cube2dall.hh:33
unsigned int size() const
number of shape functions
Definition raviartthomas0cube2dall.hh:49
Lowest order Raviart-Thomas shape functions on the reference quadrilateral.
Definition raviartthomas0cube2dall.hh:140
RT0Cube2DLocalInterpolation()
Standard constructor.
Definition raviartthomas0cube2dall.hh:144
RT0Cube2DLocalInterpolation(unsigned int s)
Make set numer s, where 0<=s<8.
Definition raviartthomas0cube2dall.hh:150
void interpolate(const F &f, std::vector< C > &out) const
Definition raviartthomas0cube2dall.hh:170
Layout map for RT0 elements on quadrilaterals.
Definition raviartthomas0cube2dall.hh:196
RT0Cube2DLocalCoefficients()
Standard constructor.
Definition raviartthomas0cube2dall.hh:199
std::size_t size() const
number of coefficients
Definition raviartthomas0cube2dall.hh:206
const LocalKey & localKey(std::size_t i) const
get i'th index
Definition raviartthomas0cube2dall.hh:212