Reference documentation for deal.II version 9.6.0
|
#include <deal.II/non_matching/mapping_info.h>
Classes | |
struct | AdditionalData |
Public Types | |
using | VectorizedArrayType |
Public Member Functions | |
MappingInfo (const Mapping< dim, spacedim > &mapping, const UpdateFlags update_flags, const AdditionalData additional_data=AdditionalData()) | |
MappingInfo (const MappingInfo &)=delete | |
MappingInfo & | operator= (const MappingInfo &)=delete |
void | reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const std::vector< Point< dim > > &unit_points) |
void | reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< dim > > &unit_points) |
void | reinit (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Quadrature< dim > &quadrature) |
template<typename ContainerType > | |
void | reinit_cells (const ContainerType &cell_iterator_range, const std::vector< std::vector< Point< dim > > > &unit_points_vector, const unsigned int n_unfiltered_cells=numbers::invalid_unsigned_int) |
template<typename ContainerType > | |
void | reinit_cells (const ContainerType &cell_iterator_range, const std::vector< Quadrature< dim > > &quadrature_vector, const unsigned int n_unfiltered_cells=numbers::invalid_unsigned_int) |
template<typename ContainerType > | |
void | reinit_surface (const ContainerType &cell_iterator_range, const std::vector< ImmersedSurfaceQuadrature< dim > > &quadrature_vector, const unsigned int n_unfiltered_cells=numbers::invalid_unsigned_int) |
template<typename ContainerType > | |
void | reinit_faces (const ContainerType &cell_iterator_range, const std::vector< std::vector< Quadrature< dim - 1 > > > &quadrature_vector, const unsigned int n_unfiltered_cells=numbers::invalid_unsigned_int) |
template<typename CellIteratorType > | |
void | reinit_faces (const std::vector< std::pair< CellIteratorType, unsigned int > > &face_iterator_range_interior, const std::vector< Quadrature< dim - 1 > > &quadrature_vector) |
bool | is_face_state () const |
unsigned int | get_face_number (const unsigned int offset, const bool is_interior) const |
const Point< dim, VectorizedArrayType > * | get_unit_point (const unsigned int offset) const |
const Point< dim - 1, VectorizedArrayType > * | get_unit_point_faces (const unsigned int offset) const |
const DerivativeForm< 1, dim, spacedim, Number > * | get_jacobian (const unsigned int offset, const bool is_interior=true) const |
const DerivativeForm< 1, spacedim, dim, Number > * | get_inverse_jacobian (const unsigned int offset, const bool is_interior=true) const |
const Tensor< 1, spacedim, Number > * | get_normal_vector (const unsigned int offset) const |
const Number * | get_JxW (const unsigned int offset) const |
const Point< spacedim, Number > * | get_real_point (const unsigned int offset) const |
const Mapping< dim, spacedim > & | get_mapping () const |
UpdateFlags | get_update_flags () const |
UpdateFlags | get_update_flags_mapping () const |
template<bool is_face> | |
unsigned int | compute_geometry_index_offset (const unsigned int cell_index, const unsigned int face_number) const |
unsigned int | compute_unit_point_index_offset (const unsigned int geometry_index) const |
unsigned int | compute_data_index_offset (const unsigned int geometry_index) const |
unsigned int | compute_compressed_data_index_offset (const unsigned int geometry_index) const |
unsigned int | get_n_q_points_unvectorized (const unsigned int geometry_index) const |
::internal::MatrixFreeFunctions::GeometryType | get_cell_type (const unsigned int geometry_index) const |
Triangulation< dim, spacedim >::cell_iterator | get_cell_iterator (const unsigned int cell_index) const |
std::size_t | memory_consumption () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
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 |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
Static Public Member Functions | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Private Types | |
enum class | State { invalid , single_cell , cell_vector , faces_on_cells_in_vector , face_vector } |
using | MappingData |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
Private Member Functions | |
template<typename NumberType > | |
unsigned int | compute_n_q_points (const unsigned int n_q_points_unvectorized) |
void | clear () |
void | resize_unit_points (const unsigned int n_unit_point_batches) |
void | resize_unit_points_faces (const unsigned int n_unit_point_batches) |
void | resize_data_fields (const unsigned int n_data_point_batches, const bool is_face_centric=false) |
void | store_unit_points (const unsigned int unit_points_index_offset, const unsigned int n_q_points, const unsigned int n_q_points_unvectorized, const std::vector< Point< dim > > &points) |
void | store_unit_points_faces (const unsigned int unit_points_index_offset, const unsigned int n_q_points, const unsigned int n_q_points_unvectorized, const std::vector< Point< dim - 1 > > &points) |
void | store_mapping_data (const unsigned int unit_points_index_offset, const unsigned int n_q_points, const unsigned int n_q_points_unvectorized, const MappingData &mapping_data, const std::vector< double > &weights, const unsigned int compressed_unit_point_index_offset, const bool affine_cell, const bool is_interior=true) |
unsigned int | compute_compressed_cell_index (const unsigned int cell_index) const |
template<typename ContainerType , typename QuadratureType > | |
void | do_reinit_cells (const ContainerType &cell_iterator_range, const std::vector< QuadratureType > &quadrature_vector, const unsigned int n_unfiltered_cells, const std::function< void(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const QuadratureType &quadrature, MappingData &mapping_data)> &compute_mapping_data) |
void | check_no_subscribers () const noexcept |
Static Private Attributes | |
static std::mutex | mutex |
This class provides the mapping information computation and mapping data storage to be used together with FEPointEvaluation.
MappingInfo is vectorized across quadrature points, which means data can be provided in vectorized format.
Two different modes are available: partially vectorized and fully vectorized across quadrature points. Partially vectorized (scalar Number template argument) means the computed mapping data is provided in scalar format and only unit points are provided vectorized (for seamless interaction with FEPointEvaluation). Fully vectorized (vectorized Number template argument, e.g. VectorizedArray<double>) provides both mapping data and unit points in vectorized format. The Number template parameter of MappingInfo and FEPointEvaluation has to be identical.
Definition at line 221 of file mapping_info.h.
using NonMatching::MappingInfo< dim, spacedim, Number >::VectorizedArrayType |
The VectorizedArray type the unit points are stored in inside this class.
Definition at line 227 of file mapping_info.h.
|
private |
Definition at line 512 of file mapping_info.h.
|
privateinherited |
The data type used in counter_map.
Definition at line 229 of file subscriptor.h.
|
privateinherited |
The iterator type used in counter_map.
Definition at line 234 of file subscriptor.h.
|
strongprivate |
Enum class for reinitialized states.
Enumerator | |
---|---|
invalid | |
single_cell | |
cell_vector | |
faces_on_cells_in_vector | |
face_vector |
Definition at line 602 of file mapping_info.h.
NonMatching::MappingInfo< dim, spacedim, Number >::MappingInfo | ( | const Mapping< dim, spacedim > & | mapping, |
const UpdateFlags | update_flags, | ||
const AdditionalData | additional_data = AdditionalData() ) |
Constructor.
mapping | The Mapping class describing the geometry of a cell. |
update_flags | Specify the quantities to be computed by the mapping during the call of reinit(). These update flags are also handed to a FEEvaluation object if you construct it with this MappingInfo object. |
additional_data | Additional data for the class to specify the behavior during reinitialization. |
Definition at line 789 of file mapping_info.h.
|
delete |
Do not allow making copies.
|
delete |
Do not allow copy assignment of this class.
void NonMatching::MappingInfo< dim, spacedim, Number >::reinit | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, |
const std::vector< Point< dim > > & | unit_points ) |
Compute the mapping information for the incoming cell and unit points. This overload is needed to resolve ambiguity.
Definition at line 838 of file mapping_info.h.
void NonMatching::MappingInfo< dim, spacedim, Number >::reinit | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, |
const ArrayView< const Point< dim > > & | unit_points ) |
Compute the mapping information for the incoming cell and unit points.
Definition at line 849 of file mapping_info.h.
void NonMatching::MappingInfo< dim, spacedim, Number >::reinit | ( | const typename Triangulation< dim, spacedim >::cell_iterator & | cell, |
const Quadrature< dim > & | quadrature ) |
Compute the mapping information for the given cell and quadrature formula. As opposed to the other reinit
function, this method allows to access a JxW
factor at the points.
Definition at line 861 of file mapping_info.h.
void NonMatching::MappingInfo< dim, spacedim, Number >::reinit_cells | ( | const ContainerType & | cell_iterator_range, |
const std::vector< std::vector< Point< dim > > > & | unit_points_vector, | ||
const unsigned int | n_unfiltered_cells = numbers::invalid_unsigned_int ) |
Compute the mapping information for the incoming iterable container of cell iterators and corresponding vector of unit points.
It is possible to give an IteratorRange<FilteredIterator> to this function and together with n_unfiltered_cells
specified this object can compress its storage while giving you access to the underlying data with the "unfiltered" index. The default argument for n_unfiltered_cells
disables this built-in compression.
Definition at line 930 of file mapping_info.h.
void NonMatching::MappingInfo< dim, spacedim, Number >::reinit_cells | ( | const ContainerType & | cell_iterator_range, |
const std::vector< Quadrature< dim > > & | quadrature_vector, | ||
const unsigned int | n_unfiltered_cells = numbers::invalid_unsigned_int ) |
Compute the mapping information for the incoming iterable container of cell iterators and corresponding vector of quadratures. As opposed to the other reinit
function, this method allows to access a JxW
factor at the points.
Definition at line 1150 of file mapping_info.h.
void NonMatching::MappingInfo< dim, spacedim, Number >::reinit_surface | ( | const ContainerType & | cell_iterator_range, |
const std::vector< ImmersedSurfaceQuadrature< dim > > & | quadrature_vector, | ||
const unsigned int | n_unfiltered_cells = numbers::invalid_unsigned_int ) |
Compute the mapping information for the incoming vector of cells and corresponding vector of ImmersedSurfaceQuadrature.
Definition at line 1183 of file mapping_info.h.
void NonMatching::MappingInfo< dim, spacedim, Number >::reinit_faces | ( | const ContainerType & | cell_iterator_range, |
const std::vector< std::vector< Quadrature< dim - 1 > > > & | quadrature_vector, | ||
const unsigned int | n_unfiltered_cells = numbers::invalid_unsigned_int ) |
Compute the mapping information for all faces of the incoming vector of cells and corresponding vector of quadratures.
Definition at line 1224 of file mapping_info.h.
void NonMatching::MappingInfo< dim, spacedim, Number >::reinit_faces | ( | const std::vector< std::pair< CellIteratorType, unsigned int > > & | face_iterator_range_interior, |
const std::vector< Quadrature< dim - 1 > > & | quadrature_vector ) |
Compute the mapping information incoming vector of faces and corresponding vector of quadratures.
Definition at line 1469 of file mapping_info.h.
bool NonMatching::MappingInfo< dim, spacedim, Number >::is_face_state | ( | ) | const |
Return if this MappingInfo object is reinitialized for faces (by reinit_faces()) or not.
Definition at line 1740 of file mapping_info.h.
|
inline |
Returns the face number of the interior/exterior face.
Definition at line 777 of file mapping_info.h.
|
inline |
Getter function for unit points. The offset can be obtained with compute_unit_point_index_offset().
Definition at line 2034 of file mapping_info.h.
|
inline |
Getter function for unit points on faces. The offset can be obtained with compute_unit_point_index_offset().
Definition at line 2046 of file mapping_info.h.
|
inline |
Getter function for Jacobians. The offset can be obtained with compute_data_index_offset().
Definition at line 2095 of file mapping_info.h.
|
inline |
Getter function for inverse Jacobians. The offset can be obtained with compute_data_index_offset().
Definition at line 2105 of file mapping_info.h.
|
inline |
Getter function for normal vectors. The offset can be obtained with compute_data_index_offset().
Definition at line 2116 of file mapping_info.h.
|
inline |
Getter function for Jacobian times quadrature weight (JxW). The offset can be obtained with compute_data_index_offset().
Definition at line 2126 of file mapping_info.h.
|
inline |
Getter function for real points. The offset can be obtained with compute_data_index_offset().
Definition at line 2056 of file mapping_info.h.
const Mapping< dim, spacedim > & NonMatching::MappingInfo< dim, spacedim, Number >::get_mapping | ( | ) | const |
Getter function for underlying mapping.
Definition at line 2135 of file mapping_info.h.
UpdateFlags NonMatching::MappingInfo< dim, spacedim, Number >::get_update_flags | ( | ) | const |
Getter function for the update flags.
Definition at line 2144 of file mapping_info.h.
UpdateFlags NonMatching::MappingInfo< dim, spacedim, Number >::get_update_flags_mapping | ( | ) | const |
Getter function for the mapping update flags.
Definition at line 2153 of file mapping_info.h.
unsigned int NonMatching::MappingInfo< dim, spacedim, Number >::compute_geometry_index_offset | ( | const unsigned int | cell_index, |
const unsigned int | face_number ) const |
Compute the geometry index offset of the current cell/face.
Definition at line 1803 of file mapping_info.h.
unsigned int NonMatching::MappingInfo< dim, spacedim, Number >::compute_unit_point_index_offset | ( | const unsigned int | geometry_index | ) | const |
Compute the unit points index offset for the current cell/face.
Definition at line 2066 of file mapping_info.h.
unsigned int NonMatching::MappingInfo< dim, spacedim, Number >::compute_data_index_offset | ( | const unsigned int | geometry_index | ) | const |
Compute the data index offset for the current cell/face.
Definition at line 2076 of file mapping_info.h.
unsigned int NonMatching::MappingInfo< dim, spacedim, Number >::compute_compressed_data_index_offset | ( | const unsigned int | geometry_index | ) | const |
Compute the data index offset for the current cell/face.
Definition at line 2085 of file mapping_info.h.
unsigned int NonMatching::MappingInfo< dim, spacedim, Number >::get_n_q_points_unvectorized | ( | const unsigned int | geometry_index | ) | const |
Get number of unvectorized quadrature points.
Definition at line 1749 of file mapping_info.h.
internal::MatrixFreeFunctions::GeometryType NonMatching::MappingInfo< dim, spacedim, Number >::get_cell_type | ( | const unsigned int | geometry_index | ) | const |
Get cell geometry type.
Definition at line 1758 of file mapping_info.h.
Triangulation< dim, spacedim >::cell_iterator NonMatching::MappingInfo< dim, spacedim, Number >::get_cell_iterator | ( | const unsigned int | cell_index | ) | const |
Return cell iterator.
Definition at line 1768 of file mapping_info.h.
std::size_t NonMatching::MappingInfo< dim, spacedim, Number >::memory_consumption | ( | ) | const |
Return the memory consumption of this class in bytes.
Definition at line 2162 of file mapping_info.h.
|
private |
Compute number of quadrature point batches depending on NumberType.
Definition at line 1785 of file mapping_info.h.
|
private |
Clears fields to make the object reusable.
Definition at line 825 of file mapping_info.h.
|
private |
Resize the unit_point data field.
Definition at line 1984 of file mapping_info.h.
|
private |
Resize the unit_point_faces data field.
Definition at line 1994 of file mapping_info.h.
|
private |
Resize the mapping data fields.
Definition at line 2004 of file mapping_info.h.
|
private |
Store the unit points.
Definition at line 1857 of file mapping_info.h.
|
private |
Store the unit points on faces.
Definition at line 1882 of file mapping_info.h.
|
private |
Store the requested mapping data.
Definition at line 1907 of file mapping_info.h.
|
private |
Compute the compressed cell index.
Definition at line 1839 of file mapping_info.h.
|
private |
Compute the mapping information for cells/surface.
Definition at line 953 of file mapping_info.h.
|
inherited |
Subscribes a user of the object by storing the pointer validity
. The subscriber may be identified by text supplied as identifier
.
Definition at line 135 of file subscriptor.cc.
|
inherited |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 155 of file subscriptor.cc.
|
inlineinherited |
Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.
Definition at line 300 of file subscriptor.h.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 317 of file subscriptor.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 203 of file subscriptor.cc.
|
inlineinherited |
Read or write the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.
This function does not actually serialize any of the member variables of this class. The reason is that what this class stores is only who subscribes to this object, but who does so at the time of storing the contents of this object does not necessarily have anything to do with who subscribes to the object when it is restored. Consequently, we do not want to overwrite the subscribers at the time of restoring, and then there is no reason to write the subscribers out in the first place.
Definition at line 309 of file subscriptor.h.
|
privatenoexceptinherited |
Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.
Definition at line 52 of file subscriptor.cc.
|
private |
Enum class that stores the currently initialized state upon the last call of reinit().
Definition at line 615 of file mapping_info.h.
|
private |
The reference points specified at reinit().
Indexed by unit_points_index
.
Definition at line 622 of file mapping_info.h.
|
private |
The reference points on faces specified at reinit().
Indexed by unit_points_index
.
Definition at line 629 of file mapping_info.h.
|
private |
Offset to point to the first unit point of a cell/face.
Definition at line 634 of file mapping_info.h.
|
private |
A pointer to the internal data of the underlying mapping.
Definition at line 640 of file mapping_info.h.
|
private |
Quadrature object that is used for the reinit(ArrayView) function.
Definition at line 645 of file mapping_info.h.
|
private |
Helper class that temporarily holds the data requested for one cell by Mapping::fill_fe_values() before it is filled into appropriate data structures for consumption by FEPointEvaluation.
Definition at line 652 of file mapping_info.h.
|
private |
A pointer to the underlying mapping.
Definition at line 657 of file mapping_info.h.
|
private |
The desired update flags for the evaluation.
Definition at line 662 of file mapping_info.h.
|
private |
The update flags for the desired mapping information.
Definition at line 667 of file mapping_info.h.
|
private |
AdditionalData for this object.
Definition at line 672 of file mapping_info.h.
|
private |
Stores whether a cell is Cartesian (cell type 0), has constant transform data (Jacobians) (cell type 1), or is general (cell type 3). Type 2 is only used for faces and no cells are assigned this value.
Definition at line 680 of file mapping_info.h.
|
private |
Stores the index offset into the arrays JxW_values
and normal_vectors
.
Definition at line 685 of file mapping_info.h.
|
private |
Stores the index offset into the arrays jacobians
and inverse_jacobians
.
Definition at line 690 of file mapping_info.h.
|
private |
The storage of the Jacobian determinant times the quadrature weight on quadrature points.
Indexed by data_index_offsets
.
Definition at line 698 of file mapping_info.h.
|
private |
Stores the normal vectors.
Indexed by data_index_offsets
.
Definition at line 705 of file mapping_info.h.
|
private |
The storage of contravariant transformation on quadrature points, i.e., the Jacobians of the transformation from the unit to the real cell.
Indexed by compressed_data_index_offsets
.
Definition at line 714 of file mapping_info.h.
|
private |
The storage of covariant transformation on quadrature points, i.e., the inverse Jacobians of the transformation from the unit to the real cell.
Indexed by compressed_data_index_offsets
.
Definition at line 724 of file mapping_info.h.
|
private |
The mapped real points.
Indexed by data_index_offsets
.
Definition at line 731 of file mapping_info.h.
|
private |
Number of unvectorized unit points per geometric entity (cell/face).
Definition at line 736 of file mapping_info.h.
|
private |
Offset to point to the first element of a cell in internal data containers.
Definition at line 742 of file mapping_info.h.
|
private |
A vector that converts the cell index to a compressed cell index for e.g. a filtered IteratorRange.
Definition at line 748 of file mapping_info.h.
|
private |
A bool that determines weather cell index compression should be done.
Definition at line 753 of file mapping_info.h.
|
private |
Reference to the triangulation passed via the cells to the reinit functions. This field is only set if AdditionalData::store_cells is enabled.
Definition at line 760 of file mapping_info.h.
|
private |
Level and indices of cells passed to the reinit functions. This vector is only filled if AdditionalData::store_cells is enabled.
Definition at line 766 of file mapping_info.h.
|
private |
Store the face number of interior or exterior faces as set up in the initialization.
Definition at line 772 of file mapping_info.h.
|
mutableprivateinherited |
Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).
The creator (and owner) of an object is counted in the map below if HE manages to supply identification.
We use the mutable
keyword in order to allow subscription to constant objects also.
This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic
class template.
Definition at line 218 of file subscriptor.h.
|
mutableprivateinherited |
In this map, we count subscriptions for each different identification string supplied to subscribe().
Definition at line 224 of file subscriptor.h.
|
mutableprivateinherited |
In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.
Definition at line 240 of file subscriptor.h.
|
mutableprivateinherited |
Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.
Definition at line 248 of file subscriptor.h.
|
staticprivateinherited |
A mutex used to ensure data consistency when accessing the mutable
members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers()
.
Definition at line 271 of file subscriptor.h.