Reference documentation for deal.II version 8.5.1
mg_transfer.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2001 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__mg_transfer_h
17 #define dealii__mg_transfer_h
18 
19 #include <deal.II/base/config.h>
20 
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>
27 
28 #include <deal.II/lac/vector_memory.h>
29 
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>
33 
34 #include <deal.II/dofs/dof_handler.h>
35 
36 #include <deal.II/base/std_cxx11/shared_ptr.h>
37 
38 
39 DEAL_II_NAMESPACE_OPEN
40 
41 
42 namespace internal
43 {
44  template <typename VectorType>
45  struct MatrixSelector
46  {
47  typedef ::SparsityPattern Sparsity;
48  typedef ::SparseMatrix<typename VectorType::value_type> Matrix;
49 
50  template <typename SparsityPatternType, typename DoFHandlerType>
51  static void reinit(Matrix &matrix, Sparsity &sparsity, int level, const SparsityPatternType &sp, const DoFHandlerType &)
52  {
53  sparsity.copy_from (sp);
54  (void)level;
55  matrix.reinit (sparsity);
56  }
57  };
58 
59 #ifdef DEAL_II_WITH_TRILINOS
60  template <typename Number>
61  struct MatrixSelector<LinearAlgebra::distributed::Vector<Number> >
62  {
63  typedef ::TrilinosWrappers::SparsityPattern Sparsity;
64  typedef ::TrilinosWrappers::SparseMatrix Matrix;
65 
66  template <typename SparsityPatternType, typename DoFHandlerType>
67  static void reinit(Matrix &matrix, Sparsity &, int level, const SparsityPatternType &sp, DoFHandlerType &dh)
68  {
69  matrix.reinit(dh.locally_owned_mg_dofs(level+1),
70  dh.locally_owned_mg_dofs(level),
71  sp, MPI_COMM_WORLD, true);
72  }
73 
74  };
75  template <>
76  struct MatrixSelector<::TrilinosWrappers::MPI::Vector>
77  {
78  typedef ::TrilinosWrappers::SparsityPattern Sparsity;
79  typedef ::TrilinosWrappers::SparseMatrix Matrix;
80 
81  template <typename SparsityPatternType, typename DoFHandlerType>
82  static void reinit(Matrix &matrix, Sparsity &, int level, const SparsityPatternType &sp, DoFHandlerType &dh)
83  {
84  matrix.reinit(dh.locally_owned_mg_dofs(level+1),
85  dh.locally_owned_mg_dofs(level),
86  sp, MPI_COMM_WORLD, true);
87  }
88 
89  };
90 
91  template <>
92  struct MatrixSelector<::TrilinosWrappers::Vector>
93  {
94  typedef ::TrilinosWrappers::SparsityPattern Sparsity;
95  typedef ::TrilinosWrappers::SparseMatrix Matrix;
96 
97  template <typename SparsityPatternType, typename DoFHandlerType>
98  static void reinit(Matrix &, Sparsity &, int /*level*/, const SparsityPatternType &, DoFHandlerType &)
99  {
100  }
101  };
102 #else
103  // ! DEAL_II_WITH_TRILINOS
104  template <typename Number>
105  struct MatrixSelector<LinearAlgebra::distributed::Vector<Number> >
106  {
107  typedef ::SparsityPattern Sparsity;
108  typedef ::SparseMatrix<Number> Matrix;
109 
110  template <typename SparsityPatternType, typename DoFHandlerType>
111  static void reinit(Matrix &, Sparsity &, int, const SparsityPatternType &, const DoFHandlerType &)
112  {
114  "ERROR: MGTransferPrebuilt with LinearAlgebra::distributed::Vector currently "
115  "needs deal.II to be configured with Trilinos."));
116  }
117  };
118 
119 #endif
120 }
121 
122 /*
123  * MGTransferBase is defined in mg_base.h
124  */
125 
128 
129 
130 
138 template <typename VectorType>
139 class MGLevelGlobalTransfer : public MGTransferBase<VectorType>
140 {
141 public:
142 
146  void clear ();
147 
152  template <int dim, class InVector, int spacedim>
153  void
154  copy_to_mg (const DoFHandler<dim,spacedim> &mg_dof,
156  const InVector &src) const;
157 
165  template <int dim, class OutVector, int spacedim>
166  void
167  copy_from_mg (const DoFHandler<dim,spacedim> &mg_dof,
168  OutVector &dst,
169  const MGLevelObject<VectorType> &src) const;
170 
176  template <int dim, class OutVector, int spacedim>
177  void
179  OutVector &dst,
180  const MGLevelObject<VectorType> &src) const;
181 
197  void
198  set_component_to_block_map (const std::vector<unsigned int> &map);
199 
203  std::size_t memory_consumption () const;
204 
208  void print_indices(std::ostream &os) const;
209 
210 protected:
211 
215  template <int dim, int spacedim>
217 
221  std::vector<types::global_dof_index> sizes;
222 
230  std::vector<std::vector<std::pair<types::global_dof_index, types::global_dof_index> > >
232 
240  std::vector<std::vector<std::pair<types::global_dof_index, types::global_dof_index> > >
242 
250  std::vector<std::vector<std::pair<types::global_dof_index, types::global_dof_index> > >
252 
260 
265  std::vector<unsigned int> component_to_block_map;
266 
271 };
272 
273 
274 
286 template <typename Number>
287 class MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number> > : public MGTransferBase<LinearAlgebra::distributed::Vector<Number> >
288 {
289 public:
290 
294  void clear ();
295 
300  template <int dim, typename Number2, int spacedim>
301  void
302  copy_to_mg (const DoFHandler<dim,spacedim> &mg_dof,
305 
313  template <int dim, typename Number2, int spacedim>
314  void
315  copy_from_mg (const DoFHandler<dim,spacedim> &mg_dof,
318 
324  template <int dim, typename Number2, int spacedim>
325  void
329 
345  void
346  set_component_to_block_map (const std::vector<unsigned int> &map);
347 
351  std::size_t memory_consumption () const;
352 
356  void print_indices(std::ostream &os) const;
357 
358 protected:
359 
363  template <int dim, int spacedim>
365 
369  std::vector<types::global_dof_index> sizes;
370 
378  std::vector<std::vector<std::pair<unsigned int, unsigned int> > >
380 
388  std::vector<std::vector<std::pair<unsigned int, unsigned int> > >
390 
398  std::vector<std::vector<std::pair<unsigned int, unsigned int> > >
400 
408 
413  std::vector<unsigned int> component_to_block_map;
414 
419 
426 
432 };
433 
434 
435 
449 template <typename VectorType>
450 class MGTransferPrebuilt : public MGLevelGlobalTransfer<VectorType>
451 {
452 public:
458 
464 
471  MGTransferPrebuilt (const ConstraintMatrix &constraints,
472  const MGConstrainedDoFs &mg_constrained_dofs) DEAL_II_DEPRECATED;
473 
477  virtual ~MGTransferPrebuilt ();
478 
483 
489  void initialize_constraints (const ConstraintMatrix &constraints,
490  const MGConstrainedDoFs &mg_constrained_dofs) DEAL_II_DEPRECATED;
491 
495  void clear ();
496 
500  template <int dim, int spacedim>
501  void build_matrices (const DoFHandler<dim,spacedim> &mg_dof);
502 
514  virtual void prolongate (const unsigned int to_level,
515  VectorType &dst,
516  const VectorType &src) const;
517 
533  virtual void restrict_and_add (const unsigned int from_level,
534  VectorType &dst,
535  const VectorType &src) const;
536 
541 
546 
550  std::size_t memory_consumption () const;
551 
555  void print_matrices(std::ostream &os) const;
556 
557 private:
558 
562  std::vector<std_cxx11::shared_ptr<typename internal::MatrixSelector<VectorType>::Sparsity> > prolongation_sparsities;
563 
569  std::vector<std_cxx11::shared_ptr<typename internal::MatrixSelector<VectorType>::Matrix> > prolongation_matrices;
570 
575  std::vector<std::vector<bool> > interface_dofs;
576 };
577 
578 
582 DEAL_II_NAMESPACE_CLOSE
583 
584 #endif
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
Definition: mg_transfer.h:231
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > copy_indices_global_mine
Definition: mg_transfer.h:389
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
Definition: mg_transfer.h:569
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
std::vector< unsigned int > component_to_block_map
Definition: mg_transfer.h:265
void print_indices(std::ostream &os) const
LinearAlgebra::distributed::Vector< Number > ghosted_global_vector
Definition: mg_transfer.h:425
static ::ExceptionBase & ExcNoProlongation()
std::vector< std_cxx11::shared_ptr< typename internal::MatrixSelector< VectorType >::Sparsity > > prolongation_sparsities
Definition: mg_transfer.h:562
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
Definition: mg_transfer.h:251
void fill_and_communicate_copy_indices(const DoFHandler< dim, spacedim > &mg_dof)
#define DeclException0(Exception0)
Definition: exceptions.h:541
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs
Definition: mg_transfer.h:418
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< VectorType > > mg_constrained_dofs
Definition: mg_transfer.h:270
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > copy_indices
Definition: mg_transfer.h:379
std::vector< types::global_dof_index > sizes
Definition: mg_transfer.h:221
std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > copy_indices_global_mine
Definition: mg_transfer.h:241
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
Definition: mg_transfer.h:575
void set_component_to_block_map(const std::vector< unsigned int > &map)
MGLevelObject< LinearAlgebra::distributed::Vector< Number > > ghosted_level_vector
Definition: mg_transfer.h:431
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > copy_indices_level_mine
Definition: mg_transfer.h:399
virtual void restrict_and_add(const unsigned int from_level, VectorType &dst, const VectorType &src) const