Reference documentation for deal.II version GIT relicensing-220-g40aaea9664 2024-03-28 18:00: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 | Static Public Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
MatrixFree< dim, Number, VectorizedArrayType > Class Template Reference

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

Inheritance diagram for MatrixFree< dim, Number, VectorizedArrayType >:
Inheritance graph
[legend]

Classes

struct  AdditionalData
 

Public Types

enum class  DataAccessOnFaces {
  none , values , values_all_faces , gradients ,
  gradients_all_faces , unspecified
}
 
using value_type = Number
 
using vectorized_value_type = VectorizedArrayType
 

Public Member Functions

template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Construction and initialization
 MatrixFree ()
 
 MatrixFree (const MatrixFree< dim, Number, VectorizedArrayType > &other)
 
 ~MatrixFree () override=default
 
template<typename QuadratureType , typename number2 , typename MappingType >
void reinit (const MappingType &mapping, const DoFHandler< dim > &dof_handler, const AffineConstraints< number2 > &constraint, const QuadratureType &quad, const AdditionalData &additional_data=AdditionalData())
 
template<typename QuadratureType , typename number2 , typename MappingType >
void reinit (const MappingType &mapping, const std::vector< const DoFHandler< dim > * > &dof_handler, const std::vector< const AffineConstraints< number2 > * > &constraint, const std::vector< QuadratureType > &quad, const AdditionalData &additional_data=AdditionalData())
 
template<typename QuadratureType , typename number2 , typename MappingType >
void reinit (const MappingType &mapping, const std::vector< const DoFHandler< dim > * > &dof_handler, const std::vector< const AffineConstraints< number2 > * > &constraint, const QuadratureType &quad, const AdditionalData &additional_data=AdditionalData())
 
void copy_from (const MatrixFree< dim, Number, VectorizedArrayType > &matrix_free_base)
 
void update_mapping (const Mapping< dim > &mapping)
 
void update_mapping (const std::shared_ptr< hp::MappingCollection< dim > > &mapping)
 
void clear ()
 
Matrix-free loops
template<typename OutVector , typename InVector >
void cell_loop (const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &cell_operation, OutVector &dst, const InVector &src, const bool zero_dst_vector=false) const
 
template<typename CLASS , typename OutVector , typename InVector >
void cell_loop (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const, const CLASS *owning_class, OutVector &dst, const InVector &src, const bool zero_dst_vector=false) const
 
template<typename CLASS , typename OutVector , typename InVector >
void cell_loop (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &), CLASS *owning_class, OutVector &dst, const InVector &src, const bool zero_dst_vector=false) const
 
template<typename CLASS , typename OutVector , typename InVector >
void cell_loop (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const, const CLASS *owning_class, OutVector &dst, const InVector &src, const std::function< void(const unsigned int, const unsigned int)> &operation_before_loop, const std::function< void(const unsigned int, const unsigned int)> &operation_after_loop, const unsigned int dof_handler_index_pre_post=0) const
 
template<typename CLASS , typename OutVector , typename InVector >
void cell_loop (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &), CLASS *owning_class, OutVector &dst, const InVector &src, const std::function< void(const unsigned int, const unsigned int)> &operation_before_loop, const std::function< void(const unsigned int, const unsigned int)> &operation_after_loop, const unsigned int dof_handler_index_pre_post=0) const
 
template<typename OutVector , typename InVector >
void cell_loop (const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &cell_operation, OutVector &dst, const InVector &src, const std::function< void(const unsigned int, const unsigned int)> &operation_before_loop, const std::function< void(const unsigned int, const unsigned int)> &operation_after_loop, const unsigned int dof_handler_index_pre_post=0) const
 
template<typename OutVector , typename InVector >
void loop (const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &cell_operation, const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &face_operation, const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &boundary_operation, OutVector &dst, const InVector &src, const bool zero_dst_vector=false, const DataAccessOnFaces dst_vector_face_access=DataAccessOnFaces::unspecified, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified) const
 
template<typename CLASS , typename OutVector , typename InVector >
void loop (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const, void(CLASS::*face_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const, void(CLASS::*boundary_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const, const CLASS *owning_class, OutVector &dst, const InVector &src, const bool zero_dst_vector=false, const DataAccessOnFaces dst_vector_face_access=DataAccessOnFaces::unspecified, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified) const
 
template<typename CLASS , typename OutVector , typename InVector >
void loop (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &), void(CLASS::*face_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &), void(CLASS::*boundary_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &), CLASS *owning_class, OutVector &dst, const InVector &src, const bool zero_dst_vector=false, const DataAccessOnFaces dst_vector_face_access=DataAccessOnFaces::unspecified, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified) const
 
template<typename CLASS , typename OutVector , typename InVector >
void loop (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const, void(CLASS::*face_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const, void(CLASS::*boundary_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const, const CLASS *owning_class, OutVector &dst, const InVector &src, const std::function< void(const unsigned int, const unsigned int)> &operation_before_loop, const std::function< void(const unsigned int, const unsigned int)> &operation_after_loop, const unsigned int dof_handler_index_pre_post=0, const DataAccessOnFaces dst_vector_face_access=DataAccessOnFaces::unspecified, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified) const
 
template<typename CLASS , typename OutVector , typename InVector >
void loop (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &), void(CLASS::*face_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &), void(CLASS::*boundary_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &), const CLASS *owning_class, OutVector &dst, const InVector &src, const std::function< void(const unsigned int, const unsigned int)> &operation_before_loop, const std::function< void(const unsigned int, const unsigned int)> &operation_after_loop, const unsigned int dof_handler_index_pre_post=0, const DataAccessOnFaces dst_vector_face_access=DataAccessOnFaces::unspecified, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified) const
 
template<typename OutVector , typename InVector >
void loop (const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &cell_operation, const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &face_operation, const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &boundary_operation, OutVector &dst, const InVector &src, const std::function< void(const unsigned int, const unsigned int)> &operation_before_loop, const std::function< void(const unsigned int, const unsigned int)> &operation_after_loop, const unsigned int dof_handler_index_pre_post=0, const DataAccessOnFaces dst_vector_face_access=DataAccessOnFaces::unspecified, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified) const
 
template<typename CLASS , typename OutVector , typename InVector >
void loop_cell_centric (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const, const CLASS *owning_class, OutVector &dst, const InVector &src, const bool zero_dst_vector=false, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified) const
 
template<typename CLASS , typename OutVector , typename InVector >
void loop_cell_centric (void(CLASS::*cell_operation)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &), CLASS *owning_class, OutVector &dst, const InVector &src, const bool zero_dst_vector=false, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified) const
 
template<typename OutVector , typename InVector >
void loop_cell_centric (const std::function< void(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &cell_operation, OutVector &dst, const InVector &src, const bool zero_dst_vector=false, const DataAccessOnFaces src_vector_face_access=DataAccessOnFaces::unspecified) const
 
std::pair< unsigned int, unsigned intcreate_cell_subrange_hp (const std::pair< unsigned int, unsigned int > &range, const unsigned int fe_degree, const unsigned int dof_handler_index=0) const
 
std::pair< unsigned int, unsigned intcreate_cell_subrange_hp_by_index (const std::pair< unsigned int, unsigned int > &range, const unsigned int fe_index, const unsigned int dof_handler_index=0) const
 
unsigned int n_active_fe_indices () const
 
unsigned int get_cell_active_fe_index (const std::pair< unsigned int, unsigned int > range) const
 
unsigned int get_face_active_fe_index (const std::pair< unsigned int, unsigned int > range, const bool is_interior_face=true) const
 
Initialization of vectors
template<typename T >
void initialize_cell_data_vector (AlignedVector< T > &vec) const
 
template<typename T >
void initialize_face_data_vector (AlignedVector< T > &vec) const
 
template<typename VectorType >
void initialize_dof_vector (VectorType &vec, const unsigned int dof_handler_index=0) const
 
template<typename Number2 >
void initialize_dof_vector (LinearAlgebra::distributed::Vector< Number2 > &vec, const unsigned int dof_handler_index=0) const
 
const std::shared_ptr< const Utilities::MPI::Partitioner > & get_vector_partitioner (const unsigned int dof_handler_index=0) const
 
const IndexSetget_locally_owned_set (const unsigned int dof_handler_index=0) const
 
const IndexSetget_ghost_set (const unsigned int dof_handler_index=0) const
 
const std::vector< unsigned int > & get_constrained_dofs (const unsigned int dof_handler_index=0) const
 
void renumber_dofs (std::vector< types::global_dof_index > &renumbering, const unsigned int dof_handler_index=0)
 
Access of internal data structure

Note: Expert mode, interface not stable between releases.

const internal::MatrixFreeFunctions::TaskInfoget_task_info () const
 
const internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType > & get_mapping_info () const
 
const internal::MatrixFreeFunctions::DoFInfoget_dof_info (const unsigned int dof_handler_index_component=0) const
 
unsigned int n_constraint_pool_entries () const
 
const Number * constraint_pool_begin (const unsigned int pool_index) const
 
const Number * constraint_pool_end (const unsigned int pool_index) const
 
const internal::MatrixFreeFunctions::ShapeInfo< VectorizedArrayType > & get_shape_info (const unsigned int dof_handler_index_component=0, const unsigned int quad_index=0, const unsigned int fe_base_element=0, const unsigned int hp_active_fe_index=0, const unsigned int hp_active_quad_index=0) const
 
const internal::MatrixFreeFunctions::FaceToCellTopology< VectorizedArrayType::size()> & get_face_info (const unsigned int face_batch_index) const
 
const Table< 3, unsigned int > & get_cell_and_face_to_plain_faces () const
 
AlignedVector< VectorizedArrayType > * acquire_scratch_data () const
 
void release_scratch_data (const AlignedVector< VectorizedArrayType > *memory) const
 
AlignedVector< Number > * acquire_scratch_data_non_threadsafe () const
 
void release_scratch_data_non_threadsafe (const AlignedVector< Number > *memory) const
 
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)
 

Static Public Attributes

static constexpr unsigned int dimension = dim
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

template<typename number2 , int q_dim>
void internal_reinit (const std::shared_ptr< hp::MappingCollection< dim > > &mapping, const std::vector< const DoFHandler< dim, dim > * > &dof_handlers, const std::vector< const AffineConstraints< number2 > * > &constraint, const std::vector< IndexSet > &locally_owned_set, const std::vector< hp::QCollection< q_dim > > &quad, const AdditionalData &additional_data)
 
template<typename number2 >
void initialize_indices (const std::vector< const AffineConstraints< number2 > * > &constraint, const std::vector< IndexSet > &locally_owned_set, const AdditionalData &additional_data)
 
void initialize_dof_handlers (const std::vector< const DoFHandler< dim, dim > * > &dof_handlers, const AdditionalData &additional_data)
 
void check_no_subscribers () const noexcept
 

Private Attributes

std::vector< SmartPointer< const DoFHandler< dim > > > dof_handlers
 
std::vector< SmartPointer< const AffineConstraints< Number > > > affine_constraints
 
std::vector< internal::MatrixFreeFunctions::DoFInfodof_info
 
std::vector< Number > constraint_pool_data
 
std::vector< unsigned intconstraint_pool_row_index
 
internal::MatrixFreeFunctions::MappingInfo< dim, Number, VectorizedArrayType > mapping_info
 
Table< 4, internal::MatrixFreeFunctions::ShapeInfo< VectorizedArrayType > > shape_info
 
std::vector< std::pair< unsigned int, unsigned int > > cell_level_index
 
std::vector< unsigned intmf_cell_indices
 
unsigned int cell_level_index_end_local
 
internal::MatrixFreeFunctions::TaskInfo task_info
 
internal::MatrixFreeFunctions::FaceInfo< VectorizedArrayType::size()> face_info
 
bool indices_are_initialized
 
bool mapping_is_initialized
 
Threads::ThreadLocalStorage< std::list< std::pair< bool, AlignedVector< VectorizedArrayType > > > > scratch_pad
 
std::list< std::pair< bool, AlignedVector< Number > > > scratch_pad_non_threadsafe
 
unsigned int mg_level
 
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
 

General information

unsigned int n_components () const
 
unsigned int n_base_elements (const unsigned int dof_handler_index) const
 
unsigned int n_physical_cells () const
 
unsigned int n_cell_batches () const
 
unsigned int n_ghost_cell_batches () const
 
unsigned int n_inner_face_batches () const
 
unsigned int n_boundary_face_batches () const
 
unsigned int n_ghost_inner_face_batches () const
 
types::boundary_id get_boundary_id (const unsigned int face_batch_index) const
 
std::array< types::boundary_id, VectorizedArrayType::size()> get_faces_by_cells_boundary_id (const unsigned int cell_batch_index, const unsigned int face_number) const
 
const DoFHandler< dim > & get_dof_handler (const unsigned int dof_handler_index=0) const
 
const AffineConstraints< Number > & get_affine_constraints (const unsigned int dof_handler_index=0) const
 
DoFHandler< dim >::cell_iterator get_cell_iterator (const unsigned int cell_batch_index, const unsigned int lane_index, const unsigned int dof_handler_index=0) const
 
std::pair< int, intget_cell_level_and_index (const unsigned int cell_batch_index, const unsigned int lane_index) const
 
unsigned int get_matrix_free_cell_index (const typename Triangulation< dim >::cell_iterator &cell) const
 
std::pair< typename DoFHandler< dim >::cell_iterator, unsigned intget_face_iterator (const unsigned int face_batch_index, const unsigned int lane_index, const bool interior=true, const unsigned int fe_component=0) const
 
bool at_irregular_cell (const unsigned int cell_batch_index) const
 
unsigned int n_active_entries_per_cell_batch (const unsigned int cell_batch_index) const
 
unsigned int n_active_entries_per_face_batch (const unsigned int face_batch_index) const
 
unsigned int get_dofs_per_cell (const unsigned int dof_handler_index=0, const unsigned int hp_active_fe_index=0) const
 
unsigned int get_n_q_points (const unsigned int quad_index=0, const unsigned int hp_active_fe_index=0) const
 
unsigned int get_dofs_per_face (const unsigned int dof_handler_index=0, const unsigned int hp_active_fe_index=0) const
 
unsigned int get_n_q_points_face (const unsigned int quad_index=0, const unsigned int hp_active_fe_index=0) const
 
const Quadrature< dim > & get_quadrature (const unsigned int quad_index=0, const unsigned int hp_active_fe_index=0) const
 
const Quadrature< dim - 1 > & get_face_quadrature (const unsigned int quad_index=0, const unsigned int hp_active_fe_index=0) const
 
unsigned int get_cell_range_category (const std::pair< unsigned int, unsigned int > cell_batch_range) const
 
std::pair< unsigned int, unsigned intget_face_range_category (const std::pair< unsigned int, unsigned int > face_batch_range) const
 
unsigned int get_cell_category (const unsigned int cell_batch_index) const
 
std::pair< unsigned int, unsigned intget_face_category (const unsigned int face_batch_index) const
 
bool indices_initialized () const
 
bool mapping_initialized () const
 
unsigned int get_mg_level () const
 
std::size_t memory_consumption () const
 
template<typename StreamType >
void print_memory_consumption (StreamType &out) const
 
void print (std::ostream &out) const
 
template<int spacedim>
static bool is_supported (const FiniteElement< dim, spacedim > &fe)
 

Detailed Description

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
class MatrixFree< dim, Number, VectorizedArrayType >

This class collects all the data that is stored for the matrix free implementation. The storage scheme is tailored towards several loops performed with the same data, i.e., typically doing many matrix-vector products or residual computations on the same mesh. The class is used in step-37 and step-48.

This class does not implement any operations involving finite element basis functions, i.e., regarding the operation performed on the cells. For these operations, the class FEEvaluation is designed to use the data collected in this class.

The stored data can be subdivided into three main components:

Besides the initialization routines, this class implements only a single operation, namely a loop over all cells (cell_loop()). This loop is scheduled in such a way that cells that share degrees of freedom are not worked on simultaneously, which implies that it is possible to write to vectors (or matrices) in parallel without having to explicitly synchronize access to these vectors and matrices. This class does not implement any shape values, all it does is to cache the respective data. To implement finite element operations, use the class FEEvaluation (or some of the related classes).

This class traverses the cells in a different order than the usual Triangulation class in deal.II, in order to be flexible with respect to parallelization in shared memory and vectorization.

Vectorization is implemented by merging several topological cells into one so-called cell batch. This enables the application of all cell-related operations for several cells with one CPU instruction and is one of the main features of this framework.

For details on usage of this class, see the description of FEEvaluation or the matrix-free module.

Definition at line 112 of file matrix_free.h.

Member Typedef Documentation

◆ value_type

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
using MatrixFree< dim, Number, VectorizedArrayType >::value_type = Number

An alias for the underlying number type specified by the template argument.

Definition at line 123 of file matrix_free.h.

◆ vectorized_value_type

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
using MatrixFree< dim, Number, VectorizedArrayType >::vectorized_value_type = VectorizedArrayType

Definition at line 124 of file matrix_free.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

◆ DataAccessOnFaces

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
enum class MatrixFree::DataAccessOnFaces
strong

This class defines the type of data access for face integrals in loop () that is passed on to the update_ghost_values and compress functions of the parallel vectors, with the purpose of being able to reduce the amount of data that must be exchanged. The data exchange is a real bottleneck in particular for high-degree DG methods, therefore a more restrictive way of exchange is clearly beneficial. Note that this selection applies to FEFaceEvaluation objects assigned to the exterior side of cells accessing FaceToCellTopology::exterior_cells only; all interior objects are available in any case.

Enumerator
none 

The loop does not involve any FEFaceEvaluation access into neighbors, as is the case with only boundary integrals (but no interior face integrals) or when doing mass matrices in a MatrixFree::cell_loop() like setup.

values 

The loop does only involve FEFaceEvaluation access into neighbors by function values, such as FEFaceEvaluation::gather_evaluate() with argument EvaluationFlags::values, but no access to shape function derivatives (which typically need to access more data). For FiniteElement types where only some of the shape functions have support on a face, such as an FE_DGQ element with Lagrange polynomials with nodes on the element surface, the data exchange is reduced from (k+1)^dim to (k+1)^(dim-1).

values_all_faces 

Same as above. To be used if data has to be accessed from exterior faces if FEFaceEvaluation was reinitialized by providing the cell batch number and a face number. This configuration is useful in the context of cell-centric loops.

Precondition
AdditionalData::hold_all_faces_to_owned_cells has to enabled.
gradients 

The loop does involve FEFaceEvaluation access into neighbors by function values and gradients, but no second derivatives, such as FEFaceEvaluation::gather_evaluate() with EvaluationFlags::values and EvaluationFlags::gradients set. For FiniteElement types where only some of the shape functions have non-zero value and first derivative on a face, such as an FE_DGQHermite element, the data exchange is reduced, e.g. from (k+1)^dim to 2(k+1)^(dim-1). Note that for bases that do not have this special property, the full neighboring data is sent anyway.

gradients_all_faces 

Same as above. To be used if data has to be accessed from exterior faces if FEFaceEvaluation was reinitialized by providing the cell batch number and a face number. This configuration is useful in the context of cell-centric loops.

Precondition
AdditionalData::hold_all_faces_to_owned_cells has to enabled.
unspecified 

General setup where the user does not want to make a restriction. This is typically more expensive than the other options, but also the most conservative one because the full data of elements behind the faces to be computed locally will be exchanged.

Definition at line 686 of file matrix_free.h.

Constructor & Destructor Documentation

◆ MatrixFree() [1/2]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
MatrixFree< dim, Number, VectorizedArrayType >::MatrixFree ( )

Default empty constructor. Does nothing.

◆ MatrixFree() [2/2]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
MatrixFree< dim, Number, VectorizedArrayType >::MatrixFree ( const MatrixFree< dim, Number, VectorizedArrayType > &  other)

Copy constructor, calls copy_from

◆ ~MatrixFree()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
MatrixFree< dim, Number, VectorizedArrayType >::~MatrixFree ( )
overridedefault

Destructor.

Member Function Documentation

◆ reinit() [1/3]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename QuadratureType , typename number2 , typename MappingType >
void MatrixFree< dim, Number, VectorizedArrayType >::reinit ( const MappingType &  mapping,
const DoFHandler< dim > &  dof_handler,
const AffineConstraints< number2 > &  constraint,
const QuadratureType &  quad,
const AdditionalData additional_data = AdditionalData() 
)

Extracts the information needed to perform loops over cells. The DoFHandler and AffineConstraints objects describe the layout of degrees of freedom, the DoFHandler and the mapping describe the transformations from unit to real cell, and the finite element underlying the DoFHandler together with the quadrature formula describe the local operations. Note that the finite element underlying the DoFHandler must either be scalar or contain several copies of the same element. Mixing several different elements into one FESystem is not allowed. In that case, use the initialization function with several DoFHandler arguments.

◆ reinit() [2/3]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename QuadratureType , typename number2 , typename MappingType >
void MatrixFree< dim, Number, VectorizedArrayType >::reinit ( const MappingType &  mapping,
const std::vector< const DoFHandler< dim > * > &  dof_handler,
const std::vector< const AffineConstraints< number2 > * > &  constraint,
const std::vector< QuadratureType > &  quad,
const AdditionalData additional_data = AdditionalData() 
)

Extracts the information needed to perform loops over cells. The DoFHandler and AffineConstraints objects describe the layout of degrees of freedom, the DoFHandler and the mapping describe the transformations from unit to real cell, and the finite element underlying the DoFHandler together with the quadrature formula describe the local operations. As opposed to the scalar case treated with the other initialization functions, this function allows for problems with two or more different finite elements. The DoFHandlers to each element must be passed as pointers to the initialization function. Alternatively, a system of several components may also be represented by a single DoFHandler with an FESystem element. The prerequisite for this case is that each base element of the FESystem must be compatible with the present class, such as the FE_Q or FE_DGQ classes.

This function also allows for using several quadrature formulas, e.g. when the description contains independent integrations of elements of different degrees. However, the number of different quadrature formulas can be sets independently from the number of DoFHandlers, when several elements are always integrated with the same quadrature formula.

◆ reinit() [3/3]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename QuadratureType , typename number2 , typename MappingType >
void MatrixFree< dim, Number, VectorizedArrayType >::reinit ( const MappingType &  mapping,
const std::vector< const DoFHandler< dim > * > &  dof_handler,
const std::vector< const AffineConstraints< number2 > * > &  constraint,
const QuadratureType &  quad,
const AdditionalData additional_data = AdditionalData() 
)

Initializes the data structures. Same as before, but now the index set description of the locally owned range of degrees of freedom is taken from the DoFHandler. Moreover, only a single quadrature formula is used, as might be necessary when several components in a vector-valued problem are integrated together based on the same quadrature formula.

◆ copy_from()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
void MatrixFree< dim, Number, VectorizedArrayType >::copy_from ( const MatrixFree< dim, Number, VectorizedArrayType > &  matrix_free_base)

Copy function. Creates a deep copy of all data structures. It is usually enough to keep the data for different operations once, so this function should not be needed very often.

◆ update_mapping() [1/2]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
void MatrixFree< dim, Number, VectorizedArrayType >::update_mapping ( const Mapping< dim > &  mapping)

Refreshes the geometry data stored in the MappingInfo fields when the underlying geometry has changed (e.g. by a mapping that can deform through a change in the spatial configuration like MappingFEField) whereas the topology of the mesh and unknowns have remained the same. Compared to reinit(), this operation only has to re-generate the geometry arrays and can thus be significantly cheaper (depending on the cost to evaluate the geometry).

◆ update_mapping() [2/2]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
void MatrixFree< dim, Number, VectorizedArrayType >::update_mapping ( const std::shared_ptr< hp::MappingCollection< dim > > &  mapping)

Same as above but with hp::MappingCollection.

◆ clear()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
void MatrixFree< dim, Number, VectorizedArrayType >::clear ( )

Clear all data fields and brings the class into a condition similar to after having called the default constructor.

◆ cell_loop() [1/6]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::cell_loop ( const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &  cell_operation,
OutVector &  dst,
const InVector &  src,
const bool  zero_dst_vector = false 
) const

This method runs the loop over all cells (in parallel) and performs the MPI data exchange on the source vector and destination vector.

Parameters
cell_operationstd::function with the signature cell_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) where the first argument passes the data of the calling class and the last argument defines the range of cells which should be worked on (typically more than one cell should be worked on in order to reduce overheads). One can pass a pointer to an object in this place if it has an operator() with the correct set of arguments since such a pointer can be converted to the function object.
dstDestination vector holding the result. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::compress() at the end of the call internally. For other vectors, including parallel Trilinos or PETSc vectors, no such call is issued. Note that Trilinos/Epetra or PETSc vectors do currently not work in parallel because the present class uses MPI-local index addressing, as opposed to the global addressing implied by those external libraries.
srcInput vector. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::update_ghost_values() at the start of the call internally to make sure all necessary data is locally available. Note, however, that the vector is reset to its original state at the end of the loop, i.e., if the vector was not ghosted upon entry of the loop, it will not be ghosted upon finishing the loop.
zero_dst_vectorIf this flag is set to true, the vector dst will be set to zero inside the loop. Use this case in case you perform a typical vmult() operation on a matrix object, as it will typically be faster than calling dst = 0; before the loop separately. This is because the vector entries are set to zero only on subranges of the vector, making sure that the vector entries stay in caches as much as possible.

◆ cell_loop() [2/6]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename CLASS , typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::cell_loop ( void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const  cell_operation,
const CLASS *  owning_class,
OutVector &  dst,
const InVector &  src,
const bool  zero_dst_vector = false 
) const

This is the second variant to run the loop over all cells, now providing a function pointer to a member function of class CLASS. This method obviates the need to define a lambda function or to call std::bind to bind the class into the given function in case the local function needs to access data in the class (i.e., it is a non-static member function).

Parameters
cell_operationPointer to member function of CLASS with the signature cell_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) where the first argument passes the data of the calling class and the last argument defines the range of cells which should be worked on (typically more than one cell should be worked on in order to reduce overheads).
owning_classThe object which provides the cell_operation call. To be compatible with this interface, the class must allow to call owning_class->cell_operation(...).
dstDestination vector holding the result. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::compress() at the end of the call internally. For other vectors, including parallel Trilinos or PETSc vectors, no such call is issued. Note that Trilinos/Epetra or PETSc vectors do currently not work in parallel because the present class uses MPI-local index addressing, as opposed to the global addressing implied by those external libraries.
srcInput vector. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::update_ghost_values() at the start of the call internally to make sure all necessary data is locally available. Note, however, that the vector is reset to its original state at the end of the loop, i.e., if the vector was not ghosted upon entry of the loop, it will not be ghosted upon finishing the loop.
zero_dst_vectorIf this flag is set to true, the vector dst will be set to zero inside the loop. Use this case in case you perform a typical vmult() operation on a matrix object, as it will typically be faster than calling dst = 0; before the loop separately. This is because the vector entries are set to zero only on subranges of the vector, making sure that the vector entries stay in caches as much as possible.

◆ cell_loop() [3/6]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename CLASS , typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::cell_loop ( void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)  cell_operation,
CLASS *  owning_class,
OutVector &  dst,
const InVector &  src,
const bool  zero_dst_vector = false 
) const

Same as above, but for class member functions which are non-const.

◆ cell_loop() [4/6]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename CLASS , typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::cell_loop ( void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const  cell_operation,
const CLASS *  owning_class,
OutVector &  dst,
const InVector &  src,
const std::function< void(const unsigned int, const unsigned int)> &  operation_before_loop,
const std::function< void(const unsigned int, const unsigned int)> &  operation_after_loop,
const unsigned int  dof_handler_index_pre_post = 0 
) const

This function is similar to the cell_loop with an std::function object to specify to operation to be performed on cells, but adds two additional functors to execute some additional work before and after the cell integrals are computed.

The two additional functors work on a range of degrees of freedom, expressed in terms of the degree-of-freedom numbering of the selected DoFHandler dof_handler_index_pre_post in MPI-local indices. The arguments to the functors represent a range of degrees of freedom at a granularity of internal::MatrixFreeFunctions::DoFInfo::chunk_size_zero_vector entries (except for the last chunk which is set to the number of locally owned entries) in the form [first, last). The idea of these functors is to bring operations on vectors closer to the point where they accessed in a matrix-free loop, with the goal to increase cache hits by temporal locality. This loop guarantees that the operation_before_loop hits all relevant unknowns before they are first touched in the cell_operation (including the MPI data exchange), allowing to execute some vector update that the src vector depends upon. The operation_after_loop is similar

  • it starts to execute on a range of DoFs once all DoFs in that range have been touched for the last time by the cell_operation (including the MPI data exchange), allowing e.g. to compute some vector operations that depend on the result of the current cell loop in dst or want to modify src. The efficiency of caching depends on the numbering of the degrees of freedom because of the granularity of the ranges.
Parameters
cell_operationPointer to member function of CLASS with the signature cell_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) where the first argument passes the data of the calling class and the last argument defines the range of cells which should be worked on (typically more than one cell should be worked on in order to reduce overheads).
owning_classThe object which provides the cell_operation call. To be compatible with this interface, the class must allow to call owning_class->cell_operation(...).
dstDestination vector holding the result. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::compress() at the end of the call internally. For other vectors, including parallel Trilinos or PETSc vectors, no such call is issued. Note that Trilinos/Epetra or PETSc vectors do currently not work in parallel because the present class uses MPI-local index addressing, as opposed to the global addressing implied by those external libraries.
srcInput vector. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::update_ghost_values() at the start of the call internally to make sure all necessary data is locally available. Note, however, that the vector is reset to its original state at the end of the loop, i.e., if the vector was not ghosted upon entry of the loop, it will not be ghosted upon finishing the loop.
operation_before_loopThis functor can be used to perform an operation on entries of the src and dst vectors (or other vectors) before the operation on cells first touches a particular DoF according to the general description in the text above. This function is passed a range of the locally owned degrees of freedom on the selected dof_handler_index_pre_post (in MPI-local numbering).
operation_after_loopThis functor can be used to perform an operation on entries of the src and dst vectors (or other vectors) after the operation on cells last touches a particular DoF according to the general description in the text above. This function is passed a range of the locally owned degrees of freedom on the selected dof_handler_index_pre_post (in MPI-local numbering).
dof_handler_index_pre_postSince MatrixFree can be initialized with a vector of DoFHandler objects, each of them will in general have different vector sizes and thus different ranges returned to operation_before_loop and operation_after_loop. Use this variable to specify which one of the DoFHandler objects the index range should be associated to. Defaults to the dof_handler_index 0.
Note
The close locality of the operation_before_loop and operation_after_loop is currently only implemented for the MPI-only case. In case threading is enabled, the complete operation_before_loop is scheduled before the parallel loop, and operation_after_loop is scheduled strictly afterwards, due to the complicated dependencies.

◆ cell_loop() [5/6]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename CLASS , typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::cell_loop ( void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)  cell_operation,
CLASS *  owning_class,
OutVector &  dst,
const InVector &  src,
const std::function< void(const unsigned int, const unsigned int)> &  operation_before_loop,
const std::function< void(const unsigned int, const unsigned int)> &  operation_after_loop,
const unsigned int  dof_handler_index_pre_post = 0 
) const

Same as above, but for class member functions which are non-const.

◆ cell_loop() [6/6]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::cell_loop ( const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &  cell_operation,
OutVector &  dst,
const InVector &  src,
const std::function< void(const unsigned int, const unsigned int)> &  operation_before_loop,
const std::function< void(const unsigned int, const unsigned int)> &  operation_after_loop,
const unsigned int  dof_handler_index_pre_post = 0 
) const

Same as above, but taking an std::function as the cell_operation rather than a class member function.

◆ loop() [1/6]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::loop ( const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &  cell_operation,
const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &  face_operation,
const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &  boundary_operation,
OutVector &  dst,
const InVector &  src,
const bool  zero_dst_vector = false,
const DataAccessOnFaces  dst_vector_face_access = DataAccessOnFaces::unspecified,
const DataAccessOnFaces  src_vector_face_access = DataAccessOnFaces::unspecified 
) const

This method runs a loop over all cells (in parallel) and performs the MPI data exchange on the source vector and destination vector. As opposed to the other variants that only runs a function on cells, this method also takes as arguments a function for the interior faces and for the boundary faces, respectively.

Parameters
cell_operationstd::function with the signature cell_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) where the first argument passes the data of the calling class and the last argument defines the range of cells which should be worked on (typically more than one cell should be worked on in order to reduce overheads). One can pass a pointer to an object in this place if it has an operator() with the correct set of arguments since such a pointer can be converted to the function object.
face_operationstd::function with the signature face_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) in analogy to cell_operation, but now the part associated to the work on interior faces. Note that the MatrixFree framework treats periodic faces as interior ones, so they will be assigned their correct neighbor after applying periodicity constraints within the face_operation calls.
boundary_operationstd::function with the signature boundary_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) in analogy to cell_operation and face_operation, but now the part associated to the work on boundary faces. Boundary faces are separated by their boundary_id and it is possible to query that id using MatrixFree::get_boundary_id(). Note that both interior and faces use the same numbering, and faces in the interior are assigned lower numbers than the boundary faces.
dstDestination vector holding the result. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::compress() at the end of the call internally.
srcInput vector. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::update_ghost_values() at the start of the call internally to make sure all necessary data is locally available. Note, however, that the vector is reset to its original state at the end of the loop, i.e., if the vector was not ghosted upon entry of the loop, it will not be ghosted upon finishing the loop.
zero_dst_vectorIf this flag is set to true, the vector dst will be set to zero inside the loop. Use this case in case you perform a typical vmult() operation on a matrix object, as it will typically be faster than calling dst = 0; before the loop separately. This is because the vector entries are set to zero only on subranges of the vector, making sure that the vector entries stay in caches as much as possible.
dst_vector_face_accessSet the type of access into the vector dst that will happen inside the body of the face_operation function. As explained in the description of the DataAccessOnFaces struct, the purpose of this selection is to reduce the amount of data that must be exchanged over the MPI network (or via memcpy if within the shared memory region of a node) to gain performance. Note that there is no way to communicate this setting with the FEFaceEvaluation class, therefore this selection must be made at this site in addition to what is implemented inside the face_operation function. As a consequence, there is also no way to check that the setting passed to this call is consistent with what is later done by FEFaceEvaluation, and it is the user's responsibility to ensure correctness of data.
src_vector_face_accessSet the type of access into the vector src that will happen inside the body of the face_operation function, in analogy to dst_vector_face_access.

◆ loop() [2/6]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename CLASS , typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::loop ( void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const  cell_operation,
void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const  face_operation,
void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const  boundary_operation,
const CLASS *  owning_class,
OutVector &  dst,
const InVector &  src,
const bool  zero_dst_vector = false,
const DataAccessOnFaces  dst_vector_face_access = DataAccessOnFaces::unspecified,
const DataAccessOnFaces  src_vector_face_access = DataAccessOnFaces::unspecified 
) const

This is the second variant to run the loop over all cells, interior faces, and boundary faces, now providing three function pointers to member functions of class CLASS with the signature operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int>&)const. This method obviates the need to define a lambda function or to call std::bind to bind the class into the given function in case the local function needs to access data in the class (i.e., it is a non-static member function).

Parameters
cell_operationPointer to member function of CLASS with the signature cell_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) where the first argument passes the data of the calling class and the last argument defines the range of cells which should be worked on (typically more than one cell should be worked on in order to reduce overheads). Note that the loop will typically split the cell_range into smaller pieces and work on cell_operation, face_operation, and boundary_operation alternately, in order to increase the potential reuse of vector entries in caches.
face_operationPointer to member function of CLASS with the signature face_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) in analogy to cell_operation, but now the part associated to the work on interior faces. Note that the MatrixFree framework treats periodic faces as interior ones, so they will be assigned their correct neighbor after applying periodicity constraints within the face_operation calls.
boundary_operationPointer to member function of CLASS with the signature boundary_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) in analogy to cell_operation and face_operation, but now the part associated to the work on boundary faces. Boundary faces are separated by their boundary_id and it is possible to query that id using MatrixFree::get_boundary_id(). Note that both interior and faces use the same numbering, and faces in the interior are assigned lower numbers than the boundary faces.
owning_classThe object which provides the cell_operation call. To be compatible with this interface, the class must allow to call owning_class->cell_operation(...), owning_class->face_operation(...), and owning_class->boundary_operation(...).
dstDestination vector holding the result. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::compress() at the end of the call internally.
srcInput vector. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::update_ghost_values() at the start of the call internally to make sure all necessary data is locally available. Note, however, that the vector is reset to its original state at the end of the loop, i.e., if the vector was not ghosted upon entry of the loop, it will not be ghosted upon finishing the loop.
zero_dst_vectorIf this flag is set to true, the vector dst will be set to zero inside the loop. Use this case in case you perform a typical vmult() operation on a matrix object, as it will typically be faster than calling dst = 0; before the loop separately. This is because the vector entries are set to zero only on subranges of the vector, making sure that the vector entries stay in caches as much as possible.
dst_vector_face_accessSet the type of access into the vector dst that will happen inside the body of the face_operation function. As explained in the description of the DataAccessOnFaces struct, the purpose of this selection is to reduce the amount of data that must be exchanged over the MPI network (or via memcpy if within the shared memory region of a node) to gain performance. Note that there is no way to communicate this setting with the FEFaceEvaluation class, therefore this selection must be made at this site in addition to what is implemented inside the face_operation function. As a consequence, there is also no way to check that the setting passed to this call is consistent with what is later done by FEFaceEvaluation, and it is the user's responsibility to ensure correctness of data.
src_vector_face_accessSet the type of access into the vector src that will happen inside the body of the face_operation function, in analogy to dst_vector_face_access.

◆ loop() [3/6]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename CLASS , typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::loop ( void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)  cell_operation,
void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)  face_operation,
void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)  boundary_operation,
CLASS *  owning_class,
OutVector &  dst,
const InVector &  src,
const bool  zero_dst_vector = false,
const DataAccessOnFaces  dst_vector_face_access = DataAccessOnFaces::unspecified,
const DataAccessOnFaces  src_vector_face_access = DataAccessOnFaces::unspecified 
) const

Same as above, but for class member functions which are non-const.

◆ loop() [4/6]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename CLASS , typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::loop ( void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const  cell_operation,
void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const  face_operation,
void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const  boundary_operation,
const CLASS *  owning_class,
OutVector &  dst,
const InVector &  src,
const std::function< void(const unsigned int, const unsigned int)> &  operation_before_loop,
const std::function< void(const unsigned int, const unsigned int)> &  operation_after_loop,
const unsigned int  dof_handler_index_pre_post = 0,
const DataAccessOnFaces  dst_vector_face_access = DataAccessOnFaces::unspecified,
const DataAccessOnFaces  src_vector_face_access = DataAccessOnFaces::unspecified 
) const

This function is similar to the loop method above, but adds two additional functors to execute some additional work before and after the cell, face and boundary integrals are computed.

The two additional functors work on a range of degrees of freedom, expressed in terms of the degree-of-freedom numbering of the selected DoFHandler dof_handler_index_pre_post in MPI-local indices. The arguments to the functors represent a range of degrees of freedom at a granularity of internal::MatrixFreeFunctions::DoFInfo::chunk_size_zero_vector entries (except for the last chunk which is set to the number of locally owned entries) in the form [first, last). The idea of these functors is to bring operations on vectors closer to the point where they accessed in a matrix-free loop, with the goal to increase cache hits by temporal locality. This loop guarantees that the operation_before_loop hits all relevant unknowns before they are first touched by any of the cell, face or boundary operations (including the MPI data exchange), allowing to execute some vector update that the src vector depends upon. The operation_after_loop is similar - it starts to execute on a range of DoFs once all DoFs in that range have been touched for the last time by the cell, face and boundary operations (including the MPI data exchange), allowing e.g. to compute some vector operations that depend on the result of the current cell loop in dst or want to modify src. The efficiency of caching depends on the numbering of the degrees of freedom because of the granularity of the ranges.

Parameters
cell_operationPointer to member function of CLASS with the signature cell_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) where the first argument passes the data of the calling class and the last argument defines the range of cells which should be worked on (typically more than one cell should be worked on in order to reduce overheads).
face_operationPointer to member function of CLASS with the signature face_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) in analogy to cell_operation, but now the part associated to the work on interior faces. Note that the MatrixFree framework treats periodic faces as interior ones, so they will be assigned their correct neighbor after applying periodicity constraints within the face_operation calls.
boundary_operationPointer to member function of CLASS with the signature boundary_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) in analogy to cell_operation and face_operation, but now the part associated to the work on boundary faces. Boundary faces are separated by their boundary_id and it is possible to query that id using MatrixFree::get_boundary_id(). Note that both interior and faces use the same numbering, and faces in the interior are assigned lower numbers than the boundary faces.
owning_classThe object which provides the cell_operation call. To be compatible with this interface, the class must allow to call owning_class->cell_operation(...).
dstDestination vector holding the result. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::compress() at the end of the call internally. For other vectors, including parallel Trilinos or PETSc vectors, no such call is issued. Note that Trilinos/Epetra or PETSc vectors do currently not work in parallel because the present class uses MPI-local index addressing, as opposed to the global addressing implied by those external libraries.
srcInput vector. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::update_ghost_values() at the start of the call internally to make sure all necessary data is locally available. Note, however, that the vector is reset to its original state at the end of the loop, i.e., if the vector was not ghosted upon entry of the loop, it will not be ghosted upon finishing the loop.
operation_before_loopThis functor can be used to perform an operation on entries of the src and dst vectors (or other vectors) before the operation on cells first touches a particular DoF according to the general description in the text above. This function is passed a range of the locally owned degrees of freedom on the selected dof_handler_index_pre_post (in MPI-local numbering).
operation_after_loopThis functor can be used to perform an operation on entries of the src and dst vectors (or other vectors) after the operation on cells last touches a particular DoF according to the general description in the text above. This function is passed a range of the locally owned degrees of freedom on the selected dof_handler_index_pre_post (in MPI-local numbering).
dof_handler_index_pre_postSince MatrixFree can be initialized with a vector of DoFHandler objects, each of them will in general have different vector sizes and thus different ranges returned to operation_before_loop and operation_after_loop. Use this variable to specify which one of the DoFHandler objects the index range should be associated to. Defaults to the dof_handler_index 0.
dst_vector_face_accessSet the type of access into the vector dst that will happen inside the body of the face_operation function. As explained in the description of the DataAccessOnFaces struct, the purpose of this selection is to reduce the amount of data that must be exchanged over the MPI network (or via memcpy if within the shared memory region of a node) to gain performance. Note that there is no way to communicate this setting with the FEFaceEvaluation class, therefore this selection must be made at this site in addition to what is implemented inside the face_operation function. As a consequence, there is also no way to check that the setting passed to this call is consistent with what is later done by FEFaceEvaluation, and it is the user's responsibility to ensure correctness of data.
src_vector_face_accessSet the type of access into the vector src that will happen inside the body of the face_operation function, in analogy to dst_vector_face_access.
Note
The close locality of the operation_before_loop and operation_after_loop is currently only implemented for the MPI-only case. In case threading is enabled, the complete operation_before_loop is scheduled before the parallel loop, and operation_after_loop is scheduled strictly afterwards, due to the complicated dependencies.

◆ loop() [5/6]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename CLASS , typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::loop ( void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)  cell_operation,
void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)  face_operation,
void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)  boundary_operation,
const CLASS *  owning_class,
OutVector &  dst,
const InVector &  src,
const std::function< void(const unsigned int, const unsigned int)> &  operation_before_loop,
const std::function< void(const unsigned int, const unsigned int)> &  operation_after_loop,
const unsigned int  dof_handler_index_pre_post = 0,
const DataAccessOnFaces  dst_vector_face_access = DataAccessOnFaces::unspecified,
const DataAccessOnFaces  src_vector_face_access = DataAccessOnFaces::unspecified 
) const

Same as above, but for class member functions which are non-const.

◆ loop() [6/6]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::loop ( const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &  cell_operation,
const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &  face_operation,
const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &  boundary_operation,
OutVector &  dst,
const InVector &  src,
const std::function< void(const unsigned int, const unsigned int)> &  operation_before_loop,
const std::function< void(const unsigned int, const unsigned int)> &  operation_after_loop,
const unsigned int  dof_handler_index_pre_post = 0,
const DataAccessOnFaces  dst_vector_face_access = DataAccessOnFaces::unspecified,
const DataAccessOnFaces  src_vector_face_access = DataAccessOnFaces::unspecified 
) const

Same as above, but taking an std::function as the cell_operation, face_operation and boundary_operation rather than a class member function.

◆ loop_cell_centric() [1/3]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename CLASS , typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::loop_cell_centric ( void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const  cell_operation,
const CLASS *  owning_class,
OutVector &  dst,
const InVector &  src,
const bool  zero_dst_vector = false,
const DataAccessOnFaces  src_vector_face_access = DataAccessOnFaces::unspecified 
) const

This method runs the loop over all cells (in parallel) similarly as cell_loop() does. However, this function is intended to be used for the case if face and boundary integrals should be also evaluated. In contrast to loop(), the user provides only a single function that should contain the cell integral over a cell (or batch of cells when vectorizing) and the face and boundary integrals over all its faces. This is referred to in the literature as element-centric loop or cell-centric loop.

To be able to evaluate all face integrals (with values or gradients from the neighboring cells), all ghost values from neighboring cells are updated. Use FEFaceEvaluation::reinit(cell, face_no) to access quantities on arbitrary faces of a cell and the respective neighbors.

Parameters
cell_operationPointer to member function of CLASS with the signature cell_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int> &) where the first argument passes the data of the calling class and the last argument defines the range of cells which should be worked on (typically more than one cell is passed in from the loop in order to reduce overheads).
owning_classThe object which provides the cell_operation call. To be compatible with this interface, the class must allow to call owning_class->cell_operation(...).
dstDestination vector holding the result. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::compress() at the end of the call internally.
srcInput vector. If the vector is of type LinearAlgebra::distributed::Vector (or composite objects thereof such as LinearAlgebra::distributed::BlockVector), the loop calls LinearAlgebra::distributed::Vector::update_ghost_values() at the start of the call internally to make sure all necessary data is locally available. Note, however, that the vector is reset to its original state at the end of the loop, i.e., if the vector was not ghosted upon entry of the loop, it will not be ghosted upon finishing the loop.
zero_dst_vectorIf this flag is set to true, the vector dst will be set to zero inside the loop. Use this case in case you perform a typical vmult() operation on a matrix object, as it will typically be faster than calling dst = 0; before the loop separately. This is because the vector entries are set to zero only on subranges of the vector, making sure that the vector entries stay in caches as much as possible.
src_vector_face_accessSet the type of access into the vector src that will happen inside the body of the cell_operation function during face integrals. As explained in the description of the DataAccessOnFaces struct, the purpose of this selection is to reduce the amount of data that must be exchanged over the MPI network (or via memcpy if within the shared memory region of a node) to gain performance. Note that there is no way to communicate this setting with the FEFaceEvaluation class, therefore this selection must be made at this site in addition to what is implemented inside the face_operation function. As a consequence, there is also no way to check that the setting passed to this call is consistent with what is later done by FEFaceEvaluation, and it is the user's responsibility to ensure correctness of data.

◆ loop_cell_centric() [2/3]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename CLASS , typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::loop_cell_centric ( void(CLASS::*)(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)  cell_operation,
CLASS *  owning_class,
OutVector &  dst,
const InVector &  src,
const bool  zero_dst_vector = false,
const DataAccessOnFaces  src_vector_face_access = DataAccessOnFaces::unspecified 
) const

Same as above, but for the class member function which is non-const.

◆ loop_cell_centric() [3/3]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename OutVector , typename InVector >
void MatrixFree< dim, Number, VectorizedArrayType >::loop_cell_centric ( const std::function< void(const MatrixFree< dim, Number, VectorizedArrayType > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &  cell_operation,
OutVector &  dst,
const InVector &  src,
const bool  zero_dst_vector = false,
const DataAccessOnFaces  src_vector_face_access = DataAccessOnFaces::unspecified 
) const

Same as above, but with std::function.

◆ create_cell_subrange_hp()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::pair< unsigned int, unsigned int > MatrixFree< dim, Number, VectorizedArrayType >::create_cell_subrange_hp ( const std::pair< unsigned int, unsigned int > &  range,
const unsigned int  fe_degree,
const unsigned int  dof_handler_index = 0 
) const

In the hp-adaptive case, a subrange of cells as computed during the cell loop might contain elements of different degrees. Use this function to compute what the subrange for an individual finite element degree is. The finite element degree is associated to the vector component given in the function call.

◆ create_cell_subrange_hp_by_index()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::pair< unsigned int, unsigned int > MatrixFree< dim, Number, VectorizedArrayType >::create_cell_subrange_hp_by_index ( const std::pair< unsigned int, unsigned int > &  range,
const unsigned int  fe_index,
const unsigned int  dof_handler_index = 0 
) const

In the hp-adaptive case, a subrange of cells as computed during the cell loop might contain elements of different degrees. Use this function to compute what the subrange for a given index the hp-finite element, as opposed to the finite element degree in the other function.

◆ n_active_fe_indices()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::n_active_fe_indices ( ) const

In the hp-adaptive case, return number of active FE indices.

◆ get_cell_active_fe_index()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::get_cell_active_fe_index ( const std::pair< unsigned int, unsigned int range) const

In the hp-adaptive case, return the active FE index of a cell range.

◆ get_face_active_fe_index()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::get_face_active_fe_index ( const std::pair< unsigned int, unsigned int range,
const bool  is_interior_face = true 
) const

In the hp-adaptive case, return the active FE index of a face range.

◆ initialize_cell_data_vector()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename T >
void MatrixFree< dim, Number, VectorizedArrayType >::initialize_cell_data_vector ( AlignedVector< T > &  vec) const

Initialize function for a vector with each entry associated with a cell batch (cell data). For reading and writing the vector use: FEEvaluationBase::read_cell_data() and FEEvaluationBase::set_cell_data().

◆ initialize_face_data_vector()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename T >
void MatrixFree< dim, Number, VectorizedArrayType >::initialize_face_data_vector ( AlignedVector< T > &  vec) const

Initialize function for a vector with each entry associated with a face batch (face data). For reading and writing the vector use: FEEvaluationBase::read_face_data() and FEEvaluationBase::set_face_data().

◆ initialize_dof_vector() [1/2]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename VectorType >
void MatrixFree< dim, Number, VectorizedArrayType >::initialize_dof_vector ( VectorType &  vec,
const unsigned int  dof_handler_index = 0 
) const

Initialize function for a general serial non-block vector. After a call to this function, the length of the vector is equal to the total number of degrees of freedom in the DoFHandler. Vector entries are initialized with zero.

If MatrixFree was set up with several DoFHandler objects, the parameter dof_handler_index defines which component is to be used.

Note
Serial vectors also include Trilinos and PETSc vectors; however in these cases, MatrixFree has to be used in a serial context, i.e., the size of the communicator has to be exactly one.

◆ initialize_dof_vector() [2/2]

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename Number2 >
void MatrixFree< dim, Number, VectorizedArrayType >::initialize_dof_vector ( LinearAlgebra::distributed::Vector< Number2 > &  vec,
const unsigned int  dof_handler_index = 0 
) const

Specialization of the method initialize_dof_vector() for the class LinearAlgebra::distributed::Vector<Number>. See the other function with the same name for the general descriptions.

Note
For the parallel vectors used with MatrixFree and in FEEvaluation, a vector needs to hold all locally active DoFs and also some of the locally relevant DoFs. The selection of DoFs is such that one can read all degrees of freedom on all locally relevant elements (locally active) plus the degrees of freedom that constraints expand into from the locally owned cells. However, not all locally relevant DoFs are stored because most of them would never be accessed in matrix-vector products and result in too much data sent around which impacts the performance.

◆ get_vector_partitioner()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const std::shared_ptr< const Utilities::MPI::Partitioner > & MatrixFree< dim, Number, VectorizedArrayType >::get_vector_partitioner ( const unsigned int  dof_handler_index = 0) const

Return the partitioner that represents the locally owned data and the ghost indices where access is needed to for the cell loop. The partitioner is constructed from the locally owned dofs and ghost dofs given by the respective fields. If you want to have specific information about these objects, you can query them with the respective access functions. If you just want to initialize a (parallel) vector, you should usually prefer this data structure as the data exchange information can be reused from one vector to another.

◆ get_locally_owned_set()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const IndexSet & MatrixFree< dim, Number, VectorizedArrayType >::get_locally_owned_set ( const unsigned int  dof_handler_index = 0) const

Return the set of cells that are owned by the processor.

◆ get_ghost_set()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const IndexSet & MatrixFree< dim, Number, VectorizedArrayType >::get_ghost_set ( const unsigned int  dof_handler_index = 0) const

Return the set of ghost cells needed but not owned by the processor.

◆ get_constrained_dofs()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const std::vector< unsigned int > & MatrixFree< dim, Number, VectorizedArrayType >::get_constrained_dofs ( const unsigned int  dof_handler_index = 0) const

Return a list of all degrees of freedom that are constrained. The list is returned in MPI-local index space for the locally owned range of the vector, not in global MPI index space that spans all MPI processors. To get numbers in global index space, call get_vector_partitioner()->local_to_global on an entry of the vector. In addition, it only returns the indices for degrees of freedom that are owned locally, not for ghosts.

◆ renumber_dofs()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
void MatrixFree< dim, Number, VectorizedArrayType >::renumber_dofs ( std::vector< types::global_dof_index > &  renumbering,
const unsigned int  dof_handler_index = 0 
)

Computes a renumbering of degrees of freedom that better fits with the data layout in MatrixFree according to the given layout of data. Note that this function does not re-arrange the information stored in this class, but rather creates a renumbering for consumption of DoFHandler::renumber_dofs. To have any effect a MatrixFree object must be set up again using the renumbered DoFHandler and AffineConstraints. Note that if a DoFHandler calls DoFHandler::renumber_dofs, all information in MatrixFree becomes invalid.

◆ is_supported()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<int spacedim>
static bool MatrixFree< dim, Number, VectorizedArrayType >::is_supported ( const FiniteElement< dim, spacedim > &  fe)
static

Return whether a given FiniteElement fe is supported by this class.

◆ n_components()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::n_components ( ) const

Return the number of different DoFHandlers specified at initialization.

◆ n_base_elements()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::n_base_elements ( const unsigned int  dof_handler_index) const

For the finite element underlying the DoFHandler specified by dof_handler_index, return the number of base elements.

◆ n_physical_cells()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::n_physical_cells ( ) const

Return the number of cells this structure is based on. If you are using a usual DoFHandler, it corresponds to the number of (locally owned) active cells. Note that most data structures in this class do not directly act on this number but rather on n_cell_batches() which gives the number of cells as seen when lumping several cells together with vectorization.

◆ n_cell_batches()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::n_cell_batches ( ) const

Return the number of cell batches that this structure works on. The batches are formed by application of vectorization over several cells in general. The cell range in cell_loop runs from zero to n_cell_batches() (exclusive), so this is the appropriate size if you want to store arrays of data for all cells to be worked on. This number is approximately n_physical_cells()/VectorizedArray::size() (depending on how many cell batches that do not get filled up completely).

◆ n_ghost_cell_batches()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::n_ghost_cell_batches ( ) const

Return the number of additional cell batches that this structure keeps for face integration. Note that not all cells that are ghosted in the triangulation are kept in this data structure, but only the ones which are necessary for evaluating face integrals from both sides.

◆ n_inner_face_batches()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::n_inner_face_batches ( ) const

Return the number of interior face batches that this structure works on. The batches are formed by application of vectorization over several faces in general. The face range in loop runs from zero to n_inner_face_batches() (exclusive), so this is the appropriate size if you want to store arrays of data for all interior faces to be worked on. Note that it returns 0 unless mapping_update_flags_inner_faces is set to a value different from UpdateFlags::update_default.

◆ n_boundary_face_batches()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::n_boundary_face_batches ( ) const

Return the number of boundary face batches that this structure works on. The batches are formed by application of vectorization over several faces in general. The face range in loop runs from n_inner_face_batches() to n_inner_face_batches()+n_boundary_face_batches() (exclusive), so if you need to store arrays that hold data for all boundary faces but not the interior ones, this number gives the appropriate size. Note that it returns 0 unless mapping_update_flags_boundary_faces is set to a value different from UpdateFlags::update_default.

◆ n_ghost_inner_face_batches()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::n_ghost_inner_face_batches ( ) const

Return the number of faces that are not processed locally but belong to locally owned faces.

◆ get_boundary_id()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
types::boundary_id MatrixFree< dim, Number, VectorizedArrayType >::get_boundary_id ( const unsigned int  face_batch_index) const

In order to apply different operators to different parts of the boundary, this method can be used to query the boundary id of a given face in the faces' own sorting by lanes in a VectorizedArray. Only valid for an index indicating a boundary face.

Note
Alternatively to this function, you can use FEFaceEvaluation::boundary_id() to get the same information if a FEFaceEvaluation object has been set up already.

◆ get_faces_by_cells_boundary_id()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::array< types::boundary_id, VectorizedArrayType::size()> MatrixFree< dim, Number, VectorizedArrayType >::get_faces_by_cells_boundary_id ( const unsigned int  cell_batch_index,
const unsigned int  face_number 
) const

Return the boundary ids for the faces within a cell, using the cells' sorting by lanes in the VectorizedArray.

◆ get_dof_handler()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const DoFHandler< dim > & MatrixFree< dim, Number, VectorizedArrayType >::get_dof_handler ( const unsigned int  dof_handler_index = 0) const

Return the DoFHandler with the index as given to the respective std::vector argument in the reinit() function.

◆ get_affine_constraints()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const AffineConstraints< Number > & MatrixFree< dim, Number, VectorizedArrayType >::get_affine_constraints ( const unsigned int  dof_handler_index = 0) const

Return the AffineConstraints with the index as given to the respective std::vector argument in the reinit() function. Only available if the AffineConstraints objects have the same template parameter Number as MatrixFree. Throws an exception otherwise.

◆ get_cell_iterator()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
DoFHandler< dim >::cell_iterator MatrixFree< dim, Number, VectorizedArrayType >::get_cell_iterator ( const unsigned int  cell_batch_index,
const unsigned int  lane_index,
const unsigned int  dof_handler_index = 0 
) const

Return the cell iterator in deal.II speak to a given cell batch (populating several lanes in a VectorizedArray) and the lane index within the vectorization across cells in the renumbering of this structure.

Note that the cell iterators in deal.II go through cells differently to what the cell loop of this class does. This is because several cells are processed together (vectorization across cells), and since cells with neighbors on different MPI processors need to be accessed at a certain time when accessing remote data and overlapping communication with computation.

◆ get_cell_level_and_index()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::pair< int, int > MatrixFree< dim, Number, VectorizedArrayType >::get_cell_level_and_index ( const unsigned int  cell_batch_index,
const unsigned int  lane_index 
) const

This returns the level and index for the cell that would be returned by get_cell_iterator() for the same arguments cell_batch_index and lane_index.

◆ get_matrix_free_cell_index()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::get_matrix_free_cell_index ( const typename Triangulation< dim >::cell_iterator &  cell) const

Get MatrixFree index associated to a deal.II cell. To get the actual cell batch index and lane, do the postprocessing index / VectorizedArrayType::size() and index % VectorizedArrayType::size().

◆ get_face_iterator()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::pair< typename DoFHandler< dim >::cell_iterator, unsigned int > MatrixFree< dim, Number, VectorizedArrayType >::get_face_iterator ( const unsigned int  face_batch_index,
const unsigned int  lane_index,
const bool  interior = true,
const unsigned int  fe_component = 0 
) const

Return the cell iterator in deal.II speak to an interior/exterior cell of a face in a pair of a face batch and lane index. The second element of the pair is the face number so that the face iterator can be accessed: pair.first->face(pair.second);

Note that the face iterators in deal.II go through cells differently to what the face/boundary loop of this class does. This is because several faces are worked on together (vectorization), and since faces with neighbor cells on different MPI processors need to be accessed at a certain time when accessing remote data and overlapping communication with computation.

◆ at_irregular_cell()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
bool MatrixFree< dim, Number, VectorizedArrayType >::at_irregular_cell ( const unsigned int  cell_batch_index) const

Since this class uses vectorized data types with usually more than one value in the data field, a situation might occur when some components of the vector type do not correspond to an actual cell in the mesh. When using only this class, one usually does not need to bother about that fact since the values are padded with zeros. However, when this class is mixed with deal.II access to cells, care needs to be taken. This function returns true if not all n_lanes cells for the given cell_batch_index correspond to actual cells of the mesh and some are merely present for padding reasons. To find out how many cells are actually used, use the function n_active_entries_per_cell_batch().

◆ n_active_entries_per_cell_batch()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::n_active_entries_per_cell_batch ( const unsigned int  cell_batch_index) const

This query returns how many cells among the VectorizedArrayType::size() many cells within a cell batch to actual cells in the mesh, rather than being present for padding reasons. For most given cell batches in n_cell_batches(), this number is equal to VectorizedArrayType::size(), but there might be one or a few cell batches in the mesh (where the numbers do not add up) where only some of the cells within a batch are used, indicated by the function at_irregular_cell().

◆ n_active_entries_per_face_batch()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::n_active_entries_per_face_batch ( const unsigned int  face_batch_index) const

Use this function to find out how many faces over the length of vectorization data types correspond to real faces (both interior and boundary faces, as those use the same indexing but with different ranges) in the mesh. For most given indices in n_inner_faces_batches() and n_boundary_face_batches(), this is just vectorization_length many, but there might be one or a few meshes (where the numbers do not add up) where there are less such lanes filled.

◆ get_dofs_per_cell()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::get_dofs_per_cell ( const unsigned int  dof_handler_index = 0,
const unsigned int  hp_active_fe_index = 0 
) const

Return the number of degrees of freedom per cell for a given hp-index.

◆ get_n_q_points()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::get_n_q_points ( const unsigned int  quad_index = 0,
const unsigned int  hp_active_fe_index = 0 
) const

Return the number of quadrature points per cell for a given hp-index.

◆ get_dofs_per_face()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::get_dofs_per_face ( const unsigned int  dof_handler_index = 0,
const unsigned int  hp_active_fe_index = 0 
) const

Return the number of degrees of freedom on each face of the cell for given hp-index.

◆ get_n_q_points_face()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::get_n_q_points_face ( const unsigned int  quad_index = 0,
const unsigned int  hp_active_fe_index = 0 
) const

Return the number of quadrature points on each face of the cell for given hp-index.

◆ get_quadrature()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const Quadrature< dim > & MatrixFree< dim, Number, VectorizedArrayType >::get_quadrature ( const unsigned int  quad_index = 0,
const unsigned int  hp_active_fe_index = 0 
) const

Return the quadrature rule for given hp-index.

◆ get_face_quadrature()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const Quadrature< dim - 1 > & MatrixFree< dim, Number, VectorizedArrayType >::get_face_quadrature ( const unsigned int  quad_index = 0,
const unsigned int  hp_active_fe_index = 0 
) const

Return the quadrature rule for given hp-index.

◆ get_cell_range_category()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::get_cell_range_category ( const std::pair< unsigned int, unsigned int cell_batch_range) const

Return the category the current batch range of cells was assigned to. Categories run between the given values in the field AdditionalData::cell_vectorization_category for the non-hp case and return the active FE index in the hp-adaptive case.

Note
Following the behaviour of get_cell_category(), we return the maximum category of any cell batch. In the hp case, it is guaranteed that all cells and as a consequence all cell batches in a range have the same category. Otherwise, there may be different categories in different cell batches.

◆ get_face_range_category()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::pair< unsigned int, unsigned int > MatrixFree< dim, Number, VectorizedArrayType >::get_face_range_category ( const std::pair< unsigned int, unsigned int face_batch_range) const

Return the category of the cells on the two sides of the current batch range of faces.

◆ get_cell_category()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::get_cell_category ( const unsigned int  cell_batch_index) const

Return the category the current batch of cells was assigned to. Categories run between the given values in the field AdditionalData::cell_vectorization_category for the non-hp case and return the active FE index in the hp-adaptive case.

Note
In the non-hp case, a category of a cell batch is given as the maximum category of any of its cell. In the hp case or the case that MatrixFree::AdditionalData::cell_vectorization_categories_strict was enabled, it is guaranteed that all cells have the same category.

◆ get_face_category()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::pair< unsigned int, unsigned int > MatrixFree< dim, Number, VectorizedArrayType >::get_face_category ( const unsigned int  face_batch_index) const

Return the category of the cells on the two sides of the current batch of faces.

◆ indices_initialized()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
bool MatrixFree< dim, Number, VectorizedArrayType >::indices_initialized ( ) const

Queries whether or not the indexation has been set.

◆ mapping_initialized()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
bool MatrixFree< dim, Number, VectorizedArrayType >::mapping_initialized ( ) const

Queries whether or not the geometry-related information for the cells has been set.

◆ get_mg_level()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::get_mg_level ( ) const

Return the level of the mesh to be worked on. Returns numbers::invalid_unsigned_int if working on active cells.

◆ memory_consumption()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::size_t MatrixFree< dim, Number, VectorizedArrayType >::memory_consumption ( ) const

Return an approximation of the memory consumption of this class in bytes.

◆ print_memory_consumption()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename StreamType >
void MatrixFree< dim, Number, VectorizedArrayType >::print_memory_consumption ( StreamType &  out) const

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

◆ print()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
void MatrixFree< dim, Number, VectorizedArrayType >::print ( std::ostream &  out) const

Prints a summary of this class to the given output stream. It is focused on the indices, and does not print all the data stored.

◆ get_task_info()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const internal::MatrixFreeFunctions::TaskInfo & MatrixFree< dim, Number, VectorizedArrayType >::get_task_info ( ) const

Return information on task graph.

◆ get_mapping_info()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const internal::MatrixFreeFunctions:: MappingInfo< dim, Number, VectorizedArrayType > & MatrixFree< dim, Number, VectorizedArrayType >::get_mapping_info ( ) const

◆ get_dof_info()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const internal::MatrixFreeFunctions::DoFInfo & MatrixFree< dim, Number, VectorizedArrayType >::get_dof_info ( const unsigned int  dof_handler_index_component = 0) const

Return information on indexation degrees of freedom.

◆ n_constraint_pool_entries()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::n_constraint_pool_entries ( ) const

Return the number of weights in the constraint pool.

◆ constraint_pool_begin()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const Number * MatrixFree< dim, Number, VectorizedArrayType >::constraint_pool_begin ( const unsigned int  pool_index) const

Return a pointer to the first number in the constraint pool data with index pool_index (to be used together with constraint_pool_end()).

◆ constraint_pool_end()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const Number * MatrixFree< dim, Number, VectorizedArrayType >::constraint_pool_end ( const unsigned int  pool_index) const

Return a pointer to one past the last number in the constraint pool data with index pool_index (to be used together with constraint_pool_begin()).

◆ get_shape_info()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const internal::MatrixFreeFunctions::ShapeInfo< VectorizedArrayType > & MatrixFree< dim, Number, VectorizedArrayType >::get_shape_info ( const unsigned int  dof_handler_index_component = 0,
const unsigned int  quad_index = 0,
const unsigned int  fe_base_element = 0,
const unsigned int  hp_active_fe_index = 0,
const unsigned int  hp_active_quad_index = 0 
) const

Return the unit cell information for given hp-index.

◆ get_face_info()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const internal::MatrixFreeFunctions::FaceToCellTopology< VectorizedArrayType::size()> & MatrixFree< dim, Number, VectorizedArrayType >::get_face_info ( const unsigned int  face_batch_index) const

Return the connectivity information of a face.

◆ get_cell_and_face_to_plain_faces()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
const Table< 3, unsigned int > & MatrixFree< dim, Number, VectorizedArrayType >::get_cell_and_face_to_plain_faces ( ) const

Return the table that translates a triple of the cell-batch number, the index of a face within a cell and the index within the cell batch of vectorization into the index within the faces array.

◆ acquire_scratch_data()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
AlignedVector< VectorizedArrayType > * MatrixFree< dim, Number, VectorizedArrayType >::acquire_scratch_data ( ) const

Obtains a scratch data object for internal use. Make sure to release it afterwards by passing the pointer you obtain from this object to the release_scratch_data() function. This interface is used by FEEvaluation objects for storing their data structures.

The organization of the internal data structure is a thread-local storage of a list of vectors. Multiple threads will each get a separate storage field and separate vectors, ensuring thread safety. The mechanism to acquire and release objects is similar to the mechanisms used for the local contributions of WorkStream, see the WorkStream paper.

◆ release_scratch_data()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
void MatrixFree< dim, Number, VectorizedArrayType >::release_scratch_data ( const AlignedVector< VectorizedArrayType > *  memory) const

Makes the object of the scratchpad available again.

◆ acquire_scratch_data_non_threadsafe()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
AlignedVector< Number > * MatrixFree< dim, Number, VectorizedArrayType >::acquire_scratch_data_non_threadsafe ( ) const

Obtains a scratch data object for internal use. Make sure to release it afterwards by passing the pointer you obtain from this object to the release_scratch_data_non_threadsafe() function. Note that, as opposed to acquire_scratch_data(), this method can only be called by a single thread at a time, but opposed to the acquire_scratch_data() it is also possible that the thread releasing the scratch data can be different than the one that acquired it.

◆ release_scratch_data_non_threadsafe()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
void MatrixFree< dim, Number, VectorizedArrayType >::release_scratch_data_non_threadsafe ( const AlignedVector< Number > *  memory) const

Makes the object of the scratch data available again.

◆ internal_reinit()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename number2 , int q_dim>
void MatrixFree< dim, Number, VectorizedArrayType >::internal_reinit ( const std::shared_ptr< hp::MappingCollection< dim > > &  mapping,
const std::vector< const DoFHandler< dim, dim > * > &  dof_handlers,
const std::vector< const AffineConstraints< number2 > * > &  constraint,
const std::vector< IndexSet > &  locally_owned_set,
const std::vector< hp::QCollection< q_dim > > &  quad,
const AdditionalData additional_data 
)
private

This is the actual reinit function that sets up the indices for the DoFHandler case.

◆ initialize_indices()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
template<typename number2 >
void MatrixFree< dim, Number, VectorizedArrayType >::initialize_indices ( const std::vector< const AffineConstraints< number2 > * > &  constraint,
const std::vector< IndexSet > &  locally_owned_set,
const AdditionalData additional_data 
)
private

Initializes the fields in DoFInfo together with the constraint pool that holds all different weights in the constraints (not part of DoFInfo because several DoFInfo classes can have the same weights which consequently only need to be stored once).

◆ initialize_dof_handlers()

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
void MatrixFree< dim, Number, VectorizedArrayType >::initialize_dof_handlers ( const std::vector< const DoFHandler< dim, dim > * > &  dof_handlers,
const AdditionalData additional_data 
)
private

Initializes the DoFHandlers based on a DoFHandler<dim> argument.

◆ 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

◆ dimension

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
constexpr unsigned int MatrixFree< dim, Number, VectorizedArrayType >::dimension = dim
staticconstexpr

The dimension set by the template argument dim.

Definition at line 129 of file matrix_free.h.

◆ dof_handlers

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::vector<SmartPointer<const DoFHandler<dim> > > MatrixFree< dim, Number, VectorizedArrayType >::dof_handlers
private

Pointers to the DoFHandlers underlying the current problem.

Definition at line 2250 of file matrix_free.h.

◆ affine_constraints

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::vector<SmartPointer<const AffineConstraints<Number> > > MatrixFree< dim, Number, VectorizedArrayType >::affine_constraints
private

Pointers to the AffineConstraints underlying the current problem. Only filled with an AffineConstraints object if objects of the same Number template parameter as the Number template of MatrixFree is passed to reinit(). Filled with nullptr otherwise.

Definition at line 2258 of file matrix_free.h.

◆ dof_info

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::vector<internal::MatrixFreeFunctions::DoFInfo> MatrixFree< dim, Number, VectorizedArrayType >::dof_info
private

Contains the information about degrees of freedom on the individual cells and constraints.

Definition at line 2264 of file matrix_free.h.

◆ constraint_pool_data

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::vector<Number> MatrixFree< dim, Number, VectorizedArrayType >::constraint_pool_data
private

Contains the weights for constraints stored in DoFInfo. Filled into a separate field since several vector components might share similar weights, which reduces memory consumption. Moreover, it obviates template arguments on DoFInfo and keeps it a plain field of indices only.

Definition at line 2272 of file matrix_free.h.

◆ constraint_pool_row_index

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::vector<unsigned int> MatrixFree< dim, Number, VectorizedArrayType >::constraint_pool_row_index
private

Contains an indicator to the start of the ith index in the constraint pool data.

Definition at line 2278 of file matrix_free.h.

◆ mapping_info

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
internal::MatrixFreeFunctions::MappingInfo<dim, Number, VectorizedArrayType> MatrixFree< dim, Number, VectorizedArrayType >::mapping_info
private

Holds information on transformation of cells from reference cell to real cell that is needed for evaluating integrals.

Definition at line 2285 of file matrix_free.h.

◆ shape_info

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
Table<4, internal::MatrixFreeFunctions::ShapeInfo<VectorizedArrayType> > MatrixFree< dim, Number, VectorizedArrayType >::shape_info
private

Contains shape value information on the unit cell.

Definition at line 2291 of file matrix_free.h.

◆ cell_level_index

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::vector<std::pair<unsigned int, unsigned int> > MatrixFree< dim, Number, VectorizedArrayType >::cell_level_index
private

Describes how the cells are gone through. With the cell level (first index in this field) and the index within the level, one can reconstruct a deal.II cell iterator and use all the traditional things deal.II offers to do with cell iterators.

Definition at line 2299 of file matrix_free.h.

◆ mf_cell_indices

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::vector<unsigned int> MatrixFree< dim, Number, VectorizedArrayType >::mf_cell_indices
private

Conversion from deal.II index (active or level index) to MatrixFree index (inverse of cell_level_index).

Definition at line 2305 of file matrix_free.h.

◆ cell_level_index_end_local

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::cell_level_index_end_local
private

For discontinuous Galerkin, the cell_level_index includes cells that are not on the local processor but that are needed to evaluate the cell integrals. In cell_level_index_end_local, we store the number of local cells.

Definition at line 2313 of file matrix_free.h.

◆ task_info

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
internal::MatrixFreeFunctions::TaskInfo MatrixFree< dim, Number, VectorizedArrayType >::task_info
private

Stores the basic layout of the cells and faces to be treated, including the task layout for the shared memory parallelization and possible overlaps between communications and computations with MPI.

Definition at line 2320 of file matrix_free.h.

◆ face_info

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
internal::MatrixFreeFunctions::FaceInfo<VectorizedArrayType::size()> MatrixFree< dim, Number, VectorizedArrayType >::face_info
private

Vector holding face information. Only initialized if build_face_info=true.

Definition at line 2327 of file matrix_free.h.

◆ indices_are_initialized

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
bool MatrixFree< dim, Number, VectorizedArrayType >::indices_are_initialized
private

Stores whether indices have been initialized.

Definition at line 2332 of file matrix_free.h.

◆ mapping_is_initialized

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
bool MatrixFree< dim, Number, VectorizedArrayType >::mapping_is_initialized
private

Stores whether indices have been initialized.

Definition at line 2337 of file matrix_free.h.

◆ scratch_pad

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
Threads::ThreadLocalStorage< std::list<std::pair<bool, AlignedVector<VectorizedArrayType> > > > MatrixFree< dim, Number, VectorizedArrayType >::scratch_pad
mutableprivate

Scratchpad memory for use in evaluation. We allow more than one evaluation object to attach to this field (this, the outer std::vector), so we need to keep tracked of whether a certain data field is already used (first part of pair) and keep a list of objects.

Definition at line 2348 of file matrix_free.h.

◆ scratch_pad_non_threadsafe

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
std::list<std::pair<bool, AlignedVector<Number> > > MatrixFree< dim, Number, VectorizedArrayType >::scratch_pad_non_threadsafe
mutableprivate

Scratchpad memory for use in evaluation and other contexts, non-thread safe variant.

Definition at line 2355 of file matrix_free.h.

◆ mg_level

template<int dim, typename Number = double, typename VectorizedArrayType = VectorizedArray<Number>>
unsigned int MatrixFree< dim, Number, VectorizedArrayType >::mg_level
private

Stored the level of the mesh to be worked on.

Definition at line 2360 of file matrix_free.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: