16 #ifndef dealii__mg_transfer_h 17 #define dealii__mg_transfer_h 19 #include <deal.II/base/config.h> 21 #include <deal.II/lac/block_vector.h> 22 #include <deal.II/lac/constraint_matrix.h> 23 #include <deal.II/lac/sparse_matrix.h> 24 #include <deal.II/lac/block_sparsity_pattern.h> 25 #include <deal.II/lac/trilinos_sparse_matrix.h> 26 #include <deal.II/lac/la_parallel_vector.h> 28 #include <deal.II/lac/vector_memory.h> 30 #include <deal.II/multigrid/mg_base.h> 31 #include <deal.II/multigrid/mg_constrained_dofs.h> 32 #include <deal.II/base/mg_level_object.h> 34 #include <deal.II/dofs/dof_handler.h> 36 #include <deal.II/base/std_cxx11/shared_ptr.h> 39 DEAL_II_NAMESPACE_OPEN
44 template <
typename VectorType>
47 typedef ::SparsityPattern Sparsity;
48 typedef ::SparseMatrix<typename VectorType::value_type> Matrix;
50 template <
typename SparsityPatternType,
typename DoFHandlerType>
51 static void reinit(Matrix &matrix, Sparsity &sparsity,
int level,
const SparsityPatternType &sp,
const DoFHandlerType &)
53 sparsity.copy_from (sp);
55 matrix.reinit (sparsity);
59 #ifdef DEAL_II_WITH_TRILINOS 60 template <
typename Number>
61 struct MatrixSelector<
LinearAlgebra::distributed::Vector<Number> >
63 typedef ::TrilinosWrappers::SparsityPattern Sparsity;
64 typedef ::TrilinosWrappers::SparseMatrix Matrix;
66 template <
typename SparsityPatternType,
typename DoFHandlerType>
67 static void reinit(Matrix &matrix, Sparsity &,
int level,
const SparsityPatternType &sp, DoFHandlerType &dh)
69 matrix.reinit(dh.locally_owned_mg_dofs(level+1),
70 dh.locally_owned_mg_dofs(level),
71 sp, MPI_COMM_WORLD,
true);
78 typedef ::TrilinosWrappers::SparsityPattern Sparsity;
79 typedef ::TrilinosWrappers::SparseMatrix Matrix;
81 template <
typename SparsityPatternType,
typename DoFHandlerType>
82 static void reinit(Matrix &matrix, Sparsity &,
int level,
const SparsityPatternType &sp, DoFHandlerType &dh)
84 matrix.reinit(dh.locally_owned_mg_dofs(level+1),
85 dh.locally_owned_mg_dofs(level),
86 sp, MPI_COMM_WORLD,
true);
94 typedef ::TrilinosWrappers::SparsityPattern Sparsity;
95 typedef ::TrilinosWrappers::SparseMatrix Matrix;
97 template <
typename SparsityPatternType,
typename DoFHandlerType>
98 static void reinit(Matrix &, Sparsity &,
int ,
const SparsityPatternType &, DoFHandlerType &)
104 template <
typename Number>
105 struct MatrixSelector<
LinearAlgebra::distributed::Vector<Number> >
107 typedef ::SparsityPattern Sparsity;
108 typedef ::SparseMatrix<Number> Matrix;
110 template <
typename SparsityPatternType,
typename DoFHandlerType>
111 static void reinit(Matrix &, Sparsity &,
int,
const SparsityPatternType &,
const DoFHandlerType &)
114 "ERROR: MGTransferPrebuilt with LinearAlgebra::distributed::Vector currently " 115 "needs deal.II to be configured with Trilinos."));
138 template <
typename VectorType>
152 template <
int dim,
class InVector,
int spacedim>
156 const InVector &src)
const;
165 template <
int dim,
class OutVector,
int spacedim>
176 template <
int dim,
class OutVector,
int spacedim>
215 template <
int dim,
int spacedim>
221 std::vector<types::global_dof_index>
sizes;
230 std::vector<std::vector<std::pair<types::global_dof_index, types::global_dof_index> > >
240 std::vector<std::vector<std::pair<types::global_dof_index, types::global_dof_index> > >
250 std::vector<std::vector<std::pair<types::global_dof_index, types::global_dof_index> > >
286 template <
typename Number>
300 template <
int dim,
typename Number2,
int spacedim>
313 template <
int dim,
typename Number2,
int spacedim>
324 template <
int dim,
typename Number2,
int spacedim>
363 template <
int dim,
int spacedim>
369 std::vector<types::global_dof_index>
sizes;
378 std::vector<std::vector<std::pair<unsigned int, unsigned int> > >
388 std::vector<std::vector<std::pair<unsigned int, unsigned int> > >
398 std::vector<std::vector<std::pair<unsigned int, unsigned int> > >
449 template <
typename VectorType>
500 template <
int dim,
int spacedim>
514 virtual void prolongate (
const unsigned int to_level,
516 const VectorType &src)
const;
535 const VectorType &src)
const;
582 DEAL_II_NAMESPACE_CLOSE
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
void copy_from_mg(const DoFHandler< dim, spacedim > &mg_dof, OutVector &dst, const MGLevelObject< VectorType > &src) const
std::size_t memory_consumption() const
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > copy_indices
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > copy_indices_global_mine
void build_matrices(const DoFHandler< dim, spacedim > &mg_dof)
void copy_from_mg_add(const DoFHandler< dim, spacedim > &mg_dof, OutVector &dst, const MGLevelObject< VectorType > &src) const
void copy_to_mg(const DoFHandler< dim, spacedim > &mg_dof, MGLevelObject< VectorType > &dst, const InVector &src) const
std::vector< std_cxx11::shared_ptr< typename internal::MatrixSelector< VectorType >::Matrix > > prolongation_matrices
std::vector< types::global_dof_index > sizes
#define AssertThrow(cond, exc)
std::vector< unsigned int > component_to_block_map
void print_indices(std::ostream &os) const
LinearAlgebra::distributed::Vector< Number > ghosted_global_vector
static ::ExceptionBase & ExcNoProlongation()
virtual ~MGTransferPrebuilt()
std::vector< std_cxx11::shared_ptr< typename internal::MatrixSelector< VectorType >::Sparsity > > prolongation_sparsities
void print_matrices(std::ostream &os) const
std::size_t memory_consumption() const
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > copy_indices_level_mine
void fill_and_communicate_copy_indices(const DoFHandler< dim, spacedim > &mg_dof)
#define DeclException0(Exception0)
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs
std::vector< unsigned int > component_to_block_map
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< VectorType > > mg_constrained_dofs
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > copy_indices
std::vector< types::global_dof_index > sizes
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > copy_indices_global_mine
static ::ExceptionBase & ExcMatricesNotBuilt()
virtual void prolongate(const unsigned int to_level, VectorType &dst, const VectorType &src) const
static ::ExceptionBase & ExcNotImplemented()
std::vector< std::vector< bool > > interface_dofs
void set_component_to_block_map(const std::vector< unsigned int > &map)
MGLevelObject< LinearAlgebra::distributed::Vector< Number > > ghosted_level_vector
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > copy_indices_level_mine
virtual void restrict_and_add(const unsigned int from_level, VectorType &dst, const VectorType &src) const