Reference documentation for deal.II version 9.0.0
|
#include <deal.II/matrix_free/matrix_free.h>
Classes | |
struct | AdditionalData |
struct | DoFHandlers |
Public Types | |
enum | DataAccessOnFaces { DataAccessOnFaces::none, DataAccessOnFaces::values, DataAccessOnFaces::gradients, DataAccessOnFaces::unspecified } |
typedef Number | value_type |
Public Member Functions | |
1: Construction and initialization | |
MatrixFree () | |
MatrixFree (const MatrixFree< dim, Number > &other) | |
~MatrixFree ()=default | |
template<typename DoFHandlerType , typename QuadratureType > | |
void | reinit (const Mapping< dim > &mapping, const DoFHandlerType &dof_handler, const ConstraintMatrix &constraint, const QuadratureType &quad, const AdditionalData additional_data=AdditionalData()) |
template<typename DoFHandlerType , typename QuadratureType > | |
void | reinit (const DoFHandlerType &dof_handler, const ConstraintMatrix &constraint, const QuadratureType &quad, const AdditionalData additional_data=AdditionalData()) |
template<typename DoFHandlerType , typename QuadratureType > | |
void | reinit (const Mapping< dim > &mapping, const DoFHandlerType &dof_handler, const ConstraintMatrix &constraint, const IndexSet &locally_owned_dofs, const QuadratureType &quad, const AdditionalData additional_data=AdditionalData()) |
template<typename DoFHandlerType , typename QuadratureType > | |
void | reinit (const Mapping< dim > &mapping, const std::vector< const DoFHandlerType *> &dof_handler, const std::vector< const ConstraintMatrix *> &constraint, const std::vector< QuadratureType > &quad, const AdditionalData additional_data=AdditionalData()) |
template<typename DoFHandlerType , typename QuadratureType > | |
void | reinit (const std::vector< const DoFHandlerType *> &dof_handler, const std::vector< const ConstraintMatrix *> &constraint, const std::vector< QuadratureType > &quad, const AdditionalData additional_data=AdditionalData()) |
template<typename DoFHandlerType , typename QuadratureType > | |
void | reinit (const Mapping< dim > &mapping, const std::vector< const DoFHandlerType *> &dof_handler, const std::vector< const ConstraintMatrix *> &constraint, const std::vector< IndexSet > &locally_owned_set, const std::vector< QuadratureType > &quad, const AdditionalData additional_data=AdditionalData()) |
template<typename DoFHandlerType , typename QuadratureType > | |
void | reinit (const Mapping< dim > &mapping, const std::vector< const DoFHandlerType *> &dof_handler, const std::vector< const ConstraintMatrix *> &constraint, const QuadratureType &quad, const AdditionalData additional_data=AdditionalData()) |
template<typename DoFHandlerType , typename QuadratureType > | |
void | reinit (const std::vector< const DoFHandlerType *> &dof_handler, const std::vector< const ConstraintMatrix *> &constraint, const QuadratureType &quad, const AdditionalData additional_data=AdditionalData()) |
void | copy_from (const MatrixFree< dim, Number > &matrix_free_base) |
void | clear () |
2: Loop over cells | |
template<typename OutVector , typename InVector > | |
void | cell_loop (const std::function< void(const MatrixFree< dim, Number > &, 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 OutVector , typename InVector > | |
void | loop (const std::function< void(const MatrixFree< dim, Number > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &cell_operation, const std::function< void(const MatrixFree< dim, Number > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &face_operation, const std::function< void(const MatrixFree< dim, Number > &, 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 |
std::pair< unsigned int, unsigned int > | 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 |
std::pair< unsigned int, unsigned int > | 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 |
3: Initialization of vectors | |
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 IndexSet & | get_locally_owned_set (const unsigned int dof_handler_index=0) const |
const IndexSet & | get_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) |
5: Access of internal data structure (expert mode, interface not stable between releases) | |
const internal::MatrixFreeFunctions::TaskInfo & | get_task_info () const |
const internal::MatrixFreeFunctions::TaskInfo & | get_size_info () const |
const internal::MatrixFreeFunctions::MappingInfo< dim, Number > & | get_mapping_info () const |
const internal::MatrixFreeFunctions::DoFInfo & | get_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< VectorizedArray< Number > > & | 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< VectorizedArray< Number >::n_array_elements > & | get_face_info (const unsigned int face_batch_number) const |
AlignedVector< VectorizedArray< Number > > * | acquire_scratch_data () const |
void | release_scratch_data (const AlignedVector< VectorizedArray< Number > > *memory) const |
AlignedVector< Number > * | acquire_scratch_data_non_threadsafe () const |
void | release_scratch_data_non_threadsafe (const AlignedVector< Number > *memory) const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Attributes | |
static const unsigned int | dimension = dim |
Private Member Functions | |
void | internal_reinit (const Mapping< dim > &mapping, const std::vector< const DoFHandler< dim > *> &dof_handler, const std::vector< const ConstraintMatrix *> &constraint, const std::vector< IndexSet > &locally_owned_set, const std::vector< hp::QCollection< 1 > > &quad, const AdditionalData &additional_data) |
void | internal_reinit (const Mapping< dim > &mapping, const std::vector< const hp::DoFHandler< dim > *> &dof_handler, const std::vector< const ConstraintMatrix *> &constraint, const std::vector< IndexSet > &locally_owned_set, const std::vector< hp::QCollection< 1 > > &quad, const AdditionalData &additional_data) |
void | initialize_indices (const std::vector< const ConstraintMatrix *> &constraint, const std::vector< IndexSet > &locally_owned_set, const AdditionalData &additional_data) |
void | initialize_dof_handlers (const std::vector< const DoFHandler< dim > *> &dof_handlers, const AdditionalData &additional_data) |
void | initialize_dof_handlers (const std::vector< const hp::DoFHandler< dim > *> &dof_handlers, const AdditionalData &additional_data) |
void | make_connectivity_graph_faces (DynamicSparsityPattern &connectivity) |
Private Attributes | |
DoFHandlers | dof_handlers |
std::vector< internal::MatrixFreeFunctions::DoFInfo > | dof_info |
std::vector< Number > | constraint_pool_data |
std::vector< unsigned int > | constraint_pool_row_index |
internal::MatrixFreeFunctions::MappingInfo< dim, Number > | mapping_info |
Table< 4, internal::MatrixFreeFunctions::ShapeInfo< VectorizedArray< Number > > > | shape_info |
std::vector< std::pair< unsigned int, unsigned int > > | cell_level_index |
unsigned int | cell_level_index_end_local |
internal::MatrixFreeFunctions::TaskInfo | task_info |
internal::MatrixFreeFunctions::FaceInfo< VectorizedArray< Number >::n_array_elements > | face_info |
bool | indices_are_initialized |
bool | mapping_is_initialized |
Threads::ThreadLocalStorage< std::list< std::pair< bool, AlignedVector< VectorizedArray< Number > > > > > | scratch_pad |
std::list< std::pair< bool, AlignedVector< Number > > > | scratch_pad_non_threadsafe |
4: 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_macro_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 macro_face) const |
std::array< types::boundary_id, VectorizedArray< Number >::n_array_elements > | get_faces_by_cells_boundary_id (const unsigned int macro_cell, const unsigned int face_number) const |
const DoFHandler< dim > & | get_dof_handler (const unsigned int dof_handler_index=0) const |
DoFHandler< dim >::cell_iterator | get_cell_iterator (const unsigned int macro_cell_number, const unsigned int vector_number, const unsigned int fe_component=0) const |
hp::DoFHandler< dim >::active_cell_iterator | get_hp_cell_iterator (const unsigned int macro_cell_number, const unsigned int vector_number, const unsigned int dof_handler_index=0) const |
bool | at_irregular_cell (const unsigned int macro_cell_number) const |
unsigned int | n_components_filled (const unsigned int cell_batch_number) const |
unsigned int | n_active_entries_per_cell_batch (const unsigned int cell_batch_number) const |
unsigned int | n_active_entries_per_face_batch (const unsigned int face_batch_number) 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 fe_component=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_category (const unsigned int macro_cell) const |
std::pair< unsigned int, unsigned int > | get_face_category (const unsigned int macro_face) const |
bool | indices_initialized () const |
bool | mapping_initialized () 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) |
Additional Inherited Members | |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
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 macro cell. 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 105 of file matrix_free.h.
typedef Number MatrixFree< dim, Number >::value_type |
A typedef for the underlying number type specified by the template argument.
Definition at line 112 of file matrix_free.h.
|
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 |
gradients | The loop does involve FEFaceEvaluation access into neighbors by function values and gradients, but no second derivatives, such as |
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 614 of file matrix_free.h.
MatrixFree< dim, Number >::MatrixFree | ( | ) |
Default empty constructor. Does nothing.
MatrixFree< dim, Number >::MatrixFree | ( | const MatrixFree< dim, Number > & | other | ) |
Copy constructor, calls copy_from
|
default |
Destructor.
void MatrixFree< dim, Number >::reinit | ( | const Mapping< dim > & | mapping, |
const DoFHandlerType & | dof_handler, | ||
const ConstraintMatrix & | constraint, | ||
const QuadratureType & | quad, | ||
const AdditionalData | additional_data = AdditionalData() |
||
) |
Extracts the information needed to perform loops over cells. The DoFHandler and ConstraintMatrix 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.
void MatrixFree< dim, Number >::reinit | ( | const DoFHandlerType & | dof_handler, |
const ConstraintMatrix & | constraint, | ||
const QuadratureType & | quad, | ||
const AdditionalData | additional_data = AdditionalData() |
||
) |
Initializes the data structures. Same as above, but using a \(Q_1\) mapping.
void MatrixFree< dim, Number >::reinit | ( | const Mapping< dim > & | mapping, |
const DoFHandlerType & | dof_handler, | ||
const ConstraintMatrix & | constraint, | ||
const IndexSet & | locally_owned_dofs, | ||
const QuadratureType & | quad, | ||
const AdditionalData | additional_data = AdditionalData() |
||
) |
Same as above.
void MatrixFree< dim, Number >::reinit | ( | const Mapping< dim > & | mapping, |
const std::vector< const DoFHandlerType *> & | dof_handler, | ||
const std::vector< const ConstraintMatrix *> & | constraint, | ||
const std::vector< QuadratureType > & | quad, | ||
const AdditionalData | additional_data = AdditionalData() |
||
) |
Extracts the information needed to perform loops over cells. The DoFHandler and ConstraintMatrix 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. Note that the finite element underlying an DoFHandler must either be scalar or contain several copies of the same element. Mixing several different elements into one FE_System
is not allowed.
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.
void MatrixFree< dim, Number >::reinit | ( | const std::vector< const DoFHandlerType *> & | dof_handler, |
const std::vector< const ConstraintMatrix *> & | constraint, | ||
const std::vector< QuadratureType > & | quad, | ||
const AdditionalData | additional_data = AdditionalData() |
||
) |
Initializes the data structures. Same as above, but using a \(Q_1\) mapping.
void MatrixFree< dim, Number >::reinit | ( | const Mapping< dim > & | mapping, |
const std::vector< const DoFHandlerType *> & | dof_handler, | ||
const std::vector< const ConstraintMatrix *> & | constraint, | ||
const std::vector< IndexSet > & | locally_owned_set, | ||
const std::vector< QuadratureType > & | quad, | ||
const AdditionalData | additional_data = AdditionalData() |
||
) |
Same as above.
void MatrixFree< dim, Number >::reinit | ( | const Mapping< dim > & | mapping, |
const std::vector< const DoFHandlerType *> & | dof_handler, | ||
const std::vector< const ConstraintMatrix *> & | 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.
void MatrixFree< dim, Number >::reinit | ( | const std::vector< const DoFHandlerType *> & | dof_handler, |
const std::vector< const ConstraintMatrix *> & | constraint, | ||
const QuadratureType & | quad, | ||
const AdditionalData | additional_data = AdditionalData() |
||
) |
Initializes the data structures. Same as above, but using a \(Q_1\) mapping.
void MatrixFree< dim, Number >::copy_from | ( | const MatrixFree< dim, Number > & | 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.
void MatrixFree< dim, Number >::clear | ( | ) |
Clear all data fields and brings the class into a condition similar to after having called the default constructor.
void MatrixFree< dim, Number >::cell_loop | ( | const std::function< void(const MatrixFree< dim, Number > &, 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.
cell_operation | std::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. |
dst | Destination 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. |
src | Input 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_vector | If 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. |
void MatrixFree< dim, Number >::cell_loop | ( | void(CLASS::*)(const MatrixFree< dim, Number > &, 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 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).
cell_operation | Pointer 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_class | The object which provides the cell_operation call. To be compatible with this interface, the class must allow to call owning_class->cell_operation(...) . |
dst | Destination 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. |
src | Input 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_vector | If 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. |
void MatrixFree< dim, Number >::cell_loop | ( | void(CLASS::*)(const MatrixFree< dim, Number > &, 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.
void MatrixFree< dim, Number >::loop | ( | const std::function< void(const MatrixFree< dim, Number > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> & | cell_operation, |
const std::function< void(const MatrixFree< dim, Number > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> & | face_operation, | ||
const std::function< void(const MatrixFree< dim, Number > &, 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.
cell_operation | std::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_operation | std::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_operation | std::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. |
dst | Destination 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. |
src | Input 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_vector | If 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_access | Set 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_access | Set 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 . |
void MatrixFree< dim, Number >::loop | ( | void(CLASS::*)(const MatrixFree< dim, Number > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const | cell_operation, |
void(CLASS::*)(const MatrixFree< dim, Number > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const | face_operation, | ||
void(CLASS::*)(const MatrixFree< dim, Number > &, 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 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).
cell_operation | Pointer 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_operation | Pointer 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_operation | Pointer 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_class | The 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(...) . |
dst | Destination 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. |
src | Input 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_vector | If 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_access | Set 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_access | Set 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 . |
void MatrixFree< dim, Number >::loop | ( | void(CLASS::*)(const MatrixFree< dim, Number > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) | cell_operation, |
void(CLASS::*)(const MatrixFree< dim, Number > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) | face_operation, | ||
void(CLASS::*)(const MatrixFree< dim, Number > &, 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.
std::pair<unsigned int,unsigned int> MatrixFree< dim, Number >::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.
std::pair<unsigned int,unsigned int> MatrixFree< dim, Number >::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.
void MatrixFree< dim, Number >::initialize_dof_vector | ( | VectorType & | vec, |
const unsigned int | dof_handler_index = 0 |
||
) | const |
Initialize function for a general vector. The length of the vector is equal to the total number of degrees in the DoFHandler. If the vector is of class LinearAlgebra::distributed::Vector<Number>, the ghost entries are set accordingly. For vector-valued problems with several DoFHandlers underlying this class, the parameter vector_component
defines which component is to be used.
For the 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.
void MatrixFree< dim, Number >::initialize_dof_vector | ( | LinearAlgebra::distributed::Vector< Number2 > & | vec, |
const unsigned int | dof_handler_index = 0 |
||
) | const |
Initialize function for a distributed vector. The length of the vector is equal to the total number of degrees in the DoFHandler. If the vector is of class LinearAlgebra::distributed::Vector<Number>, the ghost entries are set accordingly. For vector-valued problems with several DoFHandlers underlying this class, the parameter vector_component
defines which component is to be used.
For the 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.
const std::shared_ptr<const Utilities::MPI::Partitioner>& MatrixFree< dim, Number >::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.
const IndexSet& MatrixFree< dim, Number >::get_locally_owned_set | ( | const unsigned int | dof_handler_index = 0 | ) | const |
Return the set of cells that are oned by the processor.
const IndexSet& MatrixFree< dim, Number >::get_ghost_set | ( | const unsigned int | dof_handler_index = 0 | ) | const |
Return the set of ghost cells needed but not owned by the processor.
const std::vector<unsigned int>& MatrixFree< dim, Number >::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.
void MatrixFree< dim, Number >::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 ConstraintMatrix. Note that if a DoFHandler calls DoFHandler::renumber_dofs, all information in MatrixFree becomes invalid.
|
static |
Return whether a given FiniteElement fe
is supported by this class.
unsigned int MatrixFree< dim, Number >::n_components | ( | ) | const |
Return the number of different DoFHandlers specified at initialization.
unsigned int MatrixFree< dim, Number >::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.
unsigned int MatrixFree< dim, Number >::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.
unsigned int MatrixFree< dim, Number >::n_macro_cells | ( | ) | 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::n_array_elements
(depending on how many cell chunks that do not get filled up completely).
unsigned int MatrixFree< dim, Number >::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::n_array_elements
(depending on how many cell chunks that do not get filled up completely).
unsigned int MatrixFree< dim, Number >::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.
unsigned int MatrixFree< dim, Number >::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.
unsigned int MatrixFree< dim, Number >::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.
unsigned int MatrixFree< dim, Number >::n_ghost_inner_face_batches | ( | ) | const |
Return the number of faces that are not processed locally but belong to locally owned faces.
types::boundary_id MatrixFree< dim, Number >::get_boundary_id | ( | const unsigned int | macro_face | ) | 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.
std::array<types::boundary_id, VectorizedArray<Number>::n_array_elements> MatrixFree< dim, Number >::get_faces_by_cells_boundary_id | ( | const unsigned int | macro_cell, |
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.
const DoFHandler<dim>& MatrixFree< dim, Number >::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.
DoFHandler<dim>::cell_iterator MatrixFree< dim, Number >::get_cell_iterator | ( | const unsigned int | macro_cell_number, |
const unsigned int | vector_number, | ||
const unsigned int | fe_component = 0 |
||
) | const |
Return the cell iterator in deal.II speak to a given cell 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 worked on together (vectorization), 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.
hp::DoFHandler<dim>::active_cell_iterator MatrixFree< dim, Number >::get_hp_cell_iterator | ( | const unsigned int | macro_cell_number, |
const unsigned int | vector_number, | ||
const unsigned int | dof_handler_index = 0 |
||
) | const |
This returns the cell iterator in deal.II speak to a given cell in the renumbering of this structure. This function returns an exception in case the structure was not constructed based on an hp::DoFHandler.
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 worked on together (vectorization), 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.
bool MatrixFree< dim, Number >::at_irregular_cell | ( | const unsigned int | macro_cell_number | ) | 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 vectorization_length
cells for the given macro_cell
are real cells. To find out how many cells are actually used, use the function n_active_entries_per_cell_batch
.
unsigned int MatrixFree< dim, Number >::n_components_filled | ( | const unsigned int | cell_batch_number | ) | const |
This query returns how many cells over the length of vectorization data types correspond to actual cells in the mesh. For most given cell_batch_number
, 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 components filled, indicated by the function at_irregular_cell
.
unsigned int MatrixFree< dim, Number >::n_active_entries_per_cell_batch | ( | const unsigned int | cell_batch_number | ) | const |
This query returns how many cells over the length of vectorization data types correspond to actual cells in the mesh. For most given cell batches in n_cell_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 components filled, indicated by the function at_irregular_cell
.
unsigned int MatrixFree< dim, Number >::n_active_entries_per_face_batch | ( | const unsigned int | face_batch_number | ) | 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.
unsigned int MatrixFree< dim, Number >::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.
unsigned int MatrixFree< dim, Number >::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.
unsigned int MatrixFree< dim, Number >::get_dofs_per_face | ( | const unsigned int | fe_component = 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.
unsigned int MatrixFree< dim, Number >::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.
const Quadrature<dim>& MatrixFree< dim, Number >::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.
const Quadrature<dim-1>& MatrixFree< dim, Number >::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.
unsigned int MatrixFree< dim, Number >::get_cell_category | ( | const unsigned int | macro_cell | ) | 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 non-hp DoFHandler types and return the active FE index in the hp-adaptive case.
std::pair<unsigned int,unsigned int> MatrixFree< dim, Number >::get_face_category | ( | const unsigned int | macro_face | ) | const |
Return the category on the cells on the two sides of the current batch of faces.
bool MatrixFree< dim, Number >::indices_initialized | ( | ) | const |
Queries whether or not the indexation has been set.
bool MatrixFree< dim, Number >::mapping_initialized | ( | ) | const |
Queries whether or not the geometry-related information for the cells has been set.
std::size_t MatrixFree< dim, Number >::memory_consumption | ( | ) | const |
Return an approximation of the memory consumption of this class in bytes.
void MatrixFree< dim, Number >::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.
void MatrixFree< dim, Number >::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.
const internal::MatrixFreeFunctions::TaskInfo& MatrixFree< dim, Number >::get_task_info | ( | ) | const |
Return information on task graph.
const internal::MatrixFreeFunctions::TaskInfo& MatrixFree< dim, Number >::get_size_info | ( | ) | const |
Return information on system size.
const internal::MatrixFreeFunctions::DoFInfo& MatrixFree< dim, Number >::get_dof_info | ( | const unsigned int | dof_handler_index_component = 0 | ) | const |
Return information on indexation degrees of freedom.
unsigned int MatrixFree< dim, Number >::n_constraint_pool_entries | ( | ) | const |
Return the number of weights in the constraint pool.
const Number* MatrixFree< dim, Number >::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()
).
const Number* MatrixFree< dim, Number >::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()
).
const internal::MatrixFreeFunctions::ShapeInfo<VectorizedArray<Number> >& MatrixFree< dim, Number >::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.
const internal::MatrixFreeFunctions::FaceToCellTopology<VectorizedArray<Number>::n_array_elements>& MatrixFree< dim, Number >::get_face_info | ( | const unsigned int | face_batch_number | ) | const |
Return the connectivity information of a face.
AlignedVector<VectorizedArray<Number> >* MatrixFree< dim, Number >::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.
void MatrixFree< dim, Number >::release_scratch_data | ( | const AlignedVector< VectorizedArray< Number > > * | memory | ) | const |
Makes the object of the scratchpad available again.
AlignedVector<Number>* MatrixFree< dim, Number >::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.
void MatrixFree< dim, Number >::release_scratch_data_non_threadsafe | ( | const AlignedVector< Number > * | memory | ) | const |
Makes the object of the scratch data available again.
|
private |
This is the actual reinit function that sets up the indices for the DoFHandler case.
|
private |
Same as before but for hp::DoFHandler instead of generic DoFHandler type.
|
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).
|
private |
Initializes the DoFHandlers based on a DoFHandler<dim> argument.
|
private |
Initializes the DoFHandlers based on a hp::DoFHandler<dim> argument.
|
private |
Setup connectivity graph with information on the dependencies between block due to shared faces.
|
static |
The dimension set by the template argument dim
.
Definition at line 117 of file matrix_free.h.
|
private |
Pointers to the DoFHandlers underlying the current problem.
Definition at line 1601 of file matrix_free.h.
|
private |
Contains the information about degrees of freedom on the individual cells and constraints.
Definition at line 1607 of file matrix_free.h.
|
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 1615 of file matrix_free.h.
|
private |
Contains an indicator to the start of the ith index in the constraint pool data.
Definition at line 1621 of file matrix_free.h.
|
private |
Holds information on transformation of cells from reference cell to real cell that is needed for evaluating integrals.
Definition at line 1627 of file matrix_free.h.
|
private |
Contains shape value information on the unit cell.
Definition at line 1632 of file matrix_free.h.
|
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 1640 of file matrix_free.h.
|
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 1649 of file matrix_free.h.
|
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 1656 of file matrix_free.h.
|
private |
Vector holding face information. Only initialized if build_face_info=true.
Definition at line 1662 of file matrix_free.h.
|
private |
Stores whether indices have been initialized.
Definition at line 1667 of file matrix_free.h.
|
private |
Stores whether indices have been initialized.
Definition at line 1672 of file matrix_free.h.
|
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 1681 of file matrix_free.h.
|
mutableprivate |
Scratchpad memory for use in evaluation and other contexts, non-thread safe variant.
Definition at line 1687 of file matrix_free.h.