36 template <
int dim,
int spacedim>
38 compute_linear_transformation(
42 for (
unsigned int j = 0; j < spacedim; ++j)
43 for (
unsigned int i = 1; i < dim + 1; ++i)
44 result[j][i - 1] = cell->
vertex(i)[j] - cell->
vertex(0)[j];
51template <
int dim,
int spacedim>
60template <
int dim,
int spacedim>
68template <
int dim,
int spacedim>
71 const UpdateFlags update_flags,
80template <
int dim,
int spacedim>
94template <
int dim,
int spacedim>
103template <
int dim,
int spacedim>
108 return reference_cell.is_simplex();
113template <
int dim,
int spacedim>
119 UpdateFlags out = in;
128template <
int dim,
int spacedim>
129std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
133 auto data_ptr = std::make_unique<InternalData>(quadrature);
150template <
int dim,
int spacedim>
151std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
153 const UpdateFlags update_flags,
156 auto data_ptr = std::make_unique<InternalData>(
167 data_ptr->update_each = update_flags;
174template <
int dim,
int spacedim>
175std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
177 const UpdateFlags update_flags,
180 auto data_ptr = std::make_unique<InternalData>(
191 data_ptr->update_each = update_flags;
198template <
int dim,
int spacedim>
212template <
int dim,
int spacedim>
221 for (
unsigned int i = 0; i < quadrature_points.size(); ++i)
222 quadrature_points[i] =
230template <
int dim,
int spacedim>
233 const unsigned int face_no,
242 normal_vector /= normal_vector.
norm();
244 std::fill(normal_vectors.begin(), normal_vectors.end(), normal_vector);
249template <
int dim,
int spacedim>
294template <
int dim,
int spacedim>
313template <
int dim,
int spacedim>
333template <
int dim,
int spacedim>
359 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
369 Assert(spacedim == dim + 1,
370 ExcMessage(
"There is no (unique) cell normal for " +
372 "-dimensional cells in " +
374 "-dimensional space. This only works if the "
375 "space dimension is one greater than the "
376 "dimensionality of the mesh cells."));
380 if constexpr (dim == 1 && spacedim == 2)
382 else if constexpr (dim == 2 && spacedim == 3)
389 normal /= normal.
norm();
399 return cell_similarity;
404template <
int dim,
int spacedim>
409 const UpdateFlags update_flags,
421 output_data.
initialize(unit_points.size(), update_flags);
439template <
int dim,
int spacedim>
443 const unsigned int face_no,
459 quadrature_collection);
473 const double J = cell->
face(face_no)->measure() /
476 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
480 for (
unsigned int i = 0; i < output_data.
boundary_forms.size(); ++i)
491template <
int dim,
int spacedim>
495 const unsigned int face_no,
496 const unsigned int subface_no,
527 cell->
face(face_no)->measure() /
528 cell->
face(face_no)->reference_cell().volume() /
533 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
537 for (
unsigned int i = 0; i < output_data.
boundary_forms.size(); ++i)
548template <
int dim,
int spacedim>
552 const MappingKind mapping_kind,
561 switch (mapping_kind)
565 for (
unsigned int i = 0; i < output.size(); ++i)
572 for (
unsigned int i = 0; i < output.size(); ++i)
580 for (
unsigned int d = 0; d < spacedim; ++d)
582 for (
unsigned int i = 0; i < output.size(); ++i)
593template <
int dim,
int spacedim>
597 const MappingKind mapping_kind,
606 switch (mapping_kind)
610 for (
unsigned int i = 0; i < output.size(); ++i)
622template <
int dim,
int spacedim>
626 const MappingKind mapping_kind,
635 switch (mapping_kind)
641 "update_covariant_transformation"));
643 for (
unsigned int i = 0; i < output.size(); ++i)
652 "update_contravariant_transformation"));
654 for (
unsigned int i = 0; i < output.size(); ++i)
663 "update_covariant_transformation"));
665 for (
unsigned int i = 0; i < output.size(); ++i)
678 "update_covariant_transformation"));
680 for (
unsigned int i = 0; i < output.size(); ++i)
695 "update_contravariant_transformation"));
698 for (
unsigned int d = 0; d < spacedim; ++d)
700 for (
unsigned int i = 0; i < output.size(); ++i)
720template <
int dim,
int spacedim>
724 const MappingKind mapping_kind,
733 switch (mapping_kind)
739 "update_covariant_transformation"));
741 for (
unsigned int i = 0; i < output.size(); ++i)
754template <
int dim,
int spacedim>
758 const MappingKind mapping_kind,
767 switch (mapping_kind)
773 "update_covariant_transformation"));
776 "update_contravariant_transformation"));
778 for (
unsigned int i = 0; i < output.size(); ++i)
791 "update_covariant_transformation"));
793 for (
unsigned int i = 0; i < output.size(); ++i)
804 "update_covariant_transformation"));
807 "update_contravariant_transformation"));
809 for (
unsigned int i = 0; i < output.size(); ++i)
825template <
int dim,
int spacedim>
832 compute_linear_transformation<dim, spacedim>(cell);
834 return cell->
vertex(0) + sheared;
839template <
int dim,
int spacedim>
846 compute_linear_transformation<dim, spacedim>(cell)
855template <
int dim,
int spacedim>
863 compute_linear_transformation<dim, spacedim>(cell)
867 for (
unsigned int i = 0; i < real_points.size(); ++i)
874template <
int dim,
int spacedim>
875std::unique_ptr<Mapping<dim, spacedim>>
878 return std::make_unique<MappingP1<dim, spacedim>>(*this);
884#include "fe/mapping_p1.inst"
internal::SubfaceCase< dim > subface_case(const unsigned int face_no) const
TriaIterator< TriaAccessor< dim - 1, dim, spacedim > > face(const unsigned int i) const
bool direction_flag() const
DerivativeForm< 1, dim, spacedim > covariant
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
DerivativeForm< 1, dim, spacedim > linear_component
Tensor< 1, spacedim > affine_component
Quadrature< dim > quadrature
InternalData(const ArrayView< const Point< dim > > &quadrature_points)
virtual std::size_t memory_consumption() const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) 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_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 bool is_compatible_with(const ReferenceCell &reference_cell) const override
void transform_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data, const typename QProjector< dim >::DataSetDescriptor &offset, std::vector< Point< spacedim > > &quadrature_points) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
virtual bool preserves_vertex_locations() 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 void transform_points_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< spacedim > > &real_points, const ArrayView< Point< dim > > &unit_points) const override
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
void update_transformation(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data) const
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< Mapping< dim, spacedim > > clone() const override
void maybe_update_inverse_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
void fill_mapping_data_for_generic_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< dim > > &unit_points, const UpdateFlags update_flags, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
void maybe_update_jacobian_derivatives(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
void maybe_update_normal_vectors(const unsigned int face_no, const InternalData &data, std::vector< Tensor< 1, spacedim > > &normal_vectors) const
void maybe_update_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
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 std::size_t memory_consumption() const
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)
double weight(const unsigned int i) const
unsigned int size() const
double face_measure(const unsigned int face_no) const
numbers::NumberTraits< Number >::real_type norm() const
Point< spacedim > & vertex(const unsigned int i) const
ReferenceCell reference_cell() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDivideByZero()
static ::ExceptionBase & ExcMessage(std::string arg1)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
@ update_jacobian_pushed_forward_2nd_derivatives
@ update_contravariant_transformation
Contravariant transformation.
@ update_jacobian_pushed_forward_grads
@ update_jacobian_3rd_derivatives
@ 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_default
No update.
@ update_jacobian_pushed_forward_3rd_derivatives
@ update_boundary_forms
Outer normal vector, not normalized.
@ update_jacobian_2nd_derivatives
@ mapping_covariant_gradient
@ mapping_contravariant_hessian
@ mapping_covariant_hessian
@ mapping_contravariant_gradient
types::geometric_orientation combined_face_orientation(const unsigned int face) const
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
constexpr const ReferenceCell & get_simplex()
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
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 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)