Reference documentation for deal.II version 8.5.1
|
#include <deal.II/multigrid/mg_transfer_matrix_free.h>
Public Member Functions | |
MGTransferMatrixFree () | |
MGTransferMatrixFree (const MGConstrainedDoFs &mg_constrained_dofs) | |
virtual | ~MGTransferMatrixFree () |
void | initialize_constraints (const MGConstrainedDoFs &mg_constrained_dofs) |
void | clear () |
void | build (const DoFHandler< dim, dim > &mg_dof) |
virtual void | prolongate (const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const |
virtual void | restrict_and_add (const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const |
std::size_t | memory_consumption () const |
Public Member Functions inherited from MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > | |
void | clear () |
template<int dim, typename Number2 , int spacedim> | |
void | copy_to_mg (const DoFHandler< dim, spacedim > &mg_dof, MGLevelObject< LinearAlgebra::distributed::Vector< Number > > &dst, const LinearAlgebra::distributed::Vector< Number2 > &src) const |
template<int dim, typename Number2 , int spacedim> | |
void | copy_from_mg (const DoFHandler< dim, spacedim > &mg_dof, LinearAlgebra::distributed::Vector< Number2 > &dst, const MGLevelObject< LinearAlgebra::distributed::Vector< Number > > &src) const |
template<int dim, typename Number2 , int spacedim> | |
void | copy_from_mg_add (const DoFHandler< dim, spacedim > &mg_dof, LinearAlgebra::distributed::Vector< Number2 > &dst, const MGLevelObject< LinearAlgebra::distributed::Vector< Number > > &src) const |
void | set_component_to_block_map (const std::vector< unsigned int > &map) |
std::size_t | memory_consumption () const |
void | print_indices (std::ostream &os) const |
Public Member Functions inherited from MGTransferBase< LinearAlgebra::distributed::Vector< Number > > | |
virtual | ~MGTransferBase () |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (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) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcNoProlongation () |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, char *arg2, std::string &arg3) |
static ::ExceptionBase & | ExcNoSubscriber (char *arg1, char *arg2) |
Private Member Functions | |
template<int degree> | |
void | do_prolongate_add (const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const |
template<int degree> | |
void | do_restrict_add (const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const |
Private Attributes | |
unsigned int | fe_degree |
bool | element_is_continuous |
unsigned int | n_components |
unsigned int | n_child_cell_dofs |
std::vector< std::vector< unsigned int > > | level_dof_indices |
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > | parent_child_connect |
std::vector< unsigned int > | n_owned_level_cells |
AlignedVector< VectorizedArray< Number > > | prolongation_matrix_1d |
AlignedVector< VectorizedArray< Number > > | evaluation_data |
std::vector< AlignedVector< VectorizedArray< Number > > > | weights_on_refined |
std::vector< std::vector< std::vector< unsigned short > > > | dirichlet_indices |
Implementation of the MGTransferBase interface for which the transfer operations is implemented in a matrix-free way based on the interpolation matrices of the underlying finite element. This requires considerably less memory than MGTransferPrebuilt and can also be considerably faster than that variant.
This class currently only works for tensor-product finite elements based on FE_Q and FE_DGQ elements, including systems involving multiple components of one of these elements. Systems with different elements or other elements are currently not implemented.
Definition at line 52 of file mg_transfer_matrix_free.h.
MGTransferMatrixFree< dim, Number >::MGTransferMatrixFree | ( | ) |
Constructor without constraint matrices. Use this constructor only with discontinuous finite elements or with no local refinement.
Definition at line 38 of file mg_transfer_matrix_free.cc.
MGTransferMatrixFree< dim, Number >::MGTransferMatrixFree | ( | const MGConstrainedDoFs & | mg_constrained_dofs | ) |
Constructor with constraints. Equivalent to the default constructor followed by initialize_constraints().
Definition at line 49 of file mg_transfer_matrix_free.cc.
|
virtual |
Destructor.
Definition at line 62 of file mg_transfer_matrix_free.cc.
void MGTransferMatrixFree< dim, Number >::initialize_constraints | ( | const MGConstrainedDoFs & | mg_constrained_dofs | ) |
Initialize the constraints to be used in build().
Definition at line 69 of file mg_transfer_matrix_free.cc.
void MGTransferMatrixFree< dim, Number >::clear | ( | ) |
Reset the object to the state it had right after the default constructor.
Definition at line 77 of file mg_transfer_matrix_free.cc.
void MGTransferMatrixFree< dim, Number >::build | ( | const DoFHandler< dim, dim > & | mg_dof | ) |
Actually build the information for the prolongation for each level.
Definition at line 96 of file mg_transfer_matrix_free.cc.
|
virtual |
Prolongate a vector from level to_level-1
to level to_level
using the embedding matrices of the underlying finite element. The previous content of dst
is overwritten.
to_level | The index of the level to prolongate to, which is the level of dst . |
src | is a vector with as many elements as there are degrees of freedom on the coarser level involved. |
dst | has as many elements as there are degrees of freedom on the finer level. |
Implements MGTransferBase< LinearAlgebra::distributed::Vector< Number > >.
Definition at line 154 of file mg_transfer_matrix_free.cc.
|
virtual |
Restrict a vector from level from_level
to level from_level-1
using the transpose operation of the prolongate() method. If the region covered by cells on level from_level
is smaller than that of level from_level-1
(local refinement), then some degrees of freedom in dst
are active and will not be altered. For the other degrees of freedom, the result of the restriction is added.
from_level | The index of the level to restrict from, which is the level of src . |
src | is a vector with as many elements as there are degrees of freedom on the finer level involved. |
dst | has as many elements as there are degrees of freedom on the coarser level. |
Implements MGTransferBase< LinearAlgebra::distributed::Vector< Number > >.
Definition at line 218 of file mg_transfer_matrix_free.cc.
std::size_t MGTransferMatrixFree< dim, Number >::memory_consumption | ( | ) | const |
Memory used by this object.
Definition at line 542 of file mg_transfer_matrix_free.cc.
|
private |
Performs templated prolongation operation
Definition at line 360 of file mg_transfer_matrix_free.cc.
|
private |
Performs templated restriction operation
Definition at line 449 of file mg_transfer_matrix_free.cc.
|
private |
Stores the degree of the finite element contained in the DoFHandler passed to build(). The selection of the computational kernel is based on this number.
Definition at line 144 of file mg_transfer_matrix_free.h.
|
private |
Stores whether the element is continuous and there is a joint degree of freedom in the center of the 1D line.
Definition at line 150 of file mg_transfer_matrix_free.h.
|
private |
Stores the number of components in the finite element contained in the DoFHandler passed to build().
Definition at line 156 of file mg_transfer_matrix_free.h.
|
private |
Stores the number of degrees of freedom on all child cells. It is 2dim*fe.dofs_per_cell
for DG elements and somewhat less for continuous elements.
Definition at line 163 of file mg_transfer_matrix_free.h.
|
private |
Holds the indices for cells on a given level, extracted from DoFHandler for fast access. All DoF indices on a given level are stored as a plain array (since this class assumes constant DoFs per cell). To index into this array, use the cell number times dofs_per_cell.
This array first is arranged such that all locally owned level cells come first (found in the variable n_owned_level_cells) and then other cells necessary for the transfer to the next level.
Definition at line 175 of file mg_transfer_matrix_free.h.
|
private |
Stores the connectivity from parent to child cell numbers for each level.
Definition at line 180 of file mg_transfer_matrix_free.h.
|
private |
Stores the number of cells owned on a given process (sets the bounds for the worker loops) for each level.
Definition at line 186 of file mg_transfer_matrix_free.h.
|
private |
Holds the one-dimensional embedding (prolongation) matrix from mother element to all the children.
Definition at line 192 of file mg_transfer_matrix_free.h.
|
mutableprivate |
Holds the temporary values for the tensor evaluation
Definition at line 197 of file mg_transfer_matrix_free.h.
|
private |
For continuous elements, restriction is not additive and we need to weight the result at the end of prolongation (and at the start of restriction) by the valence of the degrees of freedom, i.e., on how many elements they appear. We store the data in vectorized form to allow for cheap access. Moreover, we utilize the fact that we only need to store 3dim
indices.
Data is organized in terms of each level (outer vector) and the cells on each level (inner vector).
Definition at line 210 of file mg_transfer_matrix_free.h.
|
private |
Stores the local indices of Dirichlet boundary conditions on cells for all levels (outer index), the cells within the levels (second index), and the indices on the cell (inner index).
Definition at line 217 of file mg_transfer_matrix_free.h.