dune-composites (unstable)

serendipitylocalbasis.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3
4#ifndef DUNE_LOCALFUNCTIONS_SERENDIPITYLOCALBASIS_HH
5#define DUNE_LOCALFUNCTIONS_SERENDIPITYLOCALBASIS_HH
6
7#include <dune/common/fvector.hh>
8#include <dune/common/fmatrix.hh>
9
10#include <dune/geometry/type.hh>
11
12#include <dune/localfunctions/common/localbasis.hh>
13#include <dune/localfunctions/common/localfiniteelementtraits.hh>
14
15
16namespace Dune
17{
28 template<class D, class R, int k, int d>
30 {
31 typedef typename Dune::FieldVector<int,d> DuneVector;
32
33
34 // ith Lagrange polynomial of degree k in one dimension
35 static R p (int i,int deg, D x)
36 {
37 R result(1.0);
38 for (int j=0; j<=deg; j++)
39 if (j!=i) result *= (deg*x-j)/(i-j);
40 return result;
41 }
42
43 // derivative of ith Lagrange polynomial of degree k in one dimension
44 static R dp (int i, int deg, D x)
45 {
46 R result(0.0);
47
48 for (int j=0; j<=deg; j++)
49 if (j!=i)
50 {
51 R prod( (deg*1.0)/(i-j) );
52 for (int l=0; l<=deg; l++)
53 if (l!=i && l!=j)
54 prod *= (deg*x-l)/(i-l);
55 result += prod;
56 }
57 return result;
58 }
59
60 static void multiindex (DuneVector& alpha)
61 {
62 int i = 0;
63 for (int j=0; j<d; j++){
64 i+=alpha[j]*pow(k+1,j);
65 }
66 i++;
67 for (int j=0; j<d; j++)
68 {
69 alpha[j] = i % (k+1);
70 i = i/(k+1);
71 }
72 unsigned int sum = 0;
73 for (unsigned int j=0; j<d; j++){
74 if(alpha[j] > 0 && alpha[j] < k){ sum++; }
75 }
76 if(sum > 1){
77 multiindex(alpha);
78 }
79 }
80
81 bool isEdge(DuneVector& alpha) const{
82 for(int i = 0; i < d; i++){
83 if(alpha[i]>0 && alpha[i] < k)
84 return true;
85 }
86 return false;
87 }
88
89 double edgeVal(int alpha, const double& in, int der) const{
90 //add contribution from alpha
91 if(alpha > 0 && alpha < k){
92 if(der == 0)
93 return p(alpha,k,in);
94 else
95 return dp(alpha,k,in);
96 }
97 else{
98 if(der == 0)
99 return p(alpha/k,k-1,in);
100 else
101 return dp(alpha/k,k-1,in);
102 }
103 }
104
105 double nodeVal(int alpha, const double& in, int der) const{
106 // (k-1)-th order component
107 if(der == 0)
108 return p(alpha/k,k-1,in);
109 else
110 return dp(alpha/k,k-1,in);
111 }
112
113 public:
114 typedef LocalBasisTraits<D,d,Dune::FieldVector<D,d>,R,1,Dune::FieldVector<R,1>,Dune::FieldMatrix<R,1,d> > Traits;
115
117 unsigned int size () const
118 {
119 if(d == 1){
120 return 3;
121 }
122 if(d == 2){
123 return 12;
124 }
125 if(d == 3){
126 return 20;
127 }
128 }
129
131 inline void evaluateFunction (const typename Traits::DomainType& in,
132 std::vector<typename Traits::RangeType>& out) const
133 {
134 DuneVector alpha(0);
135 out.resize(size());
136 for (size_t i=0; i<size(); i++)
137 {
138 // convert index i to multiindex
139 if(i>0){ multiindex(alpha); }
140
141 out[i] = 1.0;
142 for(int j = 0; j < d; j++){
143 if(isEdge(alpha))
144 out[i] *= edgeVal(alpha[j], in[j], 0);
145 else{ //isNode
146 out[i] *= nodeVal(alpha[j], in[j], 0);
147 }
148 }
149 }
150 if(d==2){
151 out[0] -= 0.5*(out[1]+out[3]);
152 out[2] -= 0.5*(out[1]+out[4]);
153 out[5] -= 0.5*(out[3]+out[6]);
154 out[7] -= 0.5*(out[4]+out[6]);
155 }
156 if(d==3){
157 out[0] -= 0.5*(out[1] +out[3] +out[8]);
158 out[2] -= 0.5*(out[1] +out[4] +out[9]);
159 out[5] -= 0.5*(out[3] +out[6] +out[10]);
160 out[7] -= 0.5*(out[4] +out[6] +out[11]);
161 out[12] -= 0.5*(out[8] +out[13]+out[15]);
162 out[14] -= 0.5*(out[9] +out[13]+out[16]);
163 out[17] -= 0.5*(out[10]+out[15]+out[18]);
164 out[19] -= 0.5*(out[11]+out[16]+out[18]);
165 }
166
167 //========================Debug test===========================
168 /*for(size_t h =0; h < 9; h++){
169 Dune::FieldVector<double,d> inh = {(h%3)/2.0,((h/3)%3)/2.0};
170 std::cout << "At point " << inh << std::endl;
171
172 Dune::FieldVector<double,8> outh;
173 DuneVector alpha(0);
174 for (size_t i=0; i<size(); i++)
175 {
176 // convert index i to multiindex
177 if(i>0){ multiindex(alpha); }
178
179 outh[i] = 1.0;
180 for(int j = 0; j < d; j++){
181 if(isEdge(alpha))
182 outh[i] *= edgeVal(alpha[j], inh[j], 0);
183 else{ //isNode
184 outh[i] *= nodeVal(alpha[j], inh[j], 0);
185 }
186 }
187 }
188 if(d==2){
189 outh[0] = outh[0] - 0.5*(outh[1]+outh[3]);
190 outh[2] = outh[2] - 0.5*(outh[1]+outh[4]);
191 outh[5] = outh[5] - 0.5*(outh[3]+outh[6]);
192 outh[7] = outh[7] - 0.5*(outh[4]+outh[6]);
193 }
194
195 for (size_t i=0; i<size(); i++)
196 std::cout << " value " << i << " is " << outh[i] << std::endl;
197 }*/
198 //===================End======================================
199
200 }
201
206 inline void
207 evaluateJacobian (const typename Traits::DomainType& in,
208 std::vector<typename Traits::JacobianType>& out) const
209 {
210 out.resize(size());
211 DuneVector alpha(0);
212
213 // Loop over all shape functions
214 for (size_t i=0; i<size(); i++)
215 {
216 // convert index i to multiindex
217 if(i>0){ multiindex(alpha); }
218
219 out[i] =1.0;
220 // Loop over all coordinate directions
221 for (int j=0; j<d; j++)
222 {
223 // Initialize: the overall expression is a product
224 // if j-th bit of i is set to -1, else 1
225 if(isEdge(alpha))
226 out[i][0][j] *= edgeVal(alpha[j],in[j], 1);
227 else
228 out[i][0][j] *= nodeVal(alpha[j],in[j], 1);
229 // rest of the product
230 for (int l=0; l<d; l++){
231 if (l!=j){
232 if(isEdge(alpha))
233 out[i][0][j] *= edgeVal(alpha[l],in[l],0);
234 else
235 out[i][0][j] *= nodeVal(alpha[l],in[l],0);
236 }
237 }
238 }
239 }
240 if(d==2){
241 for(int j = 0; j < d; j++){
242 out[0][0][j] -= 0.5*(out[1][0][j]+out[3][0][j]);
243 out[2][0][j] -= 0.5*(out[1][0][j]+out[4][0][j]);
244 out[5][0][j] -= 0.5*(out[3][0][j]+out[6][0][j]);
245 out[7][0][j] -= 0.5*(out[4][0][j]+out[6][0][j]);
246 }
247 }
248 if(d==3){
249 for(int j = 0; j < d; j++){
250 out[0][0][j] -= 0.5*(out[1][0][j] +out[3][0][j] +out[8][0][j]);
251 out[2][0][j] -= 0.5*(out[1][0][j] +out[4][0][j] +out[9][0][j]);
252 out[5][0][j] -= 0.5*(out[3][0][j] +out[6][0][j] +out[10][0][j]);
253 out[7][0][j] -= 0.5*(out[4][0][j] +out[6][0][j] +out[11][0][j]);
254 out[12][0][j] -= 0.5*(out[8][0][j] +out[13][0][j]+out[15][0][j]);
255 out[14][0][j] -= 0.5*(out[9][0][j] +out[13][0][j]+out[16][0][j]);
256 out[17][0][j] -= 0.5*(out[10][0][j]+out[15][0][j]+out[18][0][j]);
257 out[19][0][j] -= 0.5*(out[11][0][j]+out[16][0][j]+out[18][0][j]);
258 }
259 }
260 }
261
267 template<int diffOrder>
268 inline void evaluate(
269 const std::array<int,1>& direction,
270 const typename Traits::DomainType& in,
271 std::vector<typename Traits::RangeType>& out) const
272 {
273 static_assert(diffOrder == 1, "We only can compute first derivatives");
274 out.resize(size());
275 DuneVector alpha(0);
276
277 // Loop over all shape functions
278 for (size_t i=0; i<size(); i++)
279 {
280 // convert index i to multiindex
281 if(i>0){ multiindex(alpha);}
282
283 // Loop over all coordinate directions
284 std::size_t j = direction[0];
285
286 out[i] =1.0;
287 // Initialize: the overall expression is a product
288 // if j-th bit of i is set to -1, else 1
289 if(isEdge(alpha))
290 out[i][0] *= edgeVal(alpha[j],in[j],1);
291 else
292 out[i][0] *= nodeVal(alpha[j],in[j],1);
293
294 // rest of the product
295 for (std::size_t l=0; l<d; l++){
296 if (l!=j){
297 if(isEdge(alpha))
298 out[i][0] *= edgeVal(alpha[l],in[l],0);
299 else
300 out[i][0] *= nodeVal(alpha[l],in[l],0);
301 }
302 }
303 }
304 if(d==2){
305 out[0][0] -= 0.5*(out[1][0]+out[3][0]);
306 out[2][0] -= 0.5*(out[1][0]+out[4][0]);
307 out[5][0] -= 0.5*(out[3][0]+out[6][0]);
308 out[7][0] -= 0.5*(out[4][0]+out[6][0]);
309 }
310 if(d==3){
311 out[0][0] -= 0.5*(out[1][0] +out[3][0] +out[8][0]);
312 out[2][0] -= 0.5*(out[1][0] +out[4][0] +out[9][0]);
313 out[5][0] -= 0.5*(out[3][0] +out[6][0] +out[10][0]);
314 out[7][0] -= 0.5*(out[4][0] +out[6][0] +out[11][0]);
315 out[12][0] -= 0.5*(out[8][0] +out[13][0]+out[15][0]);
316 out[14][0] -= 0.5*(out[9][0] +out[13][0]+out[16][0]);
317 out[17][0] -= 0.5*(out[10][0]+out[15][0]+out[18][0]);
318 out[19][0] -= 0.5*(out[11][0]+out[16][0]+out[18][0]);
319 }
320 }
321
323 unsigned int order () const
324 {
325 return 2;
326 }
327 };
328}
329
330#endif
Serendipity basis functions of order k on the reference cube.
Definition: serendipitylocalbasis.hh:30
unsigned int order() const
Polynomial order of the shape functions.
Definition: serendipitylocalbasis.hh:323
void evaluateJacobian(const typename Traits::DomainType &in, std::vector< typename Traits::JacobianType > &out) const
Evaluate Jacobian of all shape functions.
Definition: serendipitylocalbasis.hh:207
void evaluate(const std::array< int, 1 > &direction, const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate derivative in a given direction.
Definition: serendipitylocalbasis.hh:268
void evaluateFunction(const typename Traits::DomainType &in, std::vector< typename Traits::RangeType > &out) const
Evaluate all shape functions.
Definition: serendipitylocalbasis.hh:131
unsigned int size() const
number of shape functions
Definition: serendipitylocalbasis.hh:117
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 3, 22:46, 2025)