Reference documentation for deal.II version 9.3.2
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
mg_transfer_prebuilt.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2003 - 2021 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/base/function.h>
18 #include <deal.II/base/logstream.h>
19 
21 #include <deal.II/dofs/dof_tools.h>
22 
23 #include <deal.II/fe/fe.h>
24 
25 #include <deal.II/grid/tria.h>
27 
38 #include <deal.II/lac/vector.h>
39 
42 
43 #include <algorithm>
44 
46 
47 
48 template <typename VectorType>
50  const MGConstrainedDoFs &mg_c)
51 {
52  this->mg_constrained_dofs = &mg_c;
53 }
54 
55 
56 
57 template <typename VectorType>
58 void
60  const MGConstrainedDoFs &mg_c)
61 {
62  this->mg_constrained_dofs = &mg_c;
63 }
64 
65 
66 
67 template <typename VectorType>
68 void
70 {
72  prolongation_matrices.resize(0);
73  prolongation_sparsities.resize(0);
74  interface_dofs.resize(0);
75 }
76 
77 
78 
79 template <typename VectorType>
80 void
81 MGTransferPrebuilt<VectorType>::prolongate(const unsigned int to_level,
82  VectorType & dst,
83  const VectorType & src) const
84 {
85  Assert((to_level >= 1) && (to_level <= prolongation_matrices.size()),
86  ExcIndexRange(to_level, 1, prolongation_matrices.size() + 1));
87 
88  if (this->mg_constrained_dofs != nullptr &&
89  this->mg_constrained_dofs->get_user_constraint_matrix(to_level - 1)
90  .get_local_lines()
91  .size() > 0)
92  {
93  VectorType copy_src(src);
94  this->mg_constrained_dofs->get_user_constraint_matrix(to_level - 1)
95  .distribute(copy_src);
96  prolongation_matrices[to_level - 1]->vmult(dst, copy_src);
97  }
98  else
99  {
100  prolongation_matrices[to_level - 1]->vmult(dst, src);
101  }
102 }
103 
104 
105 
106 template <typename VectorType>
107 void
109  VectorType & dst,
110  const VectorType & src) const
111 {
112  Assert((from_level >= 1) && (from_level <= prolongation_matrices.size()),
113  ExcIndexRange(from_level, 1, prolongation_matrices.size() + 1));
114  (void)from_level;
115 
116  prolongation_matrices[from_level - 1]->Tvmult_add(dst, src);
117 }
118 
119 
120 namespace
121 {
126  void
127  replace(const MGConstrainedDoFs * mg_constrained_dofs,
128  const unsigned int level,
129  std::vector<types::global_dof_index> &dof_indices)
130  {
131  if (mg_constrained_dofs != nullptr &&
132  mg_constrained_dofs->get_level_constraints(level).n_constraints() > 0)
133  for (auto &ind : dof_indices)
134  if (mg_constrained_dofs->get_level_constraints(level)
136  {
137  Assert(mg_constrained_dofs->get_level_constraints(level)
139  ->size() == 1,
140  ExcInternalError());
141  ind = mg_constrained_dofs->get_level_constraints(level)
143  ->front()
144  .first;
145  }
146  }
147 } // namespace
148 
149 template <typename VectorType>
150 template <int dim, int spacedim>
151 void
153  const DoFHandler<dim, spacedim> &dof_handler)
154 {
155  const unsigned int n_levels =
156  dof_handler.get_triangulation().n_global_levels();
157  const unsigned int dofs_per_cell = dof_handler.get_fe().n_dofs_per_cell();
158 
159  this->sizes.resize(n_levels);
160  for (unsigned int l = 0; l < n_levels; ++l)
161  this->sizes[l] = dof_handler.n_dofs(l);
162 
163  // reset the size of the array of
164  // matrices. call resize(0) first,
165  // in order to delete all elements
166  // and clear their memory. then
167  // repopulate these arrays
168  //
169  // note that on resize(0), the
170  // shared_ptr class takes care of
171  // deleting the object it points to
172  // by itself
173  prolongation_matrices.resize(0);
174  prolongation_sparsities.resize(0);
175  prolongation_matrices.reserve(n_levels - 1);
176  prolongation_sparsities.reserve(n_levels - 1);
177 
178  for (unsigned int i = 0; i < n_levels - 1; ++i)
179  {
180  prolongation_sparsities.emplace_back(
182  prolongation_matrices.emplace_back(
184  }
185 
186  // two fields which will store the
187  // indices of the multigrid dofs
188  // for a cell and one of its children
189  std::vector<types::global_dof_index> dof_indices_parent(dofs_per_cell);
190  std::vector<types::global_dof_index> dof_indices_child(dofs_per_cell);
191  std::vector<types::global_dof_index> entries(dofs_per_cell);
192 
193  // for each level: first build the sparsity
194  // pattern of the matrices and then build the
195  // matrices themselves. note that we only
196  // need to take care of cells on the coarser
197  // level which have children
198  for (unsigned int level = 0; level < n_levels - 1; ++level)
199  {
200  // reset the dimension of the structure. note that for the number of
201  // entries per row, the number of parent dofs coupling to a child dof is
202  // necessary. this, of course, is the number of degrees of freedom per
203  // cell
204  //
205  // increment dofs_per_cell since a useless diagonal element will be
206  // stored
207  IndexSet level_p1_relevant_dofs;
209  level + 1,
210  level_p1_relevant_dofs);
211  DynamicSparsityPattern dsp(this->sizes[level + 1],
212  this->sizes[level],
213  level_p1_relevant_dofs);
214  typename DoFHandler<dim>::cell_iterator cell,
215  endc = dof_handler.end(level);
216  for (cell = dof_handler.begin(level); cell != endc; ++cell)
217  if (cell->has_children() &&
218  (dof_handler.get_triangulation().locally_owned_subdomain() ==
220  cell->level_subdomain_id() ==
222  {
223  cell->get_mg_dof_indices(dof_indices_parent);
224 
225  replace(this->mg_constrained_dofs, level, dof_indices_parent);
226 
227  Assert(cell->n_children() ==
230  for (unsigned int child = 0; child < cell->n_children(); ++child)
231  {
232  // set an alias to the prolongation matrix for this child
233  const FullMatrix<double> &prolongation =
234  dof_handler.get_fe().get_prolongation_matrix(
235  child, cell->refinement_case());
236 
237  Assert(prolongation.n() != 0, ExcNoProlongation());
238 
239  cell->child(child)->get_mg_dof_indices(dof_indices_child);
240 
241  replace(this->mg_constrained_dofs,
242  level + 1,
243  dof_indices_child);
244 
245  // now tag the entries in the
246  // matrix which will be used
247  // for this pair of parent/child
248  for (unsigned int i = 0; i < dofs_per_cell; ++i)
249  {
250  entries.resize(0);
251  for (unsigned int j = 0; j < dofs_per_cell; ++j)
252  if (prolongation(i, j) != 0)
253  entries.push_back(dof_indices_parent[j]);
254  dsp.add_entries(dof_indices_child[i],
255  entries.begin(),
256  entries.end());
257  }
258  }
259  }
260 
261 #ifdef DEAL_II_WITH_MPI
263  VectorType>::requires_distributed_sparsity_pattern)
264  {
265  // Since PETSc matrices do not offer the functionality to fill up in-
266  // complete sparsity patterns on their own, the sparsity pattern must
267  // be manually distributed.
268 
269  // Distribute sparsity pattern
271  dsp,
272  dof_handler.locally_owned_mg_dofs(level + 1),
273  dof_handler.get_communicator(),
274  dsp.row_index_set());
275  }
276 #endif
277 
279  *prolongation_matrices[level],
280  *prolongation_sparsities[level],
281  level,
282  dsp,
283  dof_handler);
284  dsp.reinit(0, 0);
285 
286  // In the end, the entries in this object will only be real valued.
287  // Nevertheless, we have to take the underlying scalar type of the
288  // vector we want to use this class with. The global matrix the entries
289  // of this matrix are copied into has to be able to perform a
290  // matrix-vector multiplication and this is in general only implemented if
291  // the scalar type for matrix and vector is the same. Furthermore,
292  // copying entries between this local object and the global matrix is only
293  // implemented if the objects have the same scalar type.
295 
296  // now actually build the matrices
297  for (cell = dof_handler.begin(level); cell != endc; ++cell)
298  if (cell->has_children() &&
299  (dof_handler.get_triangulation().locally_owned_subdomain() ==
301  cell->level_subdomain_id() ==
303  {
304  cell->get_mg_dof_indices(dof_indices_parent);
305 
306  replace(this->mg_constrained_dofs, level, dof_indices_parent);
307 
308  Assert(cell->n_children() ==
311  for (unsigned int child = 0; child < cell->n_children(); ++child)
312  {
313  // set an alias to the prolongation matrix for this child
314  prolongation = dof_handler.get_fe().get_prolongation_matrix(
315  child, cell->refinement_case());
316 
317  if (this->mg_constrained_dofs != nullptr &&
318  this->mg_constrained_dofs->have_boundary_indices())
319  for (unsigned int j = 0; j < dofs_per_cell; ++j)
320  if (this->mg_constrained_dofs->is_boundary_index(
321  level, dof_indices_parent[j]))
322  for (unsigned int i = 0; i < dofs_per_cell; ++i)
323  prolongation(i, j) = 0.;
324 
325  cell->child(child)->get_mg_dof_indices(dof_indices_child);
326 
327  replace(this->mg_constrained_dofs,
328  level + 1,
329  dof_indices_child);
330 
331  // now set the entries in the matrix
332  for (unsigned int i = 0; i < dofs_per_cell; ++i)
333  prolongation_matrices[level]->set(dof_indices_child[i],
334  dofs_per_cell,
335  dof_indices_parent.data(),
336  &prolongation(i, 0),
337  true);
338  }
339  }
340  prolongation_matrices[level]->compress(VectorOperation::insert);
341  }
342 
343  this->fill_and_communicate_copy_indices(dof_handler);
344 }
345 
346 
347 
348 template <typename VectorType>
349 template <int dim, int spacedim>
350 void
352  const DoFHandler<dim, spacedim> &dof_handler)
353 {
354  build(dof_handler);
355 }
356 
357 
358 
359 template <typename VectorType>
360 void
362 {
363  for (unsigned int level = 0; level < prolongation_matrices.size(); ++level)
364  {
365  os << "Level " << level << std::endl;
366  prolongation_matrices[level]->print(os);
367  os << std::endl;
368  }
369 }
370 
371 
372 
373 template <typename VectorType>
374 std::size_t
376 {
378  for (unsigned int i = 0; i < prolongation_matrices.size(); ++i)
379  result += prolongation_matrices[i]->memory_consumption() +
380  prolongation_sparsities[i]->memory_consumption();
381 
382  return result;
383 }
384 
385 
386 // explicit instantiation
387 #include "mg_transfer_prebuilt.inst"
388 
389 
const std::vector< std::pair< size_type, number > > * get_constraint_entries(const size_type line_n) const
bool is_identity_constrained(const size_type line_n) const
size_type n_constraints() const
cell_iterator end() const
const Triangulation< dim, spacedim > & get_triangulation() const
types::global_dof_index n_dofs() const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
cell_iterator begin(const unsigned int level=0) const
MPI_Comm get_communicator() const
unsigned int n_dofs_per_cell() const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
size_type n() const
bool is_boundary_index(const unsigned int level, const types::global_dof_index index) const
bool have_boundary_indices() const
const AffineConstraints< double > & get_level_constraints(const unsigned int level) const
std::size_t memory_consumption() const
virtual void restrict_and_add(const unsigned int from_level, VectorType &dst, const VectorType &src) const override
void build(const DoFHandler< dim, spacedim > &dof_handler)
void build_matrices(const DoFHandler< dim, spacedim > &dof_handler)
std::size_t memory_consumption() const
MGTransferPrebuilt()=default
virtual void prolongate(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
void print_matrices(std::ostream &os) const
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
virtual types::subdomain_id locally_owned_subdomain() const
virtual unsigned int n_global_levels() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:396
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:397
unsigned int level
Definition: grid_out.cc:4590
__global__ void set(Number *val, const Number s, const size_type N)
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcInternalError()
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcNotImplemented()
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:466
void add_entries(const size_type row, ForwardIterator begin, ForwardIterator end, const bool indices_are_unique_and_sorted=false)
const IndexSet & row_index_set() const
void reinit(const size_type m, const size_type n, const IndexSet &rowset=IndexSet())
void extract_locally_relevant_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level, IndexSet &dof_set)
Definition: dof_tools.cc:1252
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void distribute_sparsity_pattern(DynamicSparsityPattern &dsp, const IndexSet &locally_owned_rows, const MPI_Comm &mpi_comm, const IndexSet &locally_relevant_rows)
const types::subdomain_id invalid_subdomain_id
Definition: types.h:276
static void reinit(Matrix &matrix, Sparsity &sparsity, int level, const SparsityPatternType &sp, const DoFHandler< dim, spacedim > &)
Definition: mg_transfer.h:58