]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add python wrapper functions to Point.
authorBruno Turcksin <bruno.turcksin@gmail.com>
Tue, 11 Apr 2017 13:41:03 +0000 (09:41 -0400)
committerBruno Turcksin <bruno.turcksin@gmail.com>
Tue, 11 Apr 2017 14:24:22 +0000 (10:24 -0400)
contrib/python-bindings/include/point_wrapper.h
contrib/python-bindings/source/export_point.cc
contrib/python-bindings/source/point_wrapper.cc

index 7c266fa0073fdcaede97a337d9b57fd076984d3c..8ea77d04fb5a15248f982db6e745cb5667718b39 100644 (file)
@@ -1,6 +1,6 @@
 // ---------------------------------------------------------------------
 //
-// Copyright (C) 2016 by the deal.II authors
+// Copyright (C) 2016 - 2017 by the deal.II authors
 //
 // This file is part of the deal.II library.
 //
@@ -20,6 +20,8 @@
 
 #ifdef DEAL_II_WITH_CXX11
 
+#include <deal.II/base/point.h>
+
 #include <boost/python.hpp>
 
 DEAL_II_NAMESPACE_OPEN
@@ -40,6 +42,12 @@ namespace python
      */
     PointWrapper(boost::python::list list);
 
+    /**
+     * Constructor. Intialize PointWrapper using a Point.
+     */
+    template <int dim>
+    PointWrapper(const Point<dim> &p);
+
     /**
      * Copy constructor.
      */
@@ -50,12 +58,89 @@ namespace python
      */
     ~PointWrapper();
 
+    /**
+     * Return the Euclidean distance of this point to the point p, i.e., the
+     * l2_norm of the difference between the vectors representing the two
+     * points.
+     */
+    double distance(const PointWrapper &p) const;
+
+    /**
+     * Return the l2_norm of the vector connecting the origin to the point.
+     */
+    double norm() const;
+
+    /**
+     * Return the sum of the absolute squares of all entries.
+     */
+    double norm_square() const;
+
     /**
      * Assignment operator. The dimension of the point is changed if it is
      * different than the one of @p other.
      */
     PointWrapper &operator= (const PointWrapper &other);
 
+    /**
+     * Test for inequality of two points.
+     */
+    bool operator!= (const PointWrapper &p) const; 
+
+    /**
+     * Test for equality of two tensors.
+     */
+    bool operator== (const PointWrapper &p) const;
+
+    /**
+     * Return the scalar product of the vectors representing two points.
+     */
+    double operator* (const PointWrapper &p) const;
+
+    /**
+     * Add an offset to a point.
+     */
+    PointWrapper operator+(const PointWrapper &) const;
+
+    /**
+     * Subtract two points.
+     */
+    PointWrapper operator-(const PointWrapper &) const;
+
+    /**
+     * The opposite point.
+     */
+    PointWrapper operator-() const;
+
+    /**
+     * Divide the coordinates of the point by a factor.
+     */
+    PointWrapper operator/(const double factor) const;
+
+    /**
+     * Multiply the coordinates of the point by a factor.
+     */
+    PointWrapper operator* (const double factor) const;
+
+    /**
+     * Add another point.
+     */
+    PointWrapper& operator+= (const PointWrapper &p);
+
+    /**
+     * Subtract another point.
+     */
+    PointWrapper& operator-= (const PointWrapper &p);
+
+    /**
+     * Scale the coordinates of the point by factor.
+     */
+    PointWrapper& operator*= (const double factor);
+
+    /**
+     * Scale the coordinates of the point by 1/factor.
+     */
+    PointWrapper& operator/= (const double factor);
+
     /**
      * Return the first component of the Point.
      */
@@ -96,6 +181,11 @@ namespace python
      */
     void *get_point();
 
+    /**
+     * Return a constant pointer that can be casted to the underlying Point.
+     */
+    const void *get_point() const;
+
   private:
     /**
      * Delete the underlying Point and free the memory.
@@ -137,6 +227,14 @@ namespace python
   {
     return point;
   }
+
+
+
+  inline
+  const void *PointWrapper::get_point() const
+  {
+    return point;
+  }
 }
 
 DEAL_II_NAMESPACE_CLOSE
index 61b840fe82f6a709b29c1f88427b51242728de00..d79f6c36b156747ad167cc3d7f2186fb419add8c 100644 (file)
@@ -25,16 +25,70 @@ DEAL_II_NAMESPACE_OPEN
 
 namespace python
 {
+  const char distance_docstring [] =
+    "Return the Euclidean distance of this point to the point p             \n"
+    ;
+
+
+
+  const char norm_docstring [] =
+    "Return the L2 norm of the vector connecting the origin to the point    \n"
+    ;
+
+
+
+  const char norm_square_docstring [] =
+    "Return the sum of the absolute squares of all entries                  \n"
+    ;
+
+
+
+  const char get_x_docstring [] =
+    "Get the x component of the point                                       \n"
+    ;
+
+
+
+  const char get_y_docstring [] =
+    "Get the y component of the point                                       \n"
+    ;
+
+
+
+  const char get_z_docstring [] =
+    "Get the z component of the point                                       \n"
+    ;
+
+
+
   void export_point()
   {
     boost::python::class_<PointWrapper> ("Point",
                                          boost::python::init<boost::python::list>())
+    .def("distance", &PointWrapper::distance, distance_docstring,
+        boost::python::args("self", "p"))
+    .def("norm", &PointWrapper::norm, norm_docstring, 
+        boost::python::args("self"))
+    .def("norm_square", &PointWrapper::norm_square, norm_square_docstring,
+        boost::python::args("self"))
+    .def(boost::python::self != boost::python::self)
+    .def(boost::python::self == boost::python::self)
+    .def(boost::python::self * boost::python::self)
+    .def(boost::python::self + boost::python::self)
+    .def(boost::python::self - boost::python::self)
+    .def(- boost::python::self)
+    .def(boost::python::self / float())
+    .def(boost::python::self * float())
+    .def(boost::python::self += boost::python::self)
+    .def(boost::python::self -= boost::python::self)
+    .def(boost::python::self *= float())
+    .def(boost::python::self /= float())
     .add_property("x", &PointWrapper::get_x, &PointWrapper::set_x,
-                  "Get the x component of the point.")
+                  get_x_docstring)
     .add_property("y", &PointWrapper::get_y, &PointWrapper::set_y,
-                  "Get the y component of the point.")
+                  get_y_docstring)
     .add_property("z", &PointWrapper::get_z, &PointWrapper::set_z,
-                  "Get the z component of the point.");
+                  get_z_docstring);
   }
 }
 
index 3fb1863090a84e3166b77cfbda4ea10a2bfd4abf..d58220b048344897e40a24ba03908d0b2bf6664a 100644 (file)
 
 #ifdef DEAL_II_WITH_CXX11
 
-#include <deal.II/base/point.h>
 
 DEAL_II_NAMESPACE_OPEN
 
 namespace python
 {
+  namespace internal
+  {
+    template <int dim>
+    double distance(const Point<dim> & p1,
+                    const Point<dim> & p2) 
+    {
+      return p1.distance(p2);
+    }
+
+
+
+    template <int dim>
+    bool not_equal(const Point<dim> &p1,
+                   const Point<dim> &p2)
+    {
+      return (p1 != p2);
+    }
+
+
+
+    template <int dim>
+    bool equal(const Point<dim> &p1,
+               const Point<dim> &p2)
+    {
+      return (p1 == p2);
+    }
+
+
+
+    template <int dim>
+    double dot_product (const Point<dim> &p1,
+                        const Point<dim> &p2)
+    {
+      return p1*p2;
+    }
+
+    
+
+    template <int dim>
+    Point<dim> add_points(const Point<dim> &p1,
+        const Point<dim> &p2)
+    {
+      return p1+p2;
+    }
+
+    
+
+    template <int dim>
+    Point<dim> subtract_point(const Point<dim> &p1,
+        const Point<dim> &p2)
+    {
+      if (dim == 2)
+        return Point<dim>(p1[0]-p2[0], p1[1]-p2[1]);
+      else
+        return Point<dim>(p1[0]-p2[0], p1[1]-p2[1], p1[2]-p2[2]);
+    }
+
+
+
+    template <int dim>
+    Point<dim> opposite_point(const Point<dim> &p)
+    {
+      return -p;
+    }
+
+
+
+    template <int dim>
+    Point<dim> divide_point(const Point<dim> &p, double factor)
+    {
+      return p/factor;
+    }
+
+
+
+    template <int dim>
+    Point<dim> multiply_point(const Point<dim> &p, double factor)
+    {
+      return p*factor;
+    }
+
+
+
+    template <int dim>
+    void add_and_set(Point<dim> &p1, const Point<dim> &p2)
+    {
+      p1 += p2;
+    }
+
+
+
+    template <int dim>
+    void subtract_and_set(Point<dim> &p1, const Point<dim> &p2)
+    {
+      p1 -= p2;
+    }
+
+
+
+    template <int dim>
+    void multiply_and_set(Point<dim> &p, double factor)
+    {
+      p *= factor;
+    }
+
+
+
+    template <int dim>
+    void divide_and_set(Point<dim> &p, double factor)
+    {
+      p /= factor;
+    }
+  }
+
   PointWrapper::PointWrapper()
     :
     dim(-1),
@@ -47,6 +160,19 @@ namespace python
 
 
 
+  template <int d>
+  PointWrapper::PointWrapper(const Point<d> &p)
+  {
+    dim = d;
+    if (dim == 2)
+      point = new Point<2>(p[0],p[1]);
+                          
+    else
+      point = new Point<3>(p[0],p[1],p[2]);
+  }
+
+
+
   PointWrapper::PointWrapper(const PointWrapper &other)
   {
     copy(other);
@@ -62,6 +188,41 @@ namespace python
 
 
 
+  double PointWrapper::distance(const PointWrapper &p) const
+  {
+    AssertThrow(p.get_dim() == dim,
+        ExcMessage("The points do not have the same dimension."));
+    
+    if (dim == 2)
+      return internal::distance(*static_cast<Point<2>*>(point), 
+          *static_cast<const Point<2>*>(p.get_point()));
+    else
+      return internal::distance(*static_cast<Point<3>*>(point), 
+          *static_cast<const Point<3>*>(p.get_point()));
+  }
+
+
+
+  double PointWrapper::norm() const
+  {
+    if (dim == 2)
+      return static_cast<Point<2>*>(point)->norm();
+    else
+      return static_cast<Point<3>*>(point)->norm();
+  }
+
+
+
+  double PointWrapper::norm_square() const
+  {
+    if (dim == 2)
+      return static_cast<Point<2>*>(point)->norm_square();
+    else
+      return static_cast<Point<3>*>(point)->norm_square();
+  }
+
+
+
   PointWrapper &PointWrapper::operator= (const PointWrapper &other)
   {
     clear();
@@ -72,6 +233,183 @@ namespace python
 
 
 
+  bool PointWrapper::operator!= (const PointWrapper &p) const
+  {
+    AssertThrow(p.get_dim() == dim,
+        ExcMessage("The points do not have the same dimension."));
+
+    if (dim == 2)
+      return internal::not_equal(*static_cast<const Point<2>*>(point),
+          *static_cast<const Point<2>*>(p.get_point()));
+    else
+      return internal::not_equal(*static_cast<const Point<3>*>(point),
+          *static_cast<const Point<3>*>(p.get_point()));
+  }
+
+
+
+  bool PointWrapper::operator== (const PointWrapper &p) const
+  {
+    AssertThrow(p.get_dim() == dim,
+        ExcMessage("The points do not have the same dimension."));
+
+    if (dim == 2)
+      return internal::equal(*static_cast<const Point<2>*>(point),
+          *static_cast<const Point<2>*>(p.get_point()));
+    else
+      return internal::equal(*static_cast<const Point<3>*>(point),
+          *static_cast<const Point<3>*>(p.get_point()));
+  }
+
+
+
+  double PointWrapper::operator* (const PointWrapper &p) const
+  {
+    AssertThrow(p.get_dim() == dim,
+        ExcMessage("The points do not have the same dimension."));
+
+    if (dim == 2)
+      return internal::dot_product(*static_cast<const Point<2>*>(point),
+          *static_cast<const Point<2>*>(p.get_point()));
+    else
+      return internal::dot_product(*static_cast<const Point<3>*>(point),
+          *static_cast<const Point<3>*>(p.get_point()));
+  }
+
+
+
+  PointWrapper PointWrapper::operator+ (const PointWrapper &p) const
+  {
+    AssertThrow(p.get_dim() == dim,
+        ExcMessage("The points do not have the same dimension."));
+
+    if (dim == 2)
+      return PointWrapper(internal::add_points(
+            *static_cast<const Point<2>*>(point),
+            *static_cast<const Point<2>*>(p.get_point())));
+    else
+      return PointWrapper(internal::add_points(
+            *static_cast<const Point<3>*>(point),
+            *static_cast<const Point<3>*>(p.get_point())));
+  }
+
+
+
+  PointWrapper PointWrapper::operator- (const PointWrapper &p) const
+  {
+    AssertThrow(p.get_dim() == dim,
+        ExcMessage("The points do not have the same dimension."));
+
+    if (dim == 2)
+      return PointWrapper(internal::subtract_point(
+            *static_cast<const Point<2>*>(point),
+            *static_cast<const Point<2>*>(p.get_point())));
+    else
+      return PointWrapper(internal::subtract_point(
+            *static_cast<const Point<3>*>(point),
+            *static_cast<const Point<3>*>(p.get_point())));
+  }
+
+
+
+  PointWrapper PointWrapper::operator-() const
+  {
+    if (dim == 2)
+      return PointWrapper(internal::opposite_point(
+            *static_cast<Point<2>*>(point)));
+    else
+      return PointWrapper(internal::opposite_point(
+            *static_cast<Point<3>*>(point)));
+  }
+
+
+
+  PointWrapper PointWrapper::operator/(const double factor) const
+  {
+    AssertThrow(factor != 0., ExcMessage("Dividing by zero."));
+
+    if (dim == 2)
+      return PointWrapper(internal::divide_point(
+            *static_cast<Point<2>*>(point), factor));
+    else
+      return PointWrapper(internal::divide_point(
+            *static_cast<Point<3>*>(point), factor));
+  }
+
+
+  PointWrapper PointWrapper::operator* (const double factor) const
+  {
+
+    if (dim == 2)
+      return PointWrapper(internal::multiply_point(
+            *static_cast<Point<2>*>(point), factor));
+    else
+      return PointWrapper(internal::multiply_point(
+            *static_cast<Point<3>*>(point), factor));
+  }
+
+
+
+  PointWrapper& PointWrapper::operator+= (const PointWrapper &p)
+  {
+    AssertThrow(p.get_dim() == dim,
+        ExcMessage("The points do not have the same dimension."));
+
+    if (dim == 2)
+      internal::add_and_set(*static_cast<Point<2>*>(point),
+          *static_cast<const Point<2>*>(p.get_point()));
+    else
+      internal::add_and_set(*static_cast<Point<3>*>(point),
+          *static_cast<const Point<3>*>(p.get_point()));
+
+    return *this;
+  }
+
+
+
+  PointWrapper& PointWrapper::operator-= (const PointWrapper &p)
+  {
+    AssertThrow(p.get_dim() == dim,
+        ExcMessage("The points do not have the same dimension."));
+
+    if (dim == 2)
+      internal::subtract_and_set(*static_cast<Point<2>*>(point),
+          *static_cast<const Point<2>*>(p.get_point()));
+    else
+      internal::subtract_and_set(*static_cast<Point<3>*>(point),
+          *static_cast<const Point<3>*>(p.get_point()));
+
+    return *this;
+  }
+
+
+
+  PointWrapper& PointWrapper::operator*= (const double factor)
+  {
+    if (dim == 2)
+      internal::multiply_and_set(*static_cast<Point<2>*>(point), factor);
+    else
+      internal::multiply_and_set(*static_cast<Point<3>*>(point), factor);
+
+    return *this;
+  }
+
+
+
+  PointWrapper& PointWrapper::operator/= (const double factor)
+  {
+    AssertThrow(factor != 0., ExcMessage("Dividing by zero."));
+
+    if (dim == 2)
+      internal::divide_and_set(*static_cast<Point<2>*>(point), factor);
+    else
+      internal::divide_and_set(*static_cast<Point<3>*>(point), factor);
+
+    return *this;
+  }
+
+
+
   double PointWrapper::get_x() const
   {
     if (dim == 2)

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.