]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Manifold: add accessors for all Manifold ctor arguments.
authorDavid Wells <drwells@email.unc.edu>
Sat, 6 Jul 2024 13:32:31 +0000 (09:32 -0400)
committerDavid Wells <drwells@email.unc.edu>
Sat, 6 Jul 2024 14:20:31 +0000 (10:20 -0400)
This makes it easier to manipulate Manifolds after they are constructed, e.g.,
by applying a translation.

doc/news/changes/minor/20240706DavidWells [new file with mode: 0644]
include/deal.II/grid/manifold.h
include/deal.II/grid/manifold_lib.h
source/grid/manifold_lib.cc

diff --git a/doc/news/changes/minor/20240706DavidWells b/doc/news/changes/minor/20240706DavidWells
new file mode 100644 (file)
index 0000000..0c04cdf
--- /dev/null
@@ -0,0 +1,5 @@
+New: All the basic Manifold objects describing a coordinate system now make
+their essential geometric properties (e.g., the center) available through member
+functions.
+<br>
+(David Wells, 2024/07/06)
index 931a1c635158978748b3638b47371bb9974a259e..623288d1585ab8a064fd2a58b8d7ff929fb96fd6 100644 (file)
@@ -779,6 +779,12 @@ public:
   const Tensor<1, spacedim> &
   get_periodicity() const;
 
+  /**
+   * Return the relative tolerance set in the constructor.
+   */
+  double
+  get_tolerance() const;
+
 private:
   /**
    * The periodicity of this Manifold. Periodicity affects the way a middle
@@ -1133,6 +1139,15 @@ Manifold<3, 3>::get_new_point_on_hex(
 
 /*---Templated functions---*/
 
+template <int dim, int spacedim>
+inline double
+FlatManifold<dim, spacedim>::get_tolerance() const
+{
+  return tolerance;
+}
+
+
+
 namespace Manifolds
 {
   template <typename MeshIteratorType>
index 1cc44a64f17012a328880333ec6483b4d60680b4..156f65a84020d687abdf752918b583781905202b 100644 (file)
@@ -134,6 +134,12 @@ public:
     const typename Triangulation<dim, spacedim>::face_iterator &face,
     const Point<spacedim> &p) const override;
 
+  /**
+   * Return the center of the spherical coordinate system.
+   */
+  const Point<spacedim> &
+  get_center() const;
+
   /**
    * The center of the spherical coordinate system.
    */
@@ -323,6 +329,12 @@ public:
   get_new_point(const ArrayView<const Point<spacedim>> &vertices,
                 const ArrayView<const double>          &weights) const override;
 
+  /**
+   * Return the center of the spherical coordinate system.
+   */
+  const Point<spacedim> &
+  get_center() const;
+
   /**
    * The center of the spherical coordinate system.
    */
@@ -446,6 +458,27 @@ public:
   get_new_point(const ArrayView<const Point<spacedim>> &surrounding_points,
                 const ArrayView<const double>          &weights) const override;
 
+  /**
+   * Get the Tensor parallel to the cylinder's axis.
+   */
+  const Tensor<1, spacedim> &
+  get_direction() const;
+
+  /**
+   * Get a point on the Cylinder's axis.
+   *
+   * @note This argument, like get_direction() and get_tolerance(), just returns
+   * the arguments set in the three-argument constructor.
+   */
+  const Point<spacedim> &
+  get_point_on_axis() const;
+
+  /**
+   * Get the tolerance which determines if a point is on the Cylinder's axis.
+   */
+  double
+  get_tolerance() const;
+
 private:
   /**
    * A vector orthogonal to the normal direction.
@@ -544,16 +577,40 @@ public:
   virtual DerivativeForm<1, spacedim, spacedim>
   push_forward_gradient(const Point<spacedim> &chart_point) const override;
 
+  /**
+   * Get the Tensor parallel to the cylinder's major axis.
+   */
+  const Tensor<1, spacedim> &
+  get_major_axis_direction() const;
+
+  /**
+   * Return the center of the elliptical coordinate system.
+   */
+  const Point<spacedim> &
+  get_center() const;
+
+  /**
+   * Return the ellipse's eccentricity.
+   */
+  double
+  get_eccentricity() const;
 
 private:
   /**
    * The direction vector of the major axis.
    */
-  Tensor<1, spacedim> direction;
+  const Tensor<1, spacedim> direction;
+
   /**
    * The center of the manifold.
    */
   const Point<spacedim> center;
+
+  /**
+   * The eccentricity.
+   */
+  const double eccentricity;
+
   /**
    * Parameters deriving from the eccentricity of the manifold.
    */
@@ -816,6 +873,18 @@ public:
   virtual DerivativeForm<1, 3, 3>
   push_forward_gradient(const Point<3> &chart_point) const override;
 
+  /**
+   * Get the radius of the centerline.
+   */
+  double
+  get_centerline_radius() const;
+
+  /**
+   * Get the inner radius of the torus.
+   */
+  double
+  get_inner_radius() const;
+
 private:
   double centerline_radius;
 
@@ -1159,6 +1228,96 @@ private:
   boost::signals2::connection clear_signal;
 };
 
+/*----------------------------- inline functions -----------------------------*/
+
+template <int dim, int spacedim>
+inline const Point<spacedim> &
+PolarManifold<dim, spacedim>::get_center() const
+{
+  return center;
+}
+
+
+
+template <int dim, int spacedim>
+inline const Point<spacedim> &
+SphericalManifold<dim, spacedim>::get_center() const
+{
+  return center;
+}
+
+
+
+template <int dim, int spacedim>
+inline const Tensor<1, spacedim> &
+CylindricalManifold<dim, spacedim>::get_direction() const
+{
+  return direction;
+}
+
+
+
+template <int dim, int spacedim>
+inline const Point<spacedim> &
+CylindricalManifold<dim, spacedim>::get_point_on_axis() const
+{
+  return point_on_axis;
+}
+
+
+
+template <int dim, int spacedim>
+inline double
+CylindricalManifold<dim, spacedim>::get_tolerance() const
+{
+  return tolerance;
+}
+
+
+
+template <int dim, int spacedim>
+inline const Tensor<1, spacedim> &
+EllipticalManifold<dim, spacedim>::get_major_axis_direction() const
+{
+  return direction;
+}
+
+
+
+template <int dim, int spacedim>
+inline const Point<spacedim> &
+EllipticalManifold<dim, spacedim>::get_center() const
+{
+  return center;
+}
+
+
+
+template <int dim, int spacedim>
+inline double
+EllipticalManifold<dim, spacedim>::get_eccentricity() const
+{
+  return eccentricity;
+}
+
+
+
+template <int dim>
+inline double
+TorusManifold<dim>::get_centerline_radius() const
+{
+  return centerline_radius;
+}
+
+
+
+template <int dim>
+inline double
+TorusManifold<dim>::get_inner_radius() const
+{
+  return inner_radius;
+}
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif
index a7b4541939185fe361695d71a7d625cf2299897d..07d5e1dbff67079b27360eacca456bd589fd883c 100644 (file)
@@ -1165,6 +1165,20 @@ CylindricalManifold<dim, spacedim>::push_forward_gradient(
 
 
 
+namespace
+{
+  template <int dim>
+  Tensor<1, dim>
+  check_and_normalize(const Tensor<1, dim> &t)
+  {
+    const double norm = t.norm();
+    Assert(norm > 0.0, ExcMessage("The major axis must have a positive norm."));
+    return t / norm;
+  }
+} // namespace
+
+
+
 // ============================================================
 // EllipticalManifold
 // ============================================================
@@ -1175,8 +1189,9 @@ EllipticalManifold<dim, spacedim>::EllipticalManifold(
   const double               eccentricity)
   : ChartManifold<dim, spacedim, spacedim>(
       EllipticalManifold<dim, spacedim>::get_periodicity())
-  , direction(major_axis_direction)
+  , direction(check_and_normalize(major_axis_direction))
   , center(center)
+  , eccentricity(eccentricity)
   , cosh_u(1.0 / eccentricity)
   , sinh_u(std::sqrt(cosh_u * cosh_u - 1.0))
 {
@@ -1186,11 +1201,6 @@ EllipticalManifold<dim, spacedim>::EllipticalManifold(
   Assert(std::signbit(cosh_u * cosh_u - 1.0) == false,
          ExcMessage(
            "Invalid eccentricity: It must satisfy 0 < eccentricity < 1."));
-  const double direction_norm = direction.norm();
-  Assert(direction_norm != 0.0,
-         ExcMessage(
-           "Invalid major axis direction vector: Null vector not allowed."));
-  direction /= direction_norm;
 }
 
 

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.