--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 - 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#ifndef dealii_quadrature_wrapper_h
+#define dealii_quadrature_wrapper_h
+
+#include <deal.II/base/config.h>
+
+#include <boost/python.hpp>
+
+#include <point_wrapper.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace python
+{
+ class QuadratureWrapper
+ {
+ public:
+ /**
+ * Copy constructor.
+ */
+ QuadratureWrapper(const QuadratureWrapper &other);
+
+ /**
+ * Constructor.
+ */
+ QuadratureWrapper(const int dim);
+
+ /**
+ * Destructor.
+ */
+ ~QuadratureWrapper();
+
+ /**
+ * Create QGauss.
+ */
+ void
+ create_gauss(const unsigned int n);
+
+ /**
+ * Create QGaussLobatto.
+ */
+ void
+ create_gauss_lobatto(const unsigned int n);
+
+ /**
+ * Return list of quadrature points.
+ */
+ boost::python::list
+ get_points() const;
+
+ /**
+ * Return list of quadrature weights.
+ */
+ boost::python::list
+ get_weights() const;
+
+ /**
+ * Return pointer to an underlying Quadrature object
+ */
+ void *
+ get_quadrature() const;
+
+ /**
+ * Return the dimension of the underlying object
+ */
+ int
+ get_dim() const;
+
+ private:
+ /**
+ * Dimension of the underlying Quadrature object.
+ */
+ int dim;
+
+ /**
+ * Pointer that can be casted to the underlying Quadrature object.
+ */
+ void *quadrature_ptr;
+ };
+
+} // namespace python
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
export_triangulation.cc
export_mapping.cc
export_manifold.cc
+ export_quadrature.cc
cell_accessor_wrapper.cc
tria_accessor_wrapper.cc
point_wrapper.cc
triangulation_wrapper.cc
mapping_wrapper.cc
manifold_wrapper.cc
+ quadrature_wrapper.cc
)
FOREACH(_build ${DEAL_II_BUILD_TYPES})
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 - 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#include <boost/python.hpp>
+
+#include <quadrature_wrapper.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace python
+{
+ const char create_gauss_docstring[] =
+ "Create Gauss quadrature with n points \n"
+ "in each space direction. \n";
+
+
+ const char create_gauss_lobatto_docstring[] =
+ "Create Gauss-Lobatto quadrature with n points \n"
+ "in each space direction. \n";
+
+
+ const char get_points_docstring[] =
+ "Return the list of quadrature points. \n";
+
+
+ const char get_weights_docstring[] =
+ "Return the list of quadrature weights. \n";
+
+
+ void
+ export_quadrature()
+ {
+ boost::python::class_<QuadratureWrapper>(
+ "Quadrature", boost::python::init<const int>(boost::python::args("dim")))
+ .def("create_gauss",
+ &QuadratureWrapper::create_gauss,
+ create_gauss_docstring,
+ boost::python::args("self", "n"))
+ .def("create_gauss_lobatto",
+ &QuadratureWrapper::create_gauss_lobatto,
+ create_gauss_lobatto_docstring,
+ boost::python::args("self", "n"))
+ .def("points",
+ &QuadratureWrapper::get_points,
+ get_points_docstring,
+ boost::python::args("self"))
+ .def("weights",
+ &QuadratureWrapper::get_weights,
+ get_weights_docstring,
+ boost::python::args("self"));
+ }
+} // namespace python
+
+DEAL_II_NAMESPACE_CLOSE
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2016 - 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+#include <deal.II/base/quadrature_lib.h>
+
+#include <quadrature_wrapper.h>
+
+DEAL_II_NAMESPACE_OPEN
+
+namespace python
+{
+ namespace internal
+ {
+ template <int dim>
+ boost::python::list
+ get_points(const void *quadrature_ptr)
+ {
+ const auto quadrature =
+ static_cast<const Quadrature<dim> *>(quadrature_ptr);
+
+ const auto points = quadrature->get_points();
+
+ boost::python::list points_list;
+ for (const auto &p : points)
+ {
+ boost::python::list python_p;
+ for (int d = 0; d < dim; ++d)
+ python_p.append(p[d]);
+
+ points_list.append(python_p);
+ }
+
+ return points_list;
+ }
+
+ template <int dim>
+ boost::python::list
+ get_weights(const void *quadrature_ptr)
+ {
+ const auto quadrature =
+ static_cast<const Quadrature<dim> *>(quadrature_ptr);
+
+ const auto weights = quadrature->get_weights();
+
+ boost::python::list weights_list;
+ for (const auto &w : weights)
+ weights_list.append(w);
+
+ return weights_list;
+ }
+ } // namespace internal
+
+
+
+ QuadratureWrapper::QuadratureWrapper(const int dim)
+ : dim(dim)
+ , quadrature_ptr(nullptr)
+ {
+ AssertThrow(dim == 2 || dim == 3, ExcMessage("Unsupported dimension."));
+ }
+
+
+
+ QuadratureWrapper::QuadratureWrapper(const QuadratureWrapper &other)
+ {
+ dim = other.dim;
+
+ AssertThrow(other.quadrature_ptr != nullptr,
+ ExcMessage("Underlying quadrature does not exist."));
+
+ if (dim == 2)
+ {
+ const auto quadrature =
+ static_cast<const Quadrature<2> *>(other.quadrature_ptr);
+ quadrature_ptr = new Quadrature<2>(*quadrature);
+ }
+ else if (dim == 3)
+ {
+ const auto quadrature =
+ static_cast<const Quadrature<3> *>(other.quadrature_ptr);
+ quadrature_ptr = new Quadrature<3>(*quadrature);
+ }
+ else
+ AssertThrow(false, ExcMessage("Given dimension is not implemented."));
+ }
+
+
+
+ QuadratureWrapper::~QuadratureWrapper()
+ {
+ if (dim != -1)
+ {
+ if (dim == 2)
+ {
+ // We cannot call delete on a void pointer so cast the void pointer
+ // back first.
+ Quadrature<2> *tmp = static_cast<Quadrature<2> *>(quadrature_ptr);
+ delete tmp;
+ }
+ else if (dim == 3)
+ {
+ Quadrature<3> *tmp = static_cast<Quadrature<3> *>(quadrature_ptr);
+ delete tmp;
+ }
+
+ dim = -1;
+ quadrature_ptr = nullptr;
+ }
+ }
+
+
+
+ void
+ QuadratureWrapper::create_gauss(const unsigned int n)
+ {
+ if (dim == 2)
+ {
+ quadrature_ptr = new QGauss<2>(n);
+ }
+ else if (dim == 3)
+ {
+ quadrature_ptr = new QGauss<3>(n);
+ }
+ else
+ AssertThrow(false, ExcMessage("Given dimension is not implemented."));
+ }
+
+
+
+ void
+ QuadratureWrapper::create_gauss_lobatto(const unsigned int n)
+ {
+ if (dim == 2)
+ {
+ quadrature_ptr = new QGaussLobatto<2>(n);
+ }
+ else if (dim == 3)
+ {
+ quadrature_ptr = new QGaussLobatto<3>(n);
+ }
+ else
+ AssertThrow(false, ExcMessage("Given dimension is not implemented."));
+ }
+
+
+
+ boost::python::list
+ QuadratureWrapper::get_points() const
+ {
+ if (dim == 2)
+ return internal::get_points<2>(quadrature_ptr);
+ else if (dim == 3)
+ return internal::get_points<3>(quadrature_ptr);
+ else
+ AssertThrow(false, ExcMessage("Given dimension is not implemented."));
+ }
+
+
+
+ boost::python::list
+ QuadratureWrapper::get_weights() const
+ {
+ if (dim == 2)
+ return internal::get_weights<2>(quadrature_ptr);
+ else if (dim == 3)
+ return internal::get_weights<3>(quadrature_ptr);
+ else
+ AssertThrow(false, ExcMessage("Given dimension is not implemented."));
+ }
+
+
+
+ void *
+ QuadratureWrapper::get_quadrature() const
+ {
+ return quadrature_ptr;
+ }
+
+
+
+ int
+ QuadratureWrapper::get_dim() const
+ {
+ return dim;
+ }
+
+} // namespace python
+
+DEAL_II_NAMESPACE_CLOSE
export_mapping();
void
export_manifold();
+ void
+ export_quadrature();
} // namespace python
DEAL_II_NAMESPACE_CLOSE
dealii::python::export_triangulation();
dealii::python::export_mapping();
dealii::python::export_manifold();
+ dealii::python::export_quadrature();
}
#else
dealii::python::export_triangulation();
dealii::python::export_mapping();
dealii::python::export_manifold();
+ dealii::python::export_quadrature();
}
#endif
--- /dev/null
+# ---------------------------------------------------------------------
+#
+# Copyright (C) 2016 by the deal.II authors
+#
+# This file is part of the deal.II library.
+#
+# The deal.II library is free software; you can use it, redistribute
+# it, and/or modify it under the terms of the GNU Lesser General
+# Public License as published by the Free Software Foundation; either
+# version 2.1 of the License, or (at your option) any later version.
+# The full text of the license can be found in the file LICENSE.md at
+# the top level directory of deal.II.
+#
+# ---------------------------------------------------------------------
+
+import unittest
+from PyDealII.Debug import *
+
+class TestQuadratureWrapper(unittest.TestCase):
+
+ def test_gauss(self):
+ quadrature = Quadrature(dim = 3)
+ quadrature.create_gauss(n = 1);
+
+ w_sum = 0.
+ for weight in quadrature.weights():
+ w_sum += weight
+
+ self.assertTrue(abs(1. - w_sum) < 1e-10)
+
+ def test_gauss_lobatto(self):
+ quadrature = Quadrature(dim = 3)
+ quadrature.create_gauss_lobatto(n = 2);
+
+ w_sum = 0.
+ for weight in quadrature.weights():
+ w_sum += weight
+
+ self.assertTrue(abs(1. - w_sum) < 1e-10)
+
+if __name__ == '__main__':
+ unittest.main()
--- /dev/null
+New: Python bindings for the Quadrature class.
+<br>
+(Alexander Grayver, 2020/01/20)