Reference documentation for deal.II version 8.5.1
Classes | Private Member Functions | Private Attributes | List of all members
MatrixFree< dim, Number > Class Template Reference

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

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

Classes

struct  AdditionalData
 
struct  DoFHandlers
 

Public Member Functions

1: Construction and initialization
 MatrixFree ()
 
 MatrixFree (const MatrixFree< dim, Number > &other)
 
 ~MatrixFree ()
 
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()) 1
 
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()) 1
 
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_cxx11::function< void(const MatrixFree< dim, Number > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &)> &cell_operation, OutVector &dst, const InVector &src) const
 
template<typename CLASS , typename OutVector , typename InVector >
void cell_loop (void(CLASS::*function_pointer)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const, const CLASS *owning_class, OutVector &dst, const InVector &src) const
 
template<typename CLASS , typename OutVector , typename InVector >
void cell_loop (void(CLASS::*function_pointer)(const MatrixFree &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &), CLASS *owning_class, OutVector &dst, const InVector &src) 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 vector_component=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 vector_component=0) const
 
3: Initialization of vectors
template<typename VectorType >
void initialize_dof_vector (VectorType &vec, const unsigned int vector_component=0) const
 
template<typename Number2 >
void initialize_dof_vector (LinearAlgebra::distributed::Vector< Number2 > &vec, const unsigned int vector_component=0) const
 
const std_cxx11::shared_ptr< const Utilities::MPI::Partitioner > & get_vector_partitioner (const unsigned int vector_component=0) const
 
const IndexSetget_locally_owned_set (const unsigned int fe_component=0) const
 
const IndexSetget_ghost_set (const unsigned int fe_component=0) const
 
const std::vector< unsigned int > & get_constrained_dofs (const unsigned int fe_component=0) const
 
void renumber_dofs (std::vector< types::global_dof_index > &renumbering, const unsigned int vector_component=0)
 
5: Access of internal data structure (expert mode)
const internal::MatrixFreeFunctions::TaskInfoget_task_info () const
 
const internal::MatrixFreeFunctions::SizeInfoget_size_info () const
 
const internal::MatrixFreeFunctions::MappingInfo< dim, Number > & get_mapping_info () const
 
const internal::MatrixFreeFunctions::DoFInfo & get_dof_info (const unsigned int fe_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< Number > & get_shape_info (const unsigned int fe_component=0, const unsigned int quad_index=0, const unsigned int hp_active_fe_index=0, const unsigned int hp_active_quad_index=0) const
 
AlignedVector< VectorizedArray< Number > > * acquire_scratch_data () const
 
void release_scratch_data (const AlignedVector< VectorizedArray< Number > > *memory) const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&)
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&)
 
void subscribe (const char *identifier=0) const
 
void unsubscribe (const char *identifier=0) const
 
unsigned int n_subscriptions () const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

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)
 
void initialize_dof_handlers (const std::vector< const DoFHandler< dim > *> &dof_handlers, const unsigned int level)
 
void initialize_dof_handlers (const std::vector< const hp::DoFHandler< dim > *> &dof_handlers, const unsigned int level)
 

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< Number > > shape_info
 
std::vector< std::pair< unsigned int, unsigned int > > cell_level_index
 
internal::MatrixFreeFunctions::SizeInfo size_info
 
internal::MatrixFreeFunctions::TaskInfo task_info
 
bool indices_are_initialized
 
bool mapping_is_initialized
 
Threads::ThreadLocalStorage< std::list< std::pair< bool, AlignedVector< VectorizedArray< Number > > > > > scratch_pad
 

4: General information

unsigned int n_components () const
 
unsigned int n_physical_cells () const
 
unsigned int n_macro_cells () const
 
const DoFHandler< dim > & get_dof_handler (const unsigned int fe_component=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 fe_component=0) const
 
bool at_irregular_cell (const unsigned int macro_cell_number) const
 
unsigned int n_components_filled (const unsigned int macro_cell_number) const
 
unsigned int get_dofs_per_cell (const unsigned int fe_component=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
 
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 ::ExceptionBaseExcInUse (int arg1, char *arg2, std::string &arg3)
 
static ::ExceptionBaseExcNoSubscriber (char *arg1, char *arg2)
 

Detailed Description

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

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.

Author
Katharina Kormann, Martin Kronbichler, 2010, 2011

Definition at line 109 of file matrix_free.h.

Constructor & Destructor Documentation

◆ MatrixFree() [1/2]

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

Default empty constructor. Does nothing.

◆ MatrixFree() [2/2]

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

Copy constructor, calls copy_from

◆ ~MatrixFree()

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

Destructor.

Member Function Documentation

◆ reinit() [1/8]

template<int dim, typename Number = double>
template<typename DoFHandlerType , typename QuadratureType >
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.

◆ reinit() [2/8]

template<int dim, typename Number = double>
template<typename DoFHandlerType , typename QuadratureType >
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.

◆ reinit() [3/8]

template<int dim, typename Number = double>
template<typename DoFHandlerType , typename QuadratureType >
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.

Deprecated:
Setting the index set specifically is not supported any more. Use the reinit function without index set argument to choose the one provided by DoFHandler::locally_owned_dofs().

◆ reinit() [4/8]

template<int dim, typename Number = double>
template<typename DoFHandlerType , typename QuadratureType >
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.

◆ reinit() [5/8]

template<int dim, typename Number = double>
template<typename DoFHandlerType , typename QuadratureType >
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.

◆ reinit() [6/8]

template<int dim, typename Number = double>
template<typename DoFHandlerType , typename QuadratureType >
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.

Deprecated:
Setting the index set specifically is not supported any more. Use the reinit function without index set argument to choose the one provided by DoFHandler::locally_owned_dofs().

◆ reinit() [7/8]

template<int dim, typename Number = double>
template<typename DoFHandlerType , typename QuadratureType >
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.

◆ reinit() [8/8]

template<int dim, typename Number = double>
template<typename DoFHandlerType , typename QuadratureType >
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.

◆ copy_from()

template<int dim, typename Number = double>
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.

◆ clear()

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

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

◆ cell_loop() [1/3]

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

This method runs the loop over all cells (in parallel) and performs the MPI data exchange on the source vector and destination vector. The first argument indicates a function object that has the following 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.

◆ cell_loop() [2/3]

template<int dim, typename Number = double>
template<typename CLASS , typename OutVector , typename InVector >
void MatrixFree< dim, Number >::cell_loop ( void(CLASS::*)(const MatrixFree< dim, Number > &, OutVector &, const InVector &, const std::pair< unsigned int, unsigned int > &) const  function_pointer,
const CLASS *  owning_class,
OutVector &  dst,
const InVector &  src 
) 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 with the signature cell_operation (const MatrixFree<dim,Number> &, OutVector &, InVector &, std::pair<unsigned int,unsigned int>&)const. This method obviates the need to call std_cxx11::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_loop() [3/3]

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

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

◆ create_cell_subrange_hp()

template<int dim, typename Number = double>
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  vector_component = 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>
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  vector_component = 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.

◆ initialize_dof_vector() [1/2]

template<int dim, typename Number = double>
template<typename VectorType >
void MatrixFree< dim, Number >::initialize_dof_vector ( VectorType &  vec,
const unsigned int  vector_component = 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 contraints 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.

◆ initialize_dof_vector() [2/2]

template<int dim, typename Number = double>
template<typename Number2 >
void MatrixFree< dim, Number >::initialize_dof_vector ( LinearAlgebra::distributed::Vector< Number2 > &  vec,
const unsigned int  vector_component = 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 contraints 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>
const std_cxx11::shared_ptr<const Utilities::MPI::Partitioner>& MatrixFree< dim, Number >::get_vector_partitioner ( const unsigned int  vector_component = 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>
const IndexSet& MatrixFree< dim, Number >::get_locally_owned_set ( const unsigned int  fe_component = 0) const

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

◆ get_ghost_set()

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

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

◆ get_constrained_dofs()

template<int dim, typename Number = double>
const std::vector<unsigned int>& MatrixFree< dim, Number >::get_constrained_dofs ( const unsigned int  fe_component = 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>
void MatrixFree< dim, Number >::renumber_dofs ( std::vector< types::global_dof_index > &  renumbering,
const unsigned int  vector_component = 0 
)

Calls renumber_dofs function in dof_info which renumbers the degrees of freedom according to the ordering for parallelization.

◆ is_supported()

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

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

◆ n_components()

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

Return the number of different DoFHandlers specified at initialization.

◆ n_physical_cells()

template<int dim, typename Number = double>
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_macro_cells() which gives the number of cells as seen when lumping several cells together with vectorization.

◆ n_macro_cells()

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

Return the number of macro cells that this structure works on, i.e., the number of cell chunks that are worked on after the application of vectorization which in general works on several cells at once. The cell range in cell_loop runs from zero to n_macro_cells() (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).

◆ get_dof_handler()

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

In case this structure was built based on a DoFHandler, this returns the DoFHandler.

◆ get_cell_iterator()

template<int dim, typename Number = double>
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

This returns 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.

◆ get_hp_cell_iterator()

template<int dim, typename Number = double>
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  fe_component = 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.

◆ at_irregular_cell()

template<int dim, typename Number = double>
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_components_filled.

◆ n_components_filled()

template<int dim, typename Number = double>
unsigned int MatrixFree< dim, Number >::n_components_filled ( const unsigned int  macro_cell_number) const

Use this function to find out how many cells over the length of vectorization data types correspond to real cells in the mesh. For most given macro_cells, 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.

◆ get_dofs_per_cell()

template<int dim, typename Number = double>
unsigned int MatrixFree< dim, Number >::get_dofs_per_cell ( const unsigned int  fe_component = 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>
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.

◆ get_dofs_per_face()

template<int dim, typename Number = double>
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.

◆ get_n_q_points_face()

template<int dim, typename Number = double>
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.

◆ get_quadrature()

template<int dim, typename Number = double>
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.

◆ get_face_quadrature()

template<int dim, typename Number = double>
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.

◆ indices_initialized()

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

Queries whether or not the indexation has been set.

◆ mapping_initialized()

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

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

◆ memory_consumption()

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

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

◆ print_memory_consumption()

template<int dim, typename Number = double>
template<typename StreamType >
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.

◆ print()

template<int dim, typename Number = double>
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.

◆ get_task_info()

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

Return information on task graph.

◆ get_size_info()

template<int dim, typename Number = double>
const internal::MatrixFreeFunctions::SizeInfo& MatrixFree< dim, Number >::get_size_info ( ) const

Return information on system size.

◆ get_dof_info()

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

Return information on indexation degrees of freedom.

◆ n_constraint_pool_entries()

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

Return the number of weights in the constraint pool.

◆ constraint_pool_begin()

template<int dim, typename Number = double>
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()).

◆ constraint_pool_end()

template<int dim, typename Number = double>
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()).

◆ get_shape_info()

template<int dim, typename Number = double>
const internal::MatrixFreeFunctions::ShapeInfo<Number>& MatrixFree< dim, Number >::get_shape_info ( const unsigned int  fe_component = 0,
const unsigned int  quad_index = 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.

◆ acquire_scratch_data()

template<int dim, typename Number = double>
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

◆ release_scratch_data()

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

Makes the object of the scratchpad available again.

◆ internal_reinit() [1/2]

template<int dim, typename Number = double>
void MatrixFree< dim, Number >::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 
)
private

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

◆ internal_reinit() [2/2]

template<int dim, typename Number = double>
void MatrixFree< dim, Number >::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 
)
private

Same as before but for hp::DoFHandler instead of generic DoFHandler type.

◆ initialize_indices()

template<int dim, typename Number = double>
void MatrixFree< dim, Number >::initialize_indices ( const std::vector< const ConstraintMatrix *> &  constraint,
const std::vector< IndexSet > &  locally_owned_set 
)
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() [1/2]

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

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

◆ initialize_dof_handlers() [2/2]

template<int dim, typename Number = double>
void MatrixFree< dim, Number >::initialize_dof_handlers ( const std::vector< const hp::DoFHandler< dim > *> &  dof_handlers,
const unsigned int  level 
)
private

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

Member Data Documentation

◆ dof_handlers

template<int dim, typename Number = double>
DoFHandlers MatrixFree< dim, Number >::dof_handlers
private

Pointers to the DoFHandlers underlying the current problem.

Definition at line 993 of file matrix_free.h.

◆ dof_info

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

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

Definition at line 999 of file matrix_free.h.

◆ constraint_pool_data

template<int dim, typename Number = double>
std::vector<Number> MatrixFree< dim, Number >::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 1007 of file matrix_free.h.

◆ constraint_pool_row_index

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

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

Definition at line 1013 of file matrix_free.h.

◆ mapping_info

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

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

Definition at line 1019 of file matrix_free.h.

◆ shape_info

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

Contains shape value information on the unit cell.

Definition at line 1024 of file matrix_free.h.

◆ cell_level_index

template<int dim, typename Number = double>
std::vector<std::pair<unsigned int,unsigned int> > MatrixFree< dim, Number >::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 1032 of file matrix_free.h.

◆ size_info

template<int dim, typename Number = double>
internal::MatrixFreeFunctions::SizeInfo MatrixFree< dim, Number >::size_info
private

Stores how many cells we have, how many cells that we see after applying vectorization (i.e., the number of macro cells), and MPI-related stuff.

Definition at line 1038 of file matrix_free.h.

◆ task_info

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

Information regarding the shared memory parallelization.

Definition at line 1043 of file matrix_free.h.

◆ indices_are_initialized

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

Stores whether indices have been initialized.

Definition at line 1048 of file matrix_free.h.

◆ mapping_is_initialized

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

Stores whether indices have been initialized.

Definition at line 1053 of file matrix_free.h.

◆ scratch_pad

template<int dim, typename Number = double>
Threads::ThreadLocalStorage<std::list<std::pair<bool, AlignedVector<VectorizedArray<Number> > > > > MatrixFree< dim, Number >::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 1062 of file matrix_free.h.


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