42 template <
int spacedim>
44 get_dof_sign_change_h_div(
45 const typename ::Triangulation<1, spacedim>::cell_iterator &,
47 const std::vector<MappingKind> &,
48 std::vector<double> &)
61 template <
int spacedim>
63 get_dof_sign_change_h_div(
64 const typename ::Triangulation<2, spacedim>::cell_iterator &cell,
66 const std::vector<MappingKind> &mapping_kind,
67 std::vector<double> &face_sign)
69 const unsigned int dim = 2;
72 const CellId this_cell_id = cell->id();
74 for (
unsigned int f = 0; f < GeometryInfo<dim>::faces_per_cell; ++f)
76 typename ::Triangulation<dim, spacedim>::face_iterator face =
78 if (!face->at_boundary())
80 const unsigned int nn = cell->neighbor_face_no(f);
81 const typename ::Triangulation<dim,
82 spacedim>::cell_iterator
83 neighbor_cell_at_face = cell->neighbor(f);
84 const CellId neighbor_cell_id = neighbor_cell_at_face->id();
88 if (((nn + f) % 2 == 0) && this_cell_id < neighbor_cell_id)
95 Assert(mapping_kind.size() == 1 ||
96 cell_j < mapping_kind.size(),
101 if ((mapping_kind.size() > 1 ?
102 mapping_kind[cell_j] :
112 template <
int spacedim>
114 get_dof_sign_change_h_div(
115 const typename ::Triangulation<3, spacedim>::cell_iterator
118 const std::vector<MappingKind> & ,
119 std::vector<double> & )
125 template <
int spacedim>
127 get_dof_sign_change_nedelec(
128 const typename ::Triangulation<1, spacedim>::cell_iterator
131 const std::vector<MappingKind> & ,
132 std::vector<double> & )
137 template <
int spacedim>
139 get_dof_sign_change_nedelec(
140 const typename ::Triangulation<2, spacedim>::cell_iterator &cell,
142 const std::vector<MappingKind> &mapping_kind,
143 std::vector<double> &line_dof_sign)
168 for (
unsigned int l = 0; l < GeometryInfo<2>::lines_per_cell; ++l)
169 if (cell->line_orientation(l) !=
177 line_dof_sign[l] = -1.0;
196 for (
unsigned int local_line_dof = 0;
197 local_line_dof < (k + 1);
199 if (local_line_dof % 2 == 0)
200 line_dof_sign[local_line_dof + l * (k + 1)] = -1.0;
205 template <
int spacedim>
207 get_dof_sign_change_nedelec(
208 const typename ::Triangulation<3, spacedim>::cell_iterator &cell,
210 const std::vector<MappingKind> &mapping_kind,
211 std::vector<double> &line_dof_sign)
239 for (
unsigned int l = 0; l < GeometryInfo<3>::lines_per_cell; ++l)
240 if (cell->line_orientation(l) !=
248 line_dof_sign[l] = -1.0;
267 for (
unsigned int local_line_dof = 0;
268 local_line_dof < (k + 1);
270 if (local_line_dof % 2 == 0)
271 line_dof_sign[local_line_dof + l * (k + 1)] = -1.0;
281template <
int dim,
int spacedim>
285 const std::vector<bool> &restriction_is_additive_flags,
286 const std::vector<ComponentMask> &nonzero_components)
288 restriction_is_additive_flags,
291 , poly_space(polynomials.clone())
293 cached_point[0] = -1;
301 for (
unsigned int comp = 0; comp < this->n_components(); ++comp)
302 this->component_to_base_table[comp].first.second = comp;
306 adjust_quad_dof_sign_for_face_orientation_table.resize(
307 this->n_unique_2d_subobjects());
309 for (
unsigned int f = 0; f < this->n_unique_2d_subobjects(); ++f)
311 adjust_quad_dof_sign_for_face_orientation_table[f] =
312 Table<2, bool>(this->n_dofs_per_quad(f),
314 adjust_quad_dof_sign_for_face_orientation_table[f].fill(
false);
321template <
int dim,
int spacedim>
323 : FiniteElement<dim, spacedim>(fe)
324 , mapping_kind(fe.mapping_kind)
325 , adjust_quad_dof_sign_for_face_orientation_table(
326 fe.adjust_quad_dof_sign_for_face_orientation_table)
327 , poly_space(fe.poly_space->clone())
328 , inverse_node_matrix(fe.inverse_node_matrix)
333template <
int dim,
int spacedim>
337 return mapping_kind.size() == 1;
341template <
int dim,
int spacedim>
344 const unsigned int index,
345 const unsigned int face,
359 Assert(adjust_quad_dof_sign_for_face_orientation_table
360 [this->n_unique_2d_subobjects() == 1 ? 0 : face]
363 this->n_dofs_per_quad(face),
366 return adjust_quad_dof_sign_for_face_orientation_table
367 [this->n_unique_2d_subobjects() == 1 ? 0 : face](
index,
368 combined_orientation);
372template <
int dim,
int spacedim>
376 if (single_mapping_kind())
377 return mapping_kind[0];
380 return mapping_kind[i];
385template <
int dim,
int spacedim>
388 const Point<dim> &)
const
397template <
int dim,
int spacedim>
400 const unsigned int i,
402 const unsigned int component)
const
407 std::lock_guard<std::mutex> lock(cache_mutex);
409 if (cached_point != p || cached_values.empty())
412 cached_values.resize(poly_space->n());
414 std::vector<Tensor<4, dim>> dummy1;
415 std::vector<Tensor<5, dim>> dummy2;
416 poly_space->evaluate(
417 p, cached_values, cached_grads, cached_grad_grads, dummy1, dummy2);
421 if (inverse_node_matrix.n_cols() == 0)
422 return cached_values[i][component];
424 for (
unsigned int j = 0; j < inverse_node_matrix.n_cols(); ++j)
425 s += inverse_node_matrix(j, i) * cached_values[j][component];
431template <
int dim,
int spacedim>
434 const Point<dim> &)
const
437 return Tensor<1, dim>();
442template <
int dim,
int spacedim>
445 const unsigned int i,
447 const unsigned int component)
const
452 std::lock_guard<std::mutex> lock(cache_mutex);
454 if (cached_point != p || cached_grads.empty())
457 cached_grads.resize(poly_space->n());
459 std::vector<Tensor<4, dim>> dummy1;
460 std::vector<Tensor<5, dim>> dummy2;
461 poly_space->evaluate(
462 p, cached_values, cached_grads, cached_grad_grads, dummy1, dummy2);
466 if (inverse_node_matrix.n_cols() == 0)
467 return cached_grads[i][component];
469 for (
unsigned int j = 0; j < inverse_node_matrix.n_cols(); ++j)
470 s += inverse_node_matrix(j, i) * cached_grads[j][component];
477template <
int dim,
int spacedim>
480 const Point<dim> &)
const
483 return Tensor<2, dim>();
488template <
int dim,
int spacedim>
491 const unsigned int i,
493 const unsigned int component)
const
498 std::lock_guard<std::mutex> lock(cache_mutex);
500 if (cached_point != p || cached_grad_grads.empty())
503 cached_grad_grads.resize(poly_space->n());
505 std::vector<Tensor<4, dim>> dummy1;
506 std::vector<Tensor<5, dim>> dummy2;
507 poly_space->evaluate(
508 p, cached_values, cached_grads, cached_grad_grads, dummy1, dummy2);
512 if (inverse_node_matrix.n_cols() == 0)
513 return cached_grad_grads[i][component];
515 for (
unsigned int j = 0; j < inverse_node_matrix.n_cols(); ++j)
516 s += inverse_node_matrix(i, j) * cached_grad_grads[j][component];
526template <
int dim,
int spacedim>
531 const Quadrature<dim> &quadrature,
532 const Mapping<dim, spacedim> &
mapping,
533 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
534 const internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
537 ::internal::FEValuesImplementation::FiniteElementRelatedData<dim,
541 Assert(
dynamic_cast<const InternalData *
>(&fe_internal) !=
nullptr,
550 static_cast<const InternalData &
>(fe_internal),
556template <
int dim,
int spacedim>
560 const unsigned int face_no,
561 const hp::QCollection<dim - 1> &quadrature,
562 const Mapping<dim, spacedim> &
mapping,
563 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
564 const internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
567 ::internal::FEValuesImplementation::FiniteElementRelatedData<dim,
572 Assert(
dynamic_cast<const InternalData *
>(&fe_internal) !=
nullptr,
582 quadrature[0].size()),
583 quadrature[0].size(),
587 static_cast<const InternalData &
>(fe_internal),
593template <
int dim,
int spacedim>
597 const unsigned int face_no,
598 const unsigned int sub_no,
599 const Quadrature<dim - 1> &quadrature,
600 const Mapping<dim, spacedim> &
mapping,
601 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
602 const internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
605 ::internal::FEValuesImplementation::FiniteElementRelatedData<dim,
609 Assert(
dynamic_cast<const InternalData *
>(&fe_internal) !=
nullptr,
626 static_cast<const InternalData &
>(fe_internal),
632template <
int dim,
int spacedim>
637 const unsigned int n_q_points,
638 const Mapping<dim, spacedim> &
mapping,
639 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
640 const internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
643 internal::FEValuesImplementation::FiniteElementRelatedData<dim, spacedim>
647 fe_data.
shape_values.size()[0] == this->n_dofs_per_cell(),
649 this->n_dofs_per_cell()));
660 internal::FE_PolyTensor::get_dof_sign_change_nedelec(
668 internal::FE_PolyTensor::get_dof_sign_change_h_div(cell,
674 const unsigned int first_quad_index = this->get_first_quad_index();
676 const unsigned int n_dofs_per_quad = this->n_dofs_per_quad();
677 const unsigned int n_quad_dofs =
680 for (
unsigned int dof_index = 0; dof_index < this->n_dofs_per_cell();
687 const bool is_quad_dof =
689 (first_quad_index <= dof_index) &&
690 (dof_index < first_quad_index + n_quad_dofs));
700 double dof_sign = 1.0;
711 unsigned int face_index_from_dof_index =
712 (dof_index - first_quad_index) / (n_dofs_per_quad);
714 unsigned int local_quad_dof_index = dof_index % n_dofs_per_quad;
717 if (adjust_quad_dof_sign_for_face_orientation(
718 local_quad_dof_index,
719 face_index_from_dof_index,
724 const MappingKind mapping_kind = get_mapping_kind(dof_index);
726 const unsigned int first =
728 [dof_index * this->n_components() +
729 this->get_nonzero_components(dof_index).first_selected_component()];
733 switch (mapping_kind)
737 for (
unsigned int k = 0; k < n_q_points; ++k)
738 for (
unsigned int d = 0;
d < dim; ++
d)
747 const ArrayView<Tensor<1, spacedim>>
748 transformed_shape_values =
758 transformed_shape_values);
760 for (
unsigned int k = 0; k < n_q_points; ++k)
761 for (
unsigned int d = 0;
d < dim; ++
d)
763 transformed_shape_values[k][
d];
770 const ArrayView<Tensor<1, spacedim>>
771 transformed_shape_values =
781 transformed_shape_values);
782 for (
unsigned int k = 0; k < n_q_points; ++k)
783 for (
unsigned int d = 0;
d < dim; ++
d)
785 dof_sign * transformed_shape_values[k][
d];
791 const ArrayView<Tensor<1, spacedim>>
792 transformed_shape_values =
802 transformed_shape_values);
804 for (
unsigned int k = 0; k < n_q_points; ++k)
805 for (
unsigned int d = 0;
d < dim; ++
d)
807 dof_sign * transformed_shape_values[k][
d];
819 const ArrayView<Tensor<2, spacedim>> transformed_shape_grads =
823 switch (mapping_kind)
833 transformed_shape_grads);
834 for (
unsigned int k = 0; k < n_q_points; ++k)
835 for (
unsigned int d = 0;
d < dim; ++
d)
837 transformed_shape_grads[k][d];
849 transformed_shape_grads);
851 for (
unsigned int k = 0; k < n_q_points; ++k)
852 for (
unsigned int d = 0;
d < spacedim; ++
d)
853 for (
unsigned int n = 0; n < spacedim; ++n)
854 transformed_shape_grads[k][d] -=
858 for (
unsigned int k = 0; k < n_q_points; ++k)
859 for (
unsigned int d = 0;
d < dim; ++
d)
861 transformed_shape_grads[k][d];
867 for (
unsigned int k = 0; k < n_q_points; ++k)
876 transformed_shape_grads);
878 for (
unsigned int k = 0; k < n_q_points; ++k)
879 for (
unsigned int d = 0;
d < spacedim; ++
d)
880 for (
unsigned int n = 0; n < spacedim; ++n)
881 transformed_shape_grads[k][d] +=
885 for (
unsigned int k = 0; k < n_q_points; ++k)
886 for (
unsigned int d = 0;
d < dim; ++
d)
888 transformed_shape_grads[k][d];
896 for (
unsigned int k = 0; k < n_q_points; ++k)
905 transformed_shape_grads);
907 for (
unsigned int k = 0; k < n_q_points; ++k)
908 for (
unsigned int d = 0;
d < spacedim; ++
d)
909 for (
unsigned int n = 0; n < spacedim; ++n)
910 transformed_shape_grads[k][d] +=
917 for (
unsigned int k = 0; k < n_q_points; ++k)
918 for (
unsigned int d = 0;
d < dim; ++
d)
920 dof_sign * transformed_shape_grads[k][d];
937 for (
unsigned int k = 0; k < n_q_points; ++k)
947 transformed_shape_grads);
949 for (
unsigned int k = 0; k < n_q_points; ++k)
950 for (
unsigned int d = 0;
d < spacedim; ++
d)
951 for (
unsigned int n = 0; n < spacedim; ++n)
952 transformed_shape_grads[k][d] -=
956 for (
unsigned int k = 0; k < n_q_points; ++k)
957 for (
unsigned int d = 0;
d < dim; ++
d)
959 dof_sign * transformed_shape_grads[k][d];
971 switch (mapping_kind)
975 const ArrayView<Tensor<3, spacedim>>
976 transformed_shape_hessians =
986 transformed_shape_hessians);
988 for (
unsigned int k = 0; k < n_q_points; ++k)
989 for (
unsigned int d = 0;
d < spacedim; ++
d)
990 for (
unsigned int n = 0; n < spacedim; ++n)
991 transformed_shape_hessians[k][d] -=
995 for (
unsigned int k = 0; k < n_q_points; ++k)
996 for (
unsigned int d = 0;
d < dim; ++
d)
998 transformed_shape_hessians[k][d];
1004 for (
unsigned int k = 0; k < n_q_points; ++k)
1008 const ArrayView<Tensor<3, spacedim>>
1009 transformed_shape_hessians =
1019 transformed_shape_hessians);
1021 for (
unsigned int k = 0; k < n_q_points; ++k)
1022 for (
unsigned int d = 0;
d < spacedim; ++
d)
1023 for (
unsigned int n = 0; n < spacedim; ++n)
1024 for (
unsigned int i = 0; i < spacedim; ++i)
1025 for (
unsigned int j = 0; j < spacedim; ++j)
1027 transformed_shape_hessians[k][
d][i][j] -=
1030 .jacobian_pushed_forward_2nd_derivatives
1043 for (
unsigned int k = 0; k < n_q_points; ++k)
1044 for (
unsigned int d = 0;
d < dim; ++
d)
1046 transformed_shape_hessians[k][d];
1053 for (
unsigned int k = 0; k < n_q_points; ++k)
1057 const ArrayView<Tensor<3, spacedim>>
1058 transformed_shape_hessians =
1068 transformed_shape_hessians);
1070 for (
unsigned int k = 0; k < n_q_points; ++k)
1071 for (
unsigned int d = 0;
d < spacedim; ++
d)
1072 for (
unsigned int n = 0; n < spacedim; ++n)
1073 for (
unsigned int i = 0; i < spacedim; ++i)
1074 for (
unsigned int j = 0; j < spacedim; ++j)
1076 transformed_shape_hessians[k][
d][i][j] +=
1079 .jacobian_pushed_forward_2nd_derivatives
1090 for (
unsigned int m = 0; m < spacedim; ++m)
1091 transformed_shape_hessians[k][d][i][j] -=
1093 .jacobian_pushed_forward_grads[k][d][i]
1096 .jacobian_pushed_forward_grads[k][m][n]
1108 for (
unsigned int k = 0; k < n_q_points; ++k)
1109 for (
unsigned int d = 0;
d < dim; ++
d)
1111 transformed_shape_hessians[k][d];
1119 for (
unsigned int k = 0; k < n_q_points; ++k)
1123 const ArrayView<Tensor<3, spacedim>>
1124 transformed_shape_hessians =
1134 transformed_shape_hessians);
1136 for (
unsigned int k = 0; k < n_q_points; ++k)
1137 for (
unsigned int d = 0;
d < spacedim; ++
d)
1138 for (
unsigned int n = 0; n < spacedim; ++n)
1139 for (
unsigned int i = 0; i < spacedim; ++i)
1140 for (
unsigned int j = 0; j < spacedim; ++j)
1142 transformed_shape_hessians[k][
d][i][j] +=
1145 .jacobian_pushed_forward_2nd_derivatives
1157 transformed_shape_hessians[k][
d][i][j] -=
1160 .jacobian_pushed_forward_2nd_derivatives
1169 for (
unsigned int m = 0; m < spacedim; ++m)
1171 transformed_shape_hessians[k][
d][i][j] -=
1187 transformed_shape_hessians[k][
d][i][j] +=
1205 for (
unsigned int k = 0; k < n_q_points; ++k)
1206 for (
unsigned int d = 0;
d < dim; ++
d)
1208 dof_sign * transformed_shape_hessians[k][d];
1215 for (
unsigned int k = 0; k < n_q_points; ++k)
1219 const ArrayView<Tensor<3, spacedim>>
1220 transformed_shape_hessians =
1230 transformed_shape_hessians);
1232 for (
unsigned int k = 0; k < n_q_points; ++k)
1233 for (
unsigned int d = 0;
d < spacedim; ++
d)
1234 for (
unsigned int n = 0; n < spacedim; ++n)
1235 for (
unsigned int i = 0; i < spacedim; ++i)
1236 for (
unsigned int j = 0; j < spacedim; ++j)
1238 transformed_shape_hessians[k][
d][i][j] -=
1241 .jacobian_pushed_forward_2nd_derivatives
1254 for (
unsigned int k = 0; k < n_q_points; ++k)
1255 for (
unsigned int d = 0;
d < dim; ++
d)
1257 dof_sign * transformed_shape_hessians[k][d];
1277template <
int dim,
int spacedim>
1284 for (
unsigned int i = 0; i < this->n_dofs_per_cell(); ++i)
1286 const MappingKind mapping_kind = get_mapping_kind(i);
1288 switch (mapping_kind)
1379#include "fe/fe_poly_tensor.inst"
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
internal::SubfaceCase< dim > subface_case(const unsigned int face_no) const
std::vector< Tensor< 3, spacedim > > transformed_shape_hessians
Table< 2, DerivativeForm< 2, dim, spacedim > > shape_grad_grads
Table< 2, DerivativeForm< 1, dim, spacedim > > shape_grads
std::vector< Tensor< 2, spacedim > > transformed_shape_grads
std::vector< Tensor< 3, dim > > untransformed_shape_hessian_tensors
std::vector< Tensor< 2, dim > > untransformed_shape_grads
Table< 2, Tensor< 1, dim > > shape_values
std::vector< Tensor< 1, spacedim > > transformed_shape_values
std::vector< double > dof_sign_change
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
bool adjust_quad_dof_sign_for_face_orientation(const unsigned int index, const unsigned int face_no, const types::geometric_orientation combined_orientation) const
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const override
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const override
void compute_fill(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename QProjector< dim >::DataSetDescriptor &offset, const unsigned int n_q_points, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const typename FE_PolyTensor< dim, spacedim >::InternalData &fe_internal, internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
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
virtual double shape_value(const unsigned int i, const Point< dim > &p) const override
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 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
MappingKind get_mapping_kind(const unsigned int i) const
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)
bool single_mapping_kind() const
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const override
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int tensor_degree() const
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const types::geometric_orientation combined_orientation=numbers::default_geometric_orientation) const
friend class InternalDataBase
virtual void transform(const ArrayView< const Tensor< 1, dim > > &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 1, spacedim > > &output) const override
unsigned int size() const
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NOT_IMPLEMENTED()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcFENotPrimitive()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_contravariant_transformation
Contravariant transformation.
@ update_jacobian_pushed_forward_grads
@ update_hessians
Second derivatives of shape functions.
@ update_values
Shape function values.
@ update_3rd_derivatives
Third derivatives of shape functions.
@ update_covariant_transformation
Covariant transformation.
@ update_gradients
Shape function gradients.
@ update_default
No update.
@ update_piola
Values needed for Piola transform.
@ mapping_covariant_gradient
@ mapping_contravariant_hessian
@ mapping_covariant_hessian
@ mapping_contravariant_gradient
types::geometric_orientation combined_face_orientation(const unsigned int face) const
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
constexpr types::geometric_orientation default_geometric_orientation
unsigned char geometric_orientation
static constexpr unsigned int faces_per_cell