16 #ifndef dealii_polynomial_space_h 17 #define dealii_polynomial_space_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/exceptions.h> 22 #include <deal.II/base/tensor.h> 23 #include <deal.II/base/point.h> 24 #include <deal.II/base/polynomial.h> 25 #include <deal.II/base/smartpointer.h> 29 DEAL_II_NAMESPACE_OPEN
119 template <
class StreamType>
126 void set_numbering(
const std::vector<unsigned int> &renumber);
142 std::vector<double> &values,
193 unsigned int n ()
const;
201 unsigned int degree ()
const;
221 unsigned int (&index)[dim>0?dim:1])
const;
250 unsigned int (&index)[1])
const;
253 unsigned int (&index)[2])
const;
256 unsigned int (&index)[3])
const;
266 polynomials (pols.begin(), pols.end()),
267 n_pols (compute_n_pols(polynomials.size())),
269 index_map_inverse(n_pols)
275 for (
unsigned int i=0; i<
n_pols; ++i)
298 return polynomials.size();
303 template <
class StreamType>
307 unsigned int ix[dim];
308 for (
unsigned int i=0; i<n_pols; ++i)
312 for (
unsigned int d=0; d<dim; ++d)
324 unsigned int indices[dim];
325 compute_index (i, indices);
327 double v [dim][order+1];
329 std::vector<double> tmp (order+1);
330 for (
unsigned int d=0; d<dim; ++d)
332 polynomials[indices[d]].value (p(d), tmp);
333 for (
unsigned int j=0; j<order+1; ++j)
344 for (
unsigned int d=0; d<dim; ++d)
346 derivative_1[d] = 1.;
347 for (
unsigned int x=0; x<dim; ++x)
349 unsigned int x_order=0;
352 derivative_1[d] *= v[x][x_order];
361 for (
unsigned int d1=0; d1<dim; ++d1)
362 for (
unsigned int d2=0; d2<dim; ++d2)
364 derivative_2[d1][d2] = 1.;
365 for (
unsigned int x=0; x<dim; ++x)
367 unsigned int x_order=0;
368 if (d1==x) ++x_order;
369 if (d2==x) ++x_order;
371 derivative_2[d1][d2] *= v[x][x_order];
380 for (
unsigned int d1=0; d1<dim; ++d1)
381 for (
unsigned int d2=0; d2<dim; ++d2)
382 for (
unsigned int d3=0; d3<dim; ++d3)
384 derivative_3[d1][d2][d3] = 1.;
385 for (
unsigned int x=0; x<dim; ++x)
387 unsigned int x_order=0;
388 if (d1==x) ++x_order;
389 if (d2==x) ++x_order;
390 if (d3==x) ++x_order;
392 derivative_3[d1][d2][d3] *= v[x][x_order];
401 for (
unsigned int d1=0; d1<dim; ++d1)
402 for (
unsigned int d2=0; d2<dim; ++d2)
403 for (
unsigned int d3=0; d3<dim; ++d3)
404 for (
unsigned int d4=0; d4<dim; ++d4)
406 derivative_4[d1][d2][d3][d4] = 1.;
407 for (
unsigned int x=0; x<dim; ++x)
409 unsigned int x_order=0;
410 if (d1==x) ++x_order;
411 if (d2==x) ++x_order;
412 if (d3==x) ++x_order;
413 if (d4==x) ++x_order;
415 derivative_4[d1][d2][d3][d4] *= v[x][x_order];
431 DEAL_II_NAMESPACE_CLOSE
void set_numbering(const std::vector< unsigned int > &renumber)
static const unsigned int dimension
std::vector< unsigned int > index_map_inverse
#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
const unsigned int n_pols
Tensor< order, dim > compute_derivative(const unsigned int i, const Point< dim > &p) const
void output_indices(StreamType &out) 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
static unsigned int compute_n_pols(const unsigned int n)
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const
unsigned int degree() const
PolynomialSpace(const std::vector< Pol > &pols)
static ::ExceptionBase & ExcNotImplemented()
std::vector< unsigned int > index_map
const std::vector< Polynomials::Polynomial< double > > polynomials
void compute_index(const unsigned int n, unsigned int(&index)[dim >0?dim:1]) const