]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Reorganize (some of) the internal grid functions. 4618/head
authorDavid Wells <wellsd2@rpi.edu>
Sun, 16 Jul 2017 19:50:21 +0000 (15:50 -0400)
committerDavid Wells <wellsd2@rpi.edu>
Sun, 16 Jul 2017 19:50:21 +0000 (15:50 -0400)
This is similar to 6310b2dee3 and 52bef21ac2: these functions should be in
specific namespaces to avoid clashes if we want to support unity builds.

source/grid/grid_refinement.cc
source/grid/manifold_lib.cc

index da97212ba475bdde370efa20a79de4df31e4c840..1fe51cb306ff440cdc2f9f43a568f3daeaa79950 100644 (file)
 DEAL_II_NAMESPACE_OPEN
 
 
-namespace
+namespace internal
 {
-  namespace internal
+  namespace GridRefinement
   {
-    template <typename number>
-    inline
-    number
-    max_element (const Vector<number> &criteria)
+    namespace
     {
-      return *std::max_element(criteria.begin(), criteria.end());
-    }
+      template <typename number>
+      inline
+      number
+      max_element (const dealii::Vector<number> &criteria)
+      {
+        return *std::max_element(criteria.begin(), criteria.end());
+      }
 
 
-    template <typename number>
-    inline
-    number
-    min_element (const Vector<number> &criteria)
-    {
-      return *std::min_element(criteria.begin(), criteria.end());
-    }
+      template <typename number>
+      inline
+      number
+      min_element (const dealii::Vector<number> &criteria)
+      {
+        return *std::min_element(criteria.begin(), criteria.end());
+      }
 
-    // Silence a (bogus) warning in clang-3.6 about the following four
-    // functions being unused:
-    DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
+      // Silence a (bogus) warning in clang-3.6 about the following four
+      // functions being unused:
+      DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
 
 #ifdef DEAL_II_WITH_TRILINOS
-    inline
-    TrilinosScalar
-    max_element (const TrilinosWrappers::MPI::Vector &criteria)
-    {
-      TrilinosScalar m = 0;
-      criteria.trilinos_vector().MaxValue(&m);
-      return m;
-    }
+      inline
+      TrilinosScalar
+      max_element (const dealii::TrilinosWrappers::MPI::Vector &criteria)
+      {
+        TrilinosScalar m = 0;
+        criteria.trilinos_vector().MaxValue(&m);
+        return m;
+      }
 
 
-    inline
-    TrilinosScalar
-    min_element (const TrilinosWrappers::MPI::Vector &criteria)
-    {
-      TrilinosScalar m = 0;
-      criteria.trilinos_vector().MinValue(&m);
-      return m;
-    }
+      inline
+      TrilinosScalar
+      min_element (const dealii::TrilinosWrappers::MPI::Vector &criteria)
+      {
+        TrilinosScalar m = 0;
+        criteria.trilinos_vector().MinValue(&m);
+        return m;
+      }
 #endif
 
-    DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
+      DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
 
-  } /* namespace internal */
-
-
-  template <typename VectorType>
-  typename constraint_and_return_value<!IsBlockVector<VectorType>::value,
-           typename VectorType::value_type>::type
-           min_element (const VectorType &criteria)
-  {
-    return internal::min_element (criteria);
-  }
-
-
-  template <typename VectorType>
-  typename constraint_and_return_value<!IsBlockVector<VectorType>::value,
-           typename VectorType::value_type>::type
-           max_element (const VectorType &criteria)
-  {
-    return internal::max_element (criteria);
-  }
+      template <typename VectorType>
+      typename constraint_and_return_value<!IsBlockVector<VectorType>::value,
+               typename VectorType::value_type>::type
+               min_element (const VectorType &criteria)
+      {
+        return min_element (criteria);
+      }
 
 
-  template <typename VectorType>
-  typename constraint_and_return_value<IsBlockVector<VectorType>::value,
-           typename VectorType::value_type>::type
-           min_element (const VectorType &criteria)
-  {
-    typename VectorType::value_type t = internal::min_element(criteria.block(0));
-    for (unsigned int b=1; b<criteria.n_blocks(); ++b)
-      t = std::min (t, internal::min_element(criteria.block(b)));
+      template <typename VectorType>
+      typename constraint_and_return_value<!IsBlockVector<VectorType>::value,
+               typename VectorType::value_type>::type
+               max_element (const VectorType &criteria)
+      {
+        return max_element (criteria);
+      }
 
-    return t;
-  }
 
+      template <typename VectorType>
+      typename constraint_and_return_value<IsBlockVector<VectorType>::value,
+               typename VectorType::value_type>::type
+               min_element (const VectorType &criteria)
+      {
+        typename VectorType::value_type t = min_element(criteria.block(0));
+        for (unsigned int b=1; b<criteria.n_blocks(); ++b)
+          t = std::min (t, min_element(criteria.block(b)));
 
-  template <typename VectorType>
-  typename constraint_and_return_value<IsBlockVector<VectorType>::value,
-           typename VectorType::value_type>::type
-           max_element (const VectorType &criteria)
-  {
-    typename VectorType::value_type t = internal::max_element(criteria.block(0));
-    for (unsigned int b=1; b<criteria.n_blocks(); ++b)
-      t = std::max (t, internal::max_element(criteria.block(b)));
+        return t;
+      }
 
-    return t;
-  }
 
-}
+      template <typename VectorType>
+      typename constraint_and_return_value<IsBlockVector<VectorType>::value,
+               typename VectorType::value_type>::type
+               max_element (const VectorType &criteria)
+      {
+        typename VectorType::value_type t = max_element(criteria.block(0));
+        for (unsigned int b=1; b<criteria.n_blocks(); ++b)
+          t = std::max (t, max_element(criteria.block(b)));
 
+        return t;
+      }
+    }
+  } /* namespace GridRefinement */
+} /* namespace internal */
 
 
 template <int dim, class VectorType, int spacedim>
@@ -439,7 +438,7 @@ GridRefinement::refine_and_coarsen_fixed_fraction (Triangulation<dim,spacedim> &
   // threshold if it equals the
   // largest indicator and the
   // top_fraction!=1
-  if ((top_threshold == max_element(criteria)) &&
+  if ((top_threshold == internal::GridRefinement::max_element(criteria)) &&
       (top_fraction != 1))
     top_threshold *= 0.999;
 
@@ -447,10 +446,10 @@ GridRefinement::refine_and_coarsen_fixed_fraction (Triangulation<dim,spacedim> &
     bottom_threshold = 0.999*top_threshold;
 
   // actually flag cells
-  if (top_threshold < max_element(criteria))
+  if (top_threshold < internal::GridRefinement::max_element(criteria))
     refine (tria, criteria, top_threshold, pp - tmp.begin());
 
-  if (bottom_threshold > min_element(criteria))
+  if (bottom_threshold > internal::GridRefinement::min_element(criteria))
     coarsen (tria, criteria, bottom_threshold);
 }
 
index 6a806d0cda5a5ca7c9e21e7b163fc6d375be318b..6b1ee081458f75986a59c95336b00f734f7d766d 100644 (file)
@@ -300,37 +300,43 @@ get_new_point (const std::vector<Point<spacedim> > &vertices,
 // CylindricalManifold
 // ============================================================
 
-namespace
+namespace internal
 {
-  // helper function to compute a vector orthogonal to a given one.
-  template <int spacedim>
-  Point<spacedim>
-  compute_normal(const Tensor<1,spacedim> &vector)
+  namespace CylindricalManifold
   {
-    AssertThrow(vector.norm() != 0.,
-                ExcMessage("The direction parameter must not be zero!"));
-    Point<3> normal;
-    if (std::abs(vector[0]) >= std::abs(vector[1])
-        && std::abs(vector[0]) >= std::abs(vector[2]))
-      {
-        normal[1]=-1.;
-        normal[2]=-1.;
-        normal[0]=(vector[1]+vector[2])/vector[0];
-      }
-    else if (std::abs(vector[1]) >= std::abs(vector[0])
-             && std::abs(vector[1]) >= std::abs(vector[2]))
-      {
-        normal[0]=-1.;
-        normal[2]=-1.;
-        normal[1]=(vector[0]+vector[2])/vector[1];
-      }
-    else
+    namespace
+    {
+      // helper function to compute a vector orthogonal to a given one.
+      template <int spacedim>
+      Point<spacedim>
+      compute_normal(const Tensor<1,spacedim> &vector)
       {
-        normal[0]=-1.;
-        normal[1]=-1.;
-        normal[2]=(vector[0]+vector[1])/vector[2];
+        AssertThrow(vector.norm() != 0.,
+                    ExcMessage("The direction parameter must not be zero!"));
+        Point<3> normal;
+        if (std::abs(vector[0]) >= std::abs(vector[1])
+            && std::abs(vector[0]) >= std::abs(vector[2]))
+          {
+            normal[1]=-1.;
+            normal[2]=-1.;
+            normal[0]=(vector[1]+vector[2])/vector[0];
+          }
+        else if (std::abs(vector[1]) >= std::abs(vector[0])
+                 && std::abs(vector[1]) >= std::abs(vector[2]))
+          {
+            normal[0]=-1.;
+            normal[2]=-1.;
+            normal[1]=(vector[0]+vector[2])/vector[1];
+          }
+        else
+          {
+            normal[0]=-1.;
+            normal[1]=-1.;
+            normal[2]=(vector[0]+vector[1])/vector[2];
+          }
+        return normal;
       }
-    return normal;
+    }
   }
 }
 
@@ -351,7 +357,7 @@ CylindricalManifold<dim, spacedim>::CylindricalManifold(const Point<spacedim> &d
                                                         const Point<spacedim> &point_on_axis_,
                                                         const double tolerance) :
   ChartManifold<dim,3,3>(Tensor<1,3>({0,2.*numbers::PI,0})),
-              normal_direction(compute_normal(direction_)),
+              normal_direction(internal::CylindricalManifold::compute_normal(direction_)),
               direction (direction_/direction_.norm()),
               point_on_axis (point_on_axis_),
               tolerance(tolerance)

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.