Reference documentation for deal.II version GIT relicensing-397-g31a1263477 2024-04-16 03:30:02+00:00
\(\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_internal.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2016 - 2023 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
17
19
20#include <deal.II/fe/fe_tools.h>
21
23
25
26#include <memory>
27
29
30namespace internal
31{
32 namespace MGTransfer
33 {
34 // Internal data structure that is used in the MPI communication in
35 // fill_copy_indices(). It represents an entry in the copy_indices* map,
36 // that associates a level dof index with a global dof index.
57
58 template <int dim, int spacedim>
59 void
61 const DoFHandler<dim, spacedim> &dof_handler,
62 const MGConstrainedDoFs *mg_constrained_dofs,
63 std::vector<std::vector<
64 std::pair<types::global_dof_index, types::global_dof_index>>>
65 &copy_indices,
66 std::vector<std::vector<
67 std::pair<types::global_dof_index, types::global_dof_index>>>
68 &copy_indices_global_mine,
69 std::vector<std::vector<
70 std::pair<types::global_dof_index, types::global_dof_index>>>
71 &copy_indices_level_mine,
72 const bool skip_interface_dofs)
73 {
74 // Now we are filling the variables copy_indices*, which are essentially
75 // maps from global to mgdof for each level stored as a std::vector of
76 // pairs. We need to split this map on each level depending on the
77 // ownership of the global and mgdof, so that we later do not access
78 // non-local elements in copy_to/from_mg.
79 // We keep track in the bitfield dof_touched which global dof has been
80 // processed already (on the current level). This is the same as the
81 // multigrid running in serial.
82
83 // map cpu_index -> vector of data
84 // that will be copied into copy_indices_level_mine
85 std::vector<DoFPair> send_data_temp;
86
87 const unsigned int n_levels =
88 dof_handler.get_triangulation().n_global_levels();
89 copy_indices.resize(n_levels);
90 copy_indices_global_mine.resize(n_levels);
91 copy_indices_level_mine.resize(n_levels);
92 const IndexSet &owned_dofs = dof_handler.locally_owned_dofs();
93
94 const unsigned int dofs_per_cell = dof_handler.get_fe().n_dofs_per_cell();
95 std::vector<types::global_dof_index> global_dof_indices(dofs_per_cell);
96 std::vector<types::global_dof_index> level_dof_indices(dofs_per_cell);
97
98 for (unsigned int level = 0; level < n_levels; ++level)
99 {
100 std::vector<bool> dof_touched(owned_dofs.n_elements(), false);
101 const IndexSet &owned_level_dofs =
102 dof_handler.locally_owned_mg_dofs(level);
103
104 // for the most common case where copy_indices are locally owned
105 // both globally and on the level, we want to skip collecting pairs
106 // and later sorting them. instead, we insert these indices into a
107 // plain vector
108 std::vector<types::global_dof_index> unrolled_copy_indices;
109
110 copy_indices_level_mine[level].clear();
111 copy_indices_global_mine[level].clear();
112
113 for (const auto &level_cell :
115 {
116 if (dof_handler.get_triangulation().locally_owned_subdomain() !=
118 (level_cell->level_subdomain_id() ==
120 level_cell->subdomain_id() ==
122 continue;
123
124 unrolled_copy_indices.resize(owned_dofs.n_elements(),
126
127 // get the dof numbers of this cell for the global and the
128 // level-wise numbering
129 level_cell->get_dof_indices(global_dof_indices);
130 level_cell->get_mg_dof_indices(level_dof_indices);
131
132 for (unsigned int i = 0; i < dofs_per_cell; ++i)
133 {
134 // we need to ignore if the DoF is on a refinement edge
135 // (hanging node)
136 if (skip_interface_dofs && mg_constrained_dofs != nullptr &&
137 mg_constrained_dofs->at_refinement_edge(
138 level, level_dof_indices[i]))
139 continue;
140
141 // First check whether we own any of the active dof index
142 // and the level one. This check involves locally owned
143 // indices which often consist only of a single range, so
144 // they are cheap to look up.
145 const types::global_dof_index global_index_in_set =
146 owned_dofs.index_within_set(global_dof_indices[i]);
147 bool global_mine =
148 global_index_in_set != numbers::invalid_dof_index;
149 bool level_mine =
150 owned_level_dofs.is_element(level_dof_indices[i]);
151
152 if (global_mine && level_mine)
153 {
154 // we own both the active dof index and the level one ->
155 // set them into the vector, indexed by the local index
156 // range of the active dof
157 unrolled_copy_indices[global_index_in_set] =
158 level_dof_indices[i];
159 }
160 else if (global_mine &&
161 dof_touched[global_index_in_set] == false)
162 {
163 copy_indices_global_mine[level].emplace_back(
164 global_dof_indices[i], level_dof_indices[i]);
165
166 // send this to the owner of the level_dof:
167 send_data_temp.emplace_back(level,
168 global_dof_indices[i],
169 level_dof_indices[i]);
170 dof_touched[global_index_in_set] = true;
171 }
172 else
173 {
174 // somebody will send those to me
175 }
176 }
177 }
178
179 // we now translate the plain vector for the copy_indices field into
180 // the expected format of a pair of indices
181 if (!unrolled_copy_indices.empty())
182 {
183 copy_indices[level].clear();
184
185 // reserve the full length in case we did not hit global-mine
186 // indices, so we expect all indices to come into copy_indices
187 if (copy_indices_global_mine[level].empty())
188 copy_indices[level].reserve(unrolled_copy_indices.size());
189
190 // owned_dofs.nth_index_in_set(i) in this query is
191 // usually cheap to look up as there are few ranges in
192 // the locally owned part
193 for (unsigned int i = 0; i < unrolled_copy_indices.size(); ++i)
194 if (unrolled_copy_indices[i] != numbers::invalid_dof_index)
195 copy_indices[level].emplace_back(
196 owned_dofs.nth_index_in_set(i), unrolled_copy_indices[i]);
197 }
198 }
199
200 const ::parallel::TriangulationBase<dim, spacedim> *tria =
201 (dynamic_cast<const ::parallel::TriangulationBase<dim, spacedim>
202 *>(&dof_handler.get_triangulation()));
204 send_data_temp.empty() || tria != nullptr,
206 "We should only be sending information with a parallel Triangulation!"));
207
208#ifdef DEAL_II_WITH_MPI
209 if (tria && Utilities::MPI::sum(send_data_temp.size(),
210 tria->get_communicator()) > 0)
211 {
212 const std::set<types::subdomain_id> &neighbors =
213 tria->level_ghost_owners();
214 std::map<int, std::vector<DoFPair>> send_data;
215
216 std::sort(send_data_temp.begin(),
217 send_data_temp.end(),
218 [](const DoFPair &lhs, const DoFPair &rhs) {
219 if (lhs.level < rhs.level)
220 return true;
221 if (lhs.level > rhs.level)
222 return false;
223
224 if (lhs.level_dof_index < rhs.level_dof_index)
225 return true;
226 if (lhs.level_dof_index > rhs.level_dof_index)
227 return false;
228
229 if (lhs.global_dof_index < rhs.global_dof_index)
230 return true;
231 else
232 return false;
233 });
234 send_data_temp.erase(
235 std::unique(send_data_temp.begin(),
236 send_data_temp.end(),
237 [](const DoFPair &lhs, const DoFPair &rhs) {
238 return (lhs.level == rhs.level) &&
239 (lhs.level_dof_index == rhs.level_dof_index) &&
240 (lhs.global_dof_index == rhs.global_dof_index);
241 }),
242 send_data_temp.end());
243
244 for (unsigned int level = 0; level < n_levels; ++level)
245 {
246 const IndexSet &owned_level_dofs =
247 dof_handler.locally_owned_mg_dofs(level);
248
249 std::vector<types::global_dof_index> level_dof_indices;
250 std::vector<types::global_dof_index> global_dof_indices;
251 for (const auto &dofpair : send_data_temp)
252 if (dofpair.level == level)
253 {
254 level_dof_indices.push_back(dofpair.level_dof_index);
255 global_dof_indices.push_back(dofpair.global_dof_index);
256 }
257
258 IndexSet is_ghost(owned_level_dofs.size());
259 is_ghost.add_indices(level_dof_indices.begin(),
260 level_dof_indices.end());
261
262 AssertThrow(level_dof_indices.size() == is_ghost.n_elements(),
263 ExcMessage("Size does not match!"));
264
265 const auto index_owner =
267 is_ghost,
268 tria->get_communicator());
269
270 AssertThrow(level_dof_indices.size() == index_owner.size(),
271 ExcMessage("Size does not match!"));
272
273 for (unsigned int i = 0; i < index_owner.size(); ++i)
274 send_data[index_owner[i]].emplace_back(level,
275 global_dof_indices[i],
276 level_dof_indices[i]);
277 }
278
279
280 // Protect the send/recv logic with a mutex:
283 mutex, tria->get_communicator());
284
285 const int mpi_tag =
287
288 // * send
289 std::vector<MPI_Request> requests;
290 {
291 for (const auto dest : neighbors)
292 {
293 requests.push_back(MPI_Request());
294 std::vector<DoFPair> &data = send_data[dest];
295
296 const int ierr =
297 MPI_Isend(data.data(),
298 data.size() * sizeof(decltype(*data.data())),
299 MPI_BYTE,
300 dest,
301 mpi_tag,
302 tria->get_communicator(),
303 &*requests.rbegin());
304 AssertThrowMPI(ierr);
305 }
306 }
307
308 // * receive
309 {
310 // We should get one message from each of our neighbors
311 std::vector<DoFPair> receive_buffer;
312 for (unsigned int counter = 0; counter < neighbors.size();
313 ++counter)
314 {
315 MPI_Status status;
316 int ierr = MPI_Probe(MPI_ANY_SOURCE,
317 mpi_tag,
318 tria->get_communicator(),
319 &status);
320 AssertThrowMPI(ierr);
321 int len;
322 ierr = MPI_Get_count(&status, MPI_BYTE, &len);
323 AssertThrowMPI(ierr);
324
325 if (len == 0)
326 {
327 ierr = MPI_Recv(nullptr,
328 0,
329 MPI_BYTE,
330 status.MPI_SOURCE,
331 status.MPI_TAG,
332 tria->get_communicator(),
333 &status);
334 AssertThrowMPI(ierr);
335 continue;
336 }
337
338 int count = len / sizeof(DoFPair);
339 Assert(static_cast<int>(count * sizeof(DoFPair)) == len,
341 receive_buffer.resize(count);
342
343 void *ptr = receive_buffer.data();
344 ierr = MPI_Recv(ptr,
345 len,
346 MPI_BYTE,
347 status.MPI_SOURCE,
348 status.MPI_TAG,
349 tria->get_communicator(),
350 &status);
351 AssertThrowMPI(ierr);
352
353 for (const auto &dof_pair : receive_buffer)
354 {
355 copy_indices_level_mine[dof_pair.level].emplace_back(
356 dof_pair.global_dof_index, dof_pair.level_dof_index);
357 }
358 }
359 }
360
361 // * wait for all MPI_Isend to complete
362 if (requests.size() > 0)
363 {
364 const int ierr = MPI_Waitall(requests.size(),
365 requests.data(),
366 MPI_STATUSES_IGNORE);
367 AssertThrowMPI(ierr);
368 requests.clear();
369 }
370# ifdef DEBUG
371 // Make sure in debug mode, that everybody sent/received all packages
372 // on this level. If a deadlock occurs here, the list of expected
373 // senders is not computed correctly.
374 const int ierr = MPI_Barrier(tria->get_communicator());
375 AssertThrowMPI(ierr);
376# endif
377 }
378#endif
379
380 // Sort the indices, except the copy_indices which already are
381 // sorted. This will produce more reliable debug output for regression
382 // tests and won't hurt performance even in release mode because the
383 // non-owned indices are a small subset of all unknowns.
384 std::less<std::pair<types::global_dof_index, types::global_dof_index>>
385 compare;
386 for (auto &level_indices : copy_indices_level_mine)
387 std::sort(level_indices.begin(), level_indices.end(), compare);
388 for (auto &level_indices : copy_indices_global_mine)
389 std::sort(level_indices.begin(), level_indices.end(), compare);
390 }
391
392
393
394 // initialize the vectors needed for the transfer (and merge with the
395 // content in copy_indices_global_mine)
396 void
398 const IndexSet &locally_owned,
399 std::vector<types::global_dof_index> &ghosted_level_dofs,
400 const std::shared_ptr<const Utilities::MPI::Partitioner>
401 &external_partitioner,
402 const MPI_Comm communicator,
403 std::shared_ptr<const Utilities::MPI::Partitioner> &target_partitioner,
404 Table<2, unsigned int> &copy_indices_global_mine)
405 {
406 std::sort(ghosted_level_dofs.begin(), ghosted_level_dofs.end());
407 IndexSet ghosted_dofs(locally_owned.size());
408 ghosted_dofs.add_indices(ghosted_level_dofs.begin(),
409 ghosted_level_dofs.end());
410 ghosted_dofs.compress();
411
412 // Add possible ghosts from the previous content in the vector
413 if (target_partitioner.get() != nullptr &&
414 target_partitioner->size() == locally_owned.size())
415 {
416 ghosted_dofs.add_indices(target_partitioner->ghost_indices());
417 }
418
419 // check if the given partitioner's ghosts represent a superset of the
420 // ghosts we require in this function
421 const int ghosts_locally_contained =
422 (external_partitioner.get() != nullptr &&
423 (external_partitioner->ghost_indices() & ghosted_dofs) ==
424 ghosted_dofs) ?
425 1 :
426 0;
427 if (external_partitioner.get() != nullptr &&
428 Utilities::MPI::min(ghosts_locally_contained, communicator) == 1)
429 {
430 // shift the local number of the copy indices according to the new
431 // partitioner that we are going to use during the access to the
432 // entries
433 if (target_partitioner.get() != nullptr &&
434 target_partitioner->size() == locally_owned.size())
435 for (unsigned int i = 0; i < copy_indices_global_mine.n_cols(); ++i)
436 copy_indices_global_mine(1, i) =
437 external_partitioner->global_to_local(
438 target_partitioner->local_to_global(
439 copy_indices_global_mine(1, i)));
440 target_partitioner = external_partitioner;
441 }
442 else
443 {
444 if (target_partitioner.get() != nullptr &&
445 target_partitioner->size() == locally_owned.size())
446 for (unsigned int i = 0; i < copy_indices_global_mine.n_cols(); ++i)
447 copy_indices_global_mine(1, i) =
448 locally_owned.n_elements() +
449 ghosted_dofs.index_within_set(
450 target_partitioner->local_to_global(
451 copy_indices_global_mine(1, i)));
452 target_partitioner =
453 std::make_shared<Utilities::MPI::Partitioner>(locally_owned,
454 ghosted_dofs,
455 communicator);
456 }
457 }
458
459
460
461 // Transform the ghost indices to local index space for the vector
462 inline void
464 const Utilities::MPI::Partitioner &part,
465 const std::vector<types::global_dof_index> &mine,
466 const std::vector<types::global_dof_index> &remote,
467 std::vector<unsigned int> &localized_indices)
468 {
469 localized_indices.resize(mine.size() + remote.size(),
471 for (unsigned int i = 0; i < mine.size(); ++i)
472 if (mine[i] != numbers::invalid_dof_index)
473 localized_indices[i] = part.global_to_local(mine[i]);
474
475 for (unsigned int i = 0; i < remote.size(); ++i)
476 if (remote[i] != numbers::invalid_dof_index)
477 localized_indices[i + mine.size()] = part.global_to_local(remote[i]);
478 }
479
480
481
482 // given the collection of child cells in lexicographic ordering as seen
483 // from the parent, compute the first index of the given child
484 template <int dim>
485 unsigned int
486 compute_shift_within_children(const unsigned int child,
487 const unsigned int fe_shift_1d,
488 const unsigned int fe_degree)
489 {
490 // we put the degrees of freedom of all child cells in lexicographic
491 // ordering
492 unsigned int c_tensor_index[dim];
493 unsigned int tmp = child;
494 for (unsigned int d = 0; d < dim; ++d)
495 {
496 c_tensor_index[d] = tmp % 2;
497 tmp /= 2;
498 }
499 const unsigned int n_child_dofs_1d = fe_degree + 1 + fe_shift_1d;
500 unsigned int factor = 1;
501 unsigned int shift = fe_shift_1d * c_tensor_index[0];
502 for (unsigned int d = 1; d < dim; ++d)
503 {
504 factor *= n_child_dofs_1d;
505 shift = shift + factor * fe_shift_1d * c_tensor_index[d];
506 }
507 return shift;
508 }
509
510
511
512 // puts the indices on the given child cell in lexicographic ordering with
513 // respect to the collection of all child cells as seen from the parent
514 template <int dim>
515 void
517 const unsigned int child,
518 const unsigned int fe_shift_1d,
519 const unsigned int fe_degree,
520 const std::vector<unsigned int> &lexicographic_numbering,
521 const std::vector<types::global_dof_index> &local_dof_indices,
522 types::global_dof_index *target_indices)
523 {
524 const unsigned int n_child_dofs_1d = fe_degree + 1 + fe_shift_1d;
525 const unsigned int shift =
526 compute_shift_within_children<dim>(child, fe_shift_1d, fe_degree);
527 const unsigned int n_components =
528 local_dof_indices.size() / Utilities::fixed_power<dim>(fe_degree + 1);
529 types::global_dof_index *indices = target_indices + shift;
530 const unsigned int n_scalar_cell_dofs =
531 Utilities::fixed_power<dim>(n_child_dofs_1d);
532 for (unsigned int c = 0, m = 0; c < n_components; ++c)
533 for (unsigned int k = 0; k < (dim > 2 ? (fe_degree + 1) : 1); ++k)
534 for (unsigned int j = 0; j < (dim > 1 ? (fe_degree + 1) : 1); ++j)
535 for (unsigned int i = 0; i < (fe_degree + 1); ++i, ++m)
536 {
537 const unsigned int index =
538 c * n_scalar_cell_dofs +
539 k * n_child_dofs_1d * n_child_dofs_1d + j * n_child_dofs_1d +
540 i;
542 indices[index] ==
543 local_dof_indices[lexicographic_numbering[m]],
545 indices[index] = local_dof_indices[lexicographic_numbering[m]];
546 }
547 }
548
549
550
551 template <int dim, typename Number>
552 void
554 const FiniteElement<1> &fe,
555 const DoFHandler<dim> &dof_handler)
556 {
557 // currently, we have only FE_Q and FE_DGQ type elements implemented
558 elem_info.n_components = dof_handler.get_fe().element_multiplicity(0);
559 AssertDimension(Utilities::fixed_power<dim>(fe.n_dofs_per_cell()) *
560 elem_info.n_components,
561 dof_handler.get_fe().n_dofs_per_cell());
562 AssertDimension(fe.degree, dof_handler.get_fe().degree);
563 elem_info.fe_degree = fe.degree;
564 elem_info.element_is_continuous = fe.n_dofs_per_vertex() > 0;
566
567 // step 1.2: get renumbering of 1d basis functions to lexicographic
568 // numbers. The distinction according to fe.n_dofs_per_vertex() is to
569 // support both continuous and discontinuous bases.
570 std::vector<unsigned int> renumbering(fe.n_dofs_per_cell());
571 {
573 renumbering[0] = 0;
574 for (unsigned int i = 0; i < fe.n_dofs_per_line(); ++i)
575 renumbering[i + fe.n_dofs_per_vertex()] =
577 if (fe.n_dofs_per_vertex() > 0)
578 renumbering[fe.n_dofs_per_cell() - fe.n_dofs_per_vertex()] =
580 }
581
582 // step 1.3: create a dummy 1d quadrature formula to extract the
583 // lexicographic numbering for the elements
584 Assert(fe.n_dofs_per_vertex() == 0 || fe.n_dofs_per_vertex() == 1,
586 const unsigned int shift = fe.n_dofs_per_cell() - fe.n_dofs_per_vertex();
587 const unsigned int n_child_dofs_1d =
588 (fe.n_dofs_per_vertex() > 0 ? (2 * fe.n_dofs_per_cell() - 1) :
589 (2 * fe.n_dofs_per_cell()));
590
591 elem_info.n_child_cell_dofs =
592 elem_info.n_components * Utilities::fixed_power<dim>(n_child_dofs_1d);
593 const Quadrature<1> dummy_quadrature(
594 std::vector<Point<1>>(1, Point<1>()));
596 shape_info.reinit(dummy_quadrature, dof_handler.get_fe(), 0);
598
599 // step 1.4: get the 1d prolongation matrix and combine from both children
600 elem_info.prolongation_matrix_1d.resize(fe.n_dofs_per_cell() *
601 n_child_dofs_1d);
602
603 for (unsigned int c = 0; c < GeometryInfo<1>::max_children_per_cell; ++c)
604 for (unsigned int i = 0; i < fe.n_dofs_per_cell(); ++i)
605 for (unsigned int j = 0; j < fe.n_dofs_per_cell(); ++j)
606 elem_info
607 .prolongation_matrix_1d[i * n_child_dofs_1d + j + c * shift] =
608 fe.get_prolongation_matrix(c)(renumbering[j], renumbering[i]);
609 }
610
611
612
613 // Sets up most of the internal data structures of the MGTransferMatrixFree
614 // class
615 template <int dim, typename Number>
616 void
618 const DoFHandler<dim> &dof_handler,
619 const MGConstrainedDoFs *mg_constrained_dofs,
620 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
621 &external_partitioners,
622 ElementInfo<Number> &elem_info,
623 std::vector<std::vector<unsigned int>> &level_dof_indices,
624 std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
625 &parent_child_connect,
626 std::vector<unsigned int> &n_owned_level_cells,
627 std::vector<std::vector<std::vector<unsigned short>>> &dirichlet_indices,
628 std::vector<std::vector<Number>> &weights_on_refined,
629 std::vector<Table<2, unsigned int>> &copy_indices_global_mine,
630 MGLevelObject<std::shared_ptr<const Utilities::MPI::Partitioner>>
631 &target_partitioners)
632 {
633 level_dof_indices.clear();
634 parent_child_connect.clear();
635 n_owned_level_cells.clear();
636 dirichlet_indices.clear();
637 weights_on_refined.clear();
638
639#ifdef DEBUG
640 if (mg_constrained_dofs)
641 {
642 const unsigned int n_levels =
643 dof_handler.get_triangulation().n_global_levels();
644
645 for (unsigned int l = 0; l < n_levels; ++l)
646 {
647 const auto &constraints =
648 mg_constrained_dofs->get_user_constraint_matrix(l);
649
650 // no inhomogeneities are supported
651 AssertDimension(constraints.n_inhomogeneities(), 0);
652
653 for (const auto dof : constraints.get_local_lines())
654 {
655 const auto *entries_ptr =
656 constraints.get_constraint_entries(dof);
657
658 if (entries_ptr == nullptr)
659 continue;
660
661 // only homogeneous or identity constraints are supported
662 Assert((entries_ptr->size() == 0) ||
663 ((entries_ptr->size() == 1) &&
664 (std::abs((*entries_ptr)[0].second - 1.) <
665 100 * std::numeric_limits<double>::epsilon())),
667 }
668 }
669 }
670#endif
671
672 // we collect all child DoFs of a mother cell together. For faster
673 // tensorized operations, we align the degrees of freedom
674 // lexicographically. We distinguish FE_Q elements and FE_DGQ elements
675
676 const ::Triangulation<dim> &tria = dof_handler.get_triangulation();
677
678 // ---------------------------- 1. Extract 1d info about the finite
679 // element step 1.1: create a 1d copy of the finite element from FETools
680 // where we substitute the template argument
681 AssertDimension(dof_handler.get_fe().n_base_elements(), 1);
682 std::string fe_name = dof_handler.get_fe().base_element(0).get_name();
683 {
684 const std::size_t template_starts = fe_name.find_first_of('<');
685 Assert(fe_name[template_starts + 1] ==
686 (dim == 1 ? '1' : (dim == 2 ? '2' : '3')),
688 fe_name[template_starts + 1] = '1';
689 }
690 const std::unique_ptr<FiniteElement<1>> fe(
691 FETools::get_fe_by_name<1, 1>(fe_name));
692
693 setup_element_info(elem_info, *fe, dof_handler);
694
695
696 // ---------- 2. Extract and match dof indices between child and parent
697 const unsigned int n_levels = tria.n_global_levels();
698 level_dof_indices.resize(n_levels);
699 parent_child_connect.resize(n_levels - 1);
700 n_owned_level_cells.resize(n_levels - 1);
701 std::vector<std::vector<unsigned int>> coarse_level_indices(n_levels - 1);
702 for (unsigned int level = 0;
703 level < std::min(tria.n_levels(), n_levels - 1);
704 ++level)
705 coarse_level_indices[level].resize(tria.n_raw_cells(level),
707 std::vector<types::global_dof_index> local_dof_indices(
708 dof_handler.get_fe().n_dofs_per_cell());
709 dirichlet_indices.resize(n_levels - 1);
710
711 AssertDimension(target_partitioners.max_level(), n_levels - 1);
712 Assert(external_partitioners.empty() ||
713 external_partitioners.size() == n_levels,
714 ExcDimensionMismatch(external_partitioners.size(), n_levels));
715
716 for (unsigned int level = n_levels - 1; level > 0; --level)
717 {
718 unsigned int counter = 0;
719 std::vector<types::global_dof_index> global_level_dof_indices;
720 std::vector<types::global_dof_index> global_level_dof_indices_remote;
721 std::vector<types::global_dof_index> ghosted_level_dofs;
722 std::vector<types::global_dof_index> global_level_dof_indices_l0;
723 std::vector<types::global_dof_index> ghosted_level_dofs_l0;
724
725 // step 2.1: loop over the cells on the coarse side
726 typename DoFHandler<dim>::cell_iterator cell,
727 endc = dof_handler.end(level - 1);
728 for (cell = dof_handler.begin(level - 1); cell != endc; ++cell)
729 {
730 // need to look into a cell if it has children and it is locally
731 // owned
732 if (!cell->has_children())
733 continue;
734
735 bool consider_cell =
736 (tria.locally_owned_subdomain() ==
738 cell->level_subdomain_id() == tria.locally_owned_subdomain());
739
740 // due to the particular way we store DoF indices (via children),
741 // we also need to add the DoF indices for coarse cells where we
742 // own at least one child
743 const bool cell_is_remote = !consider_cell;
744 for (unsigned int c = 0;
745 c < GeometryInfo<dim>::max_children_per_cell;
746 ++c)
747 if (cell->child(c)->level_subdomain_id() ==
749 {
750 consider_cell = true;
751 break;
752 }
753
754 if (!consider_cell)
755 continue;
756
757 // step 2.2: loop through children and append the dof indices to
758 // the appropriate list. We need separate lists for the owned
759 // coarse cell case (which will be part of
760 // restriction/prolongation between level-1 and level) and the
761 // remote case (which needs to store DoF indices for the
762 // operations between level and level+1).
763 AssertDimension(cell->n_children(),
765 std::vector<types::global_dof_index> &next_indices =
766 cell_is_remote ? global_level_dof_indices_remote :
767 global_level_dof_indices;
768 const std::size_t start_index = next_indices.size();
769 next_indices.resize(start_index + elem_info.n_child_cell_dofs,
771 for (unsigned int c = 0;
772 c < GeometryInfo<dim>::max_children_per_cell;
773 ++c)
774 {
775 if (cell_is_remote && cell->child(c)->level_subdomain_id() !=
777 continue;
778 cell->child(c)->get_mg_dof_indices(local_dof_indices);
779
780 resolve_identity_constraints(mg_constrained_dofs,
781 level,
782 local_dof_indices);
783
784 const IndexSet &owned_level_dofs =
785 dof_handler.locally_owned_mg_dofs(level);
786 for (const auto local_dof_index : local_dof_indices)
787 if (!owned_level_dofs.is_element(local_dof_index))
788 ghosted_level_dofs.push_back(local_dof_index);
789
790 add_child_indices<dim>(c,
791 fe->n_dofs_per_cell() -
792 fe->n_dofs_per_vertex(),
793 fe->degree,
794 elem_info.lexicographic_numbering,
795 local_dof_indices,
796 &next_indices[start_index]);
797
798 // step 2.3 store the connectivity to the parent
799 if (cell->child(c)->has_children() &&
800 (tria.locally_owned_subdomain() ==
802 cell->child(c)->level_subdomain_id() ==
804 {
805 const unsigned int child_index =
806 coarse_level_indices[level][cell->child(c)->index()];
807 AssertIndexRange(child_index,
808 parent_child_connect[level].size());
809 unsigned int parent_index = counter;
810 // remote cells, i.e., cells where we work on a further
811 // level but are not treated on the current level, need to
812 // be placed at the end of the list; however, we do not
813 // yet know the exact position in the array, so shift
814 // their parent index by the number of cells so we can set
815 // the correct number after the end of this loop
816 if (cell_is_remote)
817 parent_index =
818 start_index / elem_info.n_child_cell_dofs +
819 tria.n_cells(level);
820 parent_child_connect[level][child_index] =
821 std::make_pair(parent_index, c);
823 static_cast<unsigned short>(-1));
824
825 // set Dirichlet boundary conditions (as a list of
826 // constrained DoFs) for the child
827 if (mg_constrained_dofs != nullptr)
828 for (unsigned int i = 0;
829 i < dof_handler.get_fe().n_dofs_per_cell();
830 ++i)
831 if (mg_constrained_dofs->is_boundary_index(
832 level,
833 local_dof_indices
834 [elem_info.lexicographic_numbering[i]]))
835 dirichlet_indices[level][child_index].push_back(i);
836 }
837 }
838 if (!cell_is_remote)
839 {
840 AssertIndexRange(static_cast<unsigned int>(cell->index()),
841 coarse_level_indices[level - 1].size());
842 coarse_level_indices[level - 1][cell->index()] = counter++;
843 }
844
845 // step 2.4: include indices for the coarsest cells. we still
846 // insert the indices as if they were from a child in order to use
847 // the same code (the coarsest level does not matter much in terms
848 // of memory, so we gain in code simplicity)
849 if (level == 1 && !cell_is_remote)
850 {
851 cell->get_mg_dof_indices(local_dof_indices);
852
853 resolve_identity_constraints(mg_constrained_dofs,
854 level - 1,
855 local_dof_indices);
856
857 const IndexSet &owned_level_dofs_l0 =
858 dof_handler.locally_owned_mg_dofs(0);
859 for (const auto local_dof_index : local_dof_indices)
860 if (!owned_level_dofs_l0.is_element(local_dof_index))
861 ghosted_level_dofs_l0.push_back(local_dof_index);
862
863 const std::size_t start_index =
864 global_level_dof_indices_l0.size();
865 global_level_dof_indices_l0.resize(
866 start_index + elem_info.n_child_cell_dofs,
868 add_child_indices<dim>(
869 0,
870 fe->n_dofs_per_cell() - fe->n_dofs_per_vertex(),
871 fe->degree,
872 elem_info.lexicographic_numbering,
873 local_dof_indices,
874 &global_level_dof_indices_l0[start_index]);
875
876 dirichlet_indices[0].emplace_back();
877 if (mg_constrained_dofs != nullptr)
878 for (unsigned int i = 0;
879 i < dof_handler.get_fe().n_dofs_per_cell();
880 ++i)
881 if (mg_constrained_dofs->is_boundary_index(
882 0,
883 local_dof_indices[elem_info
885 dirichlet_indices[0].back().push_back(i);
886 }
887 }
888
889 // step 2.5: store information about the current level and prepare the
890 // Dirichlet indices and parent-child relationship for the next
891 // coarser level
892 AssertDimension(counter * elem_info.n_child_cell_dofs,
893 global_level_dof_indices.size());
894 n_owned_level_cells[level - 1] = counter;
895 dirichlet_indices[level - 1].resize(counter);
896 parent_child_connect[level - 1].resize(
897 counter,
898 std::make_pair(numbers::invalid_unsigned_int,
900
901 // step 2.6: put the cells with remotely owned parent to the end of
902 // the list (these are needed for the transfer from level to level+1
903 // but not for the transfer from level-1 to level).
904 if (level < n_levels - 1)
905 for (std::vector<std::pair<unsigned int, unsigned int>>::iterator
906 i = parent_child_connect[level].begin();
907 i != parent_child_connect[level].end();
908 ++i)
909 if (i->first >= tria.n_cells(level))
910 {
911 i->first -= tria.n_cells(level);
912 i->first += counter;
913 }
914
915 // step 2.7: Initialize the partitioner for the ghosted vector
916 //
917 // We use a vector based on the target partitioner handed in also in
918 // the base class for keeping ghosted transfer indices. To avoid
919 // keeping two very similar vectors, we keep one single ghosted
920 // vector that is augmented/filled here.
922 ghosted_level_dofs,
923 external_partitioners.empty() ?
924 nullptr :
925 external_partitioners[level],
926 tria.get_communicator(),
927 target_partitioners[level],
928 copy_indices_global_mine[level]);
929
930 copy_indices_to_mpi_local_numbers(*target_partitioners[level],
931 global_level_dof_indices,
932 global_level_dof_indices_remote,
933 level_dof_indices[level]);
934 // step 2.8: Initialize the ghosted vector for level 0
935 if (level == 1)
936 {
937 for (unsigned int i = 0; i < parent_child_connect[0].size(); ++i)
938 parent_child_connect[0][i] = std::make_pair(i, 0U);
939
941 ghosted_level_dofs_l0,
942 external_partitioners.empty() ?
943 nullptr :
944 external_partitioners[0],
945 tria.get_communicator(),
946 target_partitioners[0],
947 copy_indices_global_mine[0]);
948
950 *target_partitioners[0],
951 global_level_dof_indices_l0,
952 std::vector<types::global_dof_index>(),
953 level_dof_indices[0]);
954 }
955 }
956
957 // ---------------------- 3. compute weights to make restriction additive
958
959 const unsigned int n_child_dofs_1d =
960 fe->degree + 1 + fe->n_dofs_per_cell() - fe->n_dofs_per_vertex();
961
962 // get the valence of the individual components and compute the weights as
963 // the inverse of the valence
964 weights_on_refined.resize(n_levels - 1);
965 for (unsigned int level = 1; level < n_levels; ++level)
966 {
968 target_partitioners[level]);
969 for (unsigned int c = 0; c < n_owned_level_cells[level - 1]; ++c)
970 for (unsigned int j = 0; j < elem_info.n_child_cell_dofs; ++j)
971 touch_count.local_element(
972 level_dof_indices[level][elem_info.n_child_cell_dofs * c +
973 j]) += Number(1.);
974 touch_count.compress(VectorOperation::add);
975 touch_count.update_ghost_values();
976
977 std::vector<unsigned int> degree_to_3(n_child_dofs_1d);
978 degree_to_3[0] = 0;
979 for (unsigned int i = 1; i < n_child_dofs_1d - 1; ++i)
980 degree_to_3[i] = 1;
981 degree_to_3.back() = 2;
982
983 // we only store 3^dim weights because all dofs on a line have the
984 // same valence, and all dofs on a quad have the same valence.
985 weights_on_refined[level - 1].resize(n_owned_level_cells[level - 1] *
986 Utilities::fixed_power<dim>(3));
987 for (unsigned int c = 0; c < n_owned_level_cells[level - 1]; ++c)
988 for (unsigned int k = 0, m = 0; k < (dim > 2 ? n_child_dofs_1d : 1);
989 ++k)
990 for (unsigned int j = 0; j < (dim > 1 ? n_child_dofs_1d : 1); ++j)
991 {
992 unsigned int shift = 9 * degree_to_3[k] + 3 * degree_to_3[j];
993 for (unsigned int i = 0; i < n_child_dofs_1d; ++i, ++m)
994 weights_on_refined[level -
995 1][c * Utilities::fixed_power<dim>(3) +
996 shift + degree_to_3[i]] =
997 Number(1.) /
998 touch_count.local_element(
999 level_dof_indices[level]
1000 [elem_info.n_child_cell_dofs * c + m]);
1001 }
1002 }
1003 }
1004
1005
1006
1007 void
1009 const MGConstrainedDoFs *mg_constrained_dofs,
1010 const unsigned int level,
1011 std::vector<types::global_dof_index> &dof_indices)
1012 {
1013 if (mg_constrained_dofs != nullptr &&
1014 mg_constrained_dofs->get_level_constraints(level).n_constraints() > 0)
1015 for (auto &ind : dof_indices)
1016 if (mg_constrained_dofs->get_level_constraints(level)
1018 {
1019 Assert(mg_constrained_dofs->get_level_constraints(level)
1021 ->size() == 1,
1023 ind = mg_constrained_dofs->get_level_constraints(level)
1025 ->front()
1026 .first;
1027 }
1028 }
1029
1030 } // namespace MGTransfer
1031} // namespace internal
1032
1033// Explicit instantiations
1034
1035#include "mg_transfer_internal.inst"
1036
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 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
const IndexSet & locally_owned_dofs() const
cell_iterator begin(const unsigned int level=0) const
unsigned int n_dofs_per_vertex() const
const unsigned int degree
Definition fe_data.h:452
unsigned int n_dofs_per_cell() const
unsigned int n_dofs_per_line() const
virtual std::string get_name() const =0
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
unsigned int element_multiplicity(const unsigned int index) const
unsigned int n_base_elements() const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
size_type size() const
Definition index_set.h:1765
size_type index_within_set(const size_type global_index) const
Definition index_set.h:1990
size_type n_elements() const
Definition index_set.h:1923
bool is_element(const size_type index) const
Definition index_set.h:1883
void clear()
Definition index_set.h:1741
size_type nth_index_in_set(const size_type local_index) const
Definition index_set.h:1971
void compress() const
Definition index_set.h:1773
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Definition index_set.h:1820
Number local_element(const size_type local_index) const
void compress(VectorOperation::values operation)
bool at_refinement_edge(const unsigned int level, const types::global_dof_index index) const
bool is_boundary_index(const unsigned int level, const types::global_dof_index index) const
const AffineConstraints< double > & get_user_constraint_matrix(const unsigned int level) const
const AffineConstraints< double > & get_level_constraints(const unsigned int level) const
Definition point.h:111
virtual MPI_Comm get_communicator() const
virtual types::subdomain_id locally_owned_subdomain() const
unsigned int n_levels() const
unsigned int n_raw_cells(const unsigned int level) const
virtual unsigned int n_global_levels() const
unsigned int n_cells() const
unsigned int global_to_local(const types::global_dof_index global_index) const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
unsigned int level
Definition grid_out.cc:4626
IteratorRange< active_cell_iterator > active_cell_iterators_on_level(const unsigned int level) const
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
typename ActiveSelector::cell_iterator cell_iterator
@ mg_transfer_fill_copy_indices
mg_transfer_internal.cc: fill_copy_indices()
Definition mpi_tags.h:71
T sum(const T &t, const MPI_Comm mpi_communicator)
T min(const T &t, const MPI_Comm mpi_communicator)
std::vector< unsigned int > compute_index_owner(const IndexSet &owned_indices, const IndexSet &indices_to_look_up, const MPI_Comm comm)
Definition mpi.cc:713
void setup_element_info(ElementInfo< Number > &elem_info, const FiniteElement< 1 > &fe, const DoFHandler< dim > &dof_handler)
void resolve_identity_constraints(const MGConstrainedDoFs *mg_constrained_dofs, const unsigned int level, std::vector< types::global_dof_index > &dof_indices)
unsigned int compute_shift_within_children(const unsigned int child, const unsigned int fe_shift_1d, const unsigned int fe_degree)
void add_child_indices(const unsigned int child, const unsigned int fe_shift_1d, const unsigned int fe_degree, const std::vector< unsigned int > &lexicographic_numbering, const std::vector< types::global_dof_index > &local_dof_indices, types::global_dof_index *target_indices)
void setup_transfer(const DoFHandler< dim > &dof_handler, const MGConstrainedDoFs *mg_constrained_dofs, const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &external_partitioners, ElementInfo< Number > &elem_info, std::vector< std::vector< unsigned int > > &level_dof_indices, std::vector< std::vector< std::pair< unsigned int, unsigned int > > > &parent_child_connect, std::vector< unsigned int > &n_owned_level_cells, std::vector< std::vector< std::vector< unsigned short > > > &dirichlet_indices, std::vector< std::vector< Number > > &weights_on_refined, std::vector< Table< 2, unsigned int > > &copy_indices_global_mine, MGLevelObject< std::shared_ptr< const Utilities::MPI::Partitioner > > &vector_partitioners)
void fill_copy_indices(const DoFHandler< dim, spacedim > &dof_handler, const MGConstrainedDoFs *mg_constrained_dofs, std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > &copy_indices, std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > &copy_indices_global_mine, std::vector< std::vector< std::pair< types::global_dof_index, types::global_dof_index > > > &copy_indices_level_mine, const bool skip_interface_dofs=true)
void reinit_level_partitioner(const IndexSet &locally_owned, std::vector< types::global_dof_index > &ghosted_level_dofs, const std::shared_ptr< const Utilities::MPI::Partitioner > &external_partitioner, const MPI_Comm communicator, std::shared_ptr< const Utilities::MPI::Partitioner > &target_partitioner, Table< 2, unsigned int > &copy_indices_global_mine)
void copy_indices_to_mpi_local_numbers(const Utilities::MPI::Partitioner &part, const std::vector< types::global_dof_index > &mine, const std::vector< types::global_dof_index > &remote, std::vector< unsigned int > &localized_indices)
const types::subdomain_id artificial_subdomain_id
Definition types.h:362
const types::subdomain_id invalid_subdomain_id
Definition types.h:341
static const unsigned int invalid_unsigned_int
Definition types.h:220
const types::global_dof_index invalid_dof_index
Definition types.h:252
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
types::global_dof_index global_dof_index
types::global_dof_index level_dof_index
DoFPair(const unsigned int level, const types::global_dof_index global_dof_index, const types::global_dof_index level_dof_index)
std::vector< unsigned int > lexicographic_numbering
void reinit(const Quadrature< dim_q > &quad, const FiniteElement< dim, spacedim > &fe_dim, const unsigned int base_element=0)
std::vector< unsigned int > lexicographic_numbering
Definition shape_info.h:479