Reference documentation for deal.II version 9.6.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_global_coarsening.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2020 - 2024 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#ifndef dealii_mg_transfer_global_coarsening_h
16#define dealii_mg_transfer_global_coarsening_h
17
21
23
26
30
33
35
36
37
39
40// Forward declarations
41#ifndef DOXYGEN
42namespace internal
43{
44 class MGTwoLevelTransferImplementation;
45}
46
48{
49 template <int dim, int spacedim>
50 class Base;
51}
52
53template <int dim, typename Number>
54class MGTransferMF;
55#endif
56
57
58namespace mg
59{
71 {
77 boost::signals2::signal<void(const bool before)> prolongation_cell_loop;
78
84 boost::signals2::signal<void(const bool before)> restriction_cell_loop;
85
93 boost::signals2::signal<void(const bool before)> prolongation;
94
102 boost::signals2::signal<void(const bool before)> restriction;
103 };
104} // namespace mg
105
110{
136
141 unsigned int
143 const unsigned int degree,
144 const PolynomialCoarseningSequenceType &p_sequence);
145
151 std::vector<unsigned int>
153 const unsigned int max_degree,
154 const PolynomialCoarseningSequenceType &p_sequence);
155
166 template <int dim, int spacedim>
167 std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
169 const Triangulation<dim, spacedim> &tria);
170
187 template <int dim, int spacedim>
188 std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
192 const bool preserve_fine_triangulation,
193 const bool repartition_fine_triangulation);
194
200 template <int dim, int spacedim>
201 std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
205 const bool repartition_fine_triangulation = false);
206
207} // namespace MGTransferGlobalCoarseningTools
208
209
210
218template <typename VectorType>
220{
221public:
222 static_assert(
223 std::is_same_v<
224 VectorType,
226 "This class is currently only implemented for vectors of "
227 "type LinearAlgebra::distributed::Vector.");
228
232 using Number = typename VectorType::value_type;
233
238
242 void
243 prolongate_and_add(VectorType &dst, const VectorType &src) const;
244
248 void
249 restrict_and_add(VectorType &dst, const VectorType &src) const;
250
258 virtual void
259 interpolate(VectorType &dst, const VectorType &src) const = 0;
260
266 virtual std::pair<bool, bool>
268 const std::shared_ptr<const Utilities::MPI::Partitioner>
270 const std::shared_ptr<const Utilities::MPI::Partitioner>
271 &partitioner_fine) = 0;
272
276 virtual std::size_t
278
279protected:
283 virtual void
286 const LinearAlgebra::distributed::Vector<Number> &src) const = 0;
287
291 virtual void
294 const LinearAlgebra::distributed::Vector<Number> &src) const = 0;
295
301 void
302 update_ghost_values(const VectorType &vec) const;
303
309 void
310 compress(VectorType &vec, const VectorOperation::values op) const;
311
317 void
318 zero_out_ghost_values(const VectorType &vec) const;
319
324 template <int dim, std::size_t width, typename IndexType>
325 std::pair<bool, bool>
327 const std::shared_ptr<const Utilities::MPI::Partitioner>
329 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner_fine,
332 ConstraintInfo<dim, VectorizedArray<Number, width>, IndexType>
333 &constraint_info_coarse,
334 std::vector<unsigned int> &dof_indices_fine);
335
342
343public:
347 std::shared_ptr<const Utilities::MPI::Partitioner> partitioner_coarse;
348
352 std::shared_ptr<const Utilities::MPI::Partitioner> partitioner_fine;
353
354protected:
359
368
373
378 std::shared_ptr<const Utilities::MPI::Partitioner>
380
385 std::shared_ptr<const Utilities::MPI::Partitioner> partitioner_fine_embedded;
386
392
398};
399
400
401
427template <int dim, typename VectorType>
429{
430public:
431 static_assert(
432 std::is_same_v<
433 VectorType,
435 "This class is currently only implemented for vectors of "
436 "type LinearAlgebra::distributed::Vector.");
437
441 using Number = typename VectorType::value_type;
442
448
454 void
457 const DoFHandler<dim> &dof_handler_coarse,
458 const AffineConstraints<Number> &constraint_fine =
460 const AffineConstraints<Number> &constraint_coarse =
462 const unsigned int mg_level_fine = numbers::invalid_unsigned_int,
463 const unsigned int mg_level_coarse = numbers::invalid_unsigned_int);
464
474 void
477 const DoFHandler<dim> &dof_handler_coarse,
478 const AffineConstraints<Number> &constraint_fine =
480 const AffineConstraints<Number> &constraint_coarse =
482 const unsigned int mg_level_fine = numbers::invalid_unsigned_int,
483 const unsigned int mg_level_coarse = numbers::invalid_unsigned_int);
484
498 void
500 const DoFHandler<dim> &dof_handler_coarse,
501 const AffineConstraints<Number> &constraint_fine =
503 const AffineConstraints<Number> &constraint_coarse =
505 const unsigned int mg_level_fine = numbers::invalid_unsigned_int,
506 const unsigned int mg_level_coarse = numbers::invalid_unsigned_int);
507
517 void
518 reinit(const MatrixFree<dim, Number> &matrix_free_fine,
519 const unsigned int dof_no_fine,
520 const MatrixFree<dim, Number> &matrix_free_coarse,
521 const unsigned int dof_no_coarse);
522
531 static bool
532 fast_polynomial_transfer_supported(const unsigned int fe_degree_fine,
533 const unsigned int fe_degree_coarse);
534
538 void
539 interpolate(VectorType &dst, const VectorType &src) const override;
540
545 std::pair<bool, bool>
547 const std::shared_ptr<const Utilities::MPI::Partitioner>
549 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner_fine)
550 override;
551
555 std::size_t
556 memory_consumption() const override;
557
558protected:
559 void
561 const VectorType &src) const override;
562
563 void
565 const VectorType &src) const override;
566
567private:
627
631 std::vector<MGTransferScheme> schemes;
632
637 internal::MatrixFreeFunctions::
638 ConstraintInfo<dim, VectorizedArrayType, types::global_dof_index>
640
644 internal::MatrixFreeFunctions::
645 ConstraintInfo<dim, VectorizedArrayType, types::global_dof_index>
647
680
686 std::unique_ptr<MatrixFreeRelatedData> matrix_free_data;
687
692 std::vector<unsigned int> weights_start;
693
699
704 std::vector<unsigned char> weights_are_compressed;
705
709 unsigned int n_components;
710
715
719 unsigned int mg_level_fine;
720
722
723 friend class MGTransferMF<dim, Number>;
724};
725
726
727
731template <int dim, typename VectorType>
733{
734private:
735 static_assert(
736 std::is_same_v<
737 VectorType,
739 "This class is currently only implemented for vectors of "
740 "type LinearAlgebra::distributed::Vector.");
741
742 using Number = typename VectorType::value_type;
744
746
747public:
755 {
765 AdditionalData(const double tolerance = 1e-6,
766 const unsigned int rtree_level = 0,
767 const bool enforce_all_points_found = true)
771 {}
772
777 double tolerance;
778
784 unsigned int rtree_level;
785
794 };
795
800
805 void
807 const DoFHandler<dim> &dof_handler_coarse,
808 const Mapping<dim> &mapping_fine,
809 const Mapping<dim> &mapping_coarse,
810 const AffineConstraints<Number> &constraint_fine =
812 const AffineConstraints<Number> &constraint_coarse =
814
821 void
822 interpolate(VectorType &dst, const VectorType &src) const override;
823
828 std::pair<bool, bool>
830 const std::shared_ptr<const Utilities::MPI::Partitioner>
832 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner_fine)
833 override;
834
838 std::size_t
839 memory_consumption() const override;
840
844 boost::signals2::connection
845 connect_prolongation_cell_loop(const std::function<void(const bool)> &slot);
846
850 boost::signals2::connection
851 connect_restriction_cell_loop(const std::function<void(const bool)> &slot);
852
856 boost::signals2::connection
857 connect_prolongation(const std::function<void(const bool)> &slot);
858
862 boost::signals2::connection
863 connect_restriction(const std::function<void(const bool)> &slot);
864
865protected:
870 void
872 const VectorType &src) const override;
873
877 void
879 const VectorType &src) const override;
880
881private:
885 template <int n_components>
886 void
888 const VectorType &src) const;
889
893 template <int n_components>
894 void
895 restrict_and_add_internal_comp(VectorType &dst, const VectorType &src) const;
896
901
905 unsigned int mg_level_fine;
906
912
916 mutable std::vector<Number> rpe_input_output;
917
921 mutable std::vector<Number> rpe_buffer;
922
926 std::shared_ptr<NonMatching::MappingInfo<dim, dim, Number>> mapping_info;
927
932 internal::MatrixFreeFunctions::
933 ConstraintInfo<dim, VectorizedArrayType, unsigned int>
935
939 std::unique_ptr<FiniteElement<dim>> fe_coarse;
940
945 std::vector<unsigned int> level_dof_indices_fine;
946
952 std::vector<unsigned int> level_dof_indices_fine_ptrs;
953
954 friend class MGTransferMF<dim, Number>;
955};
956
957
958
987template <int dim, typename Number>
989 LinearAlgebra::distributed::Vector<Number>>
990{
991public:
996
1003
1018 template <typename MGTwoLevelTransferObject>
1020 const std::function<void(const unsigned int, VectorType &)>
1021 &initialize_dof_vector = {});
1022
1026 template <typename MGTwoLevelTransferObject>
1027 void
1030
1039 void
1040 build(const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
1041 &external_partitioners = {});
1042
1047 void
1048 build(const std::function<void(const unsigned int, VectorType &)>
1050
1065
1071 void
1073
1090 void
1091 build(const DoFHandler<dim> &dof_handler,
1092 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
1093 &external_partitioners = {});
1094
1104 void
1105 build(const DoFHandler<dim> &dof_handler,
1106 const std::function<void(const unsigned int, VectorType &)>
1108
1119 void
1120 prolongate(const unsigned int to_level,
1121 VectorType &dst,
1122 const VectorType &src) const override;
1123
1127 void
1128 prolongate_and_add(const unsigned int to_level,
1129 VectorType &dst,
1130 const VectorType &src) const override;
1131
1135 virtual void
1136 restrict_and_add(const unsigned int from_level,
1137 VectorType &dst,
1138 const VectorType &src) const override;
1139
1149 template <class InVector>
1150 void
1151 copy_to_mg(const DoFHandler<dim> &dof_handler,
1153 const InVector &src) const;
1154
1164 template <class OutVector>
1165 void
1166 copy_from_mg(const DoFHandler<dim> &dof_handler,
1167 OutVector &dst,
1168 const MGLevelObject<VectorType> &src) const;
1169
1188 template <class InVector>
1189 void
1192 const InVector &src) const;
1193
1202 template <class InVector>
1203 void
1204 interpolate_to_mg(MGLevelObject<VectorType> &dst, const InVector &src) const;
1205
1219 std::size_t
1221
1225 unsigned int
1226 min_level() const;
1227
1231 unsigned int
1232 max_level() const;
1233
1238 void
1240
1243private:
1249 void
1251 const DoFHandler<dim> &dof_handler,
1253
1257 std::pair<const DoFHandler<dim> *, unsigned int>
1259
1265 void
1267 const DoFHandler<dim> &dof_handler);
1268
1272 template <typename MGTwoLevelTransferObject>
1273 void
1276
1280 template <class InVector>
1281 void
1282 initialize_dof_vector(const unsigned int level,
1283 VectorType &vector,
1284 const InVector &vector_reference,
1285 const bool omit_zeroing_entries) const;
1286
1292 void
1293 assert_dof_handler(const DoFHandler<dim> &dof_handler_out) const;
1294
1301
1306
1310 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
1312};
1313
1314
1315
1321template <int dim, typename Number>
1323 : public MGTransferBlockMatrixFreeBase<dim, Number, MGTransferMF<dim, Number>>
1324{
1325public:
1330
1337
1343 MGTransferBlockMF(const MGConstrainedDoFs &mg_constrained_dofs);
1344
1350 MGTransferBlockMF(const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
1351
1357 void
1358 initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs);
1359
1365 void
1367 const std::vector<MGConstrainedDoFs> &mg_constrained_dofs);
1368
1374 void
1375 build(const DoFHandler<dim> &dof_handler);
1376
1382 void
1383 build(const std::vector<const DoFHandler<dim> *> &dof_handler);
1384
1385protected:
1387 get_matrix_free_transfer(const unsigned int b) const override;
1388
1389private:
1393 std::vector<MGTransferMF<dim, Number>> transfer_operators_internal;
1394
1398 std::vector<SmartPointer<const MGTransferMF<dim, Number>>> transfer_operators;
1399};
1400
1401
1402
1403template <int dim, typename VectorType>
1406
1407template <int dim, typename VectorType>
1410
1411
1412
1413#ifndef DOXYGEN
1414
1415/* ----------------------- Inline functions --------------------------------- */
1416
1417
1418
1419template <int dim, typename Number>
1420template <typename MGTwoLevelTransferObject>
1423 const std::function<void(const unsigned int, VectorType &)>
1424 &initialize_dof_vector)
1425{
1426 this->transfer.clear();
1427 this->internal_transfer.clear();
1428
1429 this->initialize_transfer_references(transfer);
1430 this->build(initialize_dof_vector);
1431}
1432
1433
1434
1435template <int dim, typename Number>
1436template <typename MGTwoLevelTransferObject>
1437void
1440{
1441 this->initialize_transfer_references(transfer);
1442}
1443
1444
1445
1446template <int dim, typename Number>
1447template <typename MGTwoLevelTransferObject>
1448void
1451{
1452 const unsigned int min_level = transfer.min_level();
1453 const unsigned int max_level = transfer.max_level();
1454
1455 this->transfer.resize(min_level, max_level);
1456
1457 for (unsigned int l = min_level; l <= max_level; ++l)
1458 this->transfer[l] = &const_cast<MGTwoLevelTransferBase<VectorType> &>(
1459 static_cast<const MGTwoLevelTransferBase<VectorType> &>(
1460 Utilities::get_underlying_value(transfer[l])));
1461}
1462
1463
1464
1465template <int dim, typename Number>
1466template <class InVector>
1467void
1469 const unsigned int level,
1470 VectorType &vec,
1471 const InVector &vec_reference,
1472 const bool omit_zeroing_entries) const
1473{
1474 std::shared_ptr<const Utilities::MPI::Partitioner> partitioner;
1475
1476 if (external_partitioners.empty())
1477 {
1478 partitioner = vec_reference.get_partitioner();
1479 }
1480 else
1481 {
1482 Assert(transfer.min_level() <= level && level <= transfer.max_level(),
1484
1485 partitioner = external_partitioners[level - transfer.min_level()];
1486 }
1487
1488 // check if vectors are already correctly initialized
1489
1490 // yes: same partitioners are used
1491 if (vec.get_partitioner().get() == partitioner.get())
1492 {
1493 if (omit_zeroing_entries == false)
1494 vec = 0;
1495 return; // nothing to do
1496 }
1497
1498 // yes: vectors are compatible
1499 if (vec.size() == partitioner->size() &&
1500 vec.locally_owned_size() == partitioner->locally_owned_size())
1501 {
1502 if (omit_zeroing_entries == false)
1503 vec = 0;
1504 return; // nothing to do
1505 }
1506
1507 // no
1508 vec.reinit(partitioner, omit_zeroing_entries);
1509}
1510
1511
1512
1513template <int dim, typename Number>
1514template <class InVector>
1515void
1518 const InVector &src) const
1519{
1520 assert_dof_handler(dof_handler);
1521
1522 for (unsigned int level = dst.min_level(); level <= dst.max_level(); ++level)
1523 {
1524 const bool zero_out_values =
1525 (this->perform_plain_copy == false &&
1526 this->perform_renumbered_plain_copy == false) ||
1527 level != dst.max_level();
1528
1529 this->initialize_dof_vector(level, dst[level], src, !zero_out_values);
1530 }
1531
1532 if (this->perform_plain_copy)
1533 {
1534 dst[dst.max_level()].copy_locally_owned_data_from(src);
1535 }
1536 else if (this->perform_renumbered_plain_copy)
1537 {
1538 auto &dst_level = dst[dst.max_level()];
1539
1540 for (unsigned int i = 0; i < this->copy_indices.back().n_cols(); ++i)
1541 dst_level.local_element(this->copy_indices.back()(1, i)) =
1542 src.local_element(i);
1543 }
1544 else
1545 {
1546 this->ghosted_global_vector = src;
1547 this->ghosted_global_vector.update_ghost_values();
1548
1549 for (unsigned int l = dst.max_level() + 1; l != dst.min_level();)
1550 {
1551 --l;
1552
1553 auto &dst_level = dst[l];
1554
1555 const auto copy_unknowns = [&](const auto &indices) {
1556 for (unsigned int i = 0; i < indices.n_cols(); ++i)
1557 dst_level.local_element(indices(1, i)) =
1558 this->ghosted_global_vector.local_element(indices(0, i));
1559 };
1560
1561 copy_unknowns(this->copy_indices[l]);
1562 copy_unknowns(this->copy_indices_level_mine[l]);
1563
1564 dst_level.compress(VectorOperation::insert);
1565 }
1566 }
1567}
1568
1569
1570
1571template <int dim, typename Number>
1572template <class OutVector>
1573void
1575 const DoFHandler<dim> &dof_handler,
1576 OutVector &dst,
1577 const MGLevelObject<VectorType> &src) const
1578{
1579 assert_dof_handler(dof_handler);
1580
1581 if (this->perform_plain_copy)
1582 {
1583 dst.zero_out_ghost_values();
1584 dst.copy_locally_owned_data_from(src[src.max_level()]);
1585 }
1586 else if (this->perform_renumbered_plain_copy)
1587 {
1588 const auto &src_level = src[src.max_level()];
1589 dst.zero_out_ghost_values();
1590 for (unsigned int i = 0; i < this->copy_indices.back().n_cols(); ++i)
1591 dst.local_element(i) =
1592 src_level.local_element(this->copy_indices.back()(1, i));
1593 }
1594 else
1595 {
1596 dst = 0;
1597 for (unsigned int l = src.min_level(); l <= src.max_level(); ++l)
1598 {
1599 auto &ghosted_vector = this->ghosted_level_vector[l];
1600
1601 if (this->ghosted_level_vector[l].size() > 0)
1602 ghosted_vector = src[l];
1603
1604 const auto *const ghosted_vector_ptr =
1605 (this->ghosted_level_vector[l].size() > 0) ? &ghosted_vector :
1606 &src[l];
1607
1608 ghosted_vector_ptr->update_ghost_values();
1609
1610 const auto copy_unknowns = [&](const auto &indices) {
1611 for (unsigned int i = 0; i < indices.n_cols(); ++i)
1612 dst.local_element(indices(0, i)) =
1613 ghosted_vector_ptr->local_element(indices(1, i));
1614 };
1615
1616 copy_unknowns(this->copy_indices[l]);
1617 copy_unknowns(this->copy_indices_global_mine[l]);
1618 }
1619 dst.compress(VectorOperation::insert);
1620 }
1621}
1622
1623
1624
1625template <int dim, typename Number>
1626template <class InVector>
1627void
1629 const InVector &src) const
1630{
1631 DoFHandler<dim> dof_handler_dummy;
1632
1633 this->interpolate_to_mg(dof_handler_dummy, dst, src);
1634}
1635
1636
1637
1638template <int dim, typename Number>
1639template <class InVector>
1640void
1643 const InVector &src) const
1644{
1645 assert_dof_handler(dof_handler);
1646
1647 const unsigned int min_level = transfer.min_level();
1648 const unsigned int max_level = transfer.max_level();
1649
1650 AssertDimension(min_level, dst.min_level());
1651 AssertDimension(max_level, dst.max_level());
1652
1653 for (unsigned int level = min_level; level <= max_level; ++level)
1654 {
1655 const bool zero_out_values = false;
1656 this->initialize_dof_vector(level, dst[level], src, !zero_out_values);
1657 }
1658
1659 if (this->perform_plain_copy)
1660 {
1661 dst[max_level].copy_locally_owned_data_from(src);
1662
1663 for (unsigned int l = max_level; l > min_level; --l)
1664 this->transfer[l]->interpolate(dst[l - 1], dst[l]);
1665 }
1666 else if (this->perform_renumbered_plain_copy)
1667 {
1668 auto &dst_level = dst[max_level];
1669
1670 for (unsigned int i = 0; i < this->solution_copy_indices.back().n_cols();
1671 ++i)
1672 dst_level.local_element(this->solution_copy_indices.back()(1, i)) =
1673 src.local_element(i);
1674
1675 for (unsigned int l = max_level; l > min_level; --l)
1676 this->transfer[l]->interpolate(dst[l - 1], dst[l]);
1677 }
1678 else
1679 {
1680 this->solution_ghosted_global_vector = src;
1681 this->solution_ghosted_global_vector.update_ghost_values();
1682
1683 for (unsigned int l = max_level + 1; l != min_level;)
1684 {
1685 --l;
1686
1687 auto &dst_level = dst[l];
1688
1689 const auto copy_unknowns = [&](const auto &indices) {
1690 for (unsigned int i = 0; i < indices.n_cols(); ++i)
1691 dst_level.local_element(indices(1, i)) =
1692 this->solution_ghosted_global_vector.local_element(
1693 indices(0, i));
1694 };
1695
1696 copy_unknowns(this->solution_copy_indices[l]);
1697 copy_unknowns(this->solution_copy_indices_level_mine[l]);
1698
1699 dst_level.compress(VectorOperation::insert);
1700
1701 if (l != min_level)
1702 this->transfer[l]->interpolate(dst[l - 1], dst[l]);
1703 }
1704 }
1705}
1706
1707#endif
1708
1710
1711#endif
std::function< void(const unsigned int, LinearAlgebra::distributed::Vector< Number > &)> initialize_dof_vector
void resize(const unsigned int new_minlevel, const unsigned int new_maxlevel, Args &&...args)
unsigned int max_level() const
unsigned int min_level() const
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
MGTransferBlockMF(const std::vector< MGConstrainedDoFs > &mg_constrained_dofs)
MGTransferBlockMF(const MGConstrainedDoFs &mg_constrained_dofs)
void build(const std::vector< const DoFHandler< dim > * > &dof_handler)
std::vector< MGTransferMF< dim, Number > > transfer_operators_internal
MGTransferBlockMF(const MGTransferMF< dim, Number > &transfer_operator)
MGTransferBlockMF()=default
void build(const DoFHandler< dim > &dof_handler)
const MGTransferMF< dim, Number > & get_matrix_free_transfer(const unsigned int b) const override
std::vector< SmartPointer< const MGTransferMF< dim, Number > > > transfer_operators
void initialize_constraints(const std::vector< MGConstrainedDoFs > &mg_constrained_dofs)
std::pair< const DoFHandler< dim > *, unsigned int > get_dof_handler_fine() const
void copy_to_mg(const DoFHandler< dim > &dof_handler, MGLevelObject< VectorType > &dst, const InVector &src) const
unsigned int min_level() const
void build(const std::function< void(const unsigned int, VectorType &)> &initialize_dof_vector)
void assert_dof_handler(const DoFHandler< dim > &dof_handler_out) const
void interpolate_to_mg(const DoFHandler< dim > &dof_handler, MGLevelObject< VectorType > &dst, const InVector &src) const
std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > external_partitioners
void prolongate_and_add(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
unsigned int max_level() const
MGTransferMF(const MGLevelObject< MGTwoLevelTransferObject > &transfer, const std::function< void(const unsigned int, VectorType &)> &initialize_dof_vector={})
void build(const DoFHandler< dim > &dof_handler, const std::function< void(const unsigned int, VectorType &)> &initialize_dof_vector)
void build(const DoFHandler< dim > &dof_handler, const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &external_partitioners={})
void build(const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &external_partitioners={})
void initialize_dof_vector(const unsigned int level, VectorType &vector, const InVector &vector_reference, const bool omit_zeroing_entries) const
void initialize_constraints(const MGConstrainedDoFs &mg_constrained_dofs)
void initialize_internal_transfer(const DoFHandler< dim > &dof_handler, const SmartPointer< const MGConstrainedDoFs > &mg_constrained_dofs)
void initialize_two_level_transfers(const MGLevelObject< MGTwoLevelTransferObject > &transfer)
MGTransferMF(const MGConstrainedDoFs &mg_constrained_dofs)
MGLevelObject< MGTwoLevelTransfer< dim, VectorType > > internal_transfer
void fill_and_communicate_copy_indices_global_coarsening(const DoFHandler< dim > &dof_handler)
MGLevelObject< SmartPointer< MGTwoLevelTransferBase< VectorType > > > transfer
std::size_t memory_consumption() const
virtual void restrict_and_add(const unsigned int from_level, VectorType &dst, const VectorType &src) const override
void initialize_transfer_references(const MGLevelObject< MGTwoLevelTransferObject > &transfer)
void prolongate(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
void copy_from_mg(const DoFHandler< dim > &dof_handler, OutVector &dst, const MGLevelObject< VectorType > &src) const
void interpolate_to_mg(MGLevelObject< VectorType > &dst, const InVector &src) const
AlignedVector< Number > buffer_coarse_embedded
void zero_out_ghost_values(const VectorType &vec) const
std::shared_ptr< const Utilities::MPI::Partitioner > partitioner_fine_embedded
std::shared_ptr< const Utilities::MPI::Partitioner > partitioner_coarse
virtual std::size_t memory_consumption() const =0
AlignedVector< Number > buffer_fine_embedded
void prolongate_and_add(VectorType &dst, const VectorType &src) const
LinearAlgebra::distributed::Vector< Number > vec_fine
void update_ghost_values(const VectorType &vec) const
LinearAlgebra::distributed::Vector< Number > vec_coarse
typename VectorType::value_type Number
virtual void restrict_and_add_internal(LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const =0
virtual std::pair< bool, bool > enable_inplace_operations_if_possible(const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner_coarse, const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner_fine)=0
void restrict_and_add(VectorType &dst, const VectorType &src) const
std::shared_ptr< const Utilities::MPI::Partitioner > partitioner_coarse_embedded
std::pair< bool, bool > internal_enable_inplace_operations_if_possible(const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner_coarse, const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner_fine, bool &vec_fine_needs_ghost_update, internal::MatrixFreeFunctions::ConstraintInfo< dim, VectorizedArray< Number, width >, IndexType > &constraint_info_coarse, std::vector< unsigned int > &dof_indices_fine)
std::shared_ptr< const Utilities::MPI::Partitioner > partitioner_fine
virtual void interpolate(VectorType &dst, const VectorType &src) const =0
void compress(VectorType &vec, const VectorOperation::values op) const
virtual void prolongate_and_add_internal(LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const =0
std::shared_ptr< NonMatching::MappingInfo< dim, dim, Number > > mapping_info
boost::signals2::connection connect_restriction(const std::function< void(const bool)> &slot)
void reinit(const DoFHandler< dim > &dof_handler_fine, const DoFHandler< dim > &dof_handler_coarse, const Mapping< dim > &mapping_fine, const Mapping< dim > &mapping_coarse, const AffineConstraints< Number > &constraint_fine=AffineConstraints< Number >(), const AffineConstraints< Number > &constraint_coarse=AffineConstraints< Number >())
std::unique_ptr< FiniteElement< dim > > fe_coarse
void restrict_and_add_internal_comp(VectorType &dst, const VectorType &src) const
boost::signals2::connection connect_prolongation_cell_loop(const std::function< void(const bool)> &slot)
std::size_t memory_consumption() const override
std::pair< bool, bool > enable_inplace_operations_if_possible(const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner_coarse, const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner_fine) override
Utilities::MPI::RemotePointEvaluation< dim > rpe
void restrict_and_add_internal(VectorType &dst, const VectorType &src) const override
std::vector< unsigned int > level_dof_indices_fine_ptrs
internal::MatrixFreeFunctions::ConstraintInfo< dim, VectorizedArrayType, unsigned int > constraint_info
boost::signals2::connection connect_prolongation(const std::function< void(const bool)> &slot)
boost::signals2::connection connect_restriction_cell_loop(const std::function< void(const bool)> &slot)
void prolongate_and_add_internal_comp(VectorType &dst, const VectorType &src) const
MGTwoLevelTransferNonNested(const AdditionalData &data=AdditionalData())
void interpolate(VectorType &dst, const VectorType &src) const override
std::vector< unsigned int > level_dof_indices_fine
typename VectorType::value_type Number
void prolongate_and_add_internal(VectorType &dst, const VectorType &src) const override
SmartPointer< const DoFHandler< dim > > dof_handler_fine
SmartPointer< const DoFHandler< dim > > dof_handler_fine
void reinit_geometric_transfer(const DoFHandler< dim > &dof_handler_fine, const DoFHandler< dim > &dof_handler_coarse, const AffineConstraints< Number > &constraint_fine=AffineConstraints< Number >(), const AffineConstraints< Number > &constraint_coarse=AffineConstraints< Number >(), const unsigned int mg_level_fine=numbers::invalid_unsigned_int, const unsigned int mg_level_coarse=numbers::invalid_unsigned_int)
internal::MatrixFreeFunctions::ConstraintInfo< dim, VectorizedArrayType, types::global_dof_index > constraint_info_fine
friend class internal::MGTwoLevelTransferImplementation
std::vector< MGTransferScheme > schemes
typename VectorType::value_type Number
std::vector< unsigned char > weights_are_compressed
void interpolate(VectorType &dst, const VectorType &src) const override
internal::MatrixFreeFunctions::ConstraintInfo< dim, VectorizedArrayType, types::global_dof_index > constraint_info_coarse
std::unique_ptr< MatrixFreeRelatedData > matrix_free_data
static bool fast_polynomial_transfer_supported(const unsigned int fe_degree_fine, const unsigned int fe_degree_coarse)
AlignedVector< VectorizedArrayType > weights
void reinit_polynomial_transfer(const DoFHandler< dim > &dof_handler_fine, const DoFHandler< dim > &dof_handler_coarse, const AffineConstraints< Number > &constraint_fine=AffineConstraints< Number >(), const AffineConstraints< Number > &constraint_coarse=AffineConstraints< Number >(), const unsigned int mg_level_fine=numbers::invalid_unsigned_int, const unsigned int mg_level_coarse=numbers::invalid_unsigned_int)
std::pair< bool, bool > enable_inplace_operations_if_possible(const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner_coarse, const std::shared_ptr< const Utilities::MPI::Partitioner > &partitioner_fine) override
void restrict_and_add_internal(VectorType &dst, const VectorType &src) const override
std::vector< unsigned int > weights_start
void reinit(const MatrixFree< dim, Number > &matrix_free_fine, const unsigned int dof_no_fine, const MatrixFree< dim, Number > &matrix_free_coarse, const unsigned int dof_no_coarse)
void reinit(const DoFHandler< dim > &dof_handler_fine, const DoFHandler< dim > &dof_handler_coarse, const AffineConstraints< Number > &constraint_fine=AffineConstraints< Number >(), const AffineConstraints< Number > &constraint_coarse=AffineConstraints< Number >(), const unsigned int mg_level_fine=numbers::invalid_unsigned_int, const unsigned int mg_level_coarse=numbers::invalid_unsigned_int)
std::size_t memory_consumption() const override
void prolongate_and_add_internal(VectorType &dst, const VectorType &src) const override
Abstract base class for mapping classes.
Definition mapping.h:318
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
unsigned int level
Definition grid_out.cc:4626
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
void interpolate(const DoFHandler< dim, spacedim > &dof1, const InVector &u1, const DoFHandler< dim, spacedim > &dof2, OutVector &u2)
std::vector< std::shared_ptr< const Triangulation< dim, spacedim > > > create_geometric_coarsening_sequence(const Triangulation< dim, spacedim > &tria)
unsigned int create_next_polynomial_coarsening_degree(const unsigned int degree, const PolynomialCoarseningSequenceType &p_sequence)
std::vector< unsigned int > create_polynomial_coarsening_sequence(const unsigned int max_degree, const PolynomialCoarseningSequenceType &p_sequence)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
T & get_underlying_value(T &p)
Definition utilities.h:1639
Definition mg.h:81
static const unsigned int invalid_unsigned_int
Definition types.h:220
STL namespace.
AdditionalData(const double tolerance=1e-6, const unsigned int rtree_level=0, const bool enforce_all_points_found=true)
internal::MatrixFreeFunctions::ShapeInfo< double > shape_info_coarse
std::vector< std::array< unsigned int, VectorizedArrayType::size()> > cell_list_fine_to_coarse
SmartPointer< const MatrixFree< dim, Number > > matrix_free_fine
SmartPointer< const MatrixFree< dim, Number > > matrix_free_coarse
boost::signals2::signal< void(const bool before)> prolongation_cell_loop
boost::signals2::signal< void(const bool before)> restriction
boost::signals2::signal< void(const bool before)> restriction_cell_loop
boost::signals2::signal< void(const bool before)> prolongation