Reference documentation for deal.II version 9.6.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
tensor_product_polynomials_bubbles.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2013 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_tensor_product_polynomials_bubbles_h
16#define dealii_tensor_product_polynomials_bubbles_h
17
18
19#include <deal.II/base/config.h>
20
22#include <deal.II/base/point.h>
24#include <deal.II/base/tensor.h>
27
28#include <vector>
29
31
51template <int dim>
53{
54public:
59 static constexpr unsigned int dimension = dim;
60
66 template <class Pol>
67 TensorProductPolynomialsBubbles(const std::vector<Pol> &pols);
68
72 void
73 output_indices(std::ostream &out) const;
74
80 void
81 set_numbering(const std::vector<unsigned int> &renumber);
82
86 const std::vector<unsigned int> &
88
92 const std::vector<unsigned int> &
94
107 void
108 evaluate(const Point<dim> &unit_point,
109 std::vector<double> &values,
110 std::vector<Tensor<1, dim>> &grads,
111 std::vector<Tensor<2, dim>> &grad_grads,
112 std::vector<Tensor<3, dim>> &third_derivatives,
113 std::vector<Tensor<4, dim>> &fourth_derivatives) const override;
114
127 double
128 compute_value(const unsigned int i, const Point<dim> &p) const override;
129
142 template <int order>
144 compute_derivative(const unsigned int i, const Point<dim> &p) const;
145
149 virtual Tensor<1, dim>
150 compute_1st_derivative(const unsigned int i,
151 const Point<dim> &p) const override;
152
156 virtual Tensor<2, dim>
157 compute_2nd_derivative(const unsigned int i,
158 const Point<dim> &p) const override;
159
163 virtual Tensor<3, dim>
164 compute_3rd_derivative(const unsigned int i,
165 const Point<dim> &p) const override;
166
170 virtual Tensor<4, dim>
171 compute_4th_derivative(const unsigned int i,
172 const Point<dim> &p) const override;
173
187 compute_grad(const unsigned int i, const Point<dim> &p) const override;
188
202 compute_grad_grad(const unsigned int i, const Point<dim> &p) const override;
203
210 unsigned int
211 n() const;
212
217 std::string
218 name() const override;
219
223 virtual std::unique_ptr<ScalarPolynomialsBase<dim>>
224 clone() const override;
225
226private:
231
235 std::vector<unsigned int> index_map;
236
240 std::vector<unsigned int> index_map_inverse;
241};
242
246/* ---------------- template and inline functions ---------- */
247
248#ifndef DOXYGEN
249
250template <int dim>
251template <class Pol>
253 const std::vector<Pol> &pols)
254 : ScalarPolynomialsBase<dim>(1,
255 Utilities::fixed_power<dim>(pols.size()) + dim)
256 , tensor_polys(pols)
257 , index_map(tensor_polys.n() +
258 ((tensor_polys.polynomials.size() <= 2) ? 1 : dim))
259 , index_map_inverse(tensor_polys.n() +
260 ((tensor_polys.polynomials.size() <= 2) ? 1 : dim))
261{
262 const unsigned int q_degree = tensor_polys.polynomials.size() - 1;
263 const unsigned int n_bubbles = ((q_degree <= 1) ? 1 : dim);
264 // append index for renumbering
265 for (unsigned int i = 0; i < tensor_polys.n() + n_bubbles; ++i)
266 {
267 index_map[i] = i;
268 index_map_inverse[i] = i;
269 }
270}
271
272
273template <int dim>
274inline unsigned int
276{
277 return tensor_polys.n() + dim;
278}
279
280
281template <>
282inline unsigned int
284{
286}
287
288
289template <int dim>
290inline const std::vector<unsigned int> &
292{
293 return index_map;
294}
295
296
297template <int dim>
298inline const std::vector<unsigned int> &
300{
301 return index_map_inverse;
302}
303
304
305template <int dim>
306inline std::string
308{
309 return "TensorProductPolynomialsBubbles";
310}
311
312
313template <int dim>
314template <int order>
317 const unsigned int i,
318 const Point<dim> &p) const
319{
320 const unsigned int q_degree = tensor_polys.polynomials.size() - 1;
321 const unsigned int max_q_indices = tensor_polys.n();
322 Assert(i < max_q_indices + /* n_bubbles= */ ((q_degree <= 1) ? 1 : dim),
324
325 // treat the regular basis functions
326 if (i < max_q_indices)
327 return tensor_polys.template compute_derivative<order>(i, p);
328
329 const unsigned int comp = i - tensor_polys.n();
330
331 Tensor<order, dim> derivative;
332 switch (order)
333 {
334 case 1:
335 {
336 Tensor<1, dim> &derivative_1 =
337 *reinterpret_cast<Tensor<1, dim> *>(&derivative);
338
339 for (unsigned int d = 0; d < dim; ++d)
340 {
341 derivative_1[d] = 1.;
342 // compute grad(4*\prod_{i=1}^d (x_i(1-x_i)))(p)
343 for (unsigned j = 0; j < dim; ++j)
344 derivative_1[d] *=
345 (d == j ? 4 * (1 - 2 * p[j]) : 4 * p[j] * (1 - p[j]));
346 // and multiply with (2*x_i-1)^{r-1}
347 for (unsigned int i = 0; i < q_degree - 1; ++i)
348 derivative_1[d] *= 2 * p[comp] - 1;
349 }
350
351 if (q_degree >= 2)
352 {
353 // add \prod_{i=1}^d 4*(x_i(1-x_i))(p)
354 double value = 1.;
355 for (unsigned int j = 0; j < dim; ++j)
356 value *= 4 * p[j] * (1 - p[j]);
357 // and multiply with grad(2*x_i-1)^{r-1}
358 double tmp = value * 2 * (q_degree - 1);
359 for (unsigned int i = 0; i < q_degree - 2; ++i)
360 tmp *= 2 * p[comp] - 1;
361 derivative_1[comp] += tmp;
362 }
363
364 return derivative;
365 }
366 case 2:
367 {
368 Tensor<2, dim> &derivative_2 =
369 *reinterpret_cast<Tensor<2, dim> *>(&derivative);
370
371 double v[dim + 1][3];
372 {
373 for (unsigned int c = 0; c < dim; ++c)
374 {
375 v[c][0] = 4 * p[c] * (1 - p[c]);
376 v[c][1] = 4 * (1 - 2 * p[c]);
377 v[c][2] = -8;
378 }
379
380 double tmp = 1.;
381 for (unsigned int i = 0; i < q_degree - 1; ++i)
382 tmp *= 2 * p[comp] - 1;
383 v[dim][0] = tmp;
384
385 if (q_degree >= 2)
386 {
387 double tmp = 2 * (q_degree - 1);
388 for (unsigned int i = 0; i < q_degree - 2; ++i)
389 tmp *= 2 * p[comp] - 1;
390 v[dim][1] = tmp;
391 }
392 else
393 v[dim][1] = 0.;
394
395 if (q_degree >= 3)
396 {
397 double tmp = 4 * (q_degree - 2) * (q_degree - 1);
398 for (unsigned int i = 0; i < q_degree - 3; ++i)
399 tmp *= 2 * p[comp] - 1;
400 v[dim][2] = tmp;
401 }
402 else
403 v[dim][2] = 0.;
404 }
405
406 // calculate (\partial_j \partial_k \psi) * monomial
407 Tensor<2, dim> grad_grad_1;
408 for (unsigned int d1 = 0; d1 < dim; ++d1)
409 for (unsigned int d2 = 0; d2 < dim; ++d2)
410 {
411 grad_grad_1[d1][d2] = v[dim][0];
412 for (unsigned int x = 0; x < dim; ++x)
413 {
414 unsigned int derivative = 0;
415 if (d1 == x || d2 == x)
416 {
417 if (d1 == d2)
418 derivative = 2;
419 else
420 derivative = 1;
421 }
422 grad_grad_1[d1][d2] *= v[x][derivative];
423 }
424 }
425
426 // calculate (\partial_j \psi) *(\partial_k monomial)
427 // and (\partial_k \psi) *(\partial_j monomial)
428 Tensor<2, dim> grad_grad_2;
429 Tensor<2, dim> grad_grad_3;
430 for (unsigned int d = 0; d < dim; ++d)
431 {
432 grad_grad_2[d][comp] = v[dim][1];
433 grad_grad_3[comp][d] = v[dim][1];
434 for (unsigned int x = 0; x < dim; ++x)
435 {
436 grad_grad_2[d][comp] *= v[x][d == x];
437 grad_grad_3[comp][d] *= v[x][d == x];
438 }
439 }
440
441 // calculate \psi *(\partial j \partial_k monomial) and sum
442 double psi_value = 1.;
443 for (unsigned int x = 0; x < dim; ++x)
444 psi_value *= v[x][0];
445
446 for (unsigned int d1 = 0; d1 < dim; ++d1)
447 for (unsigned int d2 = 0; d2 < dim; ++d2)
448 derivative_2[d1][d2] =
449 grad_grad_1[d1][d2] + grad_grad_2[d1][d2] + grad_grad_3[d1][d2];
450 derivative_2[comp][comp] += psi_value * v[dim][2];
451
452 return derivative;
453 }
454 default:
455 {
457 return derivative;
458 }
459 }
460}
461
462
463
464template <int dim>
465inline Tensor<1, dim>
467 const unsigned int i,
468 const Point<dim> &p) const
469{
470 return compute_derivative<1>(i, p);
471}
472
473
474
475template <int dim>
476inline Tensor<2, dim>
478 const unsigned int i,
479 const Point<dim> &p) const
480{
481 return compute_derivative<2>(i, p);
482}
483
484
485
486template <int dim>
487inline Tensor<3, dim>
489 const unsigned int i,
490 const Point<dim> &p) const
491{
492 return compute_derivative<3>(i, p);
493}
494
495
496
497template <int dim>
498inline Tensor<4, dim>
500 const unsigned int i,
501 const Point<dim> &p) const
502{
503 return compute_derivative<4>(i, p);
504}
505
506#endif // DOXYGEN
508
509#endif
Definition point.h:111
const std::vector< unsigned int > & get_numbering_inverse() const
virtual Tensor< 2, dim > compute_2nd_derivative(const unsigned int i, const Point< dim > &p) const override
std::string name() const override
void evaluate(const Point< dim > &unit_point, std::vector< double > &values, std::vector< Tensor< 1, dim > > &grads, std::vector< Tensor< 2, dim > > &grad_grads, std::vector< Tensor< 3, dim > > &third_derivatives, std::vector< Tensor< 4, dim > > &fourth_derivatives) const override
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const override
Tensor< order, dim > compute_derivative(const unsigned int i, const Point< dim > &p) const
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const override
TensorProductPolynomialsBubbles(const std::vector< Pol > &pols)
void set_numbering(const std::vector< unsigned int > &renumber)
virtual Tensor< 1, dim > compute_1st_derivative(const unsigned int i, const Point< dim > &p) const override
virtual Tensor< 3, dim > compute_3rd_derivative(const unsigned int i, const Point< dim > &p) const override
const std::vector< unsigned int > & get_numbering() const
virtual Tensor< 4, dim > compute_4th_derivative(const unsigned int i, const Point< dim > &p) const override
double compute_value(const unsigned int i, const Point< dim > &p) const override
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone() const override
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
#define DEAL_II_NOT_IMPLEMENTED()
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static const unsigned int invalid_unsigned_int
Definition types.h:220