Reference documentation for deal.II version 9.2.0
\(\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\}}\)
Public Member Functions | Static Public Member Functions | Public Attributes | List of all members
internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType > Struct Template Reference

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

Inheritance diagram for internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >:
[legend]

Public Member Functions

void initialize (const ::Triangulation< dim > &tria, const std::vector< std::pair< unsigned int, unsigned int >> &cells, const FaceInfo< VectorizedArrayType::size()> &faces, const std::vector< unsigned int > &active_fe_index, const Mapping< dim > &mapping, const std::vector<::hp::QCollection< 1 >> &quad, const UpdateFlags update_flags_cells, const UpdateFlags update_flags_boundary_faces, const UpdateFlags update_flags_inner_faces, const UpdateFlags update_flags_faces_by_cells)
 
void update_mapping (const ::Triangulation< dim > &tria, const std::vector< std::pair< unsigned int, unsigned int >> &cells, const FaceInfo< VectorizedArrayType::size()> &faces, const std::vector< unsigned int > &active_fe_index, const Mapping< dim > &mapping)
 
GeometryType get_cell_type (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 TaskInfo &task_info) const
 
void compute_mapping_q (const ::Triangulation< dim > &tria, const std::vector< std::pair< unsigned int, unsigned int >> &cells, const std::vector< FaceToCellTopology< VectorizedArrayType::size()>> &faces)
 
void initialize_cells (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)
 
void initialize_faces (const ::Triangulation< dim > &tria, const std::vector< std::pair< unsigned int, unsigned int >> &cells, const std::vector< FaceToCellTopology< VectorizedArrayType::size()>> &faces, const Mapping< dim > &mapping)
 
void initialize_faces_by_cells (const ::Triangulation< dim > &tria, const std::vector< std::pair< unsigned int, unsigned int >> &cells, const Mapping< dim > &mapping)
 

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

UpdateFlags update_flags_cells
 
UpdateFlags update_flags_boundary_faces
 
UpdateFlags update_flags_inner_faces
 
UpdateFlags update_flags_faces_by_cells
 
std::vector< GeometryTypecell_type
 
std::vector< GeometryTypeface_type
 
std::vector< MappingInfoStorage< dim, dim, Number, VectorizedArrayType > > cell_data
 
std::vector< MappingInfoStorage< dim - 1, dim, Number, VectorizedArrayType > > face_data
 
std::vector< MappingInfoStorage< dim - 1, dim, Number, VectorizedArrayType > > face_data_by_cells
 
SmartPointer< const Mapping< dim > > mapping
 

Detailed Description

template<int dim, typename Number, typename VectorizedArrayType>
struct internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >

The class that stores all geometry-dependent data related with cell interiors for use in the matrix-free class.

Author
Katharina Kormann and Martin Kronbichler, 2010, 2011, 2017

Definition at line 316 of file mapping_info.h.

Member Function Documentation

◆ initialize()

template<int dim, typename Number , typename VectorizedArrayType >
void internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::initialize ( const ::Triangulation< dim > &  tria,
const std::vector< std::pair< unsigned int, unsigned int >> &  cells,
const FaceInfo< VectorizedArrayType::size()> &  faces,
const std::vector< unsigned int > &  active_fe_index,
const Mapping< dim > &  mapping,
const std::vector<::hp::QCollection< 1 >> &  quad,
const UpdateFlags  update_flags_cells,
const UpdateFlags  update_flags_boundary_faces,
const UpdateFlags  update_flags_inner_faces,
const UpdateFlags  update_flags_faces_by_cells 
)

Compute the information in the given cells and faces. 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.

◆ update_mapping()

template<int dim, typename Number , typename VectorizedArrayType >
void internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::update_mapping ( const ::Triangulation< dim > &  tria,
const std::vector< std::pair< unsigned int, unsigned int >> &  cells,
const FaceInfo< VectorizedArrayType::size()> &  faces,
const std::vector< unsigned int > &  active_fe_index,
const Mapping< dim > &  mapping 
)

Update the information in the given cells and faces that is the result of a change in the given mapping class, keeping the cells, quadrature formulas and other unknowns unchanged. This call is only valid if MappingInfo::initialize() has been called before.

◆ get_cell_type()

template<int dim, typename Number , typename VectorizedArrayType >
GeometryType internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::get_cell_type ( const unsigned int  cell_chunk_no) const
inline

Return the type of a given cell as detected during initialization.

Definition at line 643 of file mapping_info.h.

◆ clear()

template<int dim, typename Number , typename VectorizedArrayType >
void internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::clear ( )

Clear all data fields in this class.

◆ memory_consumption()

template<int dim, typename Number , typename VectorizedArrayType >
std::size_t internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::memory_consumption ( ) const

Return the memory consumption of this class in bytes.

◆ print_memory_consumption()

template<int dim, typename Number , typename VectorizedArrayType >
template<typename StreamType >
void internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::print_memory_consumption ( StreamType &  out,
const TaskInfo task_info 
) const

Prints a detailed summary of memory consumption in the different structures of this class to the given output stream.

◆ compute_mapping_q()

template<int dim, typename Number , typename VectorizedArrayType >
void internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::compute_mapping_q ( const ::Triangulation< dim > &  tria,
const std::vector< std::pair< unsigned int, unsigned int >> &  cells,
const std::vector< FaceToCellTopology< VectorizedArrayType::size()>> &  faces 
)

Internal function to compute the geometry for the case the mapping is a MappingQ and a single quadrature formula per slot (non-hp case) is used. This method computes all data from the underlying cell quadrature points using the fast operator evaluation techniques from the matrix-free framework itself, i.e., it uses a polynomial description of the cell geometry (that is computed in a first step) and then computes all Jacobians and normal vectors based on this information. This optimized approach is much faster than going through FEValues and FEFaceValues, especially when several different quadrature formulas are involved, and consumes less memory.

Parameters
triaThe triangulation to be used for setup
cellsThe actual cells of the triangulation to be worked on, given as a tuple of the level and index within the level as used in the main initialization of the class
facesThe description of the connectivity from faces to cells as filled in the MatrixFree class

◆ initialize_cells()

template<int dim, typename Number , typename VectorizedArrayType >
void internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::initialize_cells ( 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 
)

Computes the information in the given cells, called within initialize.

◆ initialize_faces()

template<int dim, typename Number , typename VectorizedArrayType >
void internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::initialize_faces ( const ::Triangulation< dim > &  tria,
const std::vector< std::pair< unsigned int, unsigned int >> &  cells,
const std::vector< FaceToCellTopology< VectorizedArrayType::size()>> &  faces,
const Mapping< dim > &  mapping 
)

Computes the information in the given faces, called within initialize.

◆ initialize_faces_by_cells()

template<int dim, typename Number , typename VectorizedArrayType >
void internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::initialize_faces_by_cells ( const ::Triangulation< dim > &  tria,
const std::vector< std::pair< unsigned int, unsigned int >> &  cells,
const Mapping< dim > &  mapping 
)

Computes the information in the given faces, called within initialize.

◆ compute_update_flags()

template<int dim, typename Number , typename VectorizedArrayType >
static UpdateFlags internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::compute_update_flags ( const UpdateFlags  update_flags,
const std::vector<::hp::QCollection< 1 >> &  quad = std::vector<::hp::QCollection< 1 >>() 
)
static

Helper function to determine which update flags must be set in the internal functions to initialize all data as requested by the user.

Member Data Documentation

◆ update_flags_cells

template<int dim, typename Number , typename VectorizedArrayType >
UpdateFlags internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::update_flags_cells

The given update flags for computing the geometry on the cells.

Definition at line 383 of file mapping_info.h.

◆ update_flags_boundary_faces

template<int dim, typename Number , typename VectorizedArrayType >
UpdateFlags internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::update_flags_boundary_faces

The given update flags for computing the geometry on the boundary faces.

Definition at line 389 of file mapping_info.h.

◆ update_flags_inner_faces

template<int dim, typename Number , typename VectorizedArrayType >
UpdateFlags internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::update_flags_inner_faces

The given update flags for computing the geometry on the interior faces.

Definition at line 395 of file mapping_info.h.

◆ update_flags_faces_by_cells

template<int dim, typename Number , typename VectorizedArrayType >
UpdateFlags internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::update_flags_faces_by_cells

The given update flags for computing the geometry on the faces for cell-centric loops.

Definition at line 401 of file mapping_info.h.

◆ cell_type

template<int dim, typename Number , typename VectorizedArrayType >
std::vector<GeometryType> internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::cell_type

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

◆ face_type

template<int dim, typename Number , typename VectorizedArrayType >
std::vector<GeometryType> internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::face_type

Stores whether a face (and both cells adjacent to the face) is Cartesian (face type 0), whether it represents an affine situation (face type 1), whether it is a flat face where the normal vector is the same throughout the face (face type 2), or is general (face type 3).

Definition at line 418 of file mapping_info.h.

◆ cell_data

template<int dim, typename Number , typename VectorizedArrayType >
std::vector<MappingInfoStorage<dim, dim, Number, VectorizedArrayType> > internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::cell_data

The data cache for the cells.

Definition at line 424 of file mapping_info.h.

◆ face_data

template<int dim, typename Number , typename VectorizedArrayType >
std::vector<MappingInfoStorage<dim - 1, dim, Number, VectorizedArrayType> > internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::face_data

The data cache for the faces.

Definition at line 430 of file mapping_info.h.

◆ face_data_by_cells

template<int dim, typename Number , typename VectorizedArrayType >
std::vector<MappingInfoStorage<dim - 1, dim, Number, VectorizedArrayType> > internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::face_data_by_cells

The data cache for the face-associated-with-cell topology, following the cell_type variable for the cell types.

Definition at line 437 of file mapping_info.h.

◆ mapping

template<int dim, typename Number , typename VectorizedArrayType >
SmartPointer<const Mapping<dim> > internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType >::mapping

The pointer to the underlying Mapping object.

Definition at line 442 of file mapping_info.h.


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