42 "You are using MappingCartesian, but the incoming cell is not Cartesian.");
51template <
class CellType>
55 if (!cell->reference_cell().is_hyper_cube())
61 const double abs_tol = 1e-30;
62 const double rel_tol = 1e-28;
63 const auto bounding_box = cell->bounding_box();
64 const auto & bounding_vertices = bounding_box.get_boundary_points();
65 const auto bb_diagonal_length_squared =
66 bounding_vertices.first.distance_square(bounding_vertices.second);
68 for (
const unsigned int v : cell->vertex_indices())
79 const double tolerance =
std::max(abs_tol * cell->vertex(v).norm_square(),
80 rel_tol * bb_diagonal_length_squared);
82 if (cell->vertex(v).distance_square(bounding_box.vertex(v)) > tolerance)
91template <
int dim,
int spacedim>
95 , volume_element(
numbers::signaling_nan<double>())
96 , quadrature_points(q.get_points())
101template <
int dim,
int spacedim>
113template <
int dim,
int spacedim>
122template <
int dim,
int spacedim>
127 Assert(dim == reference_cell.get_dimension(),
128 ExcMessage(
"The dimension of your mapping (" +
130 ") and the reference cell cell_type (" +
132 " ) do not agree."));
134 return reference_cell.is_hyper_cube();
139template <
int dim,
int spacedim>
158template <
int dim,
int spacedim>
159std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
163 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
164 std::make_unique<InternalData>(q);
177template <
int dim,
int spacedim>
178std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
185 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
187 ReferenceCells::get_hypercube<dim>(), quadrature[0]));
197 data.update_each = update_flags;
204template <
int dim,
int spacedim>
205std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
210 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
212 ReferenceCells::get_hypercube<dim>(), quadrature));
222 data.update_each = update_flags;
229template <
int dim,
int spacedim>
265template <
int dim,
int spacedim>
270 std::vector<
Point<dim>> &quadrature_points)
const
282template <
int dim,
int spacedim>
286 const unsigned int face_no,
288 std::vector<
Point<dim>> &quadrature_points)
const
295 ReferenceCells::get_hypercube<dim>(),
297 cell->face_orientation(face_no),
298 cell->face_flip(face_no),
299 cell->face_rotation(face_no),
300 quadrature_points.size());
309template <
int dim,
int spacedim>
313 const unsigned int face_no,
314 const unsigned int sub_no,
316 std::vector<
Point<dim>> &quadrature_points)
const
320 if (cell->face(face_no)->has_children())
328 ReferenceCells::get_hypercube<dim>(),
331 cell->face_orientation(face_no),
332 cell->face_flip(face_no),
333 cell->face_rotation(face_no),
334 quadrature_points.size(),
335 cell->subface_case(face_no));
343template <
int dim,
int spacedim>
349 std::vector<
Point<dim>> &quadrature_points)
const
355 for (
unsigned int i = 0; i < quadrature_points.size(); ++i)
357 quadrature_points[i] = start;
358 for (
unsigned int d = 0; d < dim; ++d)
359 quadrature_points[i](d) +=
366template <
int dim,
int spacedim>
369 const unsigned int face_no,
377 std::fill(normal_vectors.begin(),
378 normal_vectors.end(),
385template <
int dim,
int spacedim>
396 for (
unsigned int i = 0; i < output_data.
jacobian_grads.size(); ++i)
400 for (
unsigned int i = 0;
406 for (
unsigned int i = 0;
413 for (
unsigned int i = 0;
420 for (
unsigned int i = 0;
427 for (
unsigned int i = 0;
437template <
int dim,
int spacedim>
445 for (
unsigned int d = 1; d < dim; ++d)
453template <
int dim,
int spacedim>
465 for (
unsigned int i = 0; i < output_data.
jacobians.size(); ++i)
468 for (
unsigned int j = 0; j < dim; ++j)
475template <
int dim,
int spacedim>
490 for (
unsigned int j = 0; j < dim; ++j)
497template <
int dim,
int spacedim>
528 for (
unsigned int d = 1; d < dim; ++d)
532 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
541 return cell_similarity;
546template <
int dim,
int spacedim>
565 output_data.
initialize(unit_points.size(), update_flags);
568 data.update_each = update_flags;
570 std::vector<Point<dim>>(unit_points.begin(), unit_points.end());
584template <
int dim,
int spacedim>
588 const unsigned int face_no,
617 for (
unsigned int d = 0; d < dim; ++d)
622 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
623 output_data.
JxW_values[i] = J * quadrature[0].weight(i);
626 for (
unsigned int i = 0; i < output_data.
boundary_forms.size(); ++i)
637template <
int dim,
int spacedim>
641 const unsigned int face_no,
642 const unsigned int subface_no,
666 for (
unsigned int d = 0; d < dim; ++d)
676 const unsigned int n_subfaces =
677 cell->face(face_no)->has_children() ?
678 cell->face(face_no)->n_children() :
680 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
685 for (
unsigned int i = 0; i < output_data.
boundary_forms.size(); ++i)
696template <
int dim,
int spacedim>
722 for (
unsigned int i = 0; i < output_data.
normal_vectors.size(); ++i)
727 for (
unsigned int d = 0; d < dim; ++d)
731 normal /= normal.
norm();
736 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
742 double det_jacobian = 1.;
743 for (
unsigned int d = 0; d < dim; ++d)
746 invJTxNormal[d] = ref_space_normal[d] / data.
cell_extents[d];
749 det_jacobian * invJTxNormal.
norm() * quadrature.
weight(i);
760template <
int dim,
int spacedim>
773 switch (mapping_kind)
779 "update_covariant_transformation"));
781 for (
unsigned int i = 0; i < output.size(); ++i)
782 for (
unsigned int d = 0; d < dim; ++d)
791 "update_contravariant_transformation"));
793 for (
unsigned int i = 0; i < output.size(); ++i)
794 for (
unsigned int d = 0; d < dim; ++d)
802 "update_contravariant_transformation"));
805 "update_volume_elements"));
807 for (
unsigned int i = 0; i < output.size(); ++i)
808 for (
unsigned int d = 0; d < dim; ++d)
820template <
int dim,
int spacedim>
833 switch (mapping_kind)
839 "update_covariant_transformation"));
841 for (
unsigned int i = 0; i < output.size(); ++i)
842 for (
unsigned int d1 = 0; d1 < dim; ++d1)
843 for (
unsigned int d2 = 0; d2 < dim; ++d2)
844 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2];
852 "update_contravariant_transformation"));
854 for (
unsigned int i = 0; i < output.size(); ++i)
855 for (
unsigned int d1 = 0; d1 < dim; ++d1)
856 for (
unsigned int d2 = 0; d2 < dim; ++d2)
857 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
865 "update_covariant_transformation"));
867 for (
unsigned int i = 0; i < output.size(); ++i)
868 for (
unsigned int d1 = 0; d1 < dim; ++d1)
869 for (
unsigned int d2 = 0; d2 < dim; ++d2)
870 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2] /
879 "update_contravariant_transformation"));
881 for (
unsigned int i = 0; i < output.size(); ++i)
882 for (
unsigned int d1 = 0; d1 < dim; ++d1)
883 for (
unsigned int d2 = 0; d2 < dim; ++d2)
884 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
893 "update_contravariant_transformation"));
896 "update_volume_elements"));
898 for (
unsigned int i = 0; i < output.size(); ++i)
899 for (
unsigned int d1 = 0; d1 < dim; ++d1)
900 for (
unsigned int d2 = 0; d2 < dim; ++d2)
901 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
910 "update_contravariant_transformation"));
913 "update_volume_elements"));
915 for (
unsigned int i = 0; i < output.size(); ++i)
916 for (
unsigned int d1 = 0; d1 < dim; ++d1)
917 for (
unsigned int d2 = 0; d2 < dim; ++d2)
918 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
930template <
int dim,
int spacedim>
943 switch (mapping_kind)
949 "update_covariant_transformation"));
951 for (
unsigned int i = 0; i < output.size(); ++i)
952 for (
unsigned int d1 = 0; d1 < dim; ++d1)
953 for (
unsigned int d2 = 0; d2 < dim; ++d2)
954 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2];
962 "update_contravariant_transformation"));
964 for (
unsigned int i = 0; i < output.size(); ++i)
965 for (
unsigned int d1 = 0; d1 < dim; ++d1)
966 for (
unsigned int d2 = 0; d2 < dim; ++d2)
967 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
975 "update_covariant_transformation"));
977 for (
unsigned int i = 0; i < output.size(); ++i)
978 for (
unsigned int d1 = 0; d1 < dim; ++d1)
979 for (
unsigned int d2 = 0; d2 < dim; ++d2)
980 output[i][d1][d2] = input[i][d1][d2] / data.
cell_extents[d2] /
989 "update_contravariant_transformation"));
991 for (
unsigned int i = 0; i < output.size(); ++i)
992 for (
unsigned int d1 = 0; d1 < dim; ++d1)
993 for (
unsigned int d2 = 0; d2 < dim; ++d2)
994 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
1003 "update_contravariant_transformation"));
1006 "update_volume_elements"));
1008 for (
unsigned int i = 0; i < output.size(); ++i)
1009 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1010 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1011 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
1020 "update_contravariant_transformation"));
1023 "update_volume_elements"));
1025 for (
unsigned int i = 0; i < output.size(); ++i)
1026 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1027 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1028 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
1040template <
int dim,
int spacedim>
1053 switch (mapping_kind)
1059 "update_covariant_transformation"));
1061 for (
unsigned int q = 0; q < output.size(); ++q)
1062 for (
unsigned int i = 0; i < spacedim; ++i)
1063 for (
unsigned int j = 0; j < spacedim; ++j)
1064 for (
unsigned int k = 0; k < spacedim; ++k)
1066 output[q][i][j][k] = input[q][i][j][k] /
1079template <
int dim,
int spacedim>
1092 switch (mapping_kind)
1098 "update_covariant_transformation"));
1101 "update_contravariant_transformation"));
1103 for (
unsigned int q = 0; q < output.size(); ++q)
1104 for (
unsigned int i = 0; i < spacedim; ++i)
1105 for (
unsigned int j = 0; j < spacedim; ++j)
1106 for (
unsigned int k = 0; k < spacedim; ++k)
1108 output[q][i][j][k] =
1119 "update_covariant_transformation"));
1121 for (
unsigned int q = 0; q < output.size(); ++q)
1122 for (
unsigned int i = 0; i < spacedim; ++i)
1123 for (
unsigned int j = 0; j < spacedim; ++j)
1124 for (
unsigned int k = 0; k < spacedim; ++k)
1126 output[q][i][j][k] =
1138 "update_covariant_transformation"));
1141 "update_contravariant_transformation"));
1144 "update_volume_elements"));
1146 for (
unsigned int q = 0; q < output.size(); ++q)
1147 for (
unsigned int i = 0; i < spacedim; ++i)
1148 for (
unsigned int j = 0; j < spacedim; ++j)
1149 for (
unsigned int k = 0; k < spacedim; ++k)
1151 output[q][i][j][k] =
1167template <
int dim,
int spacedim>
1180 length[0] = cell->vertex(1)(0) - start(0);
1183 length[0] = cell->vertex(1)(0) - start(0);
1184 length[1] = cell->vertex(2)(1) - start(1);
1187 length[0] = cell->vertex(1)(0) - start(0);
1188 length[1] = cell->vertex(2)(1) - start(1);
1189 length[2] = cell->vertex(4)(2) - start(2);
1196 for (
unsigned int d = 0; d < dim; ++d)
1197 p_real(d) += length[d] * p(d);
1204template <
int dim,
int spacedim>
1212 if (dim != spacedim)
1221 real(0) /= cell->vertex(1)(0) - start(0);
1224 real(0) /= cell->vertex(1)(0) - start(0);
1225 real(1) /= cell->vertex(2)(1) - start(1);
1228 real(0) /= cell->vertex(1)(0) - start(0);
1229 real(1) /= cell->vertex(2)(1) - start(1);
1230 real(2) /= cell->vertex(4)(2) - start(2);
1240template <
int dim,
int spacedim>
1241std::unique_ptr<Mapping<dim, spacedim>>
1244 return std::make_unique<MappingCartesian<dim, spacedim>>(*this);
1250#include "mapping_cartesian.inst"
std::vector< Point< dim > > quadrature_points
Tensor< 1, dim > cell_extents
virtual std::size_t memory_consumption() 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
void maybe_update_volume_elements(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 UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) const override
void update_cell_extents(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const InternalData &data) const
virtual void fill_fe_immersed_surface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const NonMatching::ImmersedSurfaceQuadrature< dim > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const override
void maybe_update_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
void maybe_update_subface_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const InternalData &data, std::vector< Point< dim > > &quadrature_points) const
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
void maybe_update_face_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const InternalData &data, std::vector< Point< dim > > &quadrature_points) const
virtual bool is_compatible_with(const ReferenceCell &reference_cell) const override
void maybe_update_cell_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data, std::vector< Point< dim > > &quadrature_points) const
virtual bool preserves_vertex_locations() const override
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
virtual std::unique_ptr< Mapping< dim, spacedim > > clone() const override
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
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 transform_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data, const typename QProjector< dim >::DataSetDescriptor &offset, std::vector< Point< dim > > &quadrature_points) const
void maybe_update_jacobian_derivatives(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
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
void maybe_update_normal_vectors(const unsigned int face_no, const InternalData &data, std::vector< Tensor< 1, dim > > &normal_vectors) const
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) const override
void maybe_update_inverse_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
Abstract base class for mapping classes.
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
static DataSetDescriptor cell()
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 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)
double weight(const unsigned int i) const
numbers::NumberTraits< Number >::real_type norm() const
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcCellNotCartesian()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
@ 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_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
bool is_cartesian(const CellType &cell)
std::enable_if_t< std::is_fundamental< T >::value, std::size_t > memory_consumption(const T &t)
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)