29 unsigned int n_pols = n;
30 for (
unsigned int i = 1; i < dim; ++i)
48std::array<unsigned int, 1>
52 return {{index_map[i]}};
58std::array<unsigned int, 2>
62 const unsigned int n = index_map[i];
67 const unsigned int n_1d = polynomials.size();
69 for (
unsigned int iy = 0; iy < n_1d; ++iy)
70 if (n < k + n_1d - iy)
84std::array<unsigned int, 3>
88 const unsigned int n = index_map[i];
96 const unsigned int n_1d = polynomials.size();
98 for (
unsigned int iz = 0; iz < n_1d; ++iz)
99 for (
unsigned int iy = 0; iy < n_1d - iz; ++iy)
100 if (n < k + n_1d - iy - iz)
102 return {{n - k, iy, iz}};
116 Assert(renumber.size() == index_map.size(),
119 index_map = renumber;
120 for (
unsigned int i = 0; i < index_map.size(); ++i)
121 index_map_inverse[index_map[i]] = i;
131 const auto ix = compute_index(i);
136 for (
unsigned int d = 0; d < dim; ++d)
137 result *= polynomials[ix[d]].value(p(d));
148 const auto ix = compute_index(i);
151 for (
unsigned int d = 0; d < dim; ++d)
155 std::vector<double> v(2);
156 for (
unsigned int d = 0; d < dim; ++d)
158 polynomials[ix[d]].value(p(d), v);
160 for (
unsigned int d1 = 0; d1 < dim; ++d1)
173 const auto ix = compute_index(i);
176 for (
unsigned int d = 0; d < dim; ++d)
177 for (
unsigned int d1 = 0; d1 < dim; ++d1)
182 std::vector<double> v(3);
183 for (
unsigned int d = 0; d < dim; ++d)
185 polynomials[ix[d]].value(p(d), v);
186 result[d][d] *= v[2];
187 for (
unsigned int d1 = 0; d1 < dim; ++d1)
191 result[d][d1] *= v[1];
192 result[d1][d] *= v[1];
193 for (
unsigned int d2 = 0; d2 < dim; ++d2)
195 result[d1][d2] *= v[0];
207 std::vector<double> & values,
213 const unsigned int n_1d = polynomials.size();
215 Assert(values.size() == this->n() || values.size() == 0,
217 Assert(grads.size() == this->n() || grads.size() == 0,
219 Assert(grad_grads.size() == this->n() || grad_grads.size() == 0,
221 Assert(third_derivatives.size() == this->n() || third_derivatives.size() == 0,
223 Assert(fourth_derivatives.size() == this->n() ||
224 fourth_derivatives.size() == 0,
227 unsigned int v_size = 0;
228 bool update_values =
false, update_grads =
false, update_grad_grads =
false;
230 if (values.size() == this->n())
235 if (grads.size() == this->n())
240 if (grad_grads.size() == this->n())
242 update_grad_grads =
true;
245 if (third_derivatives.size() == this->n())
250 if (fourth_derivatives.size() == this->n())
252 update_4th_derivatives =
true;
263 for (
unsigned int d = 0;
d < v.size()[0]; ++
d)
264 for (
unsigned int i = 0; i < v.size()[1]; ++i)
266 v(d, i).resize(v_size, 0.);
267 polynomials[i].value(p(d), v(d, i));
274 for (
unsigned int iz = 0; iz < ((dim > 2) ? n_1d : 1); ++iz)
275 for (
unsigned int iy = 0; iy < ((dim > 1) ? n_1d - iz : 1); ++iy)
276 for (
unsigned int ix = 0; ix < n_1d - iy - iz; ++ix)
277 values[index_map_inverse[k++]] = v[0][ix][0] *
278 ((dim > 1) ? v[1][iy][0] : 1.) *
279 ((dim > 2) ? v[2][iz][0] : 1.);
286 for (
unsigned int iz = 0; iz < ((dim > 2) ? n_1d : 1); ++iz)
287 for (
unsigned int iy = 0; iy < ((dim > 1) ? n_1d - iz : 1); ++iy)
288 for (
unsigned int ix = 0; ix < n_1d - iy - iz; ++ix)
290 const unsigned int k2 = index_map_inverse[k++];
291 for (
unsigned int d = 0;
d < dim; ++
d)
292 grads[k2][d] = v[0][ix][(d == 0) ? 1 : 0] *
293 ((dim > 1) ? v[1][iy][(d == 1) ? 1 : 0] : 1.) *
294 ((dim > 2) ? v[2][iz][(
d == 2) ? 1 : 0] : 1.);
298 if (update_grad_grads)
302 for (
unsigned int iz = 0; iz < ((dim > 2) ? n_1d : 1); ++iz)
303 for (
unsigned int iy = 0; iy < ((dim > 1) ? n_1d - iz : 1); ++iy)
304 for (
unsigned int ix = 0; ix < n_1d - iy - iz; ++ix)
306 const unsigned int k2 = index_map_inverse[k++];
307 for (
unsigned int d1 = 0; d1 < dim; ++d1)
308 for (
unsigned int d2 = 0; d2 < dim; ++d2)
313 const unsigned int j0 =
314 ((d1 == 0) ? 1 : 0) + ((d2 == 0) ? 1 : 0);
315 const unsigned int j1 =
316 ((d1 == 1) ? 1 : 0) + ((d2 == 1) ? 1 : 0);
317 const unsigned int j2 =
318 ((d1 == 2) ? 1 : 0) + ((d2 == 2) ? 1 : 0);
320 grad_grads[k2][d1][d2] = v[0][ix][j0] *
321 ((dim > 1) ? v[1][iy][j1] : 1.) *
322 ((dim > 2) ? v[2][iz][j2] : 1.);
331 for (
unsigned int iz = 0; iz < ((dim > 2) ? n_1d : 1); ++iz)
332 for (
unsigned int iy = 0; iy < ((dim > 1) ? n_1d - iz : 1); ++iy)
333 for (
unsigned int ix = 0; ix < n_1d - iy - iz; ++ix)
335 const unsigned int k2 = index_map_inverse[k++];
336 for (
unsigned int d1 = 0; d1 < dim; ++d1)
337 for (
unsigned int d2 = 0; d2 < dim; ++d2)
338 for (
unsigned int d3 = 0; d3 < dim; ++d3)
343 std::vector<unsigned int> deriv_order(dim, 0);
344 for (
unsigned int x = 0; x < dim; ++x)
354 third_derivatives[k2][d1][d2][d3] =
355 v[0][ix][deriv_order[0]] *
356 ((dim > 1) ? v[1][iy][deriv_order[1]] : 1.) *
357 ((dim > 2) ? v[2][iz][deriv_order[2]] : 1.);
362 if (update_4th_derivatives)
366 for (
unsigned int iz = 0; iz < ((dim > 2) ? n_1d : 1); ++iz)
367 for (
unsigned int iy = 0; iy < ((dim > 1) ? n_1d - iz : 1); ++iy)
368 for (
unsigned int ix = 0; ix < n_1d - iy - iz; ++ix)
370 const unsigned int k2 = index_map_inverse[k++];
371 for (
unsigned int d1 = 0; d1 < dim; ++d1)
372 for (
unsigned int d2 = 0; d2 < dim; ++d2)
373 for (
unsigned int d3 = 0; d3 < dim; ++d3)
374 for (
unsigned int d4 = 0; d4 < dim; ++d4)
379 std::vector<unsigned int> deriv_order(dim, 0);
380 for (
unsigned int x = 0; x < dim; ++x)
392 fourth_derivatives[k2][d1][d2][d3][d4] =
393 v[0][ix][deriv_order[0]] *
394 ((dim > 1) ? v[1][iy][deriv_order[1]] : 1.) *
395 ((dim > 2) ? v[2][iz][deriv_order[2]] : 1.);
404std::unique_ptr<ScalarPolynomialsBase<dim>>
407 return std::make_unique<PolynomialSpace<dim>>(*this);
double compute_value(const unsigned int i, const Point< dim > &p) const override
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const override
void set_numbering(const std::vector< unsigned int > &renumber)
static unsigned int n_polynomials(const unsigned int n)
std::array< unsigned int, dim > compute_index(const unsigned int n) const
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone() const override
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
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch2(std::size_t arg1, std::size_t arg2, std::size_t arg3)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
@ update_values
Shape function values.
@ update_3rd_derivatives
Third derivatives of shape functions.
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
static const unsigned int invalid_unsigned_int