Reference documentation for deal.II version 8.5.1
mg_transfer_prebuilt.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2016 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 
17 #include <deal.II/base/logstream.h>
18 #include <deal.II/base/function.h>
19 
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>
38 
39 #include <algorithm>
40 
41 DEAL_II_NAMESPACE_OPEN
42 
43 
44 template<typename VectorType>
46 {}
47 
48 
49 
50 template<typename VectorType>
52 {
53  this->mg_constrained_dofs = &mg_c;
54 }
55 
56 
57 
58 template<typename VectorType>
60 {
61  this->mg_constrained_dofs = &mg_c;
62 }
63 
64 
65 
66 template <typename VectorType>
68 {}
69 
70 
71 
72 template <typename VectorType>
74 (const MGConstrainedDoFs &mg_c)
75 {
76  this->mg_constrained_dofs = &mg_c;
77 }
78 
79 
80 
81 template <typename VectorType>
83 (const ConstraintMatrix &/*c*/, const MGConstrainedDoFs &mg_c)
84 {
85  initialize_constraints(mg_c);
86 }
87 
88 
89 
90 template <typename VectorType>
92 {
94  prolongation_matrices.resize(0);
95  prolongation_sparsities.resize(0);
96  interface_dofs.resize(0);
97 }
98 
99 
100 
101 template <typename VectorType>
102 void MGTransferPrebuilt<VectorType>::prolongate (const unsigned int to_level,
103  VectorType &dst,
104  const VectorType &src) const
105 {
106  Assert ((to_level >= 1) && (to_level<=prolongation_matrices.size()),
107  ExcIndexRange (to_level, 1, prolongation_matrices.size()+1));
108 
109  prolongation_matrices[to_level-1]->vmult (dst, src);
110 }
111 
112 
113 
114 template <typename VectorType>
115 void MGTransferPrebuilt<VectorType>::restrict_and_add (const unsigned int from_level,
116  VectorType &dst,
117  const VectorType &src) const
118 {
119  Assert ((from_level >= 1) && (from_level<=prolongation_matrices.size()),
120  ExcIndexRange (from_level, 1, prolongation_matrices.size()+1));
121  (void)from_level;
122 
123  prolongation_matrices[from_level-1]->Tvmult_add (dst, src);
124 }
125 
126 
127 
128 template <typename VectorType>
129 template <int dim, int spacedim>
132 {
133  const unsigned int n_levels = mg_dof.get_triangulation().n_global_levels();
134  const unsigned int dofs_per_cell = mg_dof.get_fe().dofs_per_cell;
135 
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);
139 
140  // reset the size of the array of
141  // matrices. call resize(0) first,
142  // in order to delete all elements
143  // and clear their memory. then
144  // repopulate these arrays
145  //
146  // note that on resize(0), the
147  // shared_ptr class takes care of
148  // deleting the object it points to
149  // by itself
150  prolongation_matrices.resize (0);
151  prolongation_sparsities.resize (0);
152 
153  for (unsigned int i=0; i<n_levels-1; ++i)
154  {
155  prolongation_sparsities.push_back
157  prolongation_matrices.push_back
159  }
160 
161  // two fields which will store the
162  // indices of the multigrid dofs
163  // for a cell and one of its children
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);
167 
168  // for each level: first build the sparsity
169  // pattern of the matrices and then build the
170  // matrices themselves. note that we only
171  // need to take care of cells on the coarser
172  // level which have children
173  for (unsigned int level=0; level<n_levels-1; ++level)
174  {
175  // reset the dimension of the structure. note that for the number of
176  // entries per row, the number of parent dofs coupling to a child dof is
177  // necessary. this, of course, is the number of degrees of freedom per
178  // cell
179  //
180  // increment dofs_per_cell since a useless diagonal element will be
181  // stored
182  IndexSet level_p1_relevant_dofs;
184  level_p1_relevant_dofs);
185  DynamicSparsityPattern dsp (this->sizes[level+1],
186  this->sizes[level],
187  level_p1_relevant_dofs);
188  typename DoFHandler<dim>::cell_iterator cell, endc = mg_dof.end(level);
189  for (cell=mg_dof.begin(level); cell != endc; ++cell)
190  if (cell->has_children() &&
191  ( mg_dof.get_triangulation().locally_owned_subdomain()==numbers::invalid_subdomain_id
192  || cell->level_subdomain_id()==mg_dof.get_triangulation().locally_owned_subdomain()
193  ))
194  {
195  cell->get_mg_dof_indices (dof_indices_parent);
196 
199  for (unsigned int child=0; child<cell->n_children(); ++child)
200  {
201  // set an alias to the prolongation matrix for this child
202  const FullMatrix<double> &prolongation
203  = mg_dof.get_fe().get_prolongation_matrix (child,
204  cell->refinement_case());
205 
206  Assert (prolongation.n() != 0, ExcNoProlongation());
207 
208  cell->child(child)->get_mg_dof_indices (dof_indices_child);
209 
210  // now tag the entries in the
211  // matrix which will be used
212  // for this pair of parent/child
213  for (unsigned int i=0; i<dofs_per_cell; ++i)
214  {
215  entries.resize(0);
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]);
219  dsp.add_entries (dof_indices_child[i],
220  entries.begin(), entries.end());
221  }
222  }
223  }
224 
225  internal::MatrixSelector<VectorType>::reinit(*prolongation_matrices[level],
226  *prolongation_sparsities[level],
227  level,
228  dsp,
229  mg_dof);
230  dsp.reinit(0,0);
231 
232  FullMatrix<double> prolongation;
233 
234  // now actually build the matrices
235  for (cell=mg_dof.begin(level); cell != endc; ++cell)
236  if (cell->has_children() &&
237  (mg_dof.get_triangulation().locally_owned_subdomain()==numbers::invalid_subdomain_id
238  || cell->level_subdomain_id()==mg_dof.get_triangulation().locally_owned_subdomain())
239  )
240  {
241  cell->get_mg_dof_indices (dof_indices_parent);
242 
245  for (unsigned int child=0; child<cell->n_children(); ++child)
246  {
247  // set an alias to the prolongation matrix for this child
248  prolongation
249  = mg_dof.get_fe().get_prolongation_matrix (child,
250  cell->refinement_case());
251 
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.;
258 
259  cell->child(child)->get_mg_dof_indices (dof_indices_child);
260 
261  // now set the entries in the matrix
262  for (unsigned int i=0; i<dofs_per_cell; ++i)
263  prolongation_matrices[level]->set (dof_indices_child[i],
264  dofs_per_cell,
265  &dof_indices_parent[0],
266  &prolongation(i,0),
267  true);
268  }
269  }
270  prolongation_matrices[level]->compress(VectorOperation::insert);
271  }
272 
273  this->fill_and_communicate_copy_indices(mg_dof);
274 }
275 
276 
277 
278 template <typename VectorType>
279 void
281 {
282  for (unsigned int level = 0; level<prolongation_matrices.size(); ++level)
283  {
284  os << "Level " << level << std::endl;
285  prolongation_matrices[level]->print(os);
286  os << std::endl;
287  }
288 }
289 
290 
291 
292 template <typename VectorType>
293 std::size_t
295 {
297  for (unsigned int i=0; i<prolongation_matrices.size(); ++i)
298  result += prolongation_matrices[i]->memory_consumption()
299  + prolongation_sparsities[i]->memory_consumption();
300 
301  return result;
302 }
303 
304 
305 // explicit instantiation
306 #include "mg_transfer_prebuilt.inst"
307 
308 
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
Definition: types.h:245
std::size_t memory_consumption() const
cell_iterator begin(const unsigned int level=0) const
Definition: dof_handler.cc:728
void build_matrices(const DoFHandler< dim, spacedim > &mg_dof)
cell_iterator end() const
Definition: dof_handler.cc:756
const FiniteElement< dim, spacedim > & get_fe() const
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
void extract_locally_relevant_level_dofs(const DoFHandlerType &dof_handler, const unsigned int level, IndexSet &dof_set)
Definition: dof_tools.cc:1019
void print_matrices(std::ostream &os) const
std::size_t memory_consumption() const
size_type n() const
#define Assert(cond, exc)
Definition: exceptions.h:313
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