]> https://gitweb.dealii.org/ - dealii.git/commitdiff
MGTransfer: don't repeat yourself 9593/head
authorNiklas Fehn <fehn@lnm.mw.tum.de>
Mon, 2 Mar 2020 17:06:58 +0000 (18:06 +0100)
committerNiklas Fehn <fehn@lnm.mw.tum.de>
Tue, 3 Mar 2020 07:28:02 +0000 (08:28 +0100)
include/deal.II/multigrid/mg_transfer.h
include/deal.II/multigrid/mg_transfer.templates.h

index fc2e630f228118bc12a1a1d7b790e23623da2636..89cef8c89c11928b4327d73a63bbe55ca7ea7359 100644 (file)
@@ -417,6 +417,14 @@ protected:
    */
   SmartPointer<const MGConstrainedDoFs, MGLevelGlobalTransfer<VectorType>>
     mg_constrained_dofs;
+
+private:
+  /**
+   * This function is called to make sure that build() has been invoked.
+   */
+  template <int dim, int spacedim>
+  void
+  assert_built(const DoFHandler<dim, spacedim> &dof_handler) const;
 };
 
 
@@ -634,6 +642,14 @@ protected:
    */
   mutable MGLevelObject<LinearAlgebra::distributed::Vector<Number>>
     solution_ghosted_level_vector;
+
+private:
+  /**
+   * This function is called to make sure that build() has been invoked.
+   */
+  template <int dim, int spacedim>
+  void
+  assert_built(const DoFHandler<dim, spacedim> &dof_handler) const;
 };
 
 
index 0cf2bbf44d266f558d22357490085b50cd8d824a..0158b6310b03b3a98fd9a448b3b4fa72dd38e63a 100644 (file)
@@ -210,9 +210,7 @@ MGLevelGlobalTransfer<VectorType>::copy_to_mg(
   MGLevelObject<VectorType> &      dst,
   const InVector &                 src) const
 {
-  Assert(copy_indices.size() ==
-           dof_handler.get_triangulation().n_global_levels(),
-         ExcMessage("MGTransfer::build() has not been called!"));
+  assert_built(dof_handler);
   AssertIndexRange(dst.max_level(),
                    dof_handler.get_triangulation().n_global_levels());
   AssertIndexRange(dst.min_level(), dst.max_level() + 1);
@@ -282,10 +280,7 @@ MGLevelGlobalTransfer<VectorType>::copy_from_mg(
   OutVector &                      dst,
   const MGLevelObject<VectorType> &src) const
 {
-  (void)dof_handler;
-  Assert(copy_indices.size() ==
-           dof_handler.get_triangulation().n_global_levels(),
-         ExcMessage("MGTransfer::build() has not been called!"));
+  assert_built(dof_handler);
   AssertIndexRange(src.max_level(),
                    dof_handler.get_triangulation().n_global_levels());
   AssertIndexRange(src.min_level(), src.max_level() + 1);
@@ -359,10 +354,7 @@ MGLevelGlobalTransfer<VectorType>::copy_from_mg_add(
   OutVector &                      dst,
   const MGLevelObject<VectorType> &src) const
 {
-  (void)dof_handler;
-  Assert(copy_indices.size() ==
-           dof_handler.get_triangulation().n_global_levels(),
-         ExcMessage("MGTransfer::build() has not been called!"));
+  assert_built(dof_handler);
   // For non-DG: degrees of freedom in the refinement face may need special
   // attention, since they belong to the coarse level, but have fine level
   // basis functions
@@ -399,6 +391,17 @@ MGLevelGlobalTransfer<VectorType>::set_component_to_block_map(
   component_to_block_map = map;
 }
 
+template <typename VectorType>
+template <int dim, int spacedim>
+void
+MGLevelGlobalTransfer<VectorType>::assert_built(
+  const DoFHandler<dim, spacedim> &dof_handler) const
+{
+  (void)dof_handler;
+  Assert(copy_indices.size() ==
+           dof_handler.get_triangulation().n_global_levels(),
+         ExcMessage("MGTransfer::build() has not been called!"));
+}
 
 
 /* --------- MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector> -------
@@ -412,9 +415,7 @@ MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number>>::copy_to_mg(
   MGLevelObject<LinearAlgebra::distributed::Vector<Number>> &dst,
   const LinearAlgebra::distributed::Vector<Number2> &        src) const
 {
-  Assert(copy_indices.size() ==
-           dof_handler.get_triangulation().n_global_levels(),
-         ExcMessage("MGTransfer::build() has not been called!"));
+  assert_built(dof_handler);
   copy_to_mg(dof_handler, dst, src, false);
 }
 
@@ -428,9 +429,7 @@ MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number>>::copy_to_mg(
   const LinearAlgebra::distributed::Vector<Number2> &        src,
   const bool solution_transfer) const
 {
-  Assert(copy_indices.size() ==
-           dof_handler.get_triangulation().n_global_levels(),
-         ExcMessage("MGTransfer::build() has not been called!"));
+  assert_built(dof_handler);
   LinearAlgebra::distributed::Vector<Number> &this_ghosted_global_vector =
     solution_transfer ? solution_ghosted_global_vector : ghosted_global_vector;
   const std::vector<Table<2, unsigned int>> &this_copy_indices =
@@ -439,8 +438,6 @@ MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number>>::copy_to_mg(
     solution_transfer ? solution_copy_indices_level_mine :
                         copy_indices_level_mine;
 
-  (void)dof_handler;
-
   AssertIndexRange(dst.max_level(),
                    dof_handler.get_triangulation().n_global_levels());
   AssertIndexRange(dst.min_level(), dst.max_level() + 1);
@@ -536,10 +533,7 @@ MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number>>::copy_from_mg(
   LinearAlgebra::distributed::Vector<Number2> &                    dst,
   const MGLevelObject<LinearAlgebra::distributed::Vector<Number>> &src) const
 {
-  (void)dof_handler;
-  Assert(copy_indices.size() ==
-           dof_handler.get_triangulation().n_global_levels(),
-         ExcMessage("MGTransfer::build() has not been called!"));
+  assert_built(dof_handler);
   AssertIndexRange(src.max_level(),
                    dof_handler.get_triangulation().n_global_levels());
   AssertIndexRange(src.min_level(), src.max_level() + 1);
@@ -614,10 +608,7 @@ MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number>>::
     LinearAlgebra::distributed::Vector<Number2> &dst,
     const MGLevelObject<LinearAlgebra::distributed::Vector<Number>> &src) const
 {
-  (void)dof_handler;
-  Assert(copy_indices.size() ==
-           dof_handler.get_triangulation().n_global_levels(),
-         ExcMessage("MGTransfer::build() has not been called!"));
+  assert_built(dof_handler);
   // For non-DG: degrees of freedom in the refinement face may need special
   // attention, since they belong to the coarse level, but have fine level
   // basis functions
@@ -663,6 +654,18 @@ MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number>>::
   component_to_block_map = map;
 }
 
+template <typename Number>
+template <int dim, int spacedim>
+void
+MGLevelGlobalTransfer<LinearAlgebra::distributed::Vector<Number>>::assert_built(
+  const DoFHandler<dim, spacedim> &dof_handler) const
+{
+  (void)dof_handler;
+  Assert(copy_indices.size() ==
+           dof_handler.get_triangulation().n_global_levels(),
+         ExcMessage("MGTransfer::build() has not been called!"));
+}
+
 
 DEAL_II_NAMESPACE_CLOSE
 

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.