Reference documentation for deal.II version 9.1.1
|
#include <deal.II/multigrid/mg_transfer.h>
Public Member Functions | |
MGTransferPrebuilt ()=default | |
MGTransferPrebuilt (const MGConstrainedDoFs &mg_constrained_dofs) | |
MGTransferPrebuilt (const AffineConstraints< double > &constraints, const MGConstrainedDoFs &mg_constrained_dofs) | |
virtual | ~MGTransferPrebuilt () override=default |
void | initialize_constraints (const MGConstrainedDoFs &mg_constrained_dofs) |
void | initialize_constraints (const AffineConstraints< double > &constraints, const MGConstrainedDoFs &mg_constrained_dofs) |
void | clear () |
template<int dim, int spacedim> | |
void | build_matrices (const DoFHandler< dim, spacedim > &mg_dof) |
virtual void | prolongate (const unsigned int to_level, VectorType &dst, const VectorType &src) const override |
virtual void | restrict_and_add (const unsigned int from_level, VectorType &dst, const VectorType &src) const override |
std::size_t | memory_consumption () const |
void | print_matrices (std::ostream &os) const |
Public Member Functions inherited from MGLevelGlobalTransfer< VectorType > | |
void | clear () |
template<int dim, class InVector , int spacedim> | |
void | copy_to_mg (const DoFHandler< dim, spacedim > &mg_dof, MGLevelObject< VectorType > &dst, const InVector &src) const |
template<int dim, class OutVector , int spacedim> | |
void | copy_from_mg (const DoFHandler< dim, spacedim > &mg_dof, OutVector &dst, const MGLevelObject< VectorType > &src) const |
template<int dim, class OutVector , int spacedim> | |
void | copy_from_mg_add (const DoFHandler< dim, spacedim > &mg_dof, OutVector &dst, const MGLevelObject< VectorType > &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< VectorType > | |
virtual | ~MGTransferBase () override=default |
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 (std::atomic< bool > *const validity, const std::string &identifier="") const |
void | unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const |
unsigned int | n_subscriptions () const |
template<typename StreamType > | |
void | list_subscribers (StreamType &stream) const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Static Public Member Functions | |
static ::ExceptionBase & | ExcNoProlongation () |
static ::ExceptionBase & | ExcMatricesNotBuilt () |
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) |
Private Attributes | |
std::vector< std::shared_ptr< typename internal::MatrixSelector< VectorType >::Sparsity > > | prolongation_sparsities |
std::vector< std::shared_ptr< typename internal::MatrixSelector< VectorType >::Matrix > > | prolongation_matrices |
std::vector< std::vector< bool > > | interface_dofs |
Additional Inherited Members | |
Protected Member Functions inherited from MGLevelGlobalTransfer< VectorType > | |
template<int dim, int spacedim> | |
void | fill_and_communicate_copy_indices (const DoFHandler< dim, spacedim > &mg_dof) |
Protected Attributes inherited from MGLevelGlobalTransfer< VectorType > | |
std::vector< types::global_dof_index > | sizes |
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > | copy_indices |
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > | copy_indices_global_mine |
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > | copy_indices_level_mine |
bool | perform_plain_copy |
std::vector< unsigned int > | component_to_block_map |
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< VectorType > > | mg_constrained_dofs |
Implementation of the MGTransferBase interface for which the transfer operations are prebuilt upon construction of the object of this class as matrices. This is the fast way, since it only needs to build the operation once by looping over all cells and storing the result in a matrix for each level, but requires additional memory.
See MGTransferBase to find out which of the transfer classes is best for your needs.
Definition at line 661 of file mg_transfer.h.
|
default |
Constructor without constraint matrices. Use this constructor only with discontinuous finite elements or with no local refinement.
MGTransferPrebuilt< VectorType >::MGTransferPrebuilt | ( | 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_prebuilt.cc.
MGTransferPrebuilt< VectorType >::MGTransferPrebuilt | ( | const AffineConstraints< double > & | constraints, |
const MGConstrainedDoFs & | mg_constrained_dofs | ||
) |
Constructor with constraints. Equivalent to the default constructor followed by initialize_constraints().
constraints
is unused. Definition at line 58 of file mg_transfer_prebuilt.cc.
|
overridevirtualdefault |
Destructor.
void MGTransferPrebuilt< VectorType >::initialize_constraints | ( | const MGConstrainedDoFs & | mg_constrained_dofs | ) |
Initialize the constraints to be used in build_matrices().
Definition at line 69 of file mg_transfer_prebuilt.cc.
void MGTransferPrebuilt< VectorType >::initialize_constraints | ( | const AffineConstraints< double > & | constraints, |
const MGConstrainedDoFs & | mg_constrained_dofs | ||
) |
Initialize the constraints to be used in build_matrices().
constraints
is unused. Definition at line 79 of file mg_transfer_prebuilt.cc.
void MGTransferPrebuilt< VectorType >::clear | ( | ) |
Reset the object to the state it had right after the default constructor.
Definition at line 90 of file mg_transfer_prebuilt.cc.
void MGTransferPrebuilt< VectorType >::build_matrices | ( | const DoFHandler< dim, spacedim > & | mg_dof | ) |
Actually build the prolongation matrices for each level.
Definition at line 160 of file mg_transfer_prebuilt.cc.
|
overridevirtual |
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.
Implements MGTransferBase< VectorType >.
Definition at line 102 of file mg_transfer_prebuilt.cc.
|
overridevirtual |
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.
Implements MGTransferBase< VectorType >.
Definition at line 116 of file mg_transfer_prebuilt.cc.
std::size_t MGTransferPrebuilt< VectorType >::memory_consumption | ( | ) | const |
Memory used by this object.
Definition at line 395 of file mg_transfer_prebuilt.cc.
void MGTransferPrebuilt< VectorType >::print_matrices | ( | std::ostream & | os | ) | const |
Print all the matrices for debugging purposes.
Definition at line 381 of file mg_transfer_prebuilt.cc.
|
private |
Sparsity patterns for transfer matrices.
Definition at line 784 of file mg_transfer.h.
|
private |
The actual prolongation matrix. column indices belong to the dof indices of the mother cell, i.e. the coarse level. while row indices belong to the child cell, i.e. the fine level.
Definition at line 793 of file mg_transfer.h.
|
private |
Degrees of freedom on the refinement edge excluding those on the boundary.
Definition at line 799 of file mg_transfer.h.