deal.II version 9.7.0
\(\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\}}\)
Loading...
Searching...
No Matches
mg_transfer_prebuilt.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2003 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
18
21
22#include <deal.II/fe/fe.h>
23
24#include <deal.II/grid/tria.h>
26
37#include <deal.II/lac/vector.h>
38
41
42#include <algorithm>
43
45
46
47template <typename VectorType>
53
54
55
56template <typename VectorType>
57void
63
64
65
66template <typename VectorType>
67void
75
76
77
78template <typename VectorType>
79void
81 VectorType &dst,
82 const VectorType &src) const
83{
84 Assert((to_level >= 1) && (to_level <= prolongation_matrices.size()),
85 ExcIndexRange(to_level, 1, prolongation_matrices.size() + 1));
86
87 if (this->mg_constrained_dofs != nullptr &&
88 this->mg_constrained_dofs->get_user_constraint_matrix(to_level - 1)
89 .get_local_lines()
90 .size() > 0)
91 {
92 VectorType copy_src(src);
93 this->mg_constrained_dofs->get_user_constraint_matrix(to_level - 1)
94 .distribute(copy_src);
95 prolongation_matrices[to_level - 1]->vmult(dst, copy_src);
96 }
97 else
98 {
99 prolongation_matrices[to_level - 1]->vmult(dst, src);
100 }
101}
102
103
104
105template <typename VectorType>
106void
108 VectorType &dst,
109 const VectorType &src) const
110{
111 Assert((from_level >= 1) && (from_level <= prolongation_matrices.size()),
112 ExcIndexRange(from_level, 1, prolongation_matrices.size() + 1));
113
114 prolongation_matrices[from_level - 1]->Tvmult_add(dst, src);
115}
116
117
118namespace
119{
124 void
125 replace(const MGConstrainedDoFs *mg_constrained_dofs,
126 const unsigned int level,
127 std::vector<types::global_dof_index> &dof_indices)
128 {
129 if (mg_constrained_dofs != nullptr &&
130 mg_constrained_dofs->get_level_constraints(level).n_constraints() > 0)
131 for (auto &ind : dof_indices)
132 if (mg_constrained_dofs->get_level_constraints(level)
134 {
135 Assert(mg_constrained_dofs->get_level_constraints(level)
137 ->size() == 1,
139 ind = mg_constrained_dofs->get_level_constraints(level)
141 ->front()
142 .first;
143 }
144 }
145} // namespace
146
147template <typename VectorType>
148template <int dim, int spacedim>
149void
151 const DoFHandler<dim, spacedim> &dof_handler)
152{
153 Assert(dof_handler.has_level_dofs(),
155 "The underlying DoFHandler object has not had its "
156 "distribute_mg_dofs() function called, but this is a prerequisite "
157 "for multigrid transfers. You will need to call this function, "
158 "probably close to where you already call distribute_dofs()."));
159
160 const unsigned int n_levels =
161 dof_handler.get_triangulation().n_global_levels();
162 const unsigned int dofs_per_cell = dof_handler.get_fe().n_dofs_per_cell();
163
164 this->sizes.resize(n_levels);
165 for (unsigned int l = 0; l < n_levels; ++l)
166 this->sizes[l] = dof_handler.n_dofs(l);
167
168 // reset the size of the array of
169 // matrices. call resize(0) first,
170 // in order to delete all elements
171 // and clear their memory. then
172 // repopulate these arrays
173 //
174 // note that on resize(0), the
175 // shared_ptr class takes care of
176 // deleting the object it points to
177 // by itself
178 prolongation_matrices.resize(0);
179 prolongation_sparsities.resize(0);
180 prolongation_matrices.reserve(n_levels - 1);
181 prolongation_sparsities.reserve(n_levels - 1);
182
183 for (unsigned int i = 0; i < n_levels - 1; ++i)
184 {
185 prolongation_sparsities.emplace_back(
187 prolongation_matrices.emplace_back(
189 }
190
191 // two fields which will store the
192 // indices of the multigrid dofs
193 // for a cell and one of its children
194 std::vector<types::global_dof_index> dof_indices_parent(dofs_per_cell);
195 std::vector<types::global_dof_index> dof_indices_child(dofs_per_cell);
196 std::vector<types::global_dof_index> entries(dofs_per_cell);
197
198 // for each level: first build the sparsity
199 // pattern of the matrices and then build the
200 // matrices themselves. note that we only
201 // need to take care of cells on the coarser
202 // level which have children
203 for (unsigned int level = 0; level < n_levels - 1; ++level)
204 {
205 // reset the dimension of the structure. note that for the number of
206 // entries per row, the number of parent dofs coupling to a child dof is
207 // necessary. this, of course, is the number of degrees of freedom per
208 // cell
209 //
210 // increment dofs_per_cell since a useless diagonal element will be
211 // stored
212 const IndexSet level_p1_relevant_dofs =
213 DoFTools::extract_locally_relevant_level_dofs(dof_handler, level + 1);
214 DynamicSparsityPattern dsp(this->sizes[level + 1],
215 this->sizes[level],
216 level_p1_relevant_dofs);
217 for (const auto &cell : dof_handler.cell_iterators_on_level(level))
218 if (cell->has_children() && cell->is_locally_owned_on_level())
219 {
220 cell->get_mg_dof_indices(dof_indices_parent);
221
222 replace(this->mg_constrained_dofs, level, dof_indices_parent);
223
224 Assert(cell->n_children() ==
227 for (unsigned int child = 0; child < cell->n_children(); ++child)
228 {
229 // set an alias to the prolongation matrix for this child
230 const FullMatrix<double> &prolongation =
231 dof_handler.get_fe().get_prolongation_matrix(
232 child, cell->refinement_case());
233
234 Assert(prolongation.n() != 0, ExcNoProlongation());
235
236 cell->child(child)->get_mg_dof_indices(dof_indices_child);
237
238 replace(this->mg_constrained_dofs,
239 level + 1,
240 dof_indices_child);
241
242 // now tag the entries in the
243 // matrix which will be used
244 // for this pair of parent/child
245 for (unsigned int i = 0; i < dofs_per_cell; ++i)
246 {
247 entries.resize(0);
248 for (unsigned int j = 0; j < dofs_per_cell; ++j)
249 if (prolongation(i, j) != 0)
250 entries.push_back(dof_indices_parent[j]);
251 dsp.add_entries(dof_indices_child[i],
252 entries.begin(),
253 entries.end());
254 }
255 }
256 }
257
258#ifdef DEAL_II_WITH_MPI
260 VectorType>::requires_distributed_sparsity_pattern)
261 {
262 // Since PETSc matrices do not offer the functionality to fill up in-
263 // complete sparsity patterns on their own, the sparsity pattern must
264 // be manually distributed.
265
266 // Distribute sparsity pattern
268 dsp,
269 dof_handler.locally_owned_mg_dofs(level + 1),
270 dof_handler.get_mpi_communicator(),
271 dsp.row_index_set());
272 }
273#endif
274
276 *prolongation_matrices[level],
278 level,
279 dsp,
280 dof_handler);
281 dsp.reinit(0, 0);
282
283 // In the end, the entries in this object will only be real valued.
284 // Nevertheless, we have to take the underlying scalar type of the
285 // vector we want to use this class with. The global matrix the entries
286 // of this matrix are copied into has to be able to perform a
287 // matrix-vector multiplication and this is in general only implemented if
288 // the scalar type for matrix and vector is the same. Furthermore,
289 // copying entries between this local object and the global matrix is only
290 // implemented if the objects have the same scalar type.
292
293 // now actually build the matrices
294 for (const auto &cell : dof_handler.cell_iterators_on_level(level))
295 if (cell->has_children() && cell->is_locally_owned_on_level())
296 {
297 cell->get_mg_dof_indices(dof_indices_parent);
298
299 replace(this->mg_constrained_dofs, level, dof_indices_parent);
300
301 Assert(cell->n_children() ==
304 for (unsigned int child = 0; child < cell->n_children(); ++child)
305 {
306 // set an alias to the prolongation matrix for this child
307 prolongation = dof_handler.get_fe().get_prolongation_matrix(
308 child, cell->refinement_case());
309
310 if (this->mg_constrained_dofs != nullptr &&
311 this->mg_constrained_dofs->have_boundary_indices())
312 for (unsigned int j = 0; j < dofs_per_cell; ++j)
313 if (this->mg_constrained_dofs->is_boundary_index(
314 level, dof_indices_parent[j]))
315 for (unsigned int i = 0; i < dofs_per_cell; ++i)
316 prolongation(i, j) = 0.;
317
318 cell->child(child)->get_mg_dof_indices(dof_indices_child);
319
320 replace(this->mg_constrained_dofs,
321 level + 1,
322 dof_indices_child);
323
324 // now set the entries in the matrix
325 for (unsigned int i = 0; i < dofs_per_cell; ++i)
326 prolongation_matrices[level]->set(dof_indices_child[i],
327 dofs_per_cell,
328 dof_indices_parent.data(),
329 &prolongation(i, 0),
330 true);
331 }
332 }
334 }
335
336 this->fill_and_communicate_copy_indices(dof_handler);
337}
338
339
340
341template <typename VectorType>
342void
344{
345 for (unsigned int level = 0; level < prolongation_matrices.size(); ++level)
346 {
347 os << "Level " << level << std::endl;
348 prolongation_matrices[level]->print(os);
349 os << std::endl;
350 }
351}
352
353
354
355template <typename VectorType>
356std::size_t
358{
360 for (unsigned int i = 0; i < prolongation_matrices.size(); ++i)
362 prolongation_sparsities[i]->memory_consumption();
363
364 return result;
365}
366
367
368// explicit instantiation
369#include "multigrid/mg_transfer_prebuilt.inst"
370
371
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
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
const Triangulation< dim, spacedim > & get_triangulation() const
types::global_dof_index n_dofs() const
MPI_Comm get_mpi_communicator() const
bool has_level_dofs() const
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())
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
ObserverPointer< const MGConstrainedDoFs > mg_constrained_dofs
void fill_and_communicate_copy_indices(const DoFHandler< dim, spacedim > &dof_handler)
std::vector< types::global_dof_index > sizes
std::vector< std::vector< bool > > interface_dofs
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)
std::vector< std::shared_ptr< typename internal::MatrixSelector< VectorType >::Sparsity > > prolongation_sparsities
std::size_t memory_consumption() const
MGTransferPrebuilt()=default
virtual void prolongate(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
std::vector< std::shared_ptr< typename internal::MatrixSelector< VectorType >::Matrix > > prolongation_matrices
void print_matrices(std::ostream &os) const
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
virtual unsigned int n_global_levels() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
static ::ExceptionBase & ExcNoProlongation()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
IndexSet extract_locally_relevant_level_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int level)
void distribute_sparsity_pattern(DynamicSparsityPattern &dsp, const IndexSet &locally_owned_rows, const MPI_Comm mpi_comm, const IndexSet &locally_relevant_rows)
static constexpr unsigned int max_children_per_cell
static void reinit(Matrix &matrix, Sparsity &sparsity, int level, const SparsityPatternType &sp, const DoFHandler< dim, spacedim > &)
Definition mg_transfer.h:57
::SparseMatrix< typename VectorType::value_type > Matrix
Definition mg_transfer.h:51
::SparsityPattern Sparsity
Definition mg_transfer.h:50