Reference documentation for deal.II version 9.4.1
|
#include <deal.II/multigrid/multigrid.h>
Public Member Functions | |
PreconditionMG (const DoFHandler< dim > &dof_handler, Multigrid< VectorType > &mg, const TRANSFER &transfer) | |
PreconditionMG (const std::vector< const DoFHandler< dim > * > &dof_handler, Multigrid< VectorType > &mg, const TRANSFER &transfer) | |
bool | empty () const |
template<class OtherVectorType > | |
void | vmult (OtherVectorType &dst, const OtherVectorType &src) const |
template<class OtherVectorType > | |
void | vmult_add (OtherVectorType &dst, const OtherVectorType &src) const |
template<class OtherVectorType > | |
void | Tvmult (OtherVectorType &dst, const OtherVectorType &src) const |
template<class OtherVectorType > | |
void | Tvmult_add (OtherVectorType &dst, const OtherVectorType &src) const |
IndexSet | locally_owned_range_indices (const unsigned int block=0) const |
IndexSet | locally_owned_domain_indices (const unsigned int block=0) const |
MPI_Comm | get_mpi_communicator () const |
boost::signals2::connection | connect_transfer_to_mg (const std::function< void(bool)> &slot) |
boost::signals2::connection | connect_transfer_to_global (const std::function< void(bool)> &slot) |
Multigrid< VectorType > & | get_multigrid () |
const Multigrid< VectorType > & | get_multigrid () const |
Private Attributes | |
std::vector< SmartPointer< const DoFHandler< dim >, PreconditionMG< dim, VectorType, TRANSFER > > > | dof_handler_vector |
std::vector< const DoFHandler< dim > * > | dof_handler_vector_raw |
SmartPointer< Multigrid< VectorType >, PreconditionMG< dim, VectorType, TRANSFER > > | multigrid |
SmartPointer< const TRANSFER, PreconditionMG< dim, VectorType, TRANSFER > > | transfer |
const bool | uses_dof_handler_vector |
mg::Signals | signals |
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) |
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) |
void | check_no_subscribers () const noexcept |
Multi-level preconditioner. Here, we collect all information needed for multi-level preconditioning and provide the standard interface for LAC iterative methods.
Furthermore, it needs functions void copy_to_mg(const VectorType&)
to store src
in the right hand side of the multi-level method and void copy_from_mg(VectorType&)
to store the result of the v-cycle in dst
.
If VectorType is in fact a block vector and the TRANSFER object supports use of a separate DoFHandler for each block, this class also allows to be initialized with a separate DoFHandler for each block.
Definition at line 499 of file multigrid.h.
PreconditionMG< dim, VectorType, TRANSFER >::PreconditionMG | ( | const DoFHandler< dim > & | dof_handler, |
Multigrid< VectorType > & | mg, | ||
const TRANSFER & | transfer | ||
) |
Constructor. Arguments are the multigrid object, pre-smoother, post- smoother and coarse grid solver.
PreconditionMG< dim, VectorType, TRANSFER >::PreconditionMG | ( | const std::vector< const DoFHandler< dim > * > & | dof_handler, |
Multigrid< VectorType > & | mg, | ||
const TRANSFER & | transfer | ||
) |
Same as above in case every component of a block vector uses its own DoFHandler.
bool PreconditionMG< dim, VectorType, TRANSFER >::empty | ( | ) | const |
Dummy function needed by other classes.
void PreconditionMG< dim, VectorType, TRANSFER >::vmult | ( | OtherVectorType & | dst, |
const OtherVectorType & | src | ||
) | const |
Preconditioning operator. Calls the vcycle
function of the MG
object passed to the constructor.
This is the operator used by LAC iterative solvers.
void PreconditionMG< dim, VectorType, TRANSFER >::vmult_add | ( | OtherVectorType & | dst, |
const OtherVectorType & | src | ||
) | const |
Preconditioning operator. Calls the vcycle
function of the MG
object passed to the constructor.
void PreconditionMG< dim, VectorType, TRANSFER >::Tvmult | ( | OtherVectorType & | dst, |
const OtherVectorType & | src | ||
) | const |
Transposed preconditioning operator.
Not implemented, but the definition may be needed.
void PreconditionMG< dim, VectorType, TRANSFER >::Tvmult_add | ( | OtherVectorType & | dst, |
const OtherVectorType & | src | ||
) | const |
Transposed preconditioning operator.
Not implemented, but the definition may be needed.
IndexSet PreconditionMG< dim, VectorType, TRANSFER >::locally_owned_range_indices | ( | const unsigned int | block = 0 | ) | const |
Return the partitioning of the range space of this preconditioner, i.e., the partitioning of the vectors that are result from matrix-vector products. By default, the respective information for the first DoFHandler object are returned.
IndexSet PreconditionMG< dim, VectorType, TRANSFER >::locally_owned_domain_indices | ( | const unsigned int | block = 0 | ) | const |
Return the partitioning of the domain space of this preconditioner, i.e., the partitioning of the vectors this matrix has to be multiplied with. By default, the respective information for the first DoFHandler object are returned.
MPI_Comm PreconditionMG< dim, VectorType, TRANSFER >::get_mpi_communicator | ( | ) | const |
Return the MPI communicator object in use with this preconditioner.
boost::signals2::connection PreconditionMG< dim, VectorType, TRANSFER >::connect_transfer_to_mg | ( | const std::function< void(bool)> & | slot | ) |
Connect a function to mg::Signals::transfer_to_mg.
boost::signals2::connection PreconditionMG< dim, VectorType, TRANSFER >::connect_transfer_to_global | ( | const std::function< void(bool)> & | slot | ) |
Connect a function to mg::Signals::transfer_to_global.
Multigrid< VectorType > & PreconditionMG< dim, VectorType, TRANSFER >::get_multigrid | ( | ) |
Return the Multigrid object passed to the constructor.
const Multigrid< VectorType > & PreconditionMG< dim, VectorType, TRANSFER >::get_multigrid | ( | ) | const |
Return the Multigrid object passed to the constructor.
|
private |
Associated DoFHandler
.
Definition at line 614 of file multigrid.h.
|
private |
Storage for the pointers to the DoFHandler objects without SmartPointer wrapper.
Definition at line 620 of file multigrid.h.
|
private |
The multigrid object.
Definition at line 626 of file multigrid.h.
|
private |
Object for grid transfer.
Definition at line 632 of file multigrid.h.
|
private |
Flag to indicate if the object is initialized with a single DoFHandler or with one for each block.
Definition at line 638 of file multigrid.h.
|
private |
Signals used by this object
Definition at line 643 of file multigrid.h.