Dune Core Modules (2.7.0)

raviartthomas0cube2dall.hh
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
11
12#include <dune/localfunctions/common/localbasis.hh>
13#include <dune/localfunctions/common/localkey.hh>
14#include <dune/localfunctions/common/localinterpolation.hh>
15
16namespace Dune
17{
26 template<class D, class R>
28 {
29 public:
32
35 {
36 std::fill(sign_.begin(), sign_.end(), 1.0);
37 }
38
40 RT0Cube2DLocalBasis (std::bitset<4> s)
41 {
42 for (int i=0; i<4; i++)
43 sign_[i] = s[i] ? -1.0 : 1.0;
44 }
45
47 unsigned int size () const
48 {
49 return 4;
50 }
51
53 inline void evaluateFunction (const typename Traits::DomainType& in,
54 std::vector<typename Traits::RangeType>& out) const
55 {
56 out.resize(4);
57 out[0] = {sign_[0]*(in[0]-1.0), 0.0};
58 out[1] = {sign_[1]*(in[0]), 0.0};
59 out[2] = {0.0, sign_[2]*(in[1]-1.0)};
60 out[3] = {0.0, sign_[3]*(in[1])};
61 }
62
64 inline void
65 evaluateJacobian (const typename Traits::DomainType& in, // position
66 std::vector<typename Traits::JacobianType>& out) const // return value
67 {
68 out.resize(4);
69 out[0][0] = {sign_[0], 0};
70 out[0][1] = {0, 0};
71
72 out[1][0] = {sign_[1], 0};
73 out[1][1] = {0, 0};
74
75 out[2][0] = {0, 0};
76 out[2][1] = {0, sign_[2]};
77
78 out[3][0] = {0, 0};
79 out[3][1] = {0, sign_[3]};
80 }
81
83 void partial (const std::array<unsigned int, 2>& order,
84 const typename Traits::DomainType& in, // position
85 std::vector<typename Traits::RangeType>& out) const // return value
86 {
87 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
88 if (totalOrder == 0) {
89 evaluateFunction(in, out);
90 } else if (totalOrder == 1) {
91 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
92 out.resize(size());
93
94 for (std::size_t i = 0; i < size(); ++i)
95 out[i] = {0, 0};
96
97 switch (direction) {
98 case 0:
99 out[0][0] = sign_[0];
100 out[1][0] = sign_[1];
101 break;
102 case 1:
103 out[2][1] = sign_[2];
104 out[3][1] = sign_[3];
105 break;
106 default:
107 DUNE_THROW(RangeError, "Component out of range.");
108 }
109 } else {
110 out.resize(size());
111 for (std::size_t i = 0; i < size(); ++i)
112 for (std::size_t j = 0; j < 2; ++j)
113 out[i][j] = 0;
114 }
115
116 }
117
119 unsigned int order () const
120 {
121 return 1;
122 }
123
124 private:
125 std::array<R,4> sign_;
126 };
127
128
136 template<class LB>
138 {
139 public:
140
142 RT0Cube2DLocalInterpolation (std::bitset<4> s = 0)
143 {
144 for (int i=0; i<4; i++)
145 sign_[i] = s[i] ? -1.0 : 1.0;
146
147 m0 = {0.0, 0.5};
148 m1 = {1.0, 0.5};
149 m2 = {0.5, 0.0};
150 m3 = {0.5, 1.0};
151
152 n0 = {-1.0, 0.0};
153 n1 = { 1.0, 0.0};
154 n2 = { 0.0, -1.0};
155 n3 = { 0.0, 1.0};
156 }
157
158 template<typename F, typename C>
159 void interpolate (const F& ff, std::vector<C>& out) const
160 {
161 // f gives v*outer normal at a point on the edge!
162 auto&& f = Impl::makeFunctionWithCallOperator<typename LB::Traits::DomainType>(ff);
163
164 out.resize(4);
165
166 // Evaluate the normal components at the edge midpoints
167 auto y = f(m0); out[0] = (y[0]*n0[0]+y[1]*n0[1])*sign_[0];
168 y = f(m1); out[1] = (y[0]*n1[0]+y[1]*n1[1])*sign_[1];
169 y = f(m2); out[2] = (y[0]*n2[0]+y[1]*n2[1])*sign_[2];
170 y = f(m3); out[3] = (y[0]*n3[0]+y[1]*n3[1])*sign_[3];
171 }
172
173 private:
174 std::array<typename LB::Traits::RangeFieldType,4> sign_;
175
176 // The four edge midpoints of the reference quadrilateral
177 typename LB::Traits::DomainType m0,m1,m2,m3;
178
179 // The four edge normals of the reference quadrilateral
180 typename LB::Traits::DomainType n0,n1,n2,n3;
181 };
182
190 {
191 public:
194 {
195 for (std::size_t i=0; i<4; i++)
196 li[i] = LocalKey(i,1,0);
197 }
198
200 std::size_t size () const
201 {
202 return 4;
203 }
204
206 const LocalKey& localKey (std::size_t i) const
207 {
208 return li[i];
209 }
210
211 private:
212 std::vector<LocalKey> li;
213 };
214
215}
216#endif // DUNE_LOCALFUNCTIONS_RAVIARTTHOMAS0_CUBE2D_ALL_HH
A dense n x m matrix.
Definition: fmatrix.hh:69
vector space out of a tensor product of fields.
Definition: fvector.hh:96
Describe position of one degree of freedom.
Definition: localkey.hh:21
Lowest order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas0cube2dall.hh:28
RT0Cube2DLocalBasis(std::bitset< 4 > s)
Constructor with a set of edge orientations.
Definition: raviartthomas0cube2dall.hh:40
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: raviartthomas0cube2dall.hh:53
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:83
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: raviartthomas0cube2dall.hh:65
unsigned int order() const
Polynomial order of the shape functions.
Definition: raviartthomas0cube2dall.hh:119
RT0Cube2DLocalBasis()
Standard constructor.
Definition: raviartthomas0cube2dall.hh:34
unsigned int size() const
number of shape functions
Definition: raviartthomas0cube2dall.hh:47
Layout map for RT0 elements on quadrilaterals.
Definition: raviartthomas0cube2dall.hh:190
RT0Cube2DLocalCoefficients()
Standard constructor.
Definition: raviartthomas0cube2dall.hh:193
std::size_t size() const
number of coefficients
Definition: raviartthomas0cube2dall.hh:200
const LocalKey & localKey(std::size_t i) const
get i'th index
Definition: raviartthomas0cube2dall.hh:206
Lowest order Raviart-Thomas shape functions on the reference quadrilateral.
Definition: raviartthomas0cube2dall.hh:138
RT0Cube2DLocalInterpolation(std::bitset< 4 > s=0)
Constructor with explicitly given edge orientations.
Definition: raviartthomas0cube2dall.hh:142
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 ...
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:290
Dune namespace.
Definition: alignedallocator.hh:14
Type traits for LocalBasisVirtualInterface.
Definition: localbasis.hh:32
D DomainType
domain type
Definition: localbasis.hh:43
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)