Reference documentation for deal.II version 9.4.1
\(\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// Copyright (C) 2020 - 2022 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#ifndef dealii_mg_transfer_global_coarsening_h
17#define dealii_mg_transfer_global_coarsening_h
18
21
23
26
29
32
34
35// Forward declarations
36#ifndef DOXYGEN
37namespace internal
38{
39 class MGTwoLevelTransferImplementation;
40}
41
43{
44 template <int dim, int spacedim>
45 class Base;
46}
47#endif
48
49
50
55{
64 {
69 bisect,
80 };
81
86 unsigned int
88 const unsigned int degree,
89 const PolynomialCoarseningSequenceType &p_sequence);
90
96 std::vector<unsigned int>
98 const unsigned int max_degree,
99 const PolynomialCoarseningSequenceType &p_sequence);
100
111 template <int dim, int spacedim>
112 std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
115
132 template <int dim, int spacedim>
133 std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
137 const bool preserve_fine_triangulation,
138 const bool repartition_fine_triangulation);
139
145 template <int dim, int spacedim>
146 std::vector<std::shared_ptr<const Triangulation<dim, spacedim>>>
150 const bool repartition_fine_triangulation = false);
151
152} // namespace MGTransferGlobalCoarseningTools
153
154
155
159template <int dim, typename VectorType>
161{
162public:
166 void
167 prolongate_and_add(VectorType &dst, const VectorType &src) const;
168
172 void
173 restrict_and_add(VectorType &dst, const VectorType &src) const;
174
181 void
182 interpolate(VectorType &dst, const VectorType &src) const;
183
188 void
190 const std::shared_ptr<const Utilities::MPI::Partitioner>
191 &partitioner_coarse,
192 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner_fine);
193
197 std::size_t
199};
200
201
202
207template <int dim, typename Number>
208class MGTwoLevelTransfer<dim, LinearAlgebra::distributed::Vector<Number>>
209{
210public:
216 void
218 const DoFHandler<dim> & dof_handler_fine,
219 const DoFHandler<dim> & dof_handler_coarse,
220 const AffineConstraints<Number> &constraint_fine =
222 const AffineConstraints<Number> &constraint_coarse =
224 const unsigned int mg_level_fine = numbers::invalid_unsigned_int,
225 const unsigned int mg_level_coarse = numbers::invalid_unsigned_int);
226
236 void
238 const DoFHandler<dim> & dof_handler_fine,
239 const DoFHandler<dim> & dof_handler_coarse,
240 const AffineConstraints<Number> &constraint_fine =
242 const AffineConstraints<Number> &constraint_coarse =
244 const unsigned int mg_level_fine = numbers::invalid_unsigned_int,
245 const unsigned int mg_level_coarse = numbers::invalid_unsigned_int);
246
260 void
261 reinit(const DoFHandler<dim> & dof_handler_fine,
262 const DoFHandler<dim> & dof_handler_coarse,
263 const AffineConstraints<Number> &constraint_fine =
265 const AffineConstraints<Number> &constraint_coarse =
267 const unsigned int mg_level_fine = numbers::invalid_unsigned_int,
268 const unsigned int mg_level_coarse = numbers::invalid_unsigned_int);
269
278 static bool
279 fast_polynomial_transfer_supported(const unsigned int fe_degree_fine,
280 const unsigned int fe_degree_coarse);
281
285 void
289
293 void
296
303 void
306
311 void
313 const std::shared_ptr<const Utilities::MPI::Partitioner>
314 &partitioner_coarse,
315 const std::shared_ptr<const Utilities::MPI::Partitioner> &partitioner_fine);
316
320 std::size_t
322
323private:
331 struct MGTransferScheme
332 {
336 unsigned int n_coarse_cells;
337
345
353
357 unsigned int degree_coarse;
358
364 unsigned int degree_fine;
365
370
375
380
385
392 };
393
397 std::vector<MGTransferScheme> schemes;
398
405
409 std::shared_ptr<const Utilities::MPI::Partitioner> partitioner_fine;
410
414 std::shared_ptr<const Utilities::MPI::Partitioner> partitioner_coarse;
415
424
429
436
440 std::vector<Number> weights;
441
447
452 std::vector<unsigned int> level_dof_indices_fine;
453
457 unsigned int n_components;
458
459 friend class internal::MGTwoLevelTransferImplementation;
460};
461
462
463
474template <int dim, typename VectorType>
476{
477public:
481 using Number = typename VectorType::value_type;
482
492 const std::function<void(const unsigned int, VectorType &)>
494
503 void
504 build(const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
506
511 void
512 build(const std::function<void(const unsigned int, VectorType &)>
514
518 void
519 prolongate(const unsigned int to_level,
520 VectorType & dst,
521 const VectorType & src) const override;
522
526 void
527 prolongate_and_add(const unsigned int to_level,
528 VectorType & dst,
529 const VectorType & src) const override;
530
534 virtual void
535 restrict_and_add(const unsigned int from_level,
536 VectorType & dst,
537 const VectorType & src) const override;
538
545 template <class InVector, int spacedim>
546 void
549 const InVector & src) const;
550
557 template <class OutVector, int spacedim>
558 void
560 OutVector & dst,
561 const MGLevelObject<VectorType> &src) const;
562
575 template <class InVector>
576 void
577 interpolate_to_mg(MGLevelObject<VectorType> &dst, const InVector &src) const;
578
585 template <class InVector, int spacedim>
586 void
589 const InVector & src) const;
590
597 std::size_t
599
600private:
604 void
605 initialize_dof_vector(const unsigned int level, VectorType &vector) const;
606
611
615 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
617};
618
619
620
626template <int dim, typename VectorType>
629 dim,
630 typename VectorType::value_type,
631 MGTransferGlobalCoarsening<dim, VectorType>>
632{
633public:
639
640protected:
642 get_matrix_free_transfer(const unsigned int b) const override;
643
644private:
649};
650
651
652
653#ifndef DOXYGEN
654
655/* ----------------------- Inline functions --------------------------------- */
656
657
658
659template <int dim, typename VectorType>
662 const std::function<void(const unsigned int, VectorType &)>
663 &initialize_dof_vector)
664 : transfer(transfer)
665{
666 this->build(initialize_dof_vector);
667}
668
669
670
671template <int dim, typename VectorType>
672void
674 const std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
675 &external_partitioners)
676{
677 this->external_partitioners = external_partitioners;
678
679 if (this->external_partitioners.size() > 0)
680 {
681 const unsigned int min_level = transfer.min_level();
682 const unsigned int max_level = transfer.max_level();
683
684 AssertDimension(this->external_partitioners.size(), transfer.n_levels());
685
686 for (unsigned int l = min_level + 1; l <= max_level; ++l)
687 transfer[l].enable_inplace_operations_if_possible(
688 this->external_partitioners[l - 1 - min_level],
689 this->external_partitioners[l - min_level]);
690 }
691}
692
693
694
695template <int dim, typename VectorType>
696void
698 const std::function<void(const unsigned int, VectorType &)>
699 &initialize_dof_vector)
700{
701 if (initialize_dof_vector)
702 {
703 const unsigned int min_level = transfer.min_level();
704 const unsigned int max_level = transfer.max_level();
705 const unsigned int n_levels = transfer.n_levels();
706
707 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
708 external_partitioners(n_levels);
709
710 for (unsigned int l = min_level; l <= max_level; ++l)
711 {
713 vector;
714 initialize_dof_vector(l, vector);
715 external_partitioners[l - min_level] = vector.get_partitioner();
716 }
717
718 this->build(external_partitioners);
719 }
720 else
721 {
722 this->build();
723 }
724}
725
726
727
728template <int dim, typename VectorType>
729void
731 const unsigned int level,
732 VectorType & vec) const
733{
734 AssertDimension(transfer.n_levels(), external_partitioners.size());
735 AssertIndexRange(level, transfer.max_level() + 1);
736
737 const auto &external_partitioner =
738 external_partitioners[level - transfer.min_level()];
739
740 if (vec.get_partitioner().get() != external_partitioner.get())
741 vec.reinit(external_partitioner);
742}
743
744
745
746template <int dim, typename VectorType>
747void
749 const unsigned int to_level,
750 VectorType & dst,
751 const VectorType & src) const
752{
753 dst = 0;
754 prolongate_and_add(to_level, dst, src);
755}
756
757
758
759template <int dim, typename VectorType>
760void
762 const unsigned int to_level,
763 VectorType & dst,
764 const VectorType & src) const
765{
766 this->transfer[to_level].prolongate_and_add(dst, src);
767}
768
769
770
771template <int dim, typename VectorType>
772void
774 const unsigned int from_level,
775 VectorType & dst,
776 const VectorType & src) const
777{
778 this->transfer[from_level].restrict_and_add(dst, src);
779}
780
781
782
783template <int dim, typename VectorType>
784template <class InVector, int spacedim>
785void
787 const DoFHandler<dim, spacedim> &dof_handler,
789 const InVector & src) const
790{
791 (void)dof_handler;
792
793 for (unsigned int level = dst.min_level(); level <= dst.max_level(); ++level)
794 {
795 initialize_dof_vector(level, dst[level]);
796
797 if (level == dst.max_level())
798 dst[level].copy_locally_owned_data_from(src);
799 else
800 dst[level] = 0.0;
801 }
802}
803
804
805
806template <int dim, typename VectorType>
807template <class OutVector, int spacedim>
808void
810 const DoFHandler<dim, spacedim> &dof_handler,
811 OutVector & dst,
812 const MGLevelObject<VectorType> &src) const
813{
814 (void)dof_handler;
815
816 dst.copy_locally_owned_data_from(src[src.max_level()]);
817}
818
819
820
821template <int dim, typename VectorType>
822template <class InVector>
823void
826 const InVector & src) const
827{
828 const unsigned int min_level = transfer.min_level();
829 const unsigned int max_level = transfer.max_level();
830
831 AssertDimension(min_level, dst.min_level());
832 AssertDimension(max_level, dst.max_level());
833
834 for (unsigned int level = min_level; level <= max_level; ++level)
835 initialize_dof_vector(level, dst[level]);
836
837 dst[transfer.max_level()].copy_locally_owned_data_from(src);
838
839 for (unsigned int l = max_level; l > min_level; --l)
840 this->transfer[l].interpolate(dst[l - 1], dst[l]);
841}
842
843
844
845template <int dim, typename VectorType>
846template <class InVector, int spacedim>
847void
849 const DoFHandler<dim, spacedim> &dof_handler,
851 const InVector & src) const
852{
853 (void)dof_handler;
854
855 this->interpolate_to_mg(dst, src);
856}
857
858
859
860template <int dim, typename VectorType>
861std::size_t
863{
864 std::size_t size = 0;
865
866 const unsigned int min_level = transfer.min_level();
867 const unsigned int max_level = transfer.max_level();
868
869 for (unsigned int l = min_level + 1; l <= max_level; ++l)
870 size += this->transfer[l].memory_consumption();
871
872 return size;
873}
874
875#endif
876
878
879#endif
unsigned int max_level() const
unsigned int min_level() const
const MGTransferGlobalCoarsening< dim, VectorType > & get_matrix_free_transfer(const unsigned int b) const override
MGTransferBlockGlobalCoarsening(const MGTransferGlobalCoarsening< dim, VectorType > &transfer_operator)
const MGTransferGlobalCoarsening< dim, VectorType > & transfer_operator
MGTransferGlobalCoarsening(const MGLevelObject< MGTwoLevelTransfer< dim, VectorType > > &transfer, const std::function< void(const unsigned int, VectorType &)> &initialize_dof_vector={})
void copy_from_mg(const DoFHandler< dim, spacedim > &dof_handler, OutVector &dst, const MGLevelObject< VectorType > &src) const
void interpolate_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< VectorType > &dst, const InVector &src) const
void prolongate(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
typename VectorType::value_type Number
std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > external_partitioners
void copy_to_mg(const DoFHandler< dim, spacedim > &dof_handler, MGLevelObject< VectorType > &dst, const InVector &src) const
virtual void restrict_and_add(const unsigned int from_level, VectorType &dst, const VectorType &src) const override
std::size_t memory_consumption() const
void initialize_dof_vector(const unsigned int level, VectorType &vector) const
void build(const std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > &external_partitioners={})
void build(const std::function< void(const unsigned int, VectorType &)> &initialize_dof_vector)
void interpolate_to_mg(MGLevelObject< VectorType > &dst, const InVector &src) const
MGLevelObject< MGTwoLevelTransfer< dim, VectorType > > transfer
void prolongate_and_add(const unsigned int to_level, VectorType &dst, const VectorType &src) const override
void 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)
void prolongate_and_add(LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
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)
static bool fast_polynomial_transfer_supported(const unsigned int fe_degree_fine, const unsigned int fe_degree_coarse)
void interpolate(LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
internal::MatrixFreeFunctions::ConstraintInfo< dim, VectorizedArray< Number > > constraint_info
void restrict_and_add(LinearAlgebra::distributed::Vector< Number > &dst, const LinearAlgebra::distributed::Vector< Number > &src) const
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)
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)
void restrict_and_add(VectorType &dst, const VectorType &src) const
void 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)
std::size_t memory_consumption() const
void interpolate(VectorType &dst, const VectorType &src) const
void prolongate_and_add(VectorType &dst, const VectorType &src) const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
unsigned int level
Definition: grid_out.cc:4606
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1667
#define AssertIndexRange(index, range)
Definition: exceptions.h:1732
const std::shared_ptr< const Utilities::MPI::Partitioner > & get_partitioner() const
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)
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)
static const unsigned int invalid_unsigned_int
Definition: types.h:201
internal::MatrixFreeFunctions::ShapeInfo< VectorizedArray< Number > > shape_info_coarse
const ::Triangulation< dim, spacedim > & tria