47template <
int dim,
int spacedim>
68template <
int dim,
int spacedim>
71 const UpdateFlags update_flags,
78 const unsigned int n_q_points = q.
size();
90 covariant.resize(n_q_points);
93 contravariant.resize(n_q_points);
96 volume_elements.resize(n_q_points);
101template <
int dim,
int spacedim>
104 const UpdateFlags update_flags,
106 const unsigned int n_original_q_points)
113 covariant.resize(n_original_q_points);
116 contravariant.resize(n_original_q_points);
119 volume_elements.resize(n_original_q_points);
131 unit_tangentials[i].resize(n_original_q_points);
132 std::fill(unit_tangentials[i].begin(),
133 unit_tangentials[i].end(),
138 .resize(n_original_q_points);
153template <
int dim,
int spacedim>
165template <
int dim,
int spacedim>
171 for (
unsigned int q = 0; q <
quad.size(); ++q)
183template <
int dim,
int spacedim>
190template <
int dim,
int spacedim>
191std::unique_ptr<Mapping<dim, spacedim>>
194 return std::make_unique<MappingManifold<dim, spacedim>>(*this);
199template <
int dim,
int spacedim>
211template <
int dim,
int spacedim>
218 std::array<double, GeometryInfo<dim>::vertices_per_cell> weights;
222 vertices[v] = cell->
vertex(v);
245template <
int dim,
int spacedim>
248 const UpdateFlags in)
const
256 UpdateFlags out = in;
257 for (
unsigned int i = 0; i < 5; ++i)
307template <
int dim,
int spacedim>
308std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
312 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
313 std::make_unique<InternalData>();
321template <
int dim,
int spacedim>
322std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
324 const UpdateFlags update_flags,
329 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
330 std::make_unique<InternalData>();
331 auto &data =
dynamic_cast<InternalData &
>(*data_ptr);
335 quadrature[0].size());
342template <
int dim,
int spacedim>
343std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
345 const UpdateFlags update_flags,
348 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
349 std::make_unique<InternalData>();
373 template <
int dim,
int spacedim>
375 maybe_compute_q_points(
377 const typename ::MappingManifold<dim, spacedim>::InternalData
381 const UpdateFlags update_flags = data.update_each;
385 for (
unsigned int point = 0; point < quadrature_points.size();
388 quadrature_points[point] = data.manifold->get_new_point(
391 data.cell_manifold_quadrature_weights[point + data_set]));
404 template <
int dim,
int spacedim>
406 maybe_update_Jacobians(
407 const typename ::QProjector<dim>::DataSetDescriptor data_set,
408 const typename ::MappingManifold<dim, spacedim>::InternalData
411 const UpdateFlags update_flags = data.update_each;
415 const unsigned int n_q_points = data.contravariant.size();
417 std::fill(data.contravariant.begin(),
418 data.contravariant.end(),
421 for (
unsigned int point = 0; point < n_q_points; ++point)
425 const Point<dim> &p = data.quad.point(point + data_set);
431 data.cell_manifold_quadrature_weights[point + data_set]));
440 for (
unsigned int i = 0; i < dim; ++i)
443 const double pi = p[i];
444 Assert(pi >= 0 && pi <= 1.0,
446 "Was expecting a quadrature point "
447 "inside the unit reference element."));
451 const double L = pi > .5 ? -pi : 1 - pi;
456 for (
const unsigned int j :
458 data.vertex_weights[j] =
466 data.manifold->get_tangent_vector(P, NP);
468 for (
unsigned int d = 0; d < spacedim; ++d)
469 data.contravariant[point][d][i] = T[d] / L;
475 const unsigned int n_q_points = data.contravariant.size();
476 for (
unsigned int point = 0; point < n_q_points; ++point)
478 data.covariant[point] =
479 (data.contravariant[point]).covariant_form();
485 const unsigned int n_q_points = data.contravariant.size();
486 for (
unsigned int point = 0; point < n_q_points; ++point)
487 data.volume_elements[point] =
488 data.contravariant[point].determinant();
498template <
int dim,
int spacedim>
513 const unsigned int n_q_points = quadrature.
size();
518 internal::MappingManifoldImplementation::maybe_compute_q_points<dim,
524 internal::MappingManifoldImplementation::maybe_update_Jacobians<dim,
529 const std::vector<double> &weights = quadrature.
get_weights();
544 for (
unsigned int point = 0; point < n_q_points; ++point)
558 cell->
center(), det, point)));
560 output_data.
JxW_values[point] = weights[point] * det;
568 for (
unsigned int i = 0; i < spacedim; ++i)
569 for (
unsigned int j = 0; j < dim; ++j)
573 for (
unsigned int i = 0; i < dim; ++i)
574 for (
unsigned int j = 0; j < dim; ++j)
575 G[i][j] = DX_t[i] * DX_t[j];
582 Assert(spacedim == dim + 1,
584 "There is no (unique) cell normal for " +
586 "-dimensional cells in " +
588 "-dimensional space. This only works if the "
589 "space dimension is one greater than the "
590 "dimensionality of the mesh cells."));
615 for (
unsigned int point = 0; point < n_q_points; ++point)
623 for (
unsigned int point = 0; point < n_q_points; ++point)
635 namespace MappingManifoldImplementation
648 template <
int dim,
int spacedim>
650 maybe_compute_face_data(
651 const ::MappingManifold<dim, spacedim> &
mapping,
652 const typename ::Triangulation<dim, spacedim>::cell_iterator
654 const unsigned int face_no,
655 const unsigned int subface_no,
656 const unsigned int n_q_points,
657 const std::vector<double> &weights,
658 const typename ::MappingManifold<dim, spacedim>::InternalData
676 for (
unsigned int d = 0; d != dim - 1; ++d)
679 data.unit_tangentials.size(),
682 data.aux[d].size() <=
684 .unit_tangentials[face_no +
692 .unit_tangentials[face_no +
703 for (
unsigned int i = 0; i < n_q_points; ++i)
712 (face_no == 0 ? -1 : +1);
736 for (
unsigned int point = 0;
point < n_q_points; ++
point)
744 data.contravariant[
point].transpose()[0];
746 (face_no == 0 ? -1. : +1.) *
754 const DerivativeForm<1, spacedim, dim> DX_t =
755 data.contravariant[
point].transpose();
757 Tensor<1, spacedim> cell_normal =
759 cell_normal /= cell_normal.
norm();
786 const double area_ratio =
788 cell->subface_case(face_no), subface_no);
800 for (
unsigned int point = 0;
point < n_q_points; ++
point)
801 output_data.
jacobians[point] = data.contravariant[point];
804 for (
unsigned int point = 0;
point < n_q_points; ++
point)
806 data.covariant[point].transpose();
817 template <
int dim,
int spacedim>
820 const ::MappingManifold<dim, spacedim> &
mapping,
821 const typename ::Triangulation<dim, spacedim>::cell_iterator
823 const unsigned int face_no,
824 const unsigned int subface_no,
825 const typename QProjector<dim>::DataSetDescriptor data_set,
826 const Quadrature<dim - 1> &quadrature,
827 const typename ::MappingManifold<dim, spacedim>::InternalData
829 internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
832 data.store_vertices(cell);
834 data.manifold = &cell->face(face_no)->get_manifold();
836 maybe_compute_q_points<dim, spacedim>(data_set,
839 maybe_update_Jacobians<dim, spacedim>(data_set, data);
851 template <
int dim,
int spacedim,
int rank>
854 const ArrayView<
const Tensor<rank, dim>> &input,
856 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
857 const ArrayView<Tensor<rank, spacedim>> &output)
860 Assert((
dynamic_cast<const typename ::
861 MappingManifold<dim, spacedim>::InternalData *
>(
862 &mapping_data) !=
nullptr),
864 const typename ::MappingManifold<dim, spacedim>::InternalData
866 static_cast<const typename ::MappingManifold<dim, spacedim>::
867 InternalData &
>(mapping_data);
869 switch (mapping_kind)
876 "update_contravariant_transformation"));
878 for (
unsigned int i = 0; i < output.size(); ++i)
890 "update_contravariant_transformation"));
894 "update_volume_elements"));
899 for (
unsigned int i = 0; i < output.size(); ++i)
903 output[i] /= data.volume_elements[i];
915 "update_covariant_transformation"));
917 for (
unsigned int i = 0; i < output.size(); ++i)
929 template <
int dim,
int spacedim,
int rank>
932 const ArrayView<
const Tensor<rank, dim>> &input,
934 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
935 const ArrayView<Tensor<rank, spacedim>> &output)
938 Assert((
dynamic_cast<const typename ::
939 MappingManifold<dim, spacedim>::InternalData *
>(
940 &mapping_data) !=
nullptr),
942 const typename ::MappingManifold<dim, spacedim>::InternalData
944 static_cast<const typename ::MappingManifold<dim, spacedim>::
945 InternalData &
>(mapping_data);
947 switch (mapping_kind)
954 "update_covariant_transformation"));
958 "update_contravariant_transformation"));
961 for (
unsigned int i = 0; i < output.size(); ++i)
963 DerivativeForm<1, spacedim, dim>
A =
978 "update_covariant_transformation"));
981 for (
unsigned int i = 0; i < output.size(); ++i)
983 DerivativeForm<1, spacedim, dim>
A =
998 "update_covariant_transformation"));
1002 "update_contravariant_transformation"));
1006 "update_volume_elements"));
1009 for (
unsigned int i = 0; i < output.size(); ++i)
1011 DerivativeForm<1, spacedim, dim>
A =
1013 Tensor<2, spacedim>
T =
1018 output[i] /= data.volume_elements[i];
1031 template <
int dim,
int spacedim>
1034 const ArrayView<
const Tensor<3, dim>> &input,
1036 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1037 const ArrayView<Tensor<3, spacedim>> &output)
1040 Assert((
dynamic_cast<const typename ::
1041 MappingManifold<dim, spacedim>::InternalData *
>(
1042 &mapping_data) !=
nullptr),
1044 const typename ::MappingManifold<dim, spacedim>::InternalData
1046 static_cast<const typename ::MappingManifold<dim, spacedim>::
1047 InternalData &
>(mapping_data);
1049 switch (mapping_kind)
1056 "update_covariant_transformation"));
1060 "update_contravariant_transformation"));
1062 for (
unsigned int q = 0; q < output.
size(); ++q)
1065 data.contravariant[q],
1076 "update_covariant_transformation"));
1078 for (
unsigned int q = 0; q < output.
size(); ++q)
1091 "update_covariant_transformation"));
1095 "update_contravariant_transformation"));
1099 "update_volume_elements"));
1101 for (
unsigned int q = 0; q < output.
size(); ++q)
1104 data.contravariant[q],
1105 data.volume_elements[q],
1118 template <
int dim,
int spacedim,
int rank>
1121 const ArrayView<
const DerivativeForm<rank, dim, spacedim>> &input,
1123 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_data,
1124 const ArrayView<Tensor<rank + 1, spacedim>> &output)
1127 Assert((
dynamic_cast<const typename ::
1128 MappingManifold<dim, spacedim>::InternalData *
>(
1129 &mapping_data) !=
nullptr),
1131 const typename ::MappingManifold<dim, spacedim>::InternalData
1133 static_cast<const typename ::MappingManifold<dim, spacedim>::
1134 InternalData &
>(mapping_data);
1136 switch (mapping_kind)
1143 "update_covariant_transformation"));
1145 for (
unsigned int i = 0; i < output.size(); ++i)
1160template <
int dim,
int spacedim>
1164 const unsigned int face_no,
1177 internal::MappingManifoldImplementation::do_fill_fe_face_values(
1186 quadrature[0].
size()),
1194template <
int dim,
int spacedim>
1198 const unsigned int face_no,
1199 const unsigned int subface_no,
1210 internal::MappingManifoldImplementation::do_fill_fe_face_values(
1229template <
int dim,
int spacedim>
1233 const MappingKind mapping_kind,
1237 internal::MappingManifoldImplementation::transform_fields(input,
1245template <
int dim,
int spacedim>
1249 const MappingKind mapping_kind,
1253 internal::MappingManifoldImplementation::transform_differential_forms(
1254 input, mapping_kind, mapping_data, output);
1259template <
int dim,
int spacedim>
1263 const MappingKind mapping_kind,
1267 switch (mapping_kind)
1270 internal::MappingManifoldImplementation::transform_fields(input,
1279 internal::MappingManifoldImplementation::transform_gradients(
1280 input, mapping_kind, mapping_data, output);
1289template <
int dim,
int spacedim>
1293 const MappingKind mapping_kind,
1302 switch (mapping_kind)
1308 "update_covariant_transformation"));
1310 for (
unsigned int q = 0; q < output.
size(); ++q)
1324template <
int dim,
int spacedim>
1328 const MappingKind mapping_kind,
1332 switch (mapping_kind)
1337 internal::MappingManifoldImplementation::transform_hessians(
1338 input, mapping_kind, mapping_data, output);
1346#include "fe/mapping_manifold.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
double diameter(const Mapping< dim, spacedim > &mapping) const
bool direction_flag() const
virtual Point< spacedim > get_new_point(const ArrayView< const Point< spacedim > > &surrounding_points, const ArrayView< const double > &weights) const
void store_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
std::vector< double > volume_elements
std::array< std::vector< Tensor< 1, dim > >, GeometryInfo< dim >::faces_per_cell *(dim - 1)> unit_tangentials
std::vector< DerivativeForm< 1, dim, spacedim > > contravariant
std::array< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > vertices
std::vector< std::array< double, GeometryInfo< dim >::vertices_per_cell > > cell_manifold_quadrature_weights
std::vector< std::vector< Tensor< 1, spacedim > > > aux
virtual std::size_t memory_consumption() const override
void initialize_face(const UpdateFlags update_flags, const Quadrature< dim > &quadrature, const unsigned int n_original_q_points)
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
ObserverPointer< const Manifold< dim, spacedim > > manifold
std::array< double, GeometryInfo< dim >::vertices_per_cell > vertex_weights
void compute_manifold_quadrature_weights(const Quadrature< dim > &quadrature)
std::vector< DerivativeForm< 1, dim, spacedim > > covariant
Triangulation< dim, spacedim >::cell_iterator cell
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
MappingManifold()=default
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
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
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
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 typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_face_data(const UpdateFlags flags, const hp::QCollection< dim - 1 > &quadrature) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
virtual CellSimilarity::Similarity fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
virtual std::size_t memory_consumption() const
static constexpr Point< dim, Number > unit_vector(const unsigned int i)
Class storing the offset index into a Quadrature rule created by project_to_all_faces() or project_to...
static DataSetDescriptor subface(const ReferenceCell &reference_cell, const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points, const internal::SubfaceCase< dim > ref_case=internal::SubfaceCase< dim >::case_isotropic)
static DataSetDescriptor cell()
static DataSetDescriptor face(const ReferenceCell &reference_cell, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const unsigned int n_quadrature_points)
static Quadrature< dim > project_to_all_subfaces(const ReferenceCell &reference_cell, const SubQuadrature &quadrature)
static Quadrature< dim > project_to_all_faces(const ReferenceCell &reference_cell, const hp::QCollection< dim - 1 > &quadrature)
const std::vector< double > & get_weights() const
const std::vector< Point< dim > > & get_points() const
unsigned int size() const
numbers::NumberTraits< Number >::real_type norm() const
const Manifold< dim, spacedim > & get_manifold() const
Point< spacedim > center(const bool respect_manifold=false, const bool interpolate_from_surrounding=false) const
Point< spacedim > & vertex(const unsigned int i) const
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
static ::ExceptionBase & ExcDistortedMappedCell(Point< spacedim > arg1, double arg2, int arg3)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_volume_elements
Determinant of the Jacobian.
@ update_contravariant_transformation
Contravariant transformation.
@ update_jacobian_pushed_forward_grads
@ update_jacobian_grads
Gradient of volume element.
@ update_normal_vectors
Normal vectors.
@ update_JxW_values
Transformed quadrature weights.
@ update_covariant_transformation
Covariant transformation.
@ update_jacobians
Volume element.
@ update_inverse_jacobians
Volume element.
@ update_quadrature_points
Transformed quadrature points.
@ update_jacobian_pushed_forward_3rd_derivatives
@ update_boundary_forms
Outer normal vector, not normalized.
@ 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
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
Point< spacedim > point(const gp_Pnt &p, const double tolerance=1e-10)
constexpr const ReferenceCell & get_hypercube()
constexpr T fixed_power(const T t)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
void transform_differential_forms(const ArrayView< const DerivativeForm< rank, dim, spacedim > > &input, const MappingKind mapping_kind, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, const ArrayView< Tensor< rank+1, spacedim > > &output)
void do_fill_fe_face_values(const ::MappingQ< dim, spacedim > &mapping, const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const typename QProjector< dim >::DataSetDescriptor data_set, const Quadrature< dim - 1 > &quadrature, const typename ::MappingQ< dim, spacedim >::InternalData &data, const std::vector< Polynomials::Polynomial< double > > &polynomials_1d, const std::vector< unsigned int > &renumber_lexicographic_to_hierarchic, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data)
void transform_fields(const ArrayView< const Tensor< rank, dim > > &input, const MappingKind mapping_kind, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, const ArrayView< Tensor< rank, spacedim > > &output)
void transform_gradients(const ArrayView< const Tensor< rank, dim > > &input, const MappingKind mapping_kind, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, const ArrayView< Tensor< rank, spacedim > > &output)
void transform_hessians(const ArrayView< const Tensor< 3, dim > > &input, const MappingKind mapping_kind, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_data, const ArrayView< Tensor< 3, spacedim > > &output)
void maybe_compute_face_data(const ::MappingQ< dim, spacedim > &mapping, const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int subface_no, const unsigned int n_q_points, const std::vector< double > &weights, const typename ::MappingQ< dim, spacedim >::InternalData &data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data)
Tensor< 3, spacedim, Number > apply_contravariant_hessian(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const DerivativeForm< 1, dim, spacedim, Number > &contravariant, const Tensor< 3, dim, Number > &input)
Tensor< 3, spacedim, Number > apply_piola_hessian(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const DerivativeForm< 1, dim, spacedim, Number > &contravariant, const Number &volume_element, const Tensor< 3, dim, Number > &input)
Tensor< 3, spacedim, Number > apply_covariant_gradient(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const DerivativeForm< 2, dim, spacedim, Number > &input)
Tensor< 3, spacedim, Number > apply_covariant_hessian(const DerivativeForm< 1, dim, spacedim, Number > &covariant, const Tensor< 3, dim, Number > &input)
constexpr unsigned int invalid_unsigned_int
::VectorizedArray< Number, width > sqrt(const ::VectorizedArray< Number, width > &)
static constexpr ndarray< Tensor< 1, dim >, faces_per_cell, dim - 1 > unit_tangential_vectors
static constexpr unsigned int vertices_per_cell
static constexpr unsigned int faces_per_cell
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()
static double subface_ratio(const internal::SubfaceCase< dim > &subface_case, const unsigned int subface_no)
static double d_linear_shape_function(const Point< dim > &xi, const unsigned int i)
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices()
DEAL_II_HOST constexpr Number determinant(const SymmetricTensor< 2, dim, Number > &)
constexpr Tensor< 1, dim, Number > cross_product_2d(const Tensor< 1, dim, Number > &src)
constexpr Tensor< 1, dim, typename ProductType< Number1, Number2 >::type > cross_product_3d(const Tensor< 1, dim, Number1 > &src1, const Tensor< 1, dim, Number2 > &src2)