Reference documentation for deal.II version 9.4.1
|
#include <deal.II/multigrid/mg_transfer_global_coarsening.h>
Public Types | |
using | Number = typename VectorType::value_type |
Public Member Functions | |
MGTransferGlobalCoarsening (const MGLevelObject< MGTwoLevelTransfer< dim, VectorType > > &transfer, const std::function< void(const unsigned int, VectorType &)> &initialize_dof_vector={}) | |
void | build (const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &external_partitioners={}) |
void | build (const std::function< void(const unsigned int, VectorType &)> &initialize_dof_vector) |
void | prolongate (const unsigned int to_level, VectorType &dst, const VectorType &src) const override |
void | prolongate_and_add (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 |
template<class InVector , int spacedim> | |
void | copy_to_mg (const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< VectorType > &dst, const InVector &src) const |
template<class OutVector , int spacedim> | |
void | copy_from_mg (const DoFHandler< dim, spacedim > &dof_handler, OutVector &dst, const MGLevelObject< VectorType > &src) const |
template<class InVector > | |
void | interpolate_to_mg (MGLevelObject< VectorType > &dst, const InVector &src) const |
template<class InVector , int spacedim> | |
void | interpolate_to_mg (const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< VectorType > &dst, const InVector &src) const |
std::size_t | memory_consumption () const |
Private Member Functions | |
void | initialize_dof_vector (const unsigned int level, VectorType &vector) const |
Private Attributes | |
MGLevelObject< MGTwoLevelTransfer< dim, VectorType > > | transfer |
std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > | external_partitioners |
Subscriptor functionality | |
Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class. | |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
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) |
void | check_no_subscribers () const noexcept |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
static std::mutex | mutex |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Implementation of the MGTransferBase. In contrast to other multigrid transfer operators, the user can provide separate transfer operators of type MGTwoLevelTransfer between each level.
This class currently only works for the tensor-product finite elements FE_Q and FE_DGQ and simplex elements FE_SimplexP and FE_SimplexDGP as well as for systems involving multiple components of one of these elements. Other elements are currently not implemented.
Definition at line 475 of file mg_transfer_global_coarsening.h.
using MGTransferGlobalCoarsening< dim, VectorType >::Number = typename VectorType::value_type |
Value type.
Definition at line 481 of file mg_transfer_global_coarsening.h.
MGTransferGlobalCoarsening< dim, VectorType >::MGTransferGlobalCoarsening | ( | const MGLevelObject< MGTwoLevelTransfer< dim, VectorType > > & | transfer, |
const std::function< void(const unsigned int, VectorType &)> & | initialize_dof_vector = {} |
||
) |
Constructor taking a collection of transfer operators (with the coarsest level kept empty in transfer
) and an optional function that initializes the internal level vectors within the function call copy_to_mg() if used in the context of PreconditionMG.
void MGTransferGlobalCoarsening< dim, VectorType >::build | ( | const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > & | external_partitioners = {} | ) |
Similar function to MGTransferMatrixFree::build() with the difference that the information for the prolongation for each level has already been built. So this function only tries to optimize the data structures of the two-level transfer operators, e.g., by enabling inplace vector operations, by checking if external_partitioners
and the internal ones are compatible.
void MGTransferGlobalCoarsening< dim, VectorType >::build | ( | const std::function< void(const unsigned int, VectorType &)> & | initialize_dof_vector | ) |
Same as above but taking a lambda for initializing vector instead of partitioners.
|
overridevirtual |
Perform prolongation.
Implements MGTransferBase< VectorType >.
|
overridevirtual |
Perform prolongation.
Reimplemented from MGTransferBase< VectorType >.
|
overridevirtual |
Perform restriction.
Implements MGTransferBase< VectorType >.
void MGTransferGlobalCoarsening< dim, VectorType >::copy_to_mg | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
MGLevelObject< VectorType > & | dst, | ||
const InVector & | src | ||
) | const |
Initialize internal vectors and copy src
vector to the finest multigrid level.
void MGTransferGlobalCoarsening< dim, VectorType >::copy_from_mg | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
OutVector & | dst, | ||
const MGLevelObject< VectorType > & | src | ||
) | const |
Initialize internal vectors and copy the values on the finest multigrid level to dst
vector.
void MGTransferGlobalCoarsening< dim, VectorType >::interpolate_to_mg | ( | MGLevelObject< VectorType > & | dst, |
const InVector & | src | ||
) | const |
Interpolate fine-mesh field src
to each multigrid level in dof_handler
and store the result in dst
. This function is different from restriction, where a weighted residual is transferred to a coarser level (transposition of prolongation matrix).
The argument dst
has to be initialized with the correct size according to the number of levels of the triangulation.
If an inner vector of dst
is empty or has incorrect locally owned size, it will be resized to locally relevant degrees of freedom on each level.
void MGTransferGlobalCoarsening< dim, VectorType >::interpolate_to_mg | ( | const DoFHandler< dim, spacedim > & | dof_handler, |
MGLevelObject< VectorType > & | dst, | ||
const InVector & | src | ||
) | const |
Like the above function but with a user-provided DoFHandler as additional argument. However, this DoFHandler is not used internally, but is required to be able to use MGTransferGlobalCoarsening and MGTransferMatrixFree as template argument.
std::size_t MGTransferGlobalCoarsening< dim, VectorType >::memory_consumption | ( | ) | const |
Return the memory consumption of the allocated memory in this class.
|
private |
Function to initialize internal level vectors.
|
private |
Collection of the two-level transfer operators.
Definition at line 610 of file mg_transfer_global_coarsening.h.
|
private |
External partitioners used during initialize_dof_vector().
Definition at line 616 of file mg_transfer_global_coarsening.h.