Reference documentation for deal.II version GIT relicensing-249-g48dc7357c7 2024-03-29 12:30:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Static Public Member Functions | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
NonMatching::MappingInfo< dim, spacedim, Number > Class Template Reference

#include <deal.II/non_matching/mapping_info.h>

Inheritance diagram for NonMatching::MappingInfo< dim, spacedim, Number >:
Inheritance graph
[legend]

Classes

struct  AdditionalData
 

Public Types

using VectorizedArrayType = typename ::internal::VectorizedArrayTrait< Number >::vectorized_value_type
 

Public Member Functions

 MappingInfo (const Mapping< dim, spacedim > &mapping, const UpdateFlags update_flags, const AdditionalData additional_data=AdditionalData())
 
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
 
boost::signals2::connection connect_is_reinitialized (const std::function< void()> &set_is_reinitialized)
 
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 ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (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 = ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim >
 
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
 

Private Attributes

State state
 
AlignedVector< Point< dim, VectorizedArrayType > > unit_points
 
AlignedVector< Point< dim - 1, VectorizedArrayType > > unit_points_faces
 
std::vector< unsigned intunit_points_index
 
std::shared_ptr< typename Mapping< dim, spacedim >::InternalDataBase > internal_mapping_data
 
const SmartPointer< const Mapping< dim, spacedim > > mapping
 
const UpdateFlags update_flags
 
UpdateFlags update_flags_mapping
 
const AdditionalData additional_data
 
std::vector<::internal::MatrixFreeFunctions::GeometryTypecell_type
 
std::vector< unsigned intdata_index_offsets
 
std::vector< unsigned intcompressed_data_index_offsets
 
AlignedVector< Number > JxW_values
 
AlignedVector< Tensor< 1, spacedim, Number > > normal_vectors
 
std::array< AlignedVector< DerivativeForm< 1, dim, spacedim, Number > >, 2 > jacobians
 
std::array< AlignedVector< DerivativeForm< 1, spacedim, dim, Number > >, 2 > inverse_jacobians
 
AlignedVector< Point< spacedim, Number > > real_points
 
std::vector< unsigned intn_q_points_unvectorized
 
std::vector< unsigned intcell_index_offset
 
std::vector< unsigned intcell_index_to_compressed_cell_index
 
bool do_cell_index_compression
 
boost::signals2::signal< void()> is_reinitialized
 
SmartPointer< const Triangulation< dim, spacedim > > triangulation
 
std::vector< std::pair< int, int > > cell_level_and_indices
 
std::vector< std::pair< unsigned char, unsigned char > > face_number
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Detailed Description

template<int dim, int spacedim = dim, typename Number = double>
class NonMatching::MappingInfo< dim, spacedim, Number >

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 246 of file mapping_info.h.

Member Typedef Documentation

◆ VectorizedArrayType

template<int dim, int spacedim = dim, typename Number = double>
using NonMatching::MappingInfo< dim, spacedim, Number >::VectorizedArrayType = typename ::internal::VectorizedArrayTrait< Number>::vectorized_value_type

The VectorizedArray type the unit points are stored in inside this class.

Definition at line 252 of file mapping_info.h.

◆ MappingData

template<int dim, int spacedim = dim, typename Number = double>
using NonMatching::MappingInfo< dim, spacedim, Number >::MappingData = ::internal::FEValuesImplementation::MappingRelatedData<dim, spacedim>
private

Definition at line 532 of file mapping_info.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Member Enumeration Documentation

◆ State

template<int dim, int spacedim = dim, typename Number = double>
enum class NonMatching::MappingInfo::State
strongprivate

Enum class for reinitialized states.

Enumerator
invalid 
single_cell 
cell_vector 
faces_on_cells_in_vector 
face_vector 

Definition at line 622 of file mapping_info.h.

Constructor & Destructor Documentation

◆ MappingInfo()

template<int dim, int spacedim, typename Number >
NonMatching::MappingInfo< dim, spacedim, Number >::MappingInfo ( const Mapping< dim, spacedim > &  mapping,
const UpdateFlags  update_flags,
const AdditionalData  additional_data = AdditionalData() 
)

Constructor.

Parameters
mappingThe Mapping class describing the geometry of a cell.
update_flagsSpecify 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_dataAdditional data for the class to specify the behavior during reinitialization.

Definition at line 799 of file mapping_info.h.

Member Function Documentation

◆ reinit() [1/3]

template<int dim, int spacedim, typename Number >
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 854 of file mapping_info.h.

◆ reinit() [2/3]

template<int dim, int spacedim, typename Number >
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 865 of file mapping_info.h.

◆ reinit() [3/3]

template<int dim, int spacedim, typename Number >
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 878 of file mapping_info.h.

◆ reinit_cells() [1/2]

template<int dim, int spacedim, typename Number >
template<typename ContainerType >
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 949 of file mapping_info.h.

◆ reinit_cells() [2/2]

template<int dim, int spacedim, typename Number >
template<typename ContainerType >
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 1171 of file mapping_info.h.

◆ reinit_surface()

template<int dim, int spacedim, typename Number >
template<typename ContainerType >
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 1204 of file mapping_info.h.

◆ reinit_faces() [1/2]

template<int dim, int spacedim, typename Number >
template<typename ContainerType >
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 1245 of file mapping_info.h.

◆ reinit_faces() [2/2]

template<int dim, int spacedim, typename Number >
template<typename CellIteratorType >
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 1492 of file mapping_info.h.

◆ is_face_state()

template<int dim, int spacedim, typename Number >
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 1764 of file mapping_info.h.

◆ get_face_number()

template<int dim, int spacedim, typename Number >
unsigned int NonMatching::MappingInfo< dim, spacedim, Number >::get_face_number ( const unsigned int  offset,
const bool  is_interior 
) const
inline

Returns the face number of the interior/exterior face.

Definition at line 787 of file mapping_info.h.

◆ get_unit_point()

template<int dim, int spacedim, typename Number >
const Point< dim, typename MappingInfo< dim, spacedim, Number >::VectorizedArrayType > * NonMatching::MappingInfo< dim, spacedim, Number >::get_unit_point ( const unsigned int  offset) const
inline

Getter function for unit points. The offset can be obtained with compute_unit_point_index_offset().

Definition at line 2058 of file mapping_info.h.

◆ get_unit_point_faces()

template<int dim, int spacedim, typename Number >
const Point< dim - 1, typename MappingInfo< dim, spacedim, Number >::VectorizedArrayType > * NonMatching::MappingInfo< dim, spacedim, Number >::get_unit_point_faces ( const unsigned int  offset) const
inline

Getter function for unit points on faces. The offset can be obtained with compute_unit_point_index_offset().

Definition at line 2070 of file mapping_info.h.

◆ get_jacobian()

template<int dim, int spacedim, typename Number >
const DerivativeForm< 1, dim, spacedim, Number > * NonMatching::MappingInfo< dim, spacedim, Number >::get_jacobian ( const unsigned int  offset,
const bool  is_interior = true 
) const
inline

Getter function for Jacobians. The offset can be obtained with compute_data_index_offset().

Definition at line 2119 of file mapping_info.h.

◆ get_inverse_jacobian()

template<int dim, int spacedim, typename Number >
const DerivativeForm< 1, spacedim, dim, Number > * NonMatching::MappingInfo< dim, spacedim, Number >::get_inverse_jacobian ( const unsigned int  offset,
const bool  is_interior = true 
) const
inline

Getter function for inverse Jacobians. The offset can be obtained with compute_data_index_offset().

Definition at line 2129 of file mapping_info.h.

◆ get_normal_vector()

template<int dim, int spacedim, typename Number >
const Tensor< 1, spacedim, Number > * NonMatching::MappingInfo< dim, spacedim, Number >::get_normal_vector ( const unsigned int  offset) const
inline

Getter function for normal vectors. The offset can be obtained with compute_data_index_offset().

Definition at line 2140 of file mapping_info.h.

◆ get_JxW()

template<int dim, int spacedim, typename Number >
const Number * NonMatching::MappingInfo< dim, spacedim, Number >::get_JxW ( const unsigned int  offset) const
inline

Getter function for Jacobian times quadrature weight (JxW). The offset can be obtained with compute_data_index_offset().

Definition at line 2150 of file mapping_info.h.

◆ get_real_point()

template<int dim, int spacedim, typename Number >
const Point< spacedim, Number > * NonMatching::MappingInfo< dim, spacedim, Number >::get_real_point ( const unsigned int  offset) const
inline

Getter function for real points. The offset can be obtained with compute_data_index_offset().

Definition at line 2080 of file mapping_info.h.

◆ get_mapping()

template<int dim, int spacedim, typename Number >
const Mapping< dim, spacedim > & NonMatching::MappingInfo< dim, spacedim, Number >::get_mapping ( ) const

Getter function for underlying mapping.

Definition at line 2159 of file mapping_info.h.

◆ get_update_flags()

template<int dim, int spacedim, typename Number >
UpdateFlags NonMatching::MappingInfo< dim, spacedim, Number >::get_update_flags ( ) const

Getter function for the update flags.

Definition at line 2168 of file mapping_info.h.

◆ get_update_flags_mapping()

template<int dim, int spacedim, typename Number >
UpdateFlags NonMatching::MappingInfo< dim, spacedim, Number >::get_update_flags_mapping ( ) const

Getter function for the mapping update flags.

Definition at line 2177 of file mapping_info.h.

◆ connect_is_reinitialized()

template<int dim, int spacedim, typename Number >
boost::signals2::connection NonMatching::MappingInfo< dim, spacedim, Number >::connect_is_reinitialized ( const std::function< void()> &  set_is_reinitialized)

Connects to is_reinitialized().

Definition at line 2186 of file mapping_info.h.

◆ compute_geometry_index_offset()

template<int dim, int spacedim, typename Number >
template<bool is_face>
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 1827 of file mapping_info.h.

◆ compute_unit_point_index_offset()

template<int dim, int spacedim, typename Number >
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 2090 of file mapping_info.h.

◆ compute_data_index_offset()

template<int dim, int spacedim, typename Number >
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 2100 of file mapping_info.h.

◆ compute_compressed_data_index_offset()

template<int dim, int spacedim, typename Number >
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 2109 of file mapping_info.h.

◆ get_n_q_points_unvectorized()

template<int dim, int spacedim, typename Number >
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 1773 of file mapping_info.h.

◆ get_cell_type()

template<int dim, int spacedim, typename Number >
internal::MatrixFreeFunctions::GeometryType NonMatching::MappingInfo< dim, spacedim, Number >::get_cell_type ( const unsigned int  geometry_index) const

Get cell geometry type.

Definition at line 1782 of file mapping_info.h.

◆ get_cell_iterator()

template<int dim, int spacedim, typename Number >
Triangulation< dim, spacedim >::cell_iterator NonMatching::MappingInfo< dim, spacedim, Number >::get_cell_iterator ( const unsigned int  cell_index) const

Return cell iterator.

Note
This call is only possible if AdditionalData::store_cells is enabled.

Definition at line 1792 of file mapping_info.h.

◆ memory_consumption()

template<int dim, int spacedim, typename Number >
std::size_t NonMatching::MappingInfo< dim, spacedim, Number >::memory_consumption ( ) const

Return the memory consumption of this class in bytes.

Definition at line 2196 of file mapping_info.h.

◆ compute_n_q_points()

template<int dim, int spacedim, typename Number >
template<typename NumberType >
unsigned int NonMatching::MappingInfo< dim, spacedim, Number >::compute_n_q_points ( const unsigned int  n_q_points_unvectorized)
private

Compute number of quadrature point batches depending on NumberType.

Definition at line 1809 of file mapping_info.h.

◆ clear()

template<int dim, int spacedim, typename Number >
void NonMatching::MappingInfo< dim, spacedim, Number >::clear ( )
private

Clears fields to make the object reusable.

Definition at line 841 of file mapping_info.h.

◆ resize_unit_points()

template<int dim, int spacedim, typename Number >
void NonMatching::MappingInfo< dim, spacedim, Number >::resize_unit_points ( const unsigned int  n_unit_point_batches)
private

Resize the unit_point data field.

Definition at line 2008 of file mapping_info.h.

◆ resize_unit_points_faces()

template<int dim, int spacedim, typename Number >
void NonMatching::MappingInfo< dim, spacedim, Number >::resize_unit_points_faces ( const unsigned int  n_unit_point_batches)
private

Resize the unit_point_faces data field.

Definition at line 2018 of file mapping_info.h.

◆ resize_data_fields()

template<int dim, int spacedim, typename Number >
void NonMatching::MappingInfo< dim, spacedim, Number >::resize_data_fields ( const unsigned int  n_data_point_batches,
const bool  is_face_centric = false 
)
private

Resize the mapping data fields.

Definition at line 2028 of file mapping_info.h.

◆ store_unit_points()

template<int dim, int spacedim, typename Number >
void NonMatching::MappingInfo< dim, spacedim, Number >::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 
)
private

Store the unit points.

Definition at line 1881 of file mapping_info.h.

◆ store_unit_points_faces()

template<int dim, int spacedim, typename Number >
void NonMatching::MappingInfo< dim, spacedim, Number >::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 
)
private

Store the unit points on faces.

Definition at line 1906 of file mapping_info.h.

◆ store_mapping_data()

template<int dim, int spacedim, typename Number >
void NonMatching::MappingInfo< dim, spacedim, Number >::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 
)
private

Store the requested mapping data.

Definition at line 1931 of file mapping_info.h.

◆ compute_compressed_cell_index()

template<int dim, int spacedim, typename Number >
unsigned int NonMatching::MappingInfo< dim, spacedim, Number >::compute_compressed_cell_index ( const unsigned int  cell_index) const
private

Compute the compressed cell index.

Definition at line 1863 of file mapping_info.h.

◆ do_reinit_cells()

template<int dim, int spacedim, typename Number >
template<typename ContainerType , typename QuadratureType >
void NonMatching::MappingInfo< dim, spacedim, Number >::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 
)
private

Compute the mapping information for cells/surface.

Definition at line 972 of file mapping_info.h.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
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.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
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.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ serialize()

template<class Archive >
void Subscriptor::serialize ( Archive &  ar,
const unsigned int  version 
)
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.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
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.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Member Data Documentation

◆ state

template<int dim, int spacedim = dim, typename Number = double>
State NonMatching::MappingInfo< dim, spacedim, Number >::state
private

Enum class that stores the currently initialized state upon the last call of reinit().

Definition at line 635 of file mapping_info.h.

◆ unit_points

template<int dim, int spacedim = dim, typename Number = double>
AlignedVector<Point<dim, VectorizedArrayType> > NonMatching::MappingInfo< dim, spacedim, Number >::unit_points
private

The reference points specified at reinit().

Indexed by unit_points_index.

Definition at line 642 of file mapping_info.h.

◆ unit_points_faces

template<int dim, int spacedim = dim, typename Number = double>
AlignedVector<Point<dim - 1, VectorizedArrayType> > NonMatching::MappingInfo< dim, spacedim, Number >::unit_points_faces
private

The reference points on faces specified at reinit().

Indexed by unit_points_index.

Definition at line 649 of file mapping_info.h.

◆ unit_points_index

template<int dim, int spacedim = dim, typename Number = double>
std::vector<unsigned int> NonMatching::MappingInfo< dim, spacedim, Number >::unit_points_index
private

Offset to point to the first unit point of a cell/face.

Definition at line 654 of file mapping_info.h.

◆ internal_mapping_data

template<int dim, int spacedim = dim, typename Number = double>
std::shared_ptr<typename Mapping<dim, spacedim>::InternalDataBase> NonMatching::MappingInfo< dim, spacedim, Number >::internal_mapping_data
private

A pointer to the internal data of the underlying mapping.

Definition at line 660 of file mapping_info.h.

◆ mapping

template<int dim, int spacedim = dim, typename Number = double>
const SmartPointer<const Mapping<dim, spacedim> > NonMatching::MappingInfo< dim, spacedim, Number >::mapping
private

A pointer to the underlying mapping.

Definition at line 665 of file mapping_info.h.

◆ update_flags

template<int dim, int spacedim = dim, typename Number = double>
const UpdateFlags NonMatching::MappingInfo< dim, spacedim, Number >::update_flags
private

The desired update flags for the evaluation.

Definition at line 670 of file mapping_info.h.

◆ update_flags_mapping

template<int dim, int spacedim = dim, typename Number = double>
UpdateFlags NonMatching::MappingInfo< dim, spacedim, Number >::update_flags_mapping
private

The update flags for the desired mapping information.

Definition at line 675 of file mapping_info.h.

◆ additional_data

template<int dim, int spacedim = dim, typename Number = double>
const AdditionalData NonMatching::MappingInfo< dim, spacedim, Number >::additional_data
private

AdditionalData for this object.

Definition at line 680 of file mapping_info.h.

◆ cell_type

template<int dim, int spacedim = dim, typename Number = double>
std::vector<::internal::MatrixFreeFunctions::GeometryType> NonMatching::MappingInfo< dim, spacedim, Number >::cell_type
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 688 of file mapping_info.h.

◆ data_index_offsets

template<int dim, int spacedim = dim, typename Number = double>
std::vector<unsigned int> NonMatching::MappingInfo< dim, spacedim, Number >::data_index_offsets
private

Stores the index offset into the arrays JxW_values and normal_vectors.

Definition at line 693 of file mapping_info.h.

◆ compressed_data_index_offsets

template<int dim, int spacedim = dim, typename Number = double>
std::vector<unsigned int> NonMatching::MappingInfo< dim, spacedim, Number >::compressed_data_index_offsets
private

Stores the index offset into the arrays jacobians and inverse_jacobians.

Definition at line 698 of file mapping_info.h.

◆ JxW_values

template<int dim, int spacedim = dim, typename Number = double>
AlignedVector<Number> NonMatching::MappingInfo< dim, spacedim, Number >::JxW_values
private

The storage of the Jacobian determinant times the quadrature weight on quadrature points.

Indexed by data_index_offsets.

Definition at line 706 of file mapping_info.h.

◆ normal_vectors

template<int dim, int spacedim = dim, typename Number = double>
AlignedVector<Tensor<1, spacedim, Number> > NonMatching::MappingInfo< dim, spacedim, Number >::normal_vectors
private

Stores the normal vectors.

Indexed by data_index_offsets.

Definition at line 713 of file mapping_info.h.

◆ jacobians

template<int dim, int spacedim = dim, typename Number = double>
std::array<AlignedVector<DerivativeForm<1, dim, spacedim, Number> >, 2> NonMatching::MappingInfo< dim, spacedim, Number >::jacobians
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 722 of file mapping_info.h.

◆ inverse_jacobians

template<int dim, int spacedim = dim, typename Number = double>
std::array<AlignedVector<DerivativeForm<1, spacedim, dim, Number> >, 2> NonMatching::MappingInfo< dim, spacedim, Number >::inverse_jacobians
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 732 of file mapping_info.h.

◆ real_points

template<int dim, int spacedim = dim, typename Number = double>
AlignedVector<Point<spacedim, Number> > NonMatching::MappingInfo< dim, spacedim, Number >::real_points
private

The mapped real points.

Indexed by data_index_offsets.

Definition at line 739 of file mapping_info.h.

◆ n_q_points_unvectorized

template<int dim, int spacedim = dim, typename Number = double>
std::vector<unsigned int> NonMatching::MappingInfo< dim, spacedim, Number >::n_q_points_unvectorized
private

Number of unvectorized unit points per geometric entity (cell/face).

Definition at line 744 of file mapping_info.h.

◆ cell_index_offset

template<int dim, int spacedim = dim, typename Number = double>
std::vector<unsigned int> NonMatching::MappingInfo< dim, spacedim, Number >::cell_index_offset
private

Offset to point to the first element of a cell in internal data containers.

Definition at line 750 of file mapping_info.h.

◆ cell_index_to_compressed_cell_index

template<int dim, int spacedim = dim, typename Number = double>
std::vector<unsigned int> NonMatching::MappingInfo< dim, spacedim, Number >::cell_index_to_compressed_cell_index
private

A vector that converts the cell index to a compressed cell index for e.g. a filtered IteratorRange.

Definition at line 756 of file mapping_info.h.

◆ do_cell_index_compression

template<int dim, int spacedim = dim, typename Number = double>
bool NonMatching::MappingInfo< dim, spacedim, Number >::do_cell_index_compression
private

A bool that determines weather cell index compression should be done.

Definition at line 761 of file mapping_info.h.

◆ is_reinitialized

template<int dim, int spacedim = dim, typename Number = double>
boost::signals2::signal<void()> NonMatching::MappingInfo< dim, spacedim, Number >::is_reinitialized
private

This signal is triggered right after this object is reinitialized, to let dependent objects know that they need to reinitialize as well.

Definition at line 767 of file mapping_info.h.

◆ triangulation

template<int dim, int spacedim = dim, typename Number = double>
SmartPointer<const Triangulation<dim, spacedim> > NonMatching::MappingInfo< dim, spacedim, Number >::triangulation
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 774 of file mapping_info.h.

◆ cell_level_and_indices

template<int dim, int spacedim = dim, typename Number = double>
std::vector<std::pair<int, int> > NonMatching::MappingInfo< dim, spacedim, Number >::cell_level_and_indices
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 780 of file mapping_info.h.

◆ face_number

template<int dim, int spacedim = dim, typename Number = double>
std::vector<std::pair<unsigned char, unsigned char> > NonMatching::MappingInfo< dim, spacedim, Number >::face_number
private

Definition at line 782 of file mapping_info.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
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.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
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.

◆ object_info

const std::type_info* Subscriptor::object_info
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.

◆ mutex

std::mutex Subscriptor::mutex
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.


The documentation for this class was generated from the following file: