Dune Core Modules (2.6.0)

brezzidouglasmarini1cube2dlocalbasis.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_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH
4#define DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH
5
6#include <array>
7#include <bitset>
8#include <numeric>
9#include <vector>
10
12
13#include "../../common/localbasis.hh"
14
15namespace Dune
16{
26 template<class D, class R>
28 {
29
30 public:
33
36 {
37 for (size_t i=0; i<4; i++)
38 sign_[i] = 1.0;
39 }
40
46 BDM1Cube2DLocalBasis (std::bitset<4> s)
47 {
48 for (size_t i=0; i<4; i++)
49 sign_[i] = s[i] ? -1.0 : 1.0;
50 }
51
53 unsigned int size () const
54 {
55 return 8;
56 }
57
64 inline void evaluateFunction (const typename Traits::DomainType& in,
65 std::vector<typename Traits::RangeType>& out) const
66 {
67 out.resize(8);
68
69 out[0][0] = sign_[0]*(in[0] - 1.0);
70 out[0][1] = 0.0;
71 out[1][0] = 6.0*in[0]*in[1] - 3.0*in[0]-6*in[1] + 3.0;
72 out[1][1] = -3.0*in[1]*in[1] + 3.0*in[1];
73 out[2][0] = sign_[1]*(in[0]);
74 out[2][1] = 0.0;
75 out[3][0] = -6.0*in[0]*in[1] + 3.0*in[0];
76 out[3][1] = 3.0*in[1]*in[1] - 3.0*in[1];
77 out[4][0] = 0.0;
78 out[4][1] = sign_[2]*(in[1] - 1.0);
79 out[5][0] = 3.0*in[0]*in[0] - 3.0*in[0];
80 out[5][1] = -6.0*in[0]*in[1] + 6.0*in[0] + 3.0*in[1] - 3.0;
81 out[6][0] = 0.0;
82 out[6][1] = sign_[3]*(in[1]);
83 out[7][0] = -3.0*in[0]*in[0] + 3.0*in[0];
84 out[7][1] = 6.0*in[0]*in[1] - 3.0*in[1];
85 }
86
93 inline void evaluateJacobian (const typename Traits::DomainType& in,
94 std::vector<typename Traits::JacobianType>& out) const
95 {
96 out.resize(8);
97
98 out[0][0][0] = sign_[0];
99 out[0][0][1] = 0.0;
100 out[0][1][0] = 0.0;
101 out[0][1][1] = 0.0;
102
103 out[1][0][0] = 6.0*in[1] - 3.0;
104 out[1][0][1] = 6.0*in[0] - 6.0;
105 out[1][1][0] = 0.0;
106 out[1][1][1] = -6.0*in[1] + 3.0;
107
108 out[2][0][0] = sign_[1];
109 out[2][0][1] = 0.0;
110 out[2][1][0] = 0.0;
111 out[2][1][1] = 0.0;
112
113 out[3][0][0] = -6.0*in[1] + 3.0;
114 out[3][0][1] = -6.0*in[0];
115 out[3][1][0] = 0.0;
116 out[3][1][1] = 6.0*in[1] - 3.0;
117
118 out[4][0][0] = 0.0;
119 out[4][0][1] = 0.0;
120 out[4][1][0] = 0.0;
121 out[4][1][1] = sign_[2];
122
123 out[5][0][0] = 6.0*in[0] - 3.0;
124 out[5][0][1] = 0.0;
125 out[5][1][0] = -6.0*in[1] + 6.0;
126 out[5][1][1] = -6.0*in[0] + 3.0;
127
128 out[6][0][0] = 0.0;
129 out[6][0][1] = 0.0;
130 out[6][1][0] = 0.0;
131 out[6][1][1] = sign_[3];
132
133 out[7][0][0] = -6.0*in[0] + 3.0;
134 out[7][0][1] = 0.0;
135 out[7][1][0] = 6.0*in[1];
136 out[7][1][1] = 6.0*in[0] - 3.0;
137 }
138
140 void partial (const std::array<unsigned int, 2>& order,
141 const typename Traits::DomainType& in, // position
142 std::vector<typename Traits::RangeType>& out) const // return value
143 {
144 auto totalOrder = std::accumulate(order.begin(), order.end(), 0);
145 if (totalOrder == 0) {
146 evaluateFunction(in, out);
147 } else if (totalOrder == 1) {
148 out.resize(size());
149 auto const direction = std::distance(order.begin(), std::find(order.begin(), order.end(), 1));
150
151 switch (direction) {
152 case 0:
153 out[0][0] = sign_[0];
154 out[0][1] = 0.0;
155
156 out[1][0] = 6.0*in[1] - 3.0;
157 out[1][1] = 0.0;
158
159 out[2][0] = sign_[1];
160 out[2][1] = 0.0;
161
162 out[3][0] = -6.0*in[1] + 3.0;
163 out[3][1] = 0.0;
164
165 out[4][0] = 0.0;
166 out[4][1] = 0.0;
167
168 out[5][0] = 6.0*in[0] - 3.0;
169 out[5][1] = -6.0*in[1] + 6.0;
170
171 out[6][0] = 0.0;
172 out[6][1] = 0.0;
173
174 out[7][0] = -6.0*in[0] + 3.0;
175 out[7][1] = 6.0*in[1];
176 break;
177 case 1:
178 out[0][0] = 0.0;
179 out[0][1] = 0.0;
180
181 out[1][0] = 6.0*in[0] - 6.0;
182 out[1][1] = -6.0*in[1] + 3.0;
183
184 out[2][0] = 0.0;
185 out[2][1] = 0.0;
186
187 out[3][0] = -6.0*in[0];
188 out[3][1] = 6.0*in[1] - 3.0;
189
190 out[4][0] = 0.0;
191 out[4][1] = sign_[2];
192
193 out[5][0] = 0.0;
194 out[5][1] = -6.0*in[0] + 3.0;
195
196 out[6][0] = 0.0;
197 out[6][1] = sign_[3];
198
199 out[7][0] = 0.0;
200 out[7][1] = 6.0*in[0] - 3.0;
201 break;
202 default:
203 DUNE_THROW(RangeError, "Component out of range.");
204 }
205 } else {
206 DUNE_THROW(NotImplemented, "Desired derivative order is not implemented");
207 }
208 }
209
211 unsigned int order () const
212 {
213 return 2;
214 }
215
216 private:
217 std::array<R,4> sign_;
218 };
219}
220#endif // DUNE_LOCALFUNCTIONS_BREZZIDOUGLASMARINI1_CUBE2D_LOCALBASIS_HH
First order Brezzi-Douglas-Marini shape functions on the reference quadrilateral.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:28
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:93
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: brezzidouglasmarini1cube2dlocalbasis.hh:140
BDM1Cube2DLocalBasis()
Standard constructor.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:35
unsigned int size() const
number of shape functions
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:53
BDM1Cube2DLocalBasis(std::bitset< 4 > s)
Make set number s, where 0 <= s < 16.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:46
unsigned int order() const
Polynomial order of the shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:211
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: brezzidouglasmarini1cube2dlocalbasis.hh:64
A dense n x m matrix.
Definition: fmatrix.hh:68
vector space out of a tensor product of fields.
Definition: fvector.hh:93
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 ...
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
T accumulate(Range &&range, T value, F &&f)
Accumulate values.
Definition: hybridutilities.hh:331
Dune namespace.
Definition: alignedallocator.hh:10
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 (Nov 24, 23:30, 2024)