v0.15.0
Loading...
Searching...
No Matches
BernsteinBezier.cpp
Go to the documentation of this file.
1/** \file BernsteinBezier.cpp
2\brief Bernstein-Bezier polynomials for H1 space
3
4Implementation based on \cite ainsworth2011bernstein
5
6*/
7
8template <int D, int Side>
9MoFEMErrorCode BernsteinBezier::generateIndicesVertex(const int N, int *alpha) {
11 alpha[(D + 1) * Side + Side] = N;
13}
14
15template <int D, int Side>
16MoFEMErrorCode BernsteinBezier::generateIndicesEdgeOnSimplex(const int N,
17 int *alpha) {
18 constexpr int edge_nodes[6][2] = {{0, 1}, {1, 2}, {2, 0},
19 {0, 3}, {1, 3}, {2, 3}};
21 std::fill(alpha, &alpha[(D + 1) * NBEDGE_H1(N)], 0);
22 for (int n = N - 1; n != 0; --n) {
23
24 // + o o o o +
25
26 alpha[edge_nodes[Side][0]] = n;
27 alpha[edge_nodes[Side][1]] = N - n;
28 alpha += D + 1;
29 }
31}
32
33template <int D, int Side>
34MoFEMErrorCode BernsteinBezier::generateIndicesTriOnSimplex(const int N,
35 int *alpha) {
36 constexpr int tri_nodes[4][3] = {{0, 1, 3}, {1, 2, 3}, {0, 2, 3}, {0, 1, 2}};
38 std::fill(alpha, &alpha[(D + 1) * NBFACETRI_H1(N)], 0);
39 for (int n0 = 1; n0 != N - 1; ++n0) {
40 for (int n1 = 1; n1 != N - n0; ++n1) {
41
42 // +
43 // + +
44 // + o +
45 // + o o +
46 // + o o o +
47 // + + + + + +
48
49 alpha[tri_nodes[Side][0]] = n0;
50 alpha[tri_nodes[Side][1]] = n1;
51 alpha[tri_nodes[Side][2]] = N - n0 - n1;
52
53 alpha += D + 1;
54 }
55 }
57}
58
59MoFEMErrorCode BernsteinBezier::generateIndicesTetOnSimplex(const int N,
60 int *alpha) {
62 std::fill(alpha, &alpha[4 * NBVOLUMETET_H1(N)], 0);
63 for (int n0 = 1; n0 != N - 1; ++n0) {
64 for (int n1 = 1; n1 != N - n0; ++n1) {
65 for (int n2 = 1; n2 != N - n0 - n1; ++n2) {
66 alpha[0] = n0;
67 alpha[1] = n1;
68 alpha[2] = n2;
69 alpha[3] = N - n0 - n1 - n2;
70 alpha += 4;
71 }
72 }
73 }
75}
76
77MoFEMErrorCode BernsteinBezier::generateIndicesVertexEdge(const int N,
78 int *alpha) {
80 CHKERR generateIndicesVertex<1, 0>(N, alpha);
81 CHKERR generateIndicesVertex<1, 1>(N, alpha);
83}
84
85MoFEMErrorCode BernsteinBezier::generateIndicesEdgeEdge(const int N,
86 int *alpha) {
87 return generateIndicesEdgeOnSimplex<1, 0>(N, alpha);
88}
89
90MoFEMErrorCode BernsteinBezier::generateIndicesVertexTri(const int N,
91 int *alpha) {
93 CHKERR generateIndicesVertex<2, 0>(N, alpha);
94 CHKERR generateIndicesVertex<2, 1>(N, alpha);
95 CHKERR generateIndicesVertex<2, 2>(N, alpha);
97}
98
99MoFEMErrorCode BernsteinBezier::generateIndicesEdgeTri(const int N[],
100 int *alpha[]) {
102 CHKERR generateIndicesEdgeOnSimplex<2, 0>(N[0], alpha[0]);
103 CHKERR generateIndicesEdgeOnSimplex<2, 1>(N[1], alpha[1]);
104 CHKERR generateIndicesEdgeOnSimplex<2, 2>(N[2], alpha[2]);
106}
107
108MoFEMErrorCode BernsteinBezier::generateIndicesEdgeTri(const int side,
109 const int N,
110 int *alpha) {
111 switch (side) {
112 case 0:
113 return generateIndicesEdgeOnSimplex<2, 0>(N, alpha);
114 case 1:
115 return generateIndicesEdgeOnSimplex<2, 1>(N, alpha);
116 case 2:
117 return generateIndicesEdgeOnSimplex<2, 2>(N, alpha);
118 default:
119 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
120 "Wrong side number, on triangle are edges from 0 to 2");
121 }
122}
123
124MoFEMErrorCode BernsteinBezier::generateIndicesTriTri(const int N, int *alpha) {
125 return generateIndicesTriOnSimplex<2, 3>(N, alpha);
126}
127
128MoFEMErrorCode BernsteinBezier::generateIndicesVertexTet(const int N,
129 int *alpha) {
131 CHKERR generateIndicesVertex<3, 0>(N, alpha);
132 CHKERR generateIndicesVertex<3, 1>(N, alpha);
133 CHKERR generateIndicesVertex<3, 2>(N, alpha);
134 CHKERR generateIndicesVertex<3, 3>(N, alpha);
136}
137
138MoFEMErrorCode BernsteinBezier::generateIndicesEdgeTet(const int N[],
139 int *alpha[]) {
141 CHKERR generateIndicesEdgeOnSimplex<3, 0>(N[0], alpha[0]);
142 CHKERR generateIndicesEdgeOnSimplex<3, 1>(N[1], alpha[1]);
143 CHKERR generateIndicesEdgeOnSimplex<3, 2>(N[2], alpha[2]);
144 CHKERR generateIndicesEdgeOnSimplex<3, 3>(N[3], alpha[3]);
145 CHKERR generateIndicesEdgeOnSimplex<3, 4>(N[4], alpha[4]);
146 CHKERR generateIndicesEdgeOnSimplex<3, 5>(N[5], alpha[5]);
148}
149
150MoFEMErrorCode BernsteinBezier::generateIndicesEdgeTet(const int side,
151 const int N,
152 int *alpha) {
153 switch (side) {
154 case 0:
155 return generateIndicesEdgeOnSimplex<3, 0>(N, alpha);
156 case 1:
157 return generateIndicesEdgeOnSimplex<3, 1>(N, alpha);
158 case 2:
159 return generateIndicesEdgeOnSimplex<3, 2>(N, alpha);
160 case 3:
161 return generateIndicesEdgeOnSimplex<3, 3>(N, alpha);
162 case 4:
163 return generateIndicesEdgeOnSimplex<3, 4>(N, alpha);
164 case 5:
165 return generateIndicesEdgeOnSimplex<3, 5>(N, alpha);
166 default:
168 "Wrong side number, on tetrahedron are edges from 0 to 5 "
169 "in generateIndicesEdgeTet");
170 }
171}
172
173MoFEMErrorCode BernsteinBezier::generateIndicesTriTet(const int N[],
174 int *alpha[]) {
176 CHKERR generateIndicesTriOnSimplex<3, 0>(N[0], alpha[0]);
177 CHKERR generateIndicesTriOnSimplex<3, 1>(N[1], alpha[1]);
178 CHKERR generateIndicesTriOnSimplex<3, 2>(N[2], alpha[2]);
179 CHKERR generateIndicesTriOnSimplex<3, 3>(N[3], alpha[3]);
181}
182
183MoFEMErrorCode BernsteinBezier::generateIndicesTriTet(const int side,
184 const int N, int *alpha) {
185 switch (side) {
186 case 0:
187 return generateIndicesTriOnSimplex<3, 0>(N, alpha);
188 case 1:
189 return generateIndicesTriOnSimplex<3, 1>(N, alpha);
190 case 2:
191 return generateIndicesTriOnSimplex<3, 2>(N, alpha);
192 case 3:
193 return generateIndicesTriOnSimplex<3, 3>(N, alpha);
194 default:
196 "Wrong side number, on tetrahedron are from from 0 to 3 "
197 "in generateIndicesTriTet");
198 }
199}
200
201MoFEMErrorCode BernsteinBezier::generateIndicesTetTet(const int N, int *alpha) {
202 return generateIndicesTetOnSimplex(N, alpha);
203}
204
205template <int D>
206MoFEMErrorCode BernsteinBezier::genrateDerivativeIndices(
207 const int N, const int n_alpha, const int *alpha, const int *diff,
208 const int n_alpha_diff, const int *alpha_diff, double *c) {
210
211 std::fill(c, &c[n_alpha * n_alpha_diff], 0);
212 int abs_diff = diff[0];
213 for (int d = 1; d != D + 1; ++d)
214 abs_diff += diff[d];
215
216 const double b = boost::math::factorial<double>(N) /
217 boost::math::factorial<double>(N - abs_diff);
218
219 for (int i = 0; i != n_alpha; ++i) {
220 for (int j = 0; j != n_alpha_diff; ++j) {
221 int d = 0;
222 for (; d != D + 1; ++d)
223 if (alpha_diff[(D + 1) * j + d] + diff[d] != alpha[(D + 1) * i + d])
224 break;
225 if (d == D + 1) {
226 c[i * n_alpha_diff + j] = b;
227 break;
228 }
229 }
230 }
231
233}
234
235MoFEMErrorCode BernsteinBezier::genrateDerivativeIndicesEdges(
236 const int N, const int n_alpha, const int *alpha, const int *diff,
237 const int n_alpha_diff, const int *alpha_diff, double *c) {
238 return genrateDerivativeIndices<1>(N, n_alpha, alpha, diff, n_alpha_diff,
239 alpha_diff, c);
240}
241
242MoFEMErrorCode BernsteinBezier::genrateDerivativeIndicesTri(
243 const int N, const int n_alpha, const int *alpha, const int *diff,
244 const int n_alpha_diff, const int *alpha_diff, double *c) {
245 return genrateDerivativeIndices<2>(N, n_alpha, alpha, diff, n_alpha_diff,
246 alpha_diff, c);
247}
248
249MoFEMErrorCode BernsteinBezier::genrateDerivativeIndicesTet(
250 const int N, const int n_alpha, const int *alpha, const int *diff,
251 const int n_alpha_diff, const int *alpha_diff, double *c) {
252 return genrateDerivativeIndices<3>(N, n_alpha, alpha, diff, n_alpha_diff,
253 alpha_diff, c);
254}
255
256template <int D>
258BernsteinBezier::domainPoints(const int N, const int n_x, const int n_alpha,
259 const int *alpha, const double *x_k,
260 double *x_alpha) {
261 FTensor::Index<'i', D> i;
263
265 x_alpha, &x_alpha[1], &x_alpha[2]);
266 for (int n0 = 0; n0 != n_alpha; ++n0) {
267 t_x_alpha(i) = 0;
269 x_k, &x_k[1], &x_k[2]);
270 for (int n1 = 0; n1 != n_x; ++n1) {
271 t_x_alpha(i) += static_cast<double>(*alpha) * t_x_k(i);
272 ++t_x_k;
273 ++alpha;
274 }
275 t_x_alpha(i) /= static_cast<double>(N);
276 ++t_x_alpha;
277 }
278
280}
281
282MoFEMErrorCode BernsteinBezier::domainPoints3d(const int N, const int n_x,
283 const int n_alpha,
284 const int *alpha,
285 const double *x_k,
286 double *x_alpha) {
287 return domainPoints<3>(N, n_x, n_alpha, alpha, x_k, x_alpha);
288}
289
290template <int D, bool GRAD_BASE>
292BernsteinBezier::baseFunctions(const int N, const int gdim, const int n_alpha,
293 const int *alpha, const double *lambda,
294 const double *grad_lambda, double *base,
295 double *grad_base) {
297
298 const int *const alpha0 = alpha;
299 const double fN = boost::math::factorial<double>(N);
300 constexpr int MAX_ALPHA = 12;
301 int max_alpha = *std::max_element(alpha, &alpha[(D + 1) * n_alpha]);
302 if(max_alpha > MAX_ALPHA)
303 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
304 "Is assumed maximal order not to be bigger than %d", MAX_ALPHA);
305 std::array<double, (D + 1) * (MAX_ALPHA + 1)> pow_alpha;
306 std::array<double, (MAX_ALPHA + 1)> factorial_alpha;
307 std::array<double, D + 1> terms, diff_terms;
308 const double *const grad_lambda0 = grad_lambda;
309
310 factorial_alpha[0] = 1;
311 if (max_alpha >= 1)
312 factorial_alpha[1] = 1;
313 if (max_alpha >= 2)
314 for (int a = 2; a <= max_alpha; ++a)
315 factorial_alpha[a] = factorial_alpha[a - 1] * a;
316
317 for (int g = 0; g != gdim; ++g) {
318
319 for (int n = 0; n != D + 1; ++n, ++lambda) {
320 const size_t shift = (MAX_ALPHA + 1) * n;
321 double *pow_alpha_ptr = &pow_alpha[shift];
322 *pow_alpha_ptr = 1;
323
324 if (max_alpha >= 1) {
325 ++pow_alpha_ptr;
326 *pow_alpha_ptr = *lambda;
327 }
328
329 if (max_alpha >= 2)
330 for (int a = 2; a <= max_alpha; ++a) {
331 const double p = (*pow_alpha_ptr) * (*lambda);
332 ++pow_alpha_ptr;
333 *pow_alpha_ptr = p;
334 }
335 }
336
337 for (int n0 = 0; n0 != n_alpha; ++n0) {
338
339 grad_lambda = grad_lambda0;
340
341 double f = factorial_alpha[(*alpha)];
342 double *terms_ptr = terms.data();
343 double *diff_terms_ptr = diff_terms.data();
344 *terms_ptr = pow_alpha[(*alpha)];
345 if (GRAD_BASE) {
346 if (*alpha > 0)
347 *diff_terms_ptr = (*alpha) * pow_alpha[(*alpha) - 1];
348 else
349 *diff_terms_ptr = 0;
350 }
351 *base = *terms_ptr;
352 ++alpha;
353 ++terms_ptr;
354 ++diff_terms_ptr;
355
356 for (int n1 = 1; n1 < D + 1;
357 ++n1, ++alpha, ++terms_ptr, ++diff_terms_ptr) {
358 f *= factorial_alpha[(*alpha)];
359 const size_t shift = (MAX_ALPHA + 1) * n1;
360 *terms_ptr = pow_alpha[shift + (*alpha)];
361 if (GRAD_BASE) {
362 if (*alpha > 0)
363 *diff_terms_ptr = (*alpha) * pow_alpha[shift + (*alpha) - 1];
364 else
365 *diff_terms_ptr = 0;
366 }
367 *base *= *terms_ptr;
368 }
369
370 const double b = fN / f;
371 *base *= b;
372 ++base;
373
374 if (GRAD_BASE) {
375 double *terms_ptr = terms.data();
376 double *diff_terms_ptr = diff_terms.data();
377 double z = *diff_terms_ptr;
378 ++terms_ptr;
379 ++diff_terms_ptr;
380 for (int n2 = 1; n2 != D + 1; ++n2, ++terms_ptr)
381 z *= *terms_ptr;
382
383 for (int d = 0; d != D; ++d, ++grad_lambda)
384 grad_base[d] = z * (*grad_lambda);
385
386 for (int n1 = 1; n1 < D + 1; ++n1) {
387 z = *diff_terms_ptr;
388
389 int n2 = 0;
390 for (terms_ptr = terms.data(); n2 != n1; ++n2, ++terms_ptr)
391 z *= *terms_ptr;
392
393 ++n2;
394 ++terms_ptr;
395
396 for (; n2 < D + 1; ++n2, ++terms_ptr)
397 z *= *terms_ptr;
398
399 for (int d = 0; d != D; ++d, ++grad_lambda)
400 grad_base[d] += z * (*grad_lambda);
401
402 ++diff_terms_ptr;
403 }
404 for (int d = 0; d != D; ++d)
405 grad_base[d] *= b;
406 grad_base += D;
407 }
408
409 }
410
411 alpha = alpha0;
412 }
413
415}
416
417MoFEMErrorCode BernsteinBezier::baseFunctionsEdge(
418 const int N, const int gdim, const int n_alpha, const int *alpha,
419 const double *lambda, const double *grad_lambda, double *base,
420 double *grad_base) {
421 if (grad_base)
422 return baseFunctions<1, true>(N, gdim, n_alpha, alpha, lambda, grad_lambda,
423 base, grad_base);
424 else
425 return baseFunctions<1, false>(N, gdim, n_alpha, alpha, lambda, grad_lambda,
426 base, grad_base);
427}
428
429MoFEMErrorCode BernsteinBezier::baseFunctionsTri(
430 const int N, const int gdim, const int n_alpha, const int *alpha,
431 const double *lambda, const double *grad_lambda, double *base,
432 double *grad_base) {
433 if (grad_base)
434 return baseFunctions<2, true>(N, gdim, n_alpha, alpha, lambda, grad_lambda,
435 base, grad_base);
436 else
437 return baseFunctions<2, false>(N, gdim, n_alpha, alpha, lambda, grad_lambda,
438 base, grad_base);
439}
440
441MoFEMErrorCode BernsteinBezier::baseFunctionsTet(
442 const int N, const int gdim, const int n_alpha, const int *alpha,
443 const double *lambda, const double *grad_lambda, double *base,
444 double *grad_base) {
445 if (grad_base)
446 return baseFunctions<3, true>(N, gdim, n_alpha, alpha, lambda, grad_lambda,
447 base, grad_base);
448 else
449 return baseFunctions<3, false>(N, gdim, n_alpha, alpha, lambda, grad_lambda,
450 base, grad_base);
451}
constexpr double a
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define NBVOLUMETET_H1(P)
Number of base functions on tetrahedron for H1 space.
#define NBEDGE_H1(P)
Number of base function on edge for H1 space.
#define NBFACETRI_H1(P)
Number of base function on triangle for H1 space.
FTensor::Index< 'i', SPACE_DIM > i
static double lambda
const double c
speed of light (cm/ns)
double D
const double n
refractive index of diffusive medium
FTensor::Index< 'j', 3 > j
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
constexpr double g
const int N
Definition speed_test.cpp:3