35 std::vector<std::vector<Point<1>>>
36 create_anisotropic_support_points(
const unsigned int dim,
37 const unsigned int normal_degree,
38 const unsigned int tangential_degree)
42 std::vector<std::vector<Point<1>>> points_normal_tangential(dim);
43 points_normal_tangential[0] =
46 for (
unsigned int d = 1;
d < dim; ++
d)
48 points_normal_tangential[
d] =
49 tangential_degree == 0 ?
53 return points_normal_tangential;
63 std::vector<std::vector<Polynomials::Polynomial<double>>>
64 create_aniso_polynomials(
const unsigned int dim,
65 const unsigned int normal_degree,
66 const unsigned int tangential_degree)
68 std::vector<std::vector<Point<1>>> points_aniso;
70 create_anisotropic_support_points(dim, normal_degree, tangential_degree);
71 std::vector<std::vector<Polynomials::Polynomial<double>>> pols(dim);
73 for (
unsigned int d = 1;
d < dim; ++
d)
85 const std::vector<unsigned int> &polynomial_ordering)
105 for (
unsigned int i = 0; i <
n_pols; ++i)
152 Assert(values.size() == this->n() || values.empty(),
154 Assert(grads.size() == this->n() || grads.empty(),
156 Assert(grad_grads.size() == this->n() || grad_grads.empty(),
158 Assert(third_derivatives.size() == this->n() || third_derivatives.empty(),
160 Assert(fourth_derivatives.size() == this->n() || fourth_derivatives.empty(),
163 std::vector<double> p_values;
164 std::vector<Tensor<1, dim>> p_grads;
165 std::vector<Tensor<2, dim>> p_grad_grads;
166 std::vector<Tensor<3, dim>> p_third_derivatives;
167 std::vector<Tensor<4, dim>> p_fourth_derivatives;
170 p_values.resize((values.empty()) ? 0 : n_sub);
171 p_grads.resize((grads.empty()) ? 0 : n_sub);
172 p_grad_grads.resize((grad_grads.empty()) ? 0 : n_sub);
173 p_third_derivatives.resize((third_derivatives.empty()) ? 0 : n_sub);
174 p_fourth_derivatives.resize((fourth_derivatives.empty()) ? 0 : n_sub);
176 for (
unsigned int d = 0; d < dim; ++d)
184 for (
unsigned int c = 0; c < dim; ++c)
185 p[c] = unit_point[(c + d) % dim];
192 p_fourth_derivatives);
194 for (
unsigned int i = 0; i < p_values.size(); ++i)
198 for (
unsigned int i = 0; i < p_grads.size(); ++i)
199 for (
unsigned int d1 = 0; d1 < dim; ++d1)
203 for (
unsigned int i = 0; i < p_grad_grads.size(); ++i)
204 for (
unsigned int d1 = 0; d1 < dim; ++d1)
205 for (
unsigned int d2 = 0; d2 < dim; ++d2)
207 [(d1 + d) % dim][(d2 + d) % dim] =
210 for (
unsigned int i = 0; i < p_third_derivatives.size(); ++i)
211 for (
unsigned int d1 = 0; d1 < dim; ++d1)
212 for (
unsigned int d2 = 0; d2 < dim; ++d2)
213 for (
unsigned int d3 = 0; d3 < dim; ++d3)
215 [(d1 + d) % dim][(d2 + d) % dim]
220 for (
unsigned int i = 0; i < p_fourth_derivatives.size(); ++i)
221 for (
unsigned int d1 = 0; d1 < dim; ++d1)
222 for (
unsigned int d2 = 0; d2 < dim; ++d2)
223 for (
unsigned int d3 = 0; d3 < dim; ++d3)
224 for (
unsigned int d4 = 0; d4 < dim; ++d4)
226 [d][(d1 + d) % dim][(d2 + d) % dim]
227 [(d3 + d) % dim][(d4 + d) % dim] =
239 return "VectorAnisotropic";
275std::unique_ptr<TensorPolynomialsBase<dim>>
278 return std::make_unique<PolynomialsVectorAnisotropic<dim>>(*this);
284std::vector<Point<dim>>
288 const std::vector<std::vector<Point<1>>> points_aniso =
291 std::vector<Point<dim>> points(dim * n_sub);
293 for (
unsigned int d = 0; d < dim; ++d)
294 for (
unsigned int i = 0; i < n_sub; ++i)
297 std::array<unsigned int, dim> indices_points_anisotropic;
298 indices_points_anisotropic[0] = renumbered_index % (
normal_degree + 1);
302 indices_points_anisotropic[1] =
306 indices_points_anisotropic[2] =
308 for (
unsigned int c = 0; c < dim; ++c)
311 points_aniso[c][indices_points_anisotropic[c]][0];
std::vector< unsigned int > hierarchic_to_lexicographic
const unsigned int tangential_degree
std::array< std::vector< unsigned int >, dim > renumber_aniso
PolynomialsVectorAnisotropic(const unsigned int degree_normal, const unsigned int degree_tangential, const std::vector< unsigned int > &polynomial_ordering)
const unsigned int normal_degree
const AnisotropicPolynomials< dim > polynomial_space
std::vector< Point< dim > > get_polynomial_support_points() const
virtual std::unique_ptr< TensorPolynomialsBase< dim > > clone() const override
unsigned int get_normal_degree() const
unsigned int get_tangential_degree() const
void evaluate(const Point< dim > &unit_point, std::vector< Tensor< 1, dim > > &values, std::vector< Tensor< 2, dim > > &grads, std::vector< Tensor< 3, dim > > &grad_grads, std::vector< Tensor< 4, dim > > &third_derivatives, std::vector< Tensor< 5, dim > > &fourth_derivatives) const override
std::vector< unsigned int > lexicographic_to_hierarchic
std::string name() const override
static unsigned int n_polynomials(const unsigned int normal_degree, const unsigned int tangential_degree)
const unsigned int n_pols
TensorPolynomialsBase(const unsigned int deg, const unsigned int n_polynomials)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 > > &points)
constexpr T pow(const T base, const int iexp)