40 std::size_t max_degree = 0;
41 for (
const auto &poly : polys)
45 const auto degree_0 = degrees[0];
46 std::size_t degree_d = 0;
47 for (
unsigned int d = 1; d < dim + 1; ++d)
48 degree_d =
std::max(degree_d, degrees[d]);
50 max_degree =
std::max(max_degree, degree_d + degree_0);
62 std::vector<PolyType> polys;
64 const auto reference_cell = ReferenceCells::get_simplex<dim>();
74 for (
const unsigned int v : reference_cell.vertex_indices())
81 for (
const unsigned int v : reference_cell.vertex_indices())
85 for (
const unsigned int l : reference_cell.line_indices())
87 const auto v0 = reference_cell.line_to_cell_vertices(l, 0);
88 const auto v1 = reference_cell.line_to_cell_vertices(l, 1);
106 const std::vector<PolyType> &polynomials)
117 for (std::size_t i = 0; i < polynomials.size(); ++i)
120 for (
unsigned int d = 0; d < dim; ++d)
121 poly_grads[i][d] = polynomials[i].derivative(d);
124 for (
unsigned int d0 = 0; d0 < dim; ++d0)
125 for (
unsigned int d1 = 0; d1 < dim; ++d1)
129 for (
unsigned int d0 = 0; d0 < dim; ++d0)
130 for (
unsigned int d1 = 0; d1 < dim; ++d1)
131 for (
unsigned int d2 = 0; d2 < dim; ++d2)
136 for (
unsigned int d0 = 0; d0 < dim; ++d0)
137 for (
unsigned int d1 = 0; d1 < dim; ++d1)
138 for (
unsigned int d2 = 0; d2 < dim; ++d2)
139 for (
unsigned int d3 = 0; d3 < dim; ++d3)
151 std::vector<double> & values,
157 Assert(values.size() == this->n() || values.size() == 0,
159 Assert(grads.size() == this->n() || grads.size() == 0,
161 Assert(grad_grads.size() == this->n() || grad_grads.size() == 0,
163 Assert(third_derivatives.size() == this->n() || third_derivatives.size() == 0,
165 Assert(fourth_derivatives.size() == this->n() ||
166 fourth_derivatives.size() == 0,
169 for (std::size_t i = 0; i < polys.size(); ++i)
171 if (values.size() == this->n())
172 values[i] = polys[i].value(unit_point);
175 if (grads.size() == this->n())
176 for (
unsigned int d = 0; d < dim; ++d)
177 grads[i][d] = poly_grads[i][d].value(unit_point);
180 if (grad_grads.size() == this->n())
181 for (
unsigned int d0 = 0; d0 < dim; ++d0)
182 for (
unsigned int d1 = 0; d1 < dim; ++d1)
183 grad_grads[i][d0][d1] = poly_hessians[i][d0][d1].value(unit_point);
186 if (third_derivatives.size() == this->n())
187 for (
unsigned int d0 = 0; d0 < dim; ++d0)
188 for (
unsigned int d1 = 0; d1 < dim; ++d1)
189 for (
unsigned int d2 = 0; d2 < dim; ++d2)
190 third_derivatives[i][d0][d1][d2] =
191 poly_third_derivatives[i][d0][d1][d2].value(unit_point);
194 if (fourth_derivatives.size() == this->n())
195 for (
unsigned int d0 = 0; d0 < dim; ++d0)
196 for (
unsigned int d1 = 0; d1 < dim; ++d1)
197 for (
unsigned int d2 = 0; d2 < dim; ++d2)
198 for (
unsigned int d3 = 0; d3 < dim; ++d3)
199 fourth_derivatives[i][d0][d1][d2][d3] =
200 poly_fourth_derivatives[i][d0][d1][d2][d3].value(unit_point);
212 return polys[i].value(p);
223 for (
unsigned int d = 0; d < dim; ++d)
224 result[d] = poly_grads[i][d].value(p);
236 for (
unsigned int d0 = 0; d0 < dim; ++d0)
237 for (
unsigned int d1 = 0; d1 < dim; ++d1)
238 result[d0][d1] = poly_hessians[i][d0][d1].value(p);
251 for (
unsigned int d0 = 0; d0 < dim; ++d0)
252 for (
unsigned int d1 = 0; d1 < dim; ++d1)
253 for (
unsigned int d2 = 0; d2 < dim; ++d2)
254 result[d0][d1][d2] = poly_third_derivatives[i][d0][d1][d2].value(p);
267 for (
unsigned int d0 = 0; d0 < dim; ++d0)
268 for (
unsigned int d1 = 0; d1 < dim; ++d1)
269 for (
unsigned int d2 = 0; d2 < dim; ++d2)
270 for (
unsigned int d3 = 0; d3 < dim; ++d3)
271 result[d0][d1][d2][d3] =
272 poly_fourth_derivatives[i][d0][d1][d2][d3].value(p);
284 return compute_1st_derivative(i, p);
294 return compute_2nd_derivative(i, p);
300std::unique_ptr<ScalarPolynomialsBase<dim>>
303 return std::make_unique<BarycentricPolynomials<dim>>(*this);
312 return "BarycentricPolynomials<" + std::to_string(dim) +
">";
321 std::size_t poly_memory = 0;
322 for (
const auto &poly : polys)
323 poly_memory += poly.memory_consumption();
static BarycentricPolynomial< dim, Number > monomial(const unsigned int d)
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const override
BarycentricPolynomials(const std::vector< BarycentricPolynomial< dim > > &polynomials)
virtual std::size_t memory_consumption() const override
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const override
std::vector< GradType > poly_grads
std::string name() const override
Tensor< 3, dim > compute_3rd_derivative(const unsigned int i, const Point< dim > &p) const override
Tensor< 1, dim > compute_1st_derivative(const unsigned int i, const Point< dim > &p) const override
Tensor< 2, dim > compute_2nd_derivative(const unsigned int i, const Point< dim > &p) const override
std::vector< PolyType > polys
void evaluate(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 override
double compute_value(const unsigned int i, const Point< dim > &p) const override
std::vector< ThirdDerivativesType > poly_third_derivatives
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone() const override
Tensor< 4, dim > compute_4th_derivative(const unsigned int i, const Point< dim > &p) const override
std::vector< HessianType > poly_hessians
static BarycentricPolynomials< dim > get_fe_p_basis(const unsigned int degree)
std::vector< FourthDerivativesType > poly_fourth_derivatives
virtual std::size_t memory_consumption() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch2(std::size_t arg1, std::size_t arg2, std::size_t arg3)
#define AssertIndexRange(index, range)
std::enable_if_t< std::is_fundamental< T >::value, std::size_t > memory_consumption(const T &t)
unsigned int get_degree(const std::vector< typename BarycentricPolynomials< dim >::PolyType > &polys)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)