#include <mapping_box.h>
Classes | |
class | InternalData |
Public Member Functions | |
MappingBox (const std::vector< BoundingBox< dim > > &local_boxes, const std::map< types::global_cell_index, types::global_cell_index > &polytope_translator) | |
virtual std::unique_ptr< Mapping< dim, spacedim > > | clone () const override |
virtual bool | preserves_vertex_locations () const override |
virtual bool | is_compatible_with (const ReferenceCell &reference_cell) const override |
Mapping points between reference and real cells | |
virtual Point< spacedim > | transform_unit_to_real_cell (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< dim > &p) 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 |
Functions to transform tensors from reference to real coordinates | |
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 | transform (const ArrayView< const DerivativeForm< 1, dim, spacedim > > &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 2, spacedim > > &output) const override |
virtual void | transform (const ArrayView< const Tensor< 2, dim > > &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 2, spacedim > > &output) const override |
virtual void | transform (const ArrayView< const DerivativeForm< 2, dim, spacedim > > &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 3, spacedim > > &output) const override |
virtual void | transform (const ArrayView< const Tensor< 3, dim > > &input, const MappingKind kind, const typename Mapping< dim, spacedim >::InternalDataBase &internal, const ArrayView< Tensor< 3, spacedim > > &output) const override |
Public Member Functions inherited from Mapping< dim, dim > | |
virtual | ~Mapping () override=default |
virtual boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_cell > | get_vertices (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const |
boost::container::small_vector< Point< spacedim >, GeometryInfo< dim >::vertices_per_face > | get_vertices (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no) const |
virtual Point< spacedim > | get_center (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const bool map_barycenter_of_reference_cell=true) const |
virtual BoundingBox< spacedim > | get_bounding_box (const typename Triangulation< dim, spacedim >::cell_iterator &cell) const |
void | serialize (Archive &ar, const unsigned int version) |
Point< dim - 1 > | project_real_point_to_unit_point_on_face (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Point< spacedim > &p) const |
Point< dim - 1 > | project_real_point_to_unit_point_on_face (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Point< spacedim > &p) const |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
unsigned int | n_subscriptions () const |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | serialize (Archive &ar, const unsigned int version) |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
void | subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
Private Member Functions | |
void | update_cell_extents (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const InternalData &data) const |
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_normal_vectors (const unsigned int face_no, const InternalData &data, std::vector< Tensor< 1, dim > > &normal_vectors) const |
void | maybe_update_jacobian_derivatives (const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const |
void | maybe_update_volume_elements (const InternalData &data) const |
void | maybe_update_jacobians (const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const |
void | maybe_update_inverse_jacobians (const InternalData &data, const CellSimilarity::Similarity cell_similarity, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const |
Interface with FEValues | |
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< 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 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_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 |
Private Attributes | |
std::vector< BoundingBox< dim > > | boxes |
std::map< types::global_cell_index, types::global_cell_index > | polytope_translator |
Additional Inherited Members | |
Static Public Member Functions inherited from Mapping< dim, dim > | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
static ::ExceptionBase & | ExcInvalidData () |
static ::ExceptionBase & | ExcInvalidData () |
static ::ExceptionBase & | ExcTransformationFailed () |
static ::ExceptionBase & | ExcTransformationFailed () |
static ::ExceptionBase & | ExcDistortedMappedCell (Point< spacedim > arg1, double arg2, int arg3) |
static ::ExceptionBase & | ExcDistortedMappedCell (Point< spacedim > arg1, double arg2, int arg3) |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Member Functions inherited from Mapping< dim, dim > | |
virtual std::unique_ptr< InternalDataBase > | get_face_data (const UpdateFlags update_flags, const hp::QCollection< dim - 1 > &quadrature) const |
virtual std::unique_ptr< InternalDataBase > | get_face_data (const UpdateFlags update_flags, const Quadrature< dim - 1 > &quadrature) const |
virtual std::unique_ptr< InternalDataBase > | get_face_data (const UpdateFlags update_flags, const hp::QCollection< dim - 1 > &quadrature) const |
virtual std::unique_ptr< InternalDataBase > | get_face_data (const UpdateFlags update_flags, const Quadrature< dim - 1 > &quadrature) const |
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 |
virtual void | fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const |
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 |
virtual void | fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const typename Mapping< dim, spacedim >::InternalDataBase &internal_data, internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &output_data) const |
A class providing a mapping from the reference cell to cells that are axiparallel, i.e., that have the shape of rectangles (in 2d) or boxes (in 3d) with edges parallel to the coordinate directions. The class therefore provides functionality that is equivalent to what, for example, MappingQ would provide for such cells. However, knowledge of the shape of cells allows this class to be substantially more efficient.
Specifically, the mapping is meant for cells for which the mapping from the reference to the real cell is a scaling along the coordinate directions: The transformation from reference coordinates $\hat {\mathbf x}$ to real coordinates $\mathbf x$ on each cell is of the form
\begin{align*} {\mathbf x}(\hat {\mathbf x}) = \begin{pmatrix} h_x & 0 \\ 0 & h_y \end{pmatrix} \hat{\mathbf x} + {\mathbf v}_0 \end{align*}
in 2d, and
\begin{align*} {\mathbf x}(\hat {\mathbf x}) = \begin{pmatrix} h_x & 0 & 0 \\ 0 & h_y & 0 \\ 0 & 0 & h_z \end{pmatrix} \hat{\mathbf x} + {\mathbf v}_0 \end{align*}
in 3d, where ${\mathbf v}_0$ is the bottom left vertex and $h_x,h_y,h_z$ are the extents of the cell along the axes.
The class is intended for efficiency, and it does not do a whole lot of error checking. If you apply this mapping to a cell that does not conform to the requirements above, you will get strange results.
Definition at line 79 of file mapping_box.h.
MappingBox< dim, spacedim >::MappingBox | ( | const std::vector< BoundingBox< dim > > & | local_boxes, |
const std::map< types::global_cell_index, types::global_cell_index > & | polytope_translator ) |
Definition at line 68 of file mapping_box.cc.
References Assert.
|
overridevirtual |
Implements Mapping< dim, dim >.
Definition at line 975 of file mapping_box.cc.
|
overrideprivatevirtual |
Reimplemented from Mapping< dim, dim >.
Definition at line 467 of file mapping_box.cc.
References Assert, AssertDimension, Quadrature< int dim >::get_points(), has_box(), internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::JxW_values, MappingBox< dim, spacedim >::maybe_update_cell_quadrature_points(), MappingBox< dim, spacedim >::maybe_update_inverse_jacobians(), MappingBox< dim, spacedim >::maybe_update_jacobian_derivatives(), MappingBox< dim, spacedim >::maybe_update_jacobians(), MappingBox< dim, spacedim >::maybe_update_volume_elements(), CellSimilarity::none, NonMatching::ImmersedSurfaceQuadrature< int dim, int spacedim >::normal_vector(), internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::normal_vectors, MappingBox< dim, spacedim >::polytope_translator, internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::quadrature_points, MappingBox< dim, spacedim >::update_cell_extents(), Mapping::InternalDataBase::update_each, update_JxW_values, update_normal_vectors, and Quadrature< int dim >::weight().
|
overrideprivatevirtual |
Implements Mapping< dim, dim >.
Definition at line 445 of file mapping_box.cc.
References DEAL_II_NOT_IMPLEMENTED.
|
overrideprivatevirtual |
Implements Mapping< dim, dim >.
Definition at line 395 of file mapping_box.cc.
References Assert, MappingBox< dim, spacedim >::InternalData::cell_extents, d(), Quadrature< int dim >::get_points(), has_box(), internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::JxW_values, MappingBox< dim, spacedim >::maybe_update_cell_quadrature_points(), MappingBox< dim, spacedim >::maybe_update_inverse_jacobians(), MappingBox< dim, spacedim >::maybe_update_jacobian_derivatives(), MappingBox< dim, spacedim >::maybe_update_jacobians(), MappingBox< dim, spacedim >::polytope_translator, internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::quadrature_points, CellSimilarity::translation, MappingBox< dim, spacedim >::update_cell_extents(), Mapping::InternalDataBase::update_each, update_JxW_values, update_volume_elements, MappingBox< dim, spacedim >::InternalData::volume_element, and Quadrature< int dim >::weight().
|
overrideprivatevirtual |
Implements Mapping< dim, dim >.
Definition at line 168 of file mapping_box.cc.
References MappingBox< dim, spacedim >::requires_update_flags().
|
overrideprivatevirtual |
Implements Mapping< dim, dim >.
Definition at line 182 of file mapping_box.cc.
References DEAL_II_NOT_IMPLEMENTED.
|
overridevirtual |
Implements Mapping< dim, dim >.
Definition at line 133 of file mapping_box.cc.
References Assert, reference_cell(), and Utilities::to_string().
|
private |
Compute the quadrature points if the UpdateFlags of the incoming InternalData object say that they should be updated.
Called from fill_fe_values.
Definition at line 246 of file mapping_box.cc.
References MappingBox< dim, spacedim >::boxes, MappingBox< dim, spacedim >::polytope_translator, Mapping::InternalDataBase::update_each, and update_quadrature_points.
Referenced by MappingBox< dim, spacedim >::fill_fe_immersed_surface_values(), and MappingBox< dim, spacedim >::fill_fe_values().
|
private |
Compute the inverse Jacobians if the UpdateFlags of the incoming InternalData object say that they should be updated.
Definition at line 372 of file mapping_box.cc.
References MappingBox< dim, spacedim >::InternalData::inverse_cell_extents, internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::inverse_jacobians, CellSimilarity::translation, Mapping::InternalDataBase::update_each, and update_inverse_jacobians.
Referenced by MappingBox< dim, spacedim >::fill_fe_immersed_surface_values(), and MappingBox< dim, spacedim >::fill_fe_values().
|
private |
Since the Jacobian is constant for this mapping all derivatives of the Jacobian are identically zero. Fill these quantities with zeros if the corresponding update flags say that they should be updated.
Definition at line 282 of file mapping_box.cc.
References internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::jacobian_2nd_derivatives, internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::jacobian_3rd_derivatives, internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::jacobian_grads, internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::jacobian_pushed_forward_2nd_derivatives, internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::jacobian_pushed_forward_3rd_derivatives, internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::jacobian_pushed_forward_grads, CellSimilarity::translation, Mapping::InternalDataBase::update_each, update_jacobian_2nd_derivatives, update_jacobian_3rd_derivatives, update_jacobian_grads, update_jacobian_pushed_forward_2nd_derivatives, update_jacobian_pushed_forward_3rd_derivatives, and update_jacobian_pushed_forward_grads.
Referenced by MappingBox< dim, spacedim >::fill_fe_immersed_surface_values(), and MappingBox< dim, spacedim >::fill_fe_values().
|
private |
Compute the Jacobians if the UpdateFlags of the incoming InternalData object say that they should be updated.
Definition at line 350 of file mapping_box.cc.
References MappingBox< dim, spacedim >::InternalData::cell_extents, internal::FEValuesImplementation::MappingRelatedData< int dim, int spacedim >::jacobians, CellSimilarity::translation, Mapping::InternalDataBase::update_each, and update_jacobians.
Referenced by MappingBox< dim, spacedim >::fill_fe_immersed_surface_values(), and MappingBox< dim, spacedim >::fill_fe_values().
|
private |
Compute the normal vectors if the UpdateFlags of the incoming InternalData object say that they should be updated.
Definition at line 263 of file mapping_box.cc.
References Assert, Mapping::InternalDataBase::update_each, and update_normal_vectors.
|
private |
Compute the volume elements if the UpdateFlags of the incoming InternalData object say that they should be updated.
Definition at line 334 of file mapping_box.cc.
References MappingBox< dim, spacedim >::InternalData::cell_extents, d(), Mapping::InternalDataBase::update_each, update_volume_elements, volume(), and MappingBox< dim, spacedim >::InternalData::volume_element.
Referenced by MappingBox< dim, spacedim >::fill_fe_immersed_surface_values().
|
overridevirtual |
Return true
because MappingBox preserves vertex locations.
Implements Mapping< dim, dim >.
Definition at line 124 of file mapping_box.cc.
|
overrideprivatevirtual |
Implements Mapping< dim, dim >.
Definition at line 150 of file mapping_box.cc.
References update_boundary_forms, and update_normal_vectors.
Referenced by MappingBox< dim, spacedim >::get_data().
|
overridevirtual |
Implements Mapping< dim, dim >.
Definition at line 569 of file mapping_box.cc.
References Assert, AssertDimension, MappingBox< dim, spacedim >::InternalData::cell_extents, DEAL_II_NOT_IMPLEMENTED, MappingBox< dim, spacedim >::InternalData::inverse_cell_extents, mapping_contravariant, mapping_contravariant_gradient, mapping_covariant, mapping_covariant_gradient, mapping_piola, mapping_piola_gradient, update_contravariant_transformation, update_covariant_transformation, Mapping::InternalDataBase::update_each, update_volume_elements, and MappingBox< dim, spacedim >::InternalData::volume_element.
|
overridevirtual |
Implements Mapping< dim, dim >.
Definition at line 795 of file mapping_box.cc.
References Assert, AssertDimension, DEAL_II_NOT_IMPLEMENTED, MappingBox< dim, spacedim >::InternalData::inverse_cell_extents, mapping_covariant_gradient, update_contravariant_transformation, and Mapping::InternalDataBase::update_each.
|
overridevirtual |
Implements Mapping< dim, dim >.
Definition at line 509 of file mapping_box.cc.
References Assert, AssertDimension, MappingBox< dim, spacedim >::InternalData::cell_extents, d(), DEAL_II_NOT_IMPLEMENTED, MappingBox< dim, spacedim >::InternalData::inverse_cell_extents, mapping_contravariant, mapping_covariant, mapping_piola, update_contravariant_transformation, update_covariant_transformation, Mapping::InternalDataBase::update_each, update_volume_elements, and MappingBox< dim, spacedim >::InternalData::volume_element.
|
overridevirtual |
Implements Mapping< dim, dim >.
Definition at line 682 of file mapping_box.cc.
References Assert, AssertDimension, MappingBox< dim, spacedim >::InternalData::cell_extents, DEAL_II_NOT_IMPLEMENTED, MappingBox< dim, spacedim >::InternalData::inverse_cell_extents, mapping_contravariant, mapping_contravariant_gradient, mapping_covariant, mapping_covariant_gradient, mapping_piola, mapping_piola_gradient, update_contravariant_transformation, update_covariant_transformation, Mapping::InternalDataBase::update_each, update_volume_elements, and MappingBox< dim, spacedim >::InternalData::volume_element.
|
overridevirtual |
Implements Mapping< dim, dim >.
Definition at line 834 of file mapping_box.cc.
References Assert, AssertDimension, MappingBox< dim, spacedim >::InternalData::cell_extents, DEAL_II_NOT_IMPLEMENTED, MappingBox< dim, spacedim >::InternalData::inverse_cell_extents, mapping_contravariant_hessian, mapping_covariant_hessian, mapping_piola_hessian, update_contravariant_transformation, update_covariant_transformation, Mapping::InternalDataBase::update_each, update_volume_elements, and MappingBox< dim, spacedim >::InternalData::volume_element.
|
overridevirtual |
Reimplemented from Mapping< dim, dim >.
Definition at line 955 of file mapping_box.cc.
References Assert, AssertDimension, MappingBox< dim, spacedim >::boxes, DEAL_II_NOT_IMPLEMENTED, has_box(), and MappingBox< dim, spacedim >::polytope_translator.
|
overridevirtual |
Implements Mapping< dim, dim >.
Definition at line 940 of file mapping_box.cc.
References Assert, MappingBox< dim, spacedim >::boxes, has_box(), and MappingBox< dim, spacedim >::polytope_translator.
|
overridevirtual |
Implements Mapping< dim, dim >.
Definition at line 925 of file mapping_box.cc.
References Assert, MappingBox< dim, spacedim >::boxes, has_box(), and MappingBox< dim, spacedim >::polytope_translator.
|
private |
Update the cell_extents field of the incoming InternalData object with the size of the incoming cell.
Definition at line 196 of file mapping_box.cc.
References Assert, MappingBox< dim, spacedim >::boxes, MappingBox< dim, spacedim >::InternalData::cell_extents, d(), BoundingBox< int spacedim, typename Number >::get_boundary_points(), MappingBox< dim, spacedim >::InternalData::inverse_cell_extents, MappingBox< dim, spacedim >::polytope_translator, BoundingBox< int spacedim, typename Number >::side_length(), CellSimilarity::translation, and MappingBox< dim, spacedim >::InternalData::traslation.
Referenced by MappingBox< dim, spacedim >::fill_fe_immersed_surface_values(), and MappingBox< dim, spacedim >::fill_fe_values().
|
private |
Vector of (local) bounding boxes
Definition at line 372 of file mapping_box.h.
Referenced by MappingBox< dim, spacedim >::maybe_update_cell_quadrature_points(), MappingBox< dim, spacedim >::transform_points_real_to_unit_cell(), MappingBox< dim, spacedim >::transform_real_to_unit_cell(), MappingBox< dim, spacedim >::transform_unit_to_real_cell(), and MappingBox< dim, spacedim >::update_cell_extents().
|
private |
Map from global cell index to bounding box index
Definition at line 378 of file mapping_box.h.
Referenced by MappingBox< dim, spacedim >::fill_fe_immersed_surface_values(), MappingBox< dim, spacedim >::fill_fe_values(), MappingBox< dim, spacedim >::maybe_update_cell_quadrature_points(), MappingBox< dim, spacedim >::transform_points_real_to_unit_cell(), MappingBox< dim, spacedim >::transform_real_to_unit_cell(), MappingBox< dim, spacedim >::transform_unit_to_real_cell(), and MappingBox< dim, spacedim >::update_cell_extents().