Reference documentation for deal.II version 9.0.0
|
#include <deal.II/matrix_free/mapping_info.h>
Public Member Functions | |
unsigned int | quad_index_from_n_q_points (const unsigned int n_q_points) const |
template<typename StreamType > | |
void | print_memory_consumption (StreamType &out, const SizeInfo &task_info) const |
std::size_t | memory_consumption () const |
Public Attributes | |
std::vector< QuadratureDescriptor > | descriptor |
AlignedVector< unsigned int > | data_index_offsets |
AlignedVector< VectorizedArray< Number > > | JxW_values |
AlignedVector< Tensor< 1, spacedim, VectorizedArray< Number > > > | normal_vectors |
AlignedVector< Tensor< 2, spacedim, VectorizedArray< Number > > > | jacobians [2] |
AlignedVector< Tensor< 1, spacedim *(spacedim+1)/2, Tensor< 1, spacedim, VectorizedArray< Number > > > > | jacobian_gradients [2] |
AlignedVector< Tensor< 1, spacedim, VectorizedArray< Number > > > | normals_times_jacobians [2] |
AlignedVector< unsigned int > | quadrature_point_offsets |
AlignedVector< Point< spacedim, VectorizedArray< Number > > > | quadrature_points |
Definition of a structure that stores all cached data related to the evaluated geometry from the mapping. In order to support hp-adaptivity and compressed storage (in particular for Jacobians, JxW values, and normals), storage length can be different for different rows. Thus, it allows to jump at the data of individual rows similar to compressed row storage in sparse matrices. We have two different start indices for fields with different sizes. The first category of offsets are the indices for Jacobians of the transformation from unit to real cell (we store the inverse Jacobian), second derivatives, JxW values, and normal vectors. We keep separate arrays for all these data structures because a user code might access only some of them. In such a case, one array will be gone through in a contiguous order with access to all entries, which makes it easy for the processor to prefetch data. Having all data in a single array would require some strides in the access pattern, which is much more complicated for the processor to predict (and indeed leads to prefetching of data that does not get used on Intel processors such as BroadwellEP).
The second category of indices are the offsets for the quadrature points. Quadrature points can be compressed less than the other fields and thus need longer fields. Quadrature point indices are often used in other contexts such as evaluation of right hand sides.
The third component is a descriptor of data from the unit cells, called QuadratureDescriptor, which contains the quadrature weights and permutations of how to go through quadrature points in case of face data. The latter comes in a vector for the support of hp adaptivity, with several data fields for the individual quadrature formulas.
Definition at line 106 of file mapping_info.h.
|
inline |
Returns the quadrature index for a given number of quadrature points. If not in hp mode or if the index is not found, this function always returns index 0. Hence, this function does not check whether the given degree is actually present.
Definition at line 485 of file mapping_info.h.
void internal::MatrixFreeFunctions::MappingInfoStorage< structdim, spacedim, Number >::print_memory_consumption | ( | StreamType & | out, |
const SizeInfo & | task_info | ||
) | const |
Prints a detailed summary of memory consumption in the different structures of this class to the given output stream.
std::size_t internal::MatrixFreeFunctions::MappingInfoStorage< structdim, spacedim, Number >::memory_consumption | ( | ) | const |
Returns the memory consumption in bytes.
std::vector<QuadratureDescriptor> internal::MatrixFreeFunctions::MappingInfoStorage< structdim, spacedim, Number >::descriptor |
A class describing the layout of the sections in the data_storage
field and also includes some data that depends on the number of quadrature points in the hp context such as the inner quadrature formula and re-indexing for faces that are not in the standard orientation.
Definition at line 165 of file mapping_info.h.
AlignedVector<unsigned int> internal::MatrixFreeFunctions::MappingInfoStorage< structdim, spacedim, Number >::data_index_offsets |
Stores the index offset into the arrays jxw_values
, jacobians
, normal_vectors
and the second derivatives. Note that affine cells have shorter fields of length 1, where the others have lengths equal to the numer of quadrature points of the given cell.
Definition at line 173 of file mapping_info.h.
AlignedVector<VectorizedArray<Number> > internal::MatrixFreeFunctions::MappingInfoStorage< structdim, spacedim, Number >::JxW_values |
The storage of the Jacobian determinant (times the quadrature weight in case the transformation is non-affine) on quadrature points.
Indexed by data_index_offsets
.
Definition at line 182 of file mapping_info.h.
AlignedVector<Tensor<1,spacedim,VectorizedArray<Number> > > internal::MatrixFreeFunctions::MappingInfoStorage< structdim, spacedim, Number >::normal_vectors |
Stores the normal vectors.
Indexed by data_index_offsets
.
Definition at line 189 of file mapping_info.h.
AlignedVector<Tensor<2,spacedim,VectorizedArray<Number> > > internal::MatrixFreeFunctions::MappingInfoStorage< structdim, spacedim, Number >::jacobians[2] |
The storage of covariant transformation on quadrature points, i.e., the inverse and transposed Jacobians of the transformation from the unit to the real cell.
Indexed by data_index_offsets
.
Contains two fields for access from both sides for interior faces, but the default case (cell integrals or boundary integrals) only fills the zeroth component and ignores the first one.
Definition at line 202 of file mapping_info.h.
AlignedVector<Tensor<1,spacedim *(spacedim+1)/2, Tensor<1,spacedim,VectorizedArray<Number> > > > internal::MatrixFreeFunctions::MappingInfoStorage< structdim, spacedim, Number >::jacobian_gradients[2] |
The storage of the gradients of the inverse Jacobian transformation. Because of symmetry, only the upper diagonal and diagonal part are needed. The first index runs through the derivatives, starting with the diagonal and then continuing row-wise, i.e., \(\partial^2/\partial x_1 \partial x_2\) first, then \(\partial^2/\partial x_1 \partial x_3\), and so on. The second index is the spatial coordinate.
Indexed by data_index_offsets
.
Contains two fields for access from both sides for interior faces, but the default case (cell integrals or boundary integrals) only fills the zeroth component and ignores the first one.
Definition at line 220 of file mapping_info.h.
AlignedVector<Tensor<1,spacedim,VectorizedArray<Number> > > internal::MatrixFreeFunctions::MappingInfoStorage< structdim, spacedim, Number >::normals_times_jacobians[2] |
Stores the Jacobian transformations times the normal vector (this represents a shortcut that is accessed often and can thus get higher performance).
Indexed by data_index_offsets
.
Definition at line 229 of file mapping_info.h.
AlignedVector<unsigned int> internal::MatrixFreeFunctions::MappingInfoStorage< structdim, spacedim, Number >::quadrature_point_offsets |
Stores the index offset of a particular cell into the quadrature points array in real coordinates. Note that Cartesian cells have shorter fields (length is n_q_points_1d
) than non-Cartesian cells (length is n_q_points
) or faces.
Definition at line 237 of file mapping_info.h.
AlignedVector<Point<spacedim,VectorizedArray<Number> > > internal::MatrixFreeFunctions::MappingInfoStorage< structdim, spacedim, Number >::quadrature_points |
Stores the quadrature points in real coordinates, including a compression scheme for Cartesian cells where we do not need to store the full data on all points.
Indexed by quadrature_point_offsets
.
Definition at line 246 of file mapping_info.h.