Reference documentation for deal.II version 9.1.0-pre
tensor_product_polynomials.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/exceptions.h>
17 #include <deal.II/base/polynomials_piecewise.h>
18 #include <deal.II/base/table.h>
19 #include <deal.II/base/tensor_product_polynomials.h>
20 
21 #include <boost/container/small_vector.hpp>
22 
23 #include <array>
24 
25 DEAL_II_NAMESPACE_OPEN
26 
27 
28 
29 /* ------------------- TensorProductPolynomials -------------- */
30 
31 
32 namespace internal
33 {
34  namespace
35  {
36  template <int dim>
37  inline void
38  compute_tensor_index(const unsigned int,
39  const unsigned int,
40  const unsigned int,
41  unsigned int (&)[dim])
42  {
43  Assert(false, ExcNotImplemented());
44  }
45 
46  inline void
47  compute_tensor_index(const unsigned int n,
48  const unsigned int,
49  const unsigned int,
50  unsigned int (&indices)[1])
51  {
52  indices[0] = n;
53  }
54 
55  inline void
56  compute_tensor_index(const unsigned int n,
57  const unsigned int n_pols_0,
58  const unsigned int,
59  unsigned int (&indices)[2])
60  {
61  indices[0] = n % n_pols_0;
62  indices[1] = n / n_pols_0;
63  }
64 
65  inline void
66  compute_tensor_index(const unsigned int n,
67  const unsigned int n_pols_0,
68  const unsigned int n_pols_1,
69  unsigned int (&indices)[3])
70  {
71  indices[0] = n % n_pols_0;
72  indices[1] = (n / n_pols_0) % n_pols_1;
73  indices[2] = n / (n_pols_0 * n_pols_1);
74  }
75  } // namespace
76 } // namespace internal
77 
78 
79 
80 template <int dim, typename PolynomialType>
81 inline void
83  const unsigned int i,
84  unsigned int (&indices)[(dim > 0 ? dim : 1)]) const
85 {
86  Assert(i < Utilities::fixed_power<dim>(polynomials.size()),
88  internal::compute_tensor_index(index_map[i],
89  polynomials.size(),
90  polynomials.size(),
91  indices);
92 }
93 
94 
95 
96 template <int dim, typename PolynomialType>
97 void
99  std::ostream &out) const
100 {
101  unsigned int ix[dim];
102  for (unsigned int i = 0; i < n_tensor_pols; ++i)
103  {
104  compute_index(i, ix);
105  out << i << "\t";
106  for (unsigned int d = 0; d < dim; ++d)
107  out << ix[d] << " ";
108  out << std::endl;
109  }
110 }
111 
112 
113 
114 template <int dim, typename PolynomialType>
115 void
117  const std::vector<unsigned int> &renumber)
118 {
119  Assert(renumber.size() == index_map.size(),
120  ExcDimensionMismatch(renumber.size(), index_map.size()));
121 
122  index_map = renumber;
123  for (unsigned int i = 0; i < index_map.size(); ++i)
124  index_map_inverse[index_map[i]] = i;
125 }
126 
127 
128 
129 template <>
130 double
132  const unsigned int,
133  const Point<0> &) const
134 {
135  Assert(false, ExcNotImplemented());
136  return 0;
137 }
138 
139 
140 
141 template <int dim, typename PolynomialType>
142 double
144  const unsigned int i,
145  const Point<dim> & p) const
146 {
147  Assert(dim > 0, ExcNotImplemented());
148 
149  unsigned int indices[dim];
150  compute_index(i, indices);
151 
152  double value = 1.;
153  for (unsigned int d = 0; d < dim; ++d)
154  value *= polynomials[indices[d]].value(p(d));
155 
156  return value;
157 }
158 
159 
160 
161 template <int dim, typename PolynomialType>
164  const unsigned int i,
165  const Point<dim> & p) const
166 {
167  unsigned int indices[dim];
168  compute_index(i, indices);
169 
170  // compute values and
171  // uni-directional derivatives at
172  // the given point in each
173  // co-ordinate direction
174  double v[dim][2];
175  {
176  std::vector<double> tmp(2);
177  for (unsigned int d = 0; d < dim; ++d)
178  {
179  polynomials[indices[d]].value(p(d), tmp);
180  v[d][0] = tmp[0];
181  v[d][1] = tmp[1];
182  }
183  }
184 
185  Tensor<1, dim> grad;
186  for (unsigned int d = 0; d < dim; ++d)
187  {
188  grad[d] = 1.;
189  for (unsigned int x = 0; x < dim; ++x)
190  grad[d] *= v[x][d == x];
191  }
192 
193  return grad;
194 }
195 
196 
197 
198 template <int dim, typename PolynomialType>
201  const unsigned int i,
202  const Point<dim> & p) const
203 {
204  unsigned int indices[dim];
205  compute_index(i, indices);
206 
207  double v[dim][3];
208  {
209  std::vector<double> tmp(3);
210  for (unsigned int d = 0; d < dim; ++d)
211  {
212  polynomials[indices[d]].value(p(d), tmp);
213  v[d][0] = tmp[0];
214  v[d][1] = tmp[1];
215  v[d][2] = tmp[2];
216  }
217  }
218 
219  Tensor<2, dim> grad_grad;
220  for (unsigned int d1 = 0; d1 < dim; ++d1)
221  for (unsigned int d2 = 0; d2 < dim; ++d2)
222  {
223  grad_grad[d1][d2] = 1.;
224  for (unsigned int x = 0; x < dim; ++x)
225  {
226  unsigned int derivative = 0;
227  if (d1 == x || d2 == x)
228  {
229  if (d1 == d2)
230  derivative = 2;
231  else
232  derivative = 1;
233  }
234  grad_grad[d1][d2] *= v[x][derivative];
235  }
236  }
237 
238  return grad_grad;
239 }
240 
241 
242 
243 template <int dim, typename PolynomialType>
244 void
246  const Point<dim> & p,
247  std::vector<double> & values,
248  std::vector<Tensor<1, dim>> &grads,
249  std::vector<Tensor<2, dim>> &grad_grads,
250  std::vector<Tensor<3, dim>> &third_derivatives,
251  std::vector<Tensor<4, dim>> &fourth_derivatives) const
252 {
253  Assert(dim <= 3, ExcNotImplemented());
254  Assert(values.size() == n_tensor_pols || values.size() == 0,
255  ExcDimensionMismatch2(values.size(), n_tensor_pols, 0));
256  Assert(grads.size() == n_tensor_pols || grads.size() == 0,
257  ExcDimensionMismatch2(grads.size(), n_tensor_pols, 0));
258  Assert(grad_grads.size() == n_tensor_pols || grad_grads.size() == 0,
259  ExcDimensionMismatch2(grad_grads.size(), n_tensor_pols, 0));
260  Assert(third_derivatives.size() == n_tensor_pols ||
261  third_derivatives.size() == 0,
262  ExcDimensionMismatch2(third_derivatives.size(), n_tensor_pols, 0));
263  Assert(fourth_derivatives.size() == n_tensor_pols ||
264  fourth_derivatives.size() == 0,
265  ExcDimensionMismatch2(fourth_derivatives.size(), n_tensor_pols, 0));
266 
267  const bool update_values = (values.size() == n_tensor_pols),
268  update_grads = (grads.size() == n_tensor_pols),
269  update_grad_grads = (grad_grads.size() == n_tensor_pols),
271  (third_derivatives.size() == n_tensor_pols),
272  update_4th_derivatives =
273  (fourth_derivatives.size() == n_tensor_pols);
274 
275  // check how many values/derivatives we have to compute
276  unsigned int n_values_and_derivatives = 0;
277  if (update_values)
278  n_values_and_derivatives = 1;
279  if (update_grads)
280  n_values_and_derivatives = 2;
281  if (update_grad_grads)
282  n_values_and_derivatives = 3;
284  n_values_and_derivatives = 4;
285  if (update_4th_derivatives)
286  n_values_and_derivatives = 5;
287 
288  // Compute the values (and derivatives, if necessary) of all 1D polynomials
289  // at this evaluation point. We need to compute dim*n_polynomials
290  // evaluations, involving an evaluation of each polynomial for each
291  // coordinate direction. Once we have those values, we perform the
292  // multiplications for the tensor product in the arbitrary dimension.
293  const unsigned int n_polynomials = polynomials.size();
294  boost::container::small_vector<std::array<std::array<double, 5>, dim>, 20>
295  values_1d(n_polynomials);
296  if (n_values_and_derivatives == 1)
297  for (unsigned int i = 0; i < n_polynomials; ++i)
298  for (unsigned int d = 0; d < dim; ++d)
299  values_1d[i][d][0] = polynomials[i].value(p(d));
300  else
301  for (unsigned int i = 0; i < n_polynomials; ++i)
302  for (unsigned d = 0; d < dim; ++d)
303  polynomials[i].value(p(d),
304  n_values_and_derivatives,
305  &values_1d[i][d][0]);
306 
307  unsigned int indices[3];
308  unsigned int ind = 0;
309  for (indices[2] = 0; indices[2] < (dim > 2 ? n_polynomials : 1); ++indices[2])
310  for (indices[1] = 0; indices[1] < (dim > 1 ? n_polynomials : 1);
311  ++indices[1])
312  if (n_values_and_derivatives == 1)
313  for (indices[0] = 0; indices[0] < n_polynomials; ++indices[0], ++ind)
314  {
315  double value = values_1d[indices[0]][0][0];
316  for (unsigned int d = 1; d < dim; ++d)
317  value *= values_1d[indices[d]][d][0];
318  values[index_map_inverse[ind]] = value;
319  }
320  else
321  for (indices[0] = 0; indices[0] < n_polynomials; ++indices[0], ++ind)
322  {
323  unsigned int i = index_map_inverse[ind];
324 
325  if (update_values)
326  {
327  double value = values_1d[indices[0]][0][0];
328  for (unsigned int x = 1; x < dim; ++x)
329  value *= values_1d[indices[x]][x][0];
330  values[i] = value;
331  }
332 
333  if (update_grads)
334  for (unsigned int d = 0; d < dim; ++d)
335  {
336  double grad = 1.;
337  for (unsigned int x = 0; x < dim; ++x)
338  grad *= values_1d[indices[x]][x][(d == x) ? 1 : 0];
339  grads[i][d] = grad;
340  }
341 
342  if (update_grad_grads)
343  for (unsigned int d1 = 0; d1 < dim; ++d1)
344  for (unsigned int d2 = 0; d2 < dim; ++d2)
345  {
346  double der2 = 1.;
347  for (unsigned int x = 0; x < dim; ++x)
348  {
349  unsigned int derivative = 0;
350  if (d1 == x)
351  ++derivative;
352  if (d2 == x)
353  ++derivative;
354 
355  der2 *= values_1d[indices[x]][x][derivative];
356  }
357  grad_grads[i][d1][d2] = der2;
358  }
359 
361  for (unsigned int d1 = 0; d1 < dim; ++d1)
362  for (unsigned int d2 = 0; d2 < dim; ++d2)
363  for (unsigned int d3 = 0; d3 < dim; ++d3)
364  {
365  double der3 = 1.;
366  for (unsigned int x = 0; x < dim; ++x)
367  {
368  unsigned int derivative = 0;
369  if (d1 == x)
370  ++derivative;
371  if (d2 == x)
372  ++derivative;
373  if (d3 == x)
374  ++derivative;
375 
376  der3 *= values_1d[indices[x]][x][derivative];
377  }
378  third_derivatives[i][d1][d2][d3] = der3;
379  }
380 
381  if (update_4th_derivatives)
382  for (unsigned int d1 = 0; d1 < dim; ++d1)
383  for (unsigned int d2 = 0; d2 < dim; ++d2)
384  for (unsigned int d3 = 0; d3 < dim; ++d3)
385  for (unsigned int d4 = 0; d4 < dim; ++d4)
386  {
387  double der4 = 1.;
388  for (unsigned int x = 0; x < dim; ++x)
389  {
390  unsigned int derivative = 0;
391  if (d1 == x)
392  ++derivative;
393  if (d2 == x)
394  ++derivative;
395  if (d3 == x)
396  ++derivative;
397  if (d4 == x)
398  ++derivative;
399 
400  der4 *= values_1d[indices[x]][x][derivative];
401  }
402  fourth_derivatives[i][d1][d2][d3][d4] = der4;
403  }
404  }
405 }
406 
407 
408 
409 /* ------------------- AnisotropicPolynomials -------------- */
410 
411 
412 template <int dim>
414  const std::vector<std::vector<Polynomials::Polynomial<double>>> &pols)
415  : polynomials(pols)
416  , n_tensor_pols(get_n_tensor_pols(pols))
417 {
418  Assert(pols.size() == dim, ExcDimensionMismatch(pols.size(), dim));
419  for (unsigned int d = 0; d < dim; ++d)
420  Assert(pols[d].size() > 0,
421  ExcMessage("The number of polynomials must be larger than zero "
422  "for all coordinate directions."));
423 }
424 
425 
426 
427 template <int dim>
428 void
430  unsigned int (&indices)[dim]) const
431 {
432 #ifdef DEBUG
433  unsigned int n_poly = 1;
434  for (unsigned int d = 0; d < dim; ++d)
435  n_poly *= polynomials[d].size();
436  Assert(i < n_poly, ExcInternalError());
437 #endif
438 
439  if (dim == 1)
440  internal::compute_tensor_index(i,
441  polynomials[0].size(),
442  0 /*not used*/,
443  indices);
444  else
445  internal::compute_tensor_index(i,
446  polynomials[0].size(),
447  polynomials[1].size(),
448  indices);
449 }
450 
451 
452 
453 template <int dim>
454 double
456  const Point<dim> & p) const
457 {
458  unsigned int indices[dim];
459  compute_index(i, indices);
460 
461  double value = 1.;
462  for (unsigned int d = 0; d < dim; ++d)
463  value *= polynomials[d][indices[d]].value(p(d));
464 
465  return value;
466 }
467 
468 
469 template <int dim>
472  const Point<dim> & p) const
473 {
474  unsigned int indices[dim];
475  compute_index(i, indices);
476 
477  // compute values and
478  // uni-directional derivatives at
479  // the given point in each
480  // co-ordinate direction
481  std::vector<std::vector<double>> v(dim, std::vector<double>(2));
482  for (unsigned int d = 0; d < dim; ++d)
483  polynomials[d][indices[d]].value(p(d), v[d]);
484 
485  Tensor<1, dim> grad;
486  for (unsigned int d = 0; d < dim; ++d)
487  {
488  grad[d] = 1.;
489  for (unsigned int x = 0; x < dim; ++x)
490  grad[d] *= v[x][d == x];
491  }
492 
493  return grad;
494 }
495 
496 
497 template <int dim>
500  const Point<dim> & p) const
501 {
502  unsigned int indices[dim];
503  compute_index(i, indices);
504 
505  std::vector<std::vector<double>> v(dim, std::vector<double>(3));
506  for (unsigned int d = 0; d < dim; ++d)
507  polynomials[d][indices[d]].value(p(d), v[d]);
508 
509  Tensor<2, dim> grad_grad;
510  for (unsigned int d1 = 0; d1 < dim; ++d1)
511  for (unsigned int d2 = 0; d2 < dim; ++d2)
512  {
513  grad_grad[d1][d2] = 1.;
514  for (unsigned int x = 0; x < dim; ++x)
515  {
516  unsigned int derivative = 0;
517  if (d1 == x || d2 == x)
518  {
519  if (d1 == d2)
520  derivative = 2;
521  else
522  derivative = 1;
523  }
524  grad_grad[d1][d2] *= v[x][derivative];
525  }
526  }
527 
528  return grad_grad;
529 }
530 
531 
532 
533 template <int dim>
534 void
536  const Point<dim> & p,
537  std::vector<double> & values,
538  std::vector<Tensor<1, dim>> &grads,
539  std::vector<Tensor<2, dim>> &grad_grads,
540  std::vector<Tensor<3, dim>> &third_derivatives,
541  std::vector<Tensor<4, dim>> &fourth_derivatives) const
542 {
543  Assert(values.size() == n_tensor_pols || values.size() == 0,
544  ExcDimensionMismatch2(values.size(), n_tensor_pols, 0));
545  Assert(grads.size() == n_tensor_pols || grads.size() == 0,
546  ExcDimensionMismatch2(grads.size(), n_tensor_pols, 0));
547  Assert(grad_grads.size() == n_tensor_pols || grad_grads.size() == 0,
548  ExcDimensionMismatch2(grad_grads.size(), n_tensor_pols, 0));
549  Assert(third_derivatives.size() == n_tensor_pols ||
550  third_derivatives.size() == 0,
551  ExcDimensionMismatch2(third_derivatives.size(), n_tensor_pols, 0));
552  Assert(fourth_derivatives.size() == n_tensor_pols ||
553  fourth_derivatives.size() == 0,
554  ExcDimensionMismatch2(fourth_derivatives.size(), n_tensor_pols, 0));
555 
556  const bool update_values = (values.size() == n_tensor_pols),
557  update_grads = (grads.size() == n_tensor_pols),
558  update_grad_grads = (grad_grads.size() == n_tensor_pols),
560  (third_derivatives.size() == n_tensor_pols),
561  update_4th_derivatives =
562  (fourth_derivatives.size() == n_tensor_pols);
563 
564  // check how many
565  // values/derivatives we have to
566  // compute
567  unsigned int n_values_and_derivatives = 0;
568  if (update_values)
569  n_values_and_derivatives = 1;
570  if (update_grads)
571  n_values_and_derivatives = 2;
572  if (update_grad_grads)
573  n_values_and_derivatives = 3;
575  n_values_and_derivatives = 4;
576  if (update_4th_derivatives)
577  n_values_and_derivatives = 5;
578 
579  // compute the values (and
580  // derivatives, if necessary) of
581  // all polynomials at this
582  // evaluation point
583  std::vector<std::vector<std::vector<double>>> v(dim);
584  for (unsigned int d = 0; d < dim; ++d)
585  {
586  v[d].resize(polynomials[d].size());
587  for (unsigned int i = 0; i < polynomials[d].size(); ++i)
588  {
589  v[d][i].resize(n_values_and_derivatives, 0.);
590  polynomials[d][i].value(p(d), v[d][i]);
591  };
592  }
593 
594  for (unsigned int i = 0; i < n_tensor_pols; ++i)
595  {
596  // first get the
597  // one-dimensional indices of
598  // this particular tensor
599  // product polynomial
600  unsigned int indices[dim];
601  compute_index(i, indices);
602 
603  if (update_values)
604  {
605  values[i] = 1;
606  for (unsigned int x = 0; x < dim; ++x)
607  values[i] *= v[x][indices[x]][0];
608  }
609 
610  if (update_grads)
611  for (unsigned int d = 0; d < dim; ++d)
612  {
613  grads[i][d] = 1.;
614  for (unsigned int x = 0; x < dim; ++x)
615  grads[i][d] *= v[x][indices[x]][d == x ? 1 : 0];
616  }
617 
618  if (update_grad_grads)
619  for (unsigned int d1 = 0; d1 < dim; ++d1)
620  for (unsigned int d2 = 0; d2 < dim; ++d2)
621  {
622  grad_grads[i][d1][d2] = 1.;
623  for (unsigned int x = 0; x < dim; ++x)
624  {
625  unsigned int derivative = 0;
626  if (d1 == x)
627  ++derivative;
628  if (d2 == x)
629  ++derivative;
630 
631  grad_grads[i][d1][d2] *= v[x][indices[x]][derivative];
632  }
633  }
634 
636  for (unsigned int d1 = 0; d1 < dim; ++d1)
637  for (unsigned int d2 = 0; d2 < dim; ++d2)
638  for (unsigned int d3 = 0; d3 < dim; ++d3)
639  {
640  third_derivatives[i][d1][d2][d3] = 1.;
641  for (unsigned int x = 0; x < dim; ++x)
642  {
643  unsigned int derivative = 0;
644  if (d1 == x)
645  ++derivative;
646  if (d2 == x)
647  ++derivative;
648  if (d3 == x)
649  ++derivative;
650 
651  third_derivatives[i][d1][d2][d3] *=
652  v[x][indices[x]][derivative];
653  }
654  }
655 
656  if (update_4th_derivatives)
657  for (unsigned int d1 = 0; d1 < dim; ++d1)
658  for (unsigned int d2 = 0; d2 < dim; ++d2)
659  for (unsigned int d3 = 0; d3 < dim; ++d3)
660  for (unsigned int d4 = 0; d4 < dim; ++d4)
661  {
662  fourth_derivatives[i][d1][d2][d3][d4] = 1.;
663  for (unsigned int x = 0; x < dim; ++x)
664  {
665  unsigned int derivative = 0;
666  if (d1 == x)
667  ++derivative;
668  if (d2 == x)
669  ++derivative;
670  if (d3 == x)
671  ++derivative;
672  if (d4 == x)
673  ++derivative;
674 
675  fourth_derivatives[i][d1][d2][d3][d4] *=
676  v[x][indices[x]][derivative];
677  }
678  }
679  }
680 }
681 
682 
683 
684 template <int dim>
685 unsigned int
687 {
688  return n_tensor_pols;
689 }
690 
691 
692 template <int dim>
693 unsigned int
695  const std::vector<std::vector<Polynomials::Polynomial<double>>> &pols)
696 {
697  unsigned int y = 1;
698  for (unsigned int d = 0; d < dim; ++d)
699  y *= pols[d].size();
700  return y;
701 }
702 
703 
704 
705 /* ------------------- explicit instantiations -------------- */
709 
710 template class TensorProductPolynomials<
711  1,
713 template class TensorProductPolynomials<
714  2,
715  Polynomials::PiecewisePolynomial<double>>;
716 template class TensorProductPolynomials<
717  3,
718  Polynomials::PiecewisePolynomial<double>>;
719 
720 template class AnisotropicPolynomials<1>;
721 template class AnisotropicPolynomials<2>;
722 template class AnisotropicPolynomials<3>;
723 
724 DEAL_II_NAMESPACE_CLOSE
Shape function values.
double compute_value(const unsigned int i, const Point< dim > &p) const
static::ExceptionBase & ExcDimensionMismatch2(int arg1, int arg2, int arg3)
void compute_index(const unsigned int i, unsigned int(&indices)[(dim > 0?dim:1)]) const
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
void set_numbering(const std::vector< unsigned int > &renumber)
static unsigned int get_n_tensor_pols(const std::vector< std::vector< Polynomials::Polynomial< double >>> &pols)
Definition: point.h:106
void compute(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
void compute(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
static::ExceptionBase & ExcMessage(std::string arg1)
Third derivatives of shape functions.
#define Assert(cond, exc)
Definition: exceptions.h:1227
static::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
const std::vector< std::vector< Polynomials::Polynomial< double > > > polynomials
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const
void output_indices(std::ostream &out) const
static::ExceptionBase & ExcNotImplemented()
AnisotropicPolynomials(const std::vector< std::vector< Polynomials::Polynomial< double >>> &base_polynomials)
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const
void compute_index(const unsigned int i, unsigned int(&indices)[dim]) const
double compute_value(const unsigned int i, const Point< dim > &p) const
static::ExceptionBase & ExcInternalError()