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