16 #ifndef dealii_fe_poly_tensor_h 17 #define dealii_fe_poly_tensor_h 140 template <
int dim,
int spacedim = dim>
178 const unsigned int component)
const override;
193 const unsigned int component)
const override;
208 const unsigned int component)
const override;
236 virtual std::unique_ptr<
244 &output_data)
const override 250 std::unique_ptr<typename FiniteElement<dim, spacedim>::InternalDataBase>
251 data_ptr = std::make_unique<InternalData>();
252 auto &data =
dynamic_cast<InternalData &
>(*data_ptr);
255 const unsigned int n_q_points = quadrature.
size();
258 std::vector<Tensor<1, dim>>
values(0);
259 std::vector<Tensor<2, dim>> grads(0);
260 std::vector<Tensor<3, dim>> grad_grads(0);
261 std::vector<Tensor<4, dim>> third_derivatives(0);
262 std::vector<Tensor<5, dim>> fourth_derivatives(0);
271 const bool update_transformed_shape_values =
272 std::any_of(this->mapping_kind.begin(),
273 this->mapping_kind.end(),
276 const bool update_transformed_shape_grads =
277 std::any_of(this->mapping_kind.begin(),
278 this->mapping_kind.end(),
284 const bool update_transformed_shape_hessian_tensors =
285 update_transformed_shape_values;
291 if (update_transformed_shape_values)
292 data.transformed_shape_values.resize(n_q_points);
299 data.transformed_shape_grads.resize(n_q_points);
301 if (update_transformed_shape_grads)
302 data.untransformed_shape_grads.resize(n_q_points);
309 data.transformed_shape_hessians.resize(n_q_points);
310 if (update_transformed_shape_hessian_tensors)
311 data.untransformed_shape_hessian_tensors.resize(n_q_points);
321 if (update_flags & (update_values | update_gradients))
322 for (
unsigned int k = 0; k < n_q_points; ++k)
331 if (update_flags & update_values)
335 data.shape_values[i][k] = values[i];
342 data.shape_values[i][k] = add_values;
346 if (update_flags & update_gradients)
350 data.shape_grads[i][k] = grads[i];
357 data.shape_grads[i][k] = add_grads;
361 if (update_flags & update_hessians)
365 data.shape_grad_grads[i][k] = grad_grads[i];
373 data.shape_grad_grads[i][k] = add_grad_grads;
387 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
393 &output_data)
const override;
400 const unsigned int face_no,
404 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
410 &output_data)
const override;
415 const unsigned int face_no,
416 const unsigned int sub_no,
420 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
426 &output_data)
const override;
479 const std::unique_ptr<const TensorPolynomialsBase<dim>>
poly_space;
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
std::vector< Tensor< 2, dim > > cached_grads
MappingKind get_mapping_kind(const unsigned int i) const
bool single_mapping_kind() const
std::vector< Tensor< 1, dim > > cached_values
std::vector< double > sign_change
const Point< dim > & point(const unsigned int i) const
FullMatrix< double > inverse_node_matrix
Table< 2, DerivativeForm< 1, dim, spacedim > > shape_grads
std::vector< Tensor< 2, dim > > untransformed_shape_grads
std::vector< Tensor< 3, spacedim > > transformed_shape_hessians
Table< 2, DerivativeForm< 2, dim, spacedim > > shape_grad_grads
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
std::vector< Tensor< 3, dim > > untransformed_shape_hessian_tensors
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
Abstract base class for mapping classes.
#define DEAL_II_NAMESPACE_CLOSE
Point< dim > cached_point
std::vector< Tensor< 1, spacedim > > transformed_shape_values
std::vector< Tensor< 3, dim > > cached_grad_grads
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
Second derivatives of shape functions.
unsigned int size() const
const std::unique_ptr< const TensorPolynomialsBase< dim > > poly_space
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
unsigned int n_dofs_per_cell() const
#define DEAL_II_NAMESPACE_OPEN
Shape function gradients.
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override
Table< 2, Tensor< 1, dim > > shape_values
const std::vector< ComponentMask > nonzero_components
FE_PolyTensor(const TensorPolynomialsBase< dim > &polynomials, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
std::vector< Tensor< 2, spacedim > > transformed_shape_grads
const std::vector< bool > restriction_is_additive_flags
std::vector< MappingKind > mapping_kind
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FiniteElement< dim, spacedim >::InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const override