Reference documentation for deal.II version 8.5.1
mg_transfer_matrix_free.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 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_matrix_free_h
17 #define dealii__mg_transfer_matrix_free_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/lac/la_parallel_vector.h>
22 #include <deal.II/multigrid/mg_base.h>
23 #include <deal.II/multigrid/mg_constrained_dofs.h>
24 #include <deal.II/base/mg_level_object.h>
25 #include <deal.II/multigrid/mg_transfer.h>
26 
27 #include <deal.II/dofs/dof_handler.h>
28 
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 
35 
51 template <int dim, typename Number>
52 class MGTransferMatrixFree : public MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number> >
53 {
54 public:
60 
66 
70  virtual ~MGTransferMatrixFree ();
71 
76 
80  void clear ();
81 
85  void build (const DoFHandler<dim,dim> &mg_dof);
86 
101  virtual void prolongate (const unsigned int to_level,
104 
123  virtual void restrict_and_add (const unsigned int from_level,
126 
131 
135  std::size_t memory_consumption () const;
136 
137 private:
138 
144  unsigned int fe_degree;
145 
151 
156  unsigned int n_components;
157 
163  unsigned int n_child_cell_dofs;
164 
175  std::vector<std::vector<unsigned int> > level_dof_indices;
176 
180  std::vector<std::vector<std::pair<unsigned int,unsigned int> > > parent_child_connect;
181 
186  std::vector<unsigned int> n_owned_level_cells;
187 
193 
198 
210  std::vector<AlignedVector<VectorizedArray<Number> > > weights_on_refined;
211 
217  std::vector<std::vector<std::vector<unsigned short> > > dirichlet_indices;
218 
222  template <int degree>
223  void do_prolongate_add(const unsigned int to_level,
226 
230  template <int degree>
231  void do_restrict_add(const unsigned int from_level,
234 };
235 
236 
240 DEAL_II_NAMESPACE_CLOSE
241 
242 #endif
void do_restrict_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
std::vector< unsigned int > n_owned_level_cells
AlignedVector< VectorizedArray< Number > > evaluation_data
static ::ExceptionBase & ExcNoProlongation()
virtual void restrict_and_add(const unsigned int from_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
std::size_t memory_consumption() const
std::vector< std::vector< std::vector< unsigned short > > > dirichlet_indices
AlignedVector< VectorizedArray< Number > > prolongation_matrix_1d
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > parent_child_connect
#define DeclException0(Exception0)
Definition: exceptions.h:541
SmartPointer< const MGConstrainedDoFs, MGLevelGlobalTransfer< LinearAlgebra::distributed::Vector< Number > > > mg_constrained_dofs
Definition: mg_transfer.h:418
std::vector< AlignedVector< VectorizedArray< Number > > > weights_on_refined
std::vector< std::vector< unsigned int > > level_dof_indices
void do_prolongate_add(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
virtual void prolongate(const unsigned int to_level, LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
void build(const DoFHandler< dim, dim > &mg_dof)