Reference documentation for deal.II version 8.5.1
|
#include <deal.II/matrix_free/mapping_info.h>
Classes | |
struct | CellData |
struct | MappingInfoDependent |
Public Member Functions | |
MappingInfo () | |
void | initialize (const ::Triangulation< dim > &tria, const std::vector< std::pair< unsigned int, unsigned int > > &cells, const std::vector< unsigned int > &active_fe_index, const Mapping< dim > &mapping, const std::vector<::hp::QCollection< 1 > > &quad, const UpdateFlags update_flags) |
CellType | get_cell_type (const unsigned int cell_chunk_no) const |
unsigned int | get_cell_data_index (const unsigned int cell_chunk_no) const |
void | clear () |
std::size_t | memory_consumption () const |
template<typename StreamType > | |
void | print_memory_consumption (StreamType &out, const SizeInfo &size_info) const |
void | evaluate_on_cell (const ::Triangulation< dim > &tria, const std::pair< unsigned int, unsigned int > *cells, const unsigned int cell, const unsigned int my_q, CellType(&cell_t_prev)[n_vector_elements], CellType(&cell_t)[n_vector_elements], ::FEValues< dim, dim > &fe_values, CellData &cell_data) const |
Static Public Member Functions | |
static UpdateFlags | compute_update_flags (const UpdateFlags update_flags, const std::vector<::hp::QCollection< 1 > > &quad=std::vector<::hp::QCollection< 1 > >()) |
Public Attributes | |
std::vector< unsigned int > | cell_type |
AlignedVector< std::pair< Tensor< 1, dim, VectorizedArray< Number > >, VectorizedArray< Number > > > | cartesian_data |
AlignedVector< std::pair< Tensor< 2, dim, VectorizedArray< Number > >, VectorizedArray< Number > > > | affine_data |
std::vector< MappingInfoDependent > | mapping_data_gen |
bool | JxW_values_initialized |
bool | second_derivatives_initialized |
bool | quadrature_points_initialized |
Static Public Attributes | |
static const std::size_t | n_cell_type_bits = 2 |
static const unsigned int | n_cell_types = 1U<<n_cell_type_bits |
static const unsigned int | n_vector_elements = VectorizedArray<Number>::n_array_elements |
The class that stores all geometry-dependent data related with cell interiors for use in the matrix-free class.
Definition at line 46 of file mapping_info.h.
internal::MatrixFreeFunctions::MappingInfo< dim, Number >::MappingInfo | ( | ) |
Empty constructor.
void internal::MatrixFreeFunctions::MappingInfo< dim, Number >::initialize | ( | const ::Triangulation< dim > & | tria, |
const std::vector< std::pair< unsigned int, unsigned int > > & | cells, | ||
const std::vector< unsigned int > & | active_fe_index, | ||
const Mapping< dim > & | mapping, | ||
const std::vector<::hp::QCollection< 1 > > & | quad, | ||
const UpdateFlags | update_flags | ||
) |
Compute the information in the given cells. The cells are specified by the level and the index within the level (as given by CellIterator::level() and CellIterator::index(), in order to allow for different kinds of iterators, e.g. standard DoFHandler, multigrid, etc.) on a fixed Triangulation. In addition, a mapping and several quadrature formulas are given.
|
static |
Helper function to determine which update flags must be set in the internal functions to initialize all data as requested by the user.
|
inline |
Return the type of a given cell as detected during initialization.
Definition at line 356 of file mapping_info.h.
|
inline |
Return the type of a given cell as detected during initialization.
Definition at line 369 of file mapping_info.h.
void internal::MatrixFreeFunctions::MappingInfo< dim, Number >::clear | ( | ) |
Clear all data fields in this class.
std::size_t internal::MatrixFreeFunctions::MappingInfo< dim, Number >::memory_consumption | ( | ) | const |
Return the memory consumption of this class in bytes.
void internal::MatrixFreeFunctions::MappingInfo< dim, Number >::print_memory_consumption | ( | StreamType & | out, |
const SizeInfo & | size_info | ||
) | const |
Prints a detailed summary of memory consumption in the different structures of this class to the given output stream.
void internal::MatrixFreeFunctions::MappingInfo< dim, Number >::evaluate_on_cell | ( | const ::Triangulation< dim > & | tria, |
const std::pair< unsigned int, unsigned int > * | cells, | ||
const unsigned int | cell, | ||
const unsigned int | my_q, | ||
CellType(&) | cell_t_prev[n_vector_elements], | ||
CellType(&) | cell_t[n_vector_elements], | ||
::FEValues< dim, dim > & | fe_values, | ||
CellData & | cell_data | ||
) | const |
Helper function called internally during the initialize function.
|
static |
Determines how many bits of an unsigned int are used to distinguish the cell types (Cartesian, with constant Jacobian, or general)
Definition at line 52 of file mapping_info.h.
|
static |
Determines how many types of different cells can be detected at most. Corresponds to the number of bits we reserved for it.
Definition at line 58 of file mapping_info.h.
|
static |
An abbreviation for the length of vector lines of the current data type.
Definition at line 64 of file mapping_info.h.
std::vector<unsigned int> internal::MatrixFreeFunctions::MappingInfo< dim, Number >::cell_type |
Stores whether a cell is Cartesian, has constant transform data (Jacobians) or is general. cell_type % 4 gives this information (0: Cartesian, 1: constant Jacobian throughout cell, 2: general cell), and cell_type / 4 gives the index in the data field of where to find the information in the fields Jacobian and JxW values (except for quadrature points, for which the index runs as usual).
Definition at line 131 of file mapping_info.h.
AlignedVector<std::pair<Tensor<1,dim,VectorizedArray<Number> >, VectorizedArray<Number> > > internal::MatrixFreeFunctions::MappingInfo< dim, Number >::cartesian_data |
The first field stores the inverse Jacobian for Cartesian cells: There, it is a diagonal rank-2 tensor, so we actually just store a rank-1 tensor. It is the same on all cells, therefore we only store it once per cell, and use similarities from one cell to another, too (on structured meshes, there are usually many cells with the same Jacobian).
The second field stores the Jacobian determinant for Cartesian cells (without the quadrature weight, which depends on the quadrature point, whereas the determinant is the same on each quadrature point).
Definition at line 146 of file mapping_info.h.
AlignedVector<std::pair<Tensor<2,dim,VectorizedArray<Number> >, VectorizedArray<Number> > > internal::MatrixFreeFunctions::MappingInfo< dim, Number >::affine_data |
The first field stores the Jacobian for non-Cartesian cells where all the Jacobians on the cell are the same (i.e., constant, which comes from a linear transformation from unit to real cell). Also use similarities from one cell to another (on structured meshes, there are usually many cells with the same Jacobian).
The second field stores the Jacobian determinant for non-Cartesian cells with constant Jacobian throughout the cell (without the quadrature weight, which depends on the quadrature point, whereas the determinant is the same on each quadrature point).
Definition at line 161 of file mapping_info.h.
std::vector<MappingInfoDependent> internal::MatrixFreeFunctions::MappingInfo< dim, Number >::mapping_data_gen |
Contains all the stuff that depends on the quadrature formula
Definition at line 290 of file mapping_info.h.
bool internal::MatrixFreeFunctions::MappingInfo< dim, Number >::JxW_values_initialized |
Stores whether JxW values have been initialized
Definition at line 295 of file mapping_info.h.
bool internal::MatrixFreeFunctions::MappingInfo< dim, Number >::second_derivatives_initialized |
Stores whether we computed second derivatives.
Definition at line 300 of file mapping_info.h.
bool internal::MatrixFreeFunctions::MappingInfo< dim, Number >::quadrature_points_initialized |
Stores whether we computed quadrature points.
Definition at line 305 of file mapping_info.h.