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> 21 DEAL_II_NAMESPACE_OPEN
34 void compute_tensor_index(
const unsigned int,
37 unsigned int ( &)[dim])
39 Assert(
false, ExcNotImplemented());
43 void compute_tensor_index(
const unsigned int n,
46 unsigned int (&indices)[1])
52 void compute_tensor_index(
const unsigned int n,
53 const unsigned int n_pols_0,
55 unsigned int (&indices)[2])
57 indices[0] = n % n_pols_0;
58 indices[1] = n / n_pols_0;
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])
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);
76 template <
int dim,
typename PolynomialType>
81 unsigned int (&indices)[(dim > 0 ? dim : 1)])
const 83 Assert (i<Utilities::fixed_power<dim>(polynomials.size()),ExcInternalError());
84 internal::compute_tensor_index(index_map[i], polynomials.size(),
85 polynomials.size(), indices);
90 template <
int dim,
typename PolynomialType>
95 for (
unsigned int i=0; i<n_tensor_pols; ++i)
99 for (
unsigned int d=0; d<dim; ++d)
107 template <
int dim,
typename PolynomialType>
110 (
const std::vector<unsigned int> &renumber)
112 Assert(renumber.size()==index_map.size(),
113 ExcDimensionMismatch(renumber.size(), index_map.size()));
116 for (
unsigned int i=0; i<index_map.size(); ++i)
117 index_map_inverse[index_map[i]]=i;
125 ::compute_value(
const unsigned int,
128 Assert (
false, ExcNotImplemented());
134 template <
int dim,
typename PolynomialType>
137 (
const unsigned int i,
140 Assert(dim>0, ExcNotImplemented());
142 unsigned int indices[dim];
143 compute_index (i, indices);
146 for (
unsigned int d=0; d<dim; ++d)
147 value *= polynomials[indices[d]].value(p(d));
154 template <
int dim,
typename PolynomialType>
159 unsigned int indices[dim];
160 compute_index (i, indices);
168 std::vector<double> tmp (2);
169 for (
unsigned int d=0; d<dim; ++d)
171 polynomials[indices[d]].value (p(d), tmp);
178 for (
unsigned int d=0; d<dim; ++d)
181 for (
unsigned int x=0; x<dim; ++x)
182 grad[d] *= v[x][d==x];
190 template <
int dim,
typename PolynomialType>
193 (
const unsigned int i,
196 unsigned int indices[dim];
197 compute_index (i, indices);
201 std::vector<double> tmp (3);
202 for (
unsigned int d=0; d<dim; ++d)
204 polynomials[indices[d]].value (p(d), tmp);
212 for (
unsigned int d1=0; d1<dim; ++d1)
213 for (
unsigned int d2=0; d2<dim; ++d2)
215 grad_grad[d1][d2] = 1.;
216 for (
unsigned int x=0; x<dim; ++x)
218 unsigned int derivative=0;
226 grad_grad[d1][d2] *= v[x][derivative];
236 template <
int dim,
typename PolynomialType>
240 std::vector<double> &values,
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));
258 update_grads = (grads.size()==n_tensor_pols),
259 update_grad_grads = (grad_grads.size()==n_tensor_pols),
261 update_4th_derivatives = (fourth_derivatives.size()==n_tensor_pols);
266 unsigned int n_values_and_derivatives = 0;
268 n_values_and_derivatives = 1;
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;
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)
293 polynomials[i].value(p(d), tmp);
294 for (
unsigned int e=0; e<n_values_and_derivatives; ++e)
299 for (
unsigned int i=0; i<n_tensor_pols; ++i)
305 unsigned int indices[dim];
306 compute_index (i, indices);
311 for (
unsigned int x=0; x<dim; ++x)
312 values[i] *= v(x,indices[x])[0];
316 for (
unsigned int d=0; d<dim; ++d)
319 for (
unsigned int x=0; x<dim; ++x)
320 grads[i][d] *= v(x,indices[x])[d==x];
323 if (update_grad_grads)
324 for (
unsigned int d1=0; d1<dim; ++d1)
325 for (
unsigned int d2=0; d2<dim; ++d2)
327 grad_grads[i][d1][d2] = 1.;
328 for (
unsigned int x=0; x<dim; ++x)
330 unsigned int derivative=0;
331 if (d1==x) ++derivative;
332 if (d2==x) ++derivative;
334 grad_grads[i][d1][d2]
335 *= v(x,indices[x])[derivative];
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)
344 third_derivatives[i][d1][d2][d3] = 1.;
345 for (
unsigned int x=0; x<dim; ++x)
347 unsigned int derivative=0;
348 if (d1==x) ++derivative;
349 if (d2==x) ++derivative;
350 if (d3==x) ++derivative;
352 third_derivatives[i][d1][d2][d3]
353 *= v(x,indices[x])[derivative];
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)
363 fourth_derivatives[i][d1][d2][d3][d4] = 1.;
364 for (
unsigned int x=0; x<dim; ++x)
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;
372 fourth_derivatives[i][d1][d2][d3][d4]
373 *= v(x,indices[x])[derivative];
390 n_tensor_pols(get_n_tensor_pols(pols))
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."));
406 unsigned int (&indices)[dim])
const 409 unsigned int n_poly = 1;
410 for (
unsigned int d=0; d<dim; ++d)
412 Assert (i < n_poly, ExcInternalError());
416 internal::compute_tensor_index(i,
polynomials[0].size(),
419 internal::compute_tensor_index(i,
polynomials[0].size(),
430 unsigned int indices[dim];
434 for (
unsigned int d=0; d<dim; ++d)
446 unsigned int indices[dim];
453 std::vector<std::vector<double> > v(dim, std::vector<double> (2));
454 for (
unsigned int d=0; d<dim; ++d)
458 for (
unsigned int d=0; d<dim; ++d)
461 for (
unsigned int x=0; x<dim; ++x)
462 grad[d] *= v[x][d==x];
474 unsigned int indices[dim];
477 std::vector<std::vector<double> > v(dim, std::vector<double> (3));
478 for (
unsigned int d=0; d<dim; ++d)
482 for (
unsigned int d1=0; d1<dim; ++d1)
483 for (
unsigned int d2=0; d2<dim; ++d2)
485 grad_grad[d1][d2] = 1.;
486 for (
unsigned int x=0; x<dim; ++x)
488 unsigned int derivative=0;
496 grad_grad[d1][d2] *= v[x][derivative];
510 std::vector<double> &values,
523 ExcDimensionMismatch2(third_derivatives.size(),
n_tensor_pols, 0));
525 ExcDimensionMismatch2(fourth_derivatives.size(),
n_tensor_pols, 0));
531 update_4th_derivatives = (fourth_derivatives.size()==
n_tensor_pols);
536 unsigned int n_values_and_derivatives = 0;
538 n_values_and_derivatives = 1;
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;
552 std::vector<std::vector<std::vector<double> > > v(dim);
553 for (
unsigned int d=0; d<dim; ++d)
556 for (
unsigned int i=0; i<
polynomials[d].size(); ++i)
558 v[d][i].resize (n_values_and_derivatives, 0.);
569 unsigned int indices[dim];
575 for (
unsigned int x=0; x<dim; ++x)
576 values[i] *= v[x][indices[x]][0];
580 for (
unsigned int d=0; d<dim; ++d)
583 for (
unsigned int x=0; x<dim; ++x)
584 grads[i][d] *= v[x][indices[x]][d==x ? 1 : 0];
587 if (update_grad_grads)
588 for (
unsigned int d1=0; d1<dim; ++d1)
589 for (
unsigned int d2=0; d2<dim; ++d2)
591 grad_grads[i][d1][d2] = 1.;
592 for (
unsigned int x=0; x<dim; ++x)
594 unsigned int derivative=0;
595 if (d1==x) ++derivative;
596 if (d2==x) ++derivative;
598 grad_grads[i][d1][d2]
599 *= v[x][indices[x]][derivative];
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)
608 third_derivatives[i][d1][d2][d3] = 1.;
609 for (
unsigned int x=0; x<dim; ++x)
611 unsigned int derivative=0;
612 if (d1==x) ++derivative;
613 if (d2==x) ++derivative;
614 if (d3==x) ++derivative;
616 third_derivatives[i][d1][d2][d3]
617 *= v[x][indices[x]][derivative];
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)
627 fourth_derivatives[i][d1][d2][d3][d4] = 1.;
628 for (
unsigned int x=0; x<dim; ++x)
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;
636 fourth_derivatives[i][d1][d2][d3][d4]
637 *= v[x][indices[x]][derivative];
659 for (
unsigned int d=0; d<dim; ++d)
679 DEAL_II_NAMESPACE_CLOSE
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
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
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.
unsigned int n_tensor_pols
#define Assert(cond, exc)
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
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