40 ExcCellNotAssociatedWithBox,
41 "You are using MappingBox, but the incoming element is not associated with a"
42 "Bounding Box Cartesian.");
50template <
typename CellType>
53 const std::map<types::global_cell_index, types::global_cell_index>
56 Assert((cell->reference_cell().is_hyper_cube() ||
57 cell->reference_cell().is_simplex()),
59 Assert((translator.find(cell->active_cell_index()) != translator.cend()),
60 ExcCellNotAssociatedWithBox());
67template <
int dim,
int spacedim>
70 const std::map<types::global_cell_index, types::global_cell_index>
73 Assert(input_boxes.size() > 0,
74 ExcMessage(
"Invalid number of bounding boxes."));
77 boxes.resize(input_boxes.size());
78 for (
unsigned int i = 0; i < input_boxes.size(); ++i)
79 boxes[i] = input_boxes[i];
85template <
int dim,
int spacedim>
96template <
int dim,
int spacedim>
109template <
int dim,
int spacedim>
122template <
int dim,
int spacedim>
131template <
int dim,
int spacedim>
133#if DEAL_II_VERSION_GTE(9, 8, 0)
142 ExcMessage(
"The dimension of your mapping (" +
144 ") and the reference cell cell_type (" +
146 " ) do not agree."));
153template <
int dim,
int spacedim>
171template <
int dim,
int spacedim>
172std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
176 std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase> data_ptr =
177 std::make_unique<InternalData>();
185template <
int dim,
int spacedim>
186std::unique_ptr<typename Mapping<dim, spacedim>::InternalDataBase>
199template <
int dim,
int spacedim>
212 const std::pair<Point<dim>,
Point<dim>> &bdary_points =
215 for (
unsigned int d = 0;
d < dim; ++
d)
217 const double cell_extent_d = current_box.
side_length(
d);
218 data.cell_extents[
d] = cell_extent_d;
221 .5 * (bdary_points.first[
d] +
222 bdary_points.second[
d]);
224 Assert(cell_extent_d != 0.,
225 ExcMessage(
"Cell does not appear to be Cartesian!"));
226 data.inverse_cell_extents[
d] = 1. / cell_extent_d;
237 transform_quadrature_points(
242 for (
unsigned int i = 0; i < quadrature_points.size(); ++i)
243 quadrature_points[i] = box.
unit_to_real(unit_quadrature_points[i]);
249template <
int dim,
int spacedim>
255 std::vector<
Point<dim>> &quadrature_points)
const
258 transform_quadrature_points(
260 unit_quadrature_points,
266template <
int dim,
int spacedim>
269 const unsigned int face_no,
277 std::fill(normal_vectors.begin(),
278 normal_vectors.end(),
285template <
int dim,
int spacedim>
296 for (
unsigned int i = 0; i < output_data.
jacobian_grads.size(); ++i)
300 for (
unsigned int i = 0;
306 for (
unsigned int i = 0;
313 for (
unsigned int i = 0;
320 for (
unsigned int i = 0;
327 for (
unsigned int i = 0;
337template <
int dim,
int spacedim>
345 for (
unsigned int d = 1;
d < dim; ++
d)
353template <
int dim,
int spacedim>
365 for (
unsigned int i = 0; i < output_data.
jacobians.size(); ++i)
368 for (
unsigned int j = 0; j < dim; ++j)
375template <
int dim,
int spacedim>
390 for (
unsigned int j = 0; j < dim; ++j)
392 data.inverse_cell_extents[j];
398template <
int dim,
int spacedim>
429 double J =
data.cell_extents[0];
430 for (
unsigned int d = 1;
d < dim; ++
d)
431 J *=
data.cell_extents[
d];
432 data.volume_element = J;
434 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
443 return cell_similarity;
448template <
int dim,
int spacedim>
452 const unsigned int face_no,
453 const unsigned int subface_no,
470template <
int dim,
int spacedim>
497 for (
unsigned int i = 0; i < output_data.
normal_vectors.size(); ++i)
501 for (
unsigned int i = 0; i < output_data.
JxW_values.size(); ++i)
512template <
int dim,
int spacedim>
525 switch (mapping_kind)
531 "update_covariant_transformation"));
533 for (
unsigned int i = 0; i < output.size(); ++i)
534 for (
unsigned int d = 0;
d < dim; ++
d)
535 output[i][
d] = input[i][
d] *
data.inverse_cell_extents[
d];
543 "update_contravariant_transformation"));
545 for (
unsigned int i = 0; i < output.size(); ++i)
546 for (
unsigned int d = 0;
d < dim; ++
d)
547 output[i][
d] = input[i][
d] *
data.cell_extents[
d];
554 "update_contravariant_transformation"));
557 "update_volume_elements"));
559 for (
unsigned int i = 0; i < output.size(); ++i)
560 for (
unsigned int d = 0;
d < dim; ++
d)
562 input[i][
d] *
data.cell_extents[
d] /
data.volume_element;
572template <
int dim,
int spacedim>
585 switch (mapping_kind)
591 "update_covariant_transformation"));
593 for (
unsigned int i = 0; i < output.size(); ++i)
594 for (
unsigned int d1 = 0; d1 < dim; ++d1)
595 for (
unsigned int d2 = 0; d2 < dim; ++d2)
597 input[i][d1][d2] *
data.inverse_cell_extents[d2];
605 "update_contravariant_transformation"));
607 for (
unsigned int i = 0; i < output.size(); ++i)
608 for (
unsigned int d1 = 0; d1 < dim; ++d1)
609 for (
unsigned int d2 = 0; d2 < dim; ++d2)
610 output[i][d1][d2] = input[i][d1][d2] *
data.cell_extents[d2];
618 "update_covariant_transformation"));
620 for (
unsigned int i = 0; i < output.size(); ++i)
621 for (
unsigned int d1 = 0; d1 < dim; ++d1)
622 for (
unsigned int d2 = 0; d2 < dim; ++d2)
623 output[i][d1][d2] = input[i][d1][d2] *
624 data.inverse_cell_extents[d2] *
625 data.inverse_cell_extents[d1];
633 "update_contravariant_transformation"));
635 for (
unsigned int i = 0; i < output.size(); ++i)
636 for (
unsigned int d1 = 0; d1 < dim; ++d1)
637 for (
unsigned int d2 = 0; d2 < dim; ++d2)
638 output[i][d1][d2] = input[i][d1][d2] *
data.cell_extents[d2] *
639 data.inverse_cell_extents[d1];
647 "update_contravariant_transformation"));
650 "update_volume_elements"));
652 for (
unsigned int i = 0; i < output.size(); ++i)
653 for (
unsigned int d1 = 0; d1 < dim; ++d1)
654 for (
unsigned int d2 = 0; d2 < dim; ++d2)
655 output[i][d1][d2] = input[i][d1][d2] *
data.cell_extents[d2] /
664 "update_contravariant_transformation"));
667 "update_volume_elements"));
669 for (
unsigned int i = 0; i < output.size(); ++i)
670 for (
unsigned int d1 = 0; d1 < dim; ++d1)
671 for (
unsigned int d2 = 0; d2 < dim; ++d2)
672 output[i][d1][d2] = input[i][d1][d2] *
data.cell_extents[d2] *
673 data.inverse_cell_extents[d1] /
685template <
int dim,
int spacedim>
698 switch (mapping_kind)
704 "update_covariant_transformation"));
706 for (
unsigned int i = 0; i < output.size(); ++i)
707 for (
unsigned int d1 = 0; d1 < dim; ++d1)
708 for (
unsigned int d2 = 0; d2 < dim; ++d2)
710 input[i][d1][d2] *
data.inverse_cell_extents[d2];
718 "update_contravariant_transformation"));
720 for (
unsigned int i = 0; i < output.size(); ++i)
721 for (
unsigned int d1 = 0; d1 < dim; ++d1)
722 for (
unsigned int d2 = 0; d2 < dim; ++d2)
723 output[i][d1][d2] = input[i][d1][d2] *
data.cell_extents[d2];
731 "update_covariant_transformation"));
733 for (
unsigned int i = 0; i < output.size(); ++i)
734 for (
unsigned int d1 = 0; d1 < dim; ++d1)
735 for (
unsigned int d2 = 0; d2 < dim; ++d2)
736 output[i][d1][d2] = input[i][d1][d2] *
737 data.inverse_cell_extents[d2] *
738 data.inverse_cell_extents[d1];
746 "update_contravariant_transformation"));
748 for (
unsigned int i = 0; i < output.size(); ++i)
749 for (
unsigned int d1 = 0; d1 < dim; ++d1)
750 for (
unsigned int d2 = 0; d2 < dim; ++d2)
751 output[i][d1][d2] = input[i][d1][d2] *
data.cell_extents[d2] *
752 data.inverse_cell_extents[d1];
760 "update_contravariant_transformation"));
763 "update_volume_elements"));
765 for (
unsigned int i = 0; i < output.size(); ++i)
766 for (
unsigned int d1 = 0; d1 < dim; ++d1)
767 for (
unsigned int d2 = 0; d2 < dim; ++d2)
768 output[i][d1][d2] = input[i][d1][d2] *
data.cell_extents[d2] /
777 "update_contravariant_transformation"));
780 "update_volume_elements"));
782 for (
unsigned int i = 0; i < output.size(); ++i)
783 for (
unsigned int d1 = 0; d1 < dim; ++d1)
784 for (
unsigned int d2 = 0; d2 < dim; ++d2)
785 output[i][d1][d2] = input[i][d1][d2] *
data.cell_extents[d2] *
786 data.inverse_cell_extents[d1] /
798template <
int dim,
int spacedim>
811 switch (mapping_kind)
817 "update_covariant_transformation"));
819 for (
unsigned int q = 0; q < output.size(); ++q)
820 for (
unsigned int i = 0; i < spacedim; ++i)
821 for (
unsigned int j = 0; j < spacedim; ++j)
822 for (
unsigned int k = 0; k < spacedim; ++k)
824 output[q][i][j][k] = input[q][i][j][k] *
825 data.inverse_cell_extents[j] *
826 data.inverse_cell_extents[k];
837template <
int dim,
int spacedim>
850 switch (mapping_kind)
856 "update_covariant_transformation"));
859 "update_contravariant_transformation"));
861 for (
unsigned int q = 0; q < output.size(); ++q)
862 for (
unsigned int i = 0; i < spacedim; ++i)
863 for (
unsigned int j = 0; j < spacedim; ++j)
864 for (
unsigned int k = 0; k < spacedim; ++k)
866 output[q][i][j][k] = input[q][i][j][k] *
867 data.cell_extents[i] *
868 data.inverse_cell_extents[j] *
869 data.inverse_cell_extents[k];
878 "update_covariant_transformation"));
880 for (
unsigned int q = 0; q < output.size(); ++q)
881 for (
unsigned int i = 0; i < spacedim; ++i)
882 for (
unsigned int j = 0; j < spacedim; ++j)
883 for (
unsigned int k = 0; k < spacedim; ++k)
885 output[q][i][j][k] = input[q][i][j][k] *
886 (
data.inverse_cell_extents[i] *
887 data.inverse_cell_extents[j]) *
888 data.inverse_cell_extents[k];
898 "update_covariant_transformation"));
901 "update_contravariant_transformation"));
904 "update_volume_elements"));
906 for (
unsigned int q = 0; q < output.size(); ++q)
907 for (
unsigned int i = 0; i < spacedim; ++i)
908 for (
unsigned int j = 0; j < spacedim; ++j)
909 for (
unsigned int k = 0; k < spacedim; ++k)
913 (
data.cell_extents[i] /
data.volume_element *
914 data.inverse_cell_extents[j]) *
915 data.inverse_cell_extents[k];
928template <
int dim,
int spacedim>
935 Assert(dim == spacedim, ExcNotImplemented());
943template <
int dim,
int spacedim>
950 Assert(dim == spacedim, ExcNotImplemented());
958template <
int dim,
int spacedim>
970 for (
unsigned int i = 0; i < real_points.size(); ++i)
978template <
int dim,
int spacedim>
979std::unique_ptr<Mapping<dim, spacedim>>
982 return std::make_unique<MappingBox<dim, spacedim>>(*this);
std::pair< Point< spacedim, Number >, Point< spacedim, Number > > & get_boundary_points()
Number side_length(const unsigned int direction) const
Point< spacedim, Number > unit_to_real(const Point< spacedim, Number > &point) const
virtual std::size_t memory_consumption() const
const Tensor< 1, spacedim > & normal_vector(const unsigned int i) const
double weight(const unsigned int i) const
const std::vector< Point< dim > > & get_points() const
virtual std::size_t memory_consumption() const override
Tensor< 1, dim > traslation
virtual void reinit(const UpdateFlags update_flags, const Quadrature< dim > &quadrature) override
Tensor< 1, dim > cell_extents
std::vector< Point< dim > > quadrature_points
Tensor< 1, dim > inverse_cell_extents
void maybe_update_volume_elements(const InternalData &data) const
virtual bool is_compatible_with(const ReferenceCell &reference_cell) const override
MappingBox(const std::vector< BoundingBox< dim > > &local_boxes, const std::map< types::global_cell_index, types::global_cell_index > &polytope_translator)
void maybe_update_jacobian_derivatives(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
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 Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const override
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_data(const UpdateFlags, const Quadrature< dim > &quadrature) 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 update_cell_extents(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const InternalData &data) const
virtual Point< spacedim > transform_unit_to_real_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) const override
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
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::unique_ptr< Mapping< dim, spacedim > > clone() const override
void maybe_update_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
std::vector< BoundingBox< dim > > boxes
virtual std::unique_ptr< typename Mapping< dim, spacedim >::InternalDataBase > get_subface_data(const UpdateFlags flags, const Quadrature< dim - 1 > &quadrature) 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
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const override
void maybe_update_inverse_jacobians(const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const
void maybe_update_normal_vectors(const unsigned int face_no, const InternalData &data, std::vector< Tensor< 1, dim > > &normal_vectors) 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
std::map< types::global_cell_index, types::global_cell_index > polytope_translator
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcAccessToUninitializedField(std::string arg1)
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
update_jacobian_pushed_forward_2nd_derivatives
update_contravariant_transformation
update_jacobian_pushed_forward_grads
update_jacobian_3rd_derivatives
update_covariant_transformation
update_jacobian_pushed_forward_3rd_derivatives
update_jacobian_2nd_derivatives
mapping_covariant_gradient
mapping_contravariant_hessian
mapping_covariant_hessian
mapping_contravariant_gradient
DEAL_II_NAMESPACE_OPEN DeclExceptionMsg(ExcCellNotAssociatedWithBox, "You are using MappingBox, but the incoming element is not associated with a" "Bounding Box Cartesian.")
bool has_box(const CellType &cell, const std::map< types::global_cell_index, types::global_cell_index > &translator)
std::vector< index_type > data
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell< dim > &reference_cell)
std::enable_if_t< std::is_fundamental_v< T >, std::size_t > memory_consumption(const T &t)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
static constexpr unsigned int faces_per_cell
static constexpr std::array< Tensor< 1, dim >, faces_per_cell > unit_normal_vector