Reference documentation for deal.II version 9.0.0
|
#include <deal.II/multigrid/mg_transfer_block.h>
Public Member Functions | |
MGTransferBlock () | |
virtual | ~MGTransferBlock () |
void | initialize (const std::vector< number > &factors, VectorMemory< Vector< number > > &memory) |
template<int dim, int spacedim> | |
void | build_matrices (const DoFHandler< dim, spacedim > &dof, const DoFHandler< dim, spacedim > &mg_dof, const std::vector< bool > &selected) |
virtual void | prolongate (const unsigned int to_level, BlockVector< number > &dst, const BlockVector< number > &src) const |
virtual void | restrict_and_add (const unsigned int from_level, BlockVector< number > &dst, const BlockVector< number > &src) const |
template<int dim, typename number2 , int spacedim> | |
void | copy_to_mg (const DoFHandler< dim, spacedim > &mg_dof, MGLevelObject< BlockVector< number > > &dst, const BlockVector< number2 > &src) const |
template<int dim, typename number2 , int spacedim> | |
void | copy_from_mg (const DoFHandler< dim, spacedim > &mg_dof, BlockVector< number2 > &dst, const MGLevelObject< BlockVector< number > > &src) const |
template<int dim, typename number2 , int spacedim> | |
void | copy_from_mg_add (const DoFHandler< dim, spacedim > &mg_dof, BlockVector< number2 > &dst, const MGLevelObject< BlockVector< number > > &src) const |
std::size_t | memory_consumption () const |
Public Member Functions inherited from MGTransferBase< BlockVector< number > > | |
virtual | ~MGTransferBase ()=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 (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) |
Private Attributes | |
std::vector< number > | factors |
SmartPointer< VectorMemory< Vector< number > >, MGTransferBlock< number > > | memory |
Private Attributes inherited from MGTransferBlockBase | |
std::vector< bool > | selected |
unsigned int | n_mg_blocks |
std::vector< unsigned int > | mg_block |
std::vector< std::vector< types::global_dof_index > > | sizes |
std::vector< types::global_dof_index > | block_start |
std::vector< std::vector< types::global_dof_index > > | mg_block_start |
std::vector< std::shared_ptr< BlockSparseMatrix< double > > > | prolongation_matrices |
std::vector< std::vector< std::vector< std::pair< unsigned int, unsigned int > > > > | copy_indices |
SmartPointer< const MGConstrainedDoFs, MGTransferBlockBase > | mg_constrained_dofs |
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) |
Private Member Functions inherited from MGTransferBlockBase | |
MGTransferBlockBase () | |
MGTransferBlockBase (const MGConstrainedDoFs &mg_constrained_dofs) | |
MGTransferBlockBase (const ConstraintMatrix &constraints, const MGConstrainedDoFs &mg_constrained_dofs) | |
std::size_t | memory_consumption () const |
template<int dim, int spacedim> | |
void | build_matrices (const DoFHandler< dim, spacedim > &dof, const DoFHandler< dim, spacedim > &mg_dof) |
Static Private Member Functions inherited from MGTransferBlockBase | |
static ::ExceptionBase & | ExcMatricesNotBuilt () |
Implementation of the MGTransferBase interface for block matrices and block vectors.
In addition to the functionality of MGTransferPrebuilt, the operation may be restricted to certain blocks of the vector.
If the restricted mode is chosen, block vectors used in the transfer routines may only have as many blocks as there are trues
in the selected-field.
See MGTransferBase to find out which of the transfer classes is best for your needs.
Definition at line 191 of file mg_transfer_block.h.
MGTransferBlock< number >::MGTransferBlock | ( | ) |
Default constructor.
Definition at line 62 of file multigrid.cc.
|
virtual |
Destructor.
Definition at line 69 of file multigrid.cc.
void MGTransferBlock< number >::initialize | ( | const std::vector< number > & | factors, |
VectorMemory< Vector< number > > & | memory | ||
) |
Initialize additional factors and memory if the restriction of the blocks is to be weighted differently.
Definition at line 77 of file multigrid.cc.
void MGTransferBlock< number >::build_matrices | ( | const DoFHandler< dim, spacedim > & | dof, |
const DoFHandler< dim, spacedim > & | mg_dof, | ||
const std::vector< bool > & | selected | ||
) |
Build the prolongation matrices for each level.
This function is a front-end for the same function in MGTransferBlockBase.
Definition at line 536 of file mg_transfer_block.cc.
|
virtual |
Prolongate a vector from level to_level-1
to level to_level
. The previous content of dst
is overwritten.
Implements MGTransferBase< BlockVector< number > >.
Definition at line 86 of file multigrid.cc.
|
virtual |
Restrict a vector from level from_level
to level from_level-1
and add this restriction to dst
. 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< BlockVector< number > >.
Definition at line 111 of file multigrid.cc.
void MGTransferBlock< number >::copy_to_mg | ( | const DoFHandler< dim, spacedim > & | mg_dof, |
MGLevelObject< BlockVector< number > > & | dst, | ||
const BlockVector< number2 > & | src | ||
) | const |
Transfer from a vector on the global grid to a multilevel vector for the active degrees of freedom. In particular, for a globally refined mesh only the finest level in dst
is filled as a plain copy of src
. All the other level objects are left untouched.
The action for discontinuous elements is as follows: on an active mesh cell, the global vector entries are simply copied to the corresponding entries of the level vector. Then, these values are restricted down to the coarsest level.
Definition at line 177 of file mg_transfer_block.cc.
void MGTransferBlock< number >::copy_from_mg | ( | const DoFHandler< dim, spacedim > & | mg_dof, |
BlockVector< number2 > & | dst, | ||
const MGLevelObject< BlockVector< number > > & | src | ||
) | const |
Transfer from multi-level vector to normal vector.
Copies data from active portions of a multilevel vector into the respective positions of a global vector.
void MGTransferBlock< number >::copy_from_mg_add | ( | const DoFHandler< dim, spacedim > & | mg_dof, |
BlockVector< number2 > & | dst, | ||
const MGLevelObject< BlockVector< number > > & | src | ||
) | const |
Add a multi-level vector to a normal vector.
Works as the previous function, but probably not for continuous elements.
std::size_t MGTransferBlockBase::memory_consumption |
Memory used by this object.
Definition at line 176 of file multigrid.cc.
|
private |
Optional multiplication factors for each block. Requires initialization of memory.
Definition at line 278 of file mg_transfer_block.h.
|
private |
Memory pool required if additional multiplication using factors is desired.
Definition at line 284 of file mg_transfer_block.h.