17 #include <deal.II/base/logstream.h> 18 #include <deal.II/base/function.h> 20 #include <deal.II/lac/vector.h> 21 #include <deal.II/lac/block_vector.h> 22 #include <deal.II/lac/la_parallel_vector.h> 23 #include <deal.II/lac/la_parallel_block_vector.h> 24 #include <deal.II/lac/petsc_vector.h> 25 #include <deal.II/lac/petsc_block_vector.h> 26 #include <deal.II/lac/trilinos_vector.h> 27 #include <deal.II/lac/trilinos_block_vector.h> 28 #include <deal.II/lac/sparse_matrix.h> 29 #include <deal.II/lac/block_sparse_matrix.h> 30 #include <deal.II/lac/dynamic_sparsity_pattern.h> 31 #include <deal.II/grid/tria.h> 32 #include <deal.II/grid/tria_iterator.h> 33 #include <deal.II/dofs/dof_tools.h> 34 #include <deal.II/fe/fe.h> 35 #include <deal.II/dofs/dof_accessor.h> 36 #include <deal.II/multigrid/mg_tools.h> 37 #include <deal.II/multigrid/mg_transfer.h> 41 DEAL_II_NAMESPACE_OPEN
44 template<
typename VectorType>
50 template<
typename VectorType>
53 this->mg_constrained_dofs = &mg_c;
58 template<
typename VectorType>
61 this->mg_constrained_dofs = &mg_c;
66 template <
typename VectorType>
72 template <
typename VectorType>
76 this->mg_constrained_dofs = &mg_c;
81 template <
typename VectorType>
85 initialize_constraints(mg_c);
90 template <
typename VectorType>
94 prolongation_matrices.resize(0);
95 prolongation_sparsities.resize(0);
96 interface_dofs.resize(0);
101 template <
typename VectorType>
104 const VectorType &src)
const 106 Assert ((to_level >= 1) && (to_level<=prolongation_matrices.size()),
107 ExcIndexRange (to_level, 1, prolongation_matrices.size()+1));
109 prolongation_matrices[to_level-1]->vmult (dst, src);
114 template <
typename VectorType>
117 const VectorType &src)
const 119 Assert ((from_level >= 1) && (from_level<=prolongation_matrices.size()),
120 ExcIndexRange (from_level, 1, prolongation_matrices.size()+1));
123 prolongation_matrices[from_level-1]->Tvmult_add (dst, src);
128 template <
typename VectorType>
129 template <
int dim,
int spacedim>
134 const unsigned int dofs_per_cell = mg_dof.
get_fe().dofs_per_cell;
136 this->sizes.resize(n_levels);
137 for (
unsigned int l=0; l<n_levels; ++l)
138 this->sizes[l] = mg_dof.
n_dofs(l);
150 prolongation_matrices.resize (0);
151 prolongation_sparsities.resize (0);
153 for (
unsigned int i=0; i<n_levels-1; ++i)
155 prolongation_sparsities.push_back
157 prolongation_matrices.push_back
164 std::vector<types::global_dof_index> dof_indices_parent (dofs_per_cell);
165 std::vector<types::global_dof_index> dof_indices_child (dofs_per_cell);
166 std::vector<types::global_dof_index> entries (dofs_per_cell);
173 for (
unsigned int level=0; level<n_levels-1; ++level)
184 level_p1_relevant_dofs);
187 level_p1_relevant_dofs);
189 for (cell=mg_dof.
begin(level); cell != endc; ++cell)
190 if (cell->has_children() &&
192 || cell->level_subdomain_id()==mg_dof.
get_triangulation().locally_owned_subdomain()
195 cell->get_mg_dof_indices (dof_indices_parent);
199 for (
unsigned int child=0; child<cell->n_children(); ++child)
203 = mg_dof.
get_fe().get_prolongation_matrix (child,
204 cell->refinement_case());
206 Assert (prolongation.
n() != 0, ExcNoProlongation());
208 cell->child(child)->get_mg_dof_indices (dof_indices_child);
213 for (
unsigned int i=0; i<dofs_per_cell; ++i)
216 for (
unsigned int j=0; j<dofs_per_cell; ++j)
217 if (prolongation(i,j) != 0)
218 entries.push_back (dof_indices_parent[j]);
220 entries.begin(), entries.end());
225 internal::MatrixSelector<VectorType>::reinit(*prolongation_matrices[level],
226 *prolongation_sparsities[level],
235 for (cell=mg_dof.
begin(level); cell != endc; ++cell)
236 if (cell->has_children() &&
238 || cell->level_subdomain_id()==mg_dof.
get_triangulation().locally_owned_subdomain())
241 cell->get_mg_dof_indices (dof_indices_parent);
245 for (
unsigned int child=0; child<cell->n_children(); ++child)
249 = mg_dof.
get_fe().get_prolongation_matrix (child,
250 cell->refinement_case());
252 if (this->mg_constrained_dofs != 0 &&
253 this->mg_constrained_dofs->have_boundary_indices())
254 for (
unsigned int j=0; j<dofs_per_cell; ++j)
255 if (this->mg_constrained_dofs->is_boundary_index(level, dof_indices_parent[j]))
256 for (
unsigned int i=0; i<dofs_per_cell; ++i)
257 prolongation(i,j) = 0.;
259 cell->child(child)->get_mg_dof_indices (dof_indices_child);
262 for (
unsigned int i=0; i<dofs_per_cell; ++i)
263 prolongation_matrices[level]->
set (dof_indices_child[i],
265 &dof_indices_parent[0],
273 this->fill_and_communicate_copy_indices(mg_dof);
278 template <
typename VectorType>
282 for (
unsigned int level = 0; level<prolongation_matrices.size(); ++level)
284 os <<
"Level " << level << std::endl;
285 prolongation_matrices[level]->print(os);
292 template <
typename VectorType>
297 for (
unsigned int i=0; i<prolongation_matrices.size(); ++i)
298 result += prolongation_matrices[i]->memory_consumption()
299 + prolongation_sparsities[i]->memory_consumption();
306 #include "mg_transfer_prebuilt.inst" 309 DEAL_II_NAMESPACE_CLOSE
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
void add_entries(const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_unique_and_sorted=false)
const types::subdomain_id invalid_subdomain_id
std::size_t memory_consumption() const
cell_iterator begin(const unsigned int level=0) const
void build_matrices(const DoFHandler< dim, spacedim > &mg_dof)
cell_iterator end() const
const FiniteElement< dim, spacedim > & get_fe() const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
virtual ~MGTransferPrebuilt()
void print_matrices(std::ostream &os) const
std::size_t memory_consumption() const
#define Assert(cond, exc)
void reinit(const size_type m, const size_type n, const IndexSet &rowset=IndexSet())
types::global_dof_index n_dofs() const
virtual void prolongate(const unsigned int to_level, VectorType &dst, const VectorType &src) const
static ::ExceptionBase & ExcNotImplemented()
const Triangulation< dim, spacedim > & get_triangulation() const
virtual void restrict_and_add(const unsigned int from_level, VectorType &dst, const VectorType &src) const