41 "You are using MappingCartesian, but the incoming cell is not Cartesian.");
50template <
typename CellType>
54 if (!cell->reference_cell().is_hyper_cube())
60 const double abs_tol = 1e-30;
61 const double rel_tol = 1e-28;
62 const auto bounding_box = cell->bounding_box();
63 const auto &bounding_vertices = bounding_box.get_boundary_points();
64 const auto bb_diagonal_length_squared =
65 bounding_vertices.first.distance_square(bounding_vertices.second);
67 for (
const unsigned int v : cell->vertex_indices())
78 const double tolerance =
std::max(abs_tol * cell->vertex(v).norm_square(),
79 rel_tol * bb_diagonal_length_squared);
81 if (cell->vertex(v).distance_square(bounding_box.vertex(v)) > tolerance)
90template <
int dim,
int spacedim>
94 , inverse_cell_extents(
numbers::signaling_nan<
Tensor<1, dim>>())
95 , volume_element(
numbers::signaling_nan<double>())
96 , quadrature_points(q.get_points())
101template <
int dim,
int spacedim>
110 this->update_each = update_flags;
115template <
int dim,
int spacedim>
127template <
int dim,
int spacedim>
136template <
int dim,
int spacedim>
141 Assert(dim == reference_cell.get_dimension(),
142 ExcMessage(
"The dimension of your mapping (" +
144 ") and the reference cell cell_type (" +
146 " ) do not agree."));
148 return reference_cell.is_hyper_cube();
153template <
int dim,
int spacedim>
172template <
int dim,
int spacedim>
173std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
177 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
178 std::make_unique<InternalData>();
186template <
int dim,
int spacedim>
187std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
194 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
206 data.update_each = update_flags;
213template <
int dim,
int spacedim>
214std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
219 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
231 data.update_each = update_flags;
238template <
int dim,
int spacedim>
250 for (
unsigned int d = 0; d < dim; ++d)
252 const double cell_extent_d = cell->vertex(1 << d)[d] - start[d];
254 Assert(cell_extent_d != 0.,
255 ExcMessage(
"Cell does not appear to be Cartesian!"));
267 transform_quadrature_points(
274 for (
unsigned int i = 0; i < quadrature_points.size(); ++i)
276 quadrature_points[i] = first_vertex;
277 for (
unsigned int d = 0; d < dim; ++d)
278 quadrature_points[i][d] +=
279 cell_extents[d] * unit_quadrature_points[i + offset][d];
286template <
int dim,
int spacedim>
292 std::vector<
Point<dim>> &quadrature_points)
const
298 transform_quadrature_points(cell->vertex(0),
300 unit_quadrature_points,
308template <
int dim,
int spacedim>
312 const unsigned int face_no,
314 std::vector<
Point<dim>> &quadrature_points)
const
323 cell->combined_face_orientation(face_no),
324 quadrature_points.size());
327 transform_quadrature_points(cell->vertex(0),
337template <
int dim,
int spacedim>
341 const unsigned int face_no,
342 const unsigned int sub_no,
344 std::vector<
Point<dim>> &quadrature_points)
const
348 if (cell->face(face_no)->has_children())
359 cell->combined_face_orientation(face_no),
360 quadrature_points.size(),
361 cell->subface_case(face_no));
363 transform_quadrature_points(cell->vertex(0),
373template <
int dim,
int spacedim>
376 const unsigned int face_no,
384 std::fill(normal_vectors.begin(),
385 normal_vectors.end(),
392template <
int dim,
int spacedim>
403 for (
unsigned int i = 0; i < output_data.
jacobian_grads.size(); ++i)
407 for (
unsigned int i = 0;
413 for (
unsigned int i = 0;
420 for (
unsigned int i = 0;
427 for (
unsigned int i = 0;
434 for (
unsigned int i = 0;
444template <
int dim,
int spacedim>
452 for (
unsigned int d = 1; d < dim; ++d)
460template <
int dim,
int spacedim>
472 for (
unsigned int i = 0; i < output_data.
jacobians.size(); ++i)
475 for (
unsigned int j = 0; j < dim; ++j)
482template <
int dim,
int spacedim>
497 for (
unsigned int j = 0; j < dim; ++j)
505template <
int dim,
int spacedim>
537 for (
unsigned int d = 1; d < dim; ++d)
541 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
550 return cell_similarity;
555template <
int dim,
int spacedim>
574 output_data.
initialize(unit_points.size(), update_flags);
592template <
int dim,
int spacedim>
596 const unsigned int face_no,
625 for (
unsigned int d = 0; d < dim; ++d)
630 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
631 output_data.
JxW_values[i] = J * quadrature[0].weight(i);
634 for (
unsigned int i = 0; i < output_data.
boundary_forms.size(); ++i)
645template <
int dim,
int spacedim>
649 const unsigned int face_no,
650 const unsigned int subface_no,
674 for (
unsigned int d = 0; d < dim; ++d)
684 const unsigned int n_subfaces =
685 cell->face(face_no)->has_children() ?
686 cell->face(face_no)->n_children() :
688 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
693 for (
unsigned int i = 0; i < output_data.
boundary_forms.size(); ++i)
704template <
int dim,
int spacedim>
727 quadrature.get_points(),
731 for (
unsigned int i = 0; i < output_data.
normal_vectors.size(); ++i)
736 for (
unsigned int d = 0; d < dim; ++d)
740 normal /= normal.
norm();
745 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
751 double det_jacobian = 1.;
752 for (
unsigned int d = 0; d < dim; ++d)
759 det_jacobian * invJTxNormal.
norm() * quadrature.weight(i);
770template <
int dim,
int spacedim>
783 switch (mapping_kind)
789 "update_covariant_transformation"));
791 for (
unsigned int i = 0; i < output.size(); ++i)
792 for (
unsigned int d = 0; d < dim; ++d)
801 "update_contravariant_transformation"));
803 for (
unsigned int i = 0; i < output.size(); ++i)
804 for (
unsigned int d = 0; d < dim; ++d)
812 "update_contravariant_transformation"));
815 "update_volume_elements"));
817 for (
unsigned int i = 0; i < output.size(); ++i)
818 for (
unsigned int d = 0; d < dim; ++d)
830template <
int dim,
int spacedim>
843 switch (mapping_kind)
849 "update_covariant_transformation"));
851 for (
unsigned int i = 0; i < output.size(); ++i)
852 for (
unsigned int d1 = 0; d1 < dim; ++d1)
853 for (
unsigned int d2 = 0; d2 < dim; ++d2)
863 "update_contravariant_transformation"));
865 for (
unsigned int i = 0; i < output.size(); ++i)
866 for (
unsigned int d1 = 0; d1 < dim; ++d1)
867 for (
unsigned int d2 = 0; d2 < dim; ++d2)
868 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
876 "update_covariant_transformation"));
878 for (
unsigned int i = 0; i < output.size(); ++i)
879 for (
unsigned int d1 = 0; d1 < dim; ++d1)
880 for (
unsigned int d2 = 0; d2 < dim; ++d2)
881 output[i][d1][d2] = input[i][d1][d2] *
891 "update_contravariant_transformation"));
893 for (
unsigned int i = 0; i < output.size(); ++i)
894 for (
unsigned int d1 = 0; d1 < dim; ++d1)
895 for (
unsigned int d2 = 0; d2 < dim; ++d2)
896 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] *
905 "update_contravariant_transformation"));
908 "update_volume_elements"));
910 for (
unsigned int i = 0; i < output.size(); ++i)
911 for (
unsigned int d1 = 0; d1 < dim; ++d1)
912 for (
unsigned int d2 = 0; d2 < dim; ++d2)
913 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
922 "update_contravariant_transformation"));
925 "update_volume_elements"));
927 for (
unsigned int i = 0; i < output.size(); ++i)
928 for (
unsigned int d1 = 0; d1 < dim; ++d1)
929 for (
unsigned int d2 = 0; d2 < dim; ++d2)
930 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] *
943template <
int dim,
int spacedim>
956 switch (mapping_kind)
962 "update_covariant_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)
976 "update_contravariant_transformation"));
978 for (
unsigned int i = 0; i < output.size(); ++i)
979 for (
unsigned int d1 = 0; d1 < dim; ++d1)
980 for (
unsigned int d2 = 0; d2 < dim; ++d2)
981 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2];
989 "update_covariant_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] *
1004 "update_contravariant_transformation"));
1006 for (
unsigned int i = 0; i < output.size(); ++i)
1007 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1008 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1009 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] *
1018 "update_contravariant_transformation"));
1021 "update_volume_elements"));
1023 for (
unsigned int i = 0; i < output.size(); ++i)
1024 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1025 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1026 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] /
1035 "update_contravariant_transformation"));
1038 "update_volume_elements"));
1040 for (
unsigned int i = 0; i < output.size(); ++i)
1041 for (
unsigned int d1 = 0; d1 < dim; ++d1)
1042 for (
unsigned int d2 = 0; d2 < dim; ++d2)
1043 output[i][d1][d2] = input[i][d1][d2] * data.
cell_extents[d2] *
1056template <
int dim,
int spacedim>
1069 switch (mapping_kind)
1075 "update_covariant_transformation"));
1077 for (
unsigned int q = 0; q < output.size(); ++q)
1078 for (
unsigned int i = 0; i < spacedim; ++i)
1079 for (
unsigned int j = 0; j < spacedim; ++j)
1080 for (
unsigned int k = 0; k < spacedim; ++k)
1082 output[q][i][j][k] = input[q][i][j][k] *
1095template <
int dim,
int spacedim>
1108 switch (mapping_kind)
1114 "update_covariant_transformation"));
1117 "update_contravariant_transformation"));
1119 for (
unsigned int q = 0; q < output.size(); ++q)
1120 for (
unsigned int i = 0; i < spacedim; ++i)
1121 for (
unsigned int j = 0; j < spacedim; ++j)
1122 for (
unsigned int k = 0; k < spacedim; ++k)
1124 output[q][i][j][k] = input[q][i][j][k] *
1136 "update_covariant_transformation"));
1138 for (
unsigned int q = 0; q < output.size(); ++q)
1139 for (
unsigned int i = 0; i < spacedim; ++i)
1140 for (
unsigned int j = 0; j < spacedim; ++j)
1141 for (
unsigned int k = 0; k < spacedim; ++k)
1143 output[q][i][j][k] = input[q][i][j][k] *
1156 "update_covariant_transformation"));
1159 "update_contravariant_transformation"));
1162 "update_volume_elements"));
1164 for (
unsigned int q = 0; q < output.size(); ++q)
1165 for (
unsigned int i = 0; i < spacedim; ++i)
1166 for (
unsigned int j = 0; j < spacedim; ++j)
1167 for (
unsigned int k = 0; k < spacedim; ++k)
1169 output[q][i][j][k] =
1186template <
int dim,
int spacedim>
1198 for (
unsigned int d = 0; d < dim; ++d)
1199 unit[d] += (cell->vertex(1 << d)[d] - unit[d]) * p[d];
1206template <
int dim,
int spacedim>
1219 for (
unsigned int d = 0; d < dim; ++d)
1220 real[d] = (real[d] - start[d]) / (cell->vertex(1 << d)[d] - start[d]);
1227template <
int dim,
int spacedim>
1237 if (dim != spacedim)
1243 std::array<double, dim> inverse_lengths;
1244 for (
unsigned int d = 0; d < dim; ++d)
1245 inverse_lengths[d] = 1. / (cell->vertex(1 << d)[d] - start[d]);
1247 for (
unsigned int i = 0; i < real_points.size(); ++i)
1248 for (
unsigned int d = 0; d < dim; ++d)
1249 unit_points[i][d] = (real_points[i][d] - start[d]) * inverse_lengths[d];
1254template <
int dim,
int spacedim>
1255std::unique_ptr<Mapping<dim, spacedim>>
1258 return std::make_unique<MappingCartesian<dim, spacedim>>(*this);
1264#include "mapping_cartesian.inst"
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
std::vector< Point< dim > > quadrature_points
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
Tensor< 1, dim > cell_extents
Tensor< 1, dim > inverse_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_cell_quadrature_points(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const InternalData &data, const ArrayView< const Point< dim > > &unit_quadrature_points, std::vector< Point< dim > > &quadrature_points) const
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
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
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
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 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 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)
double weight(const unsigned int i) const
const std::vector< Point< dim > > & get_points() 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
#define DEAL_II_NOT_IMPLEMENTED()
bool is_cartesian(const CellType &cell)
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
constexpr const ReferenceCell & get_hypercube()
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 > &)