]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Implement the ScalarPolynomialsBase class 8651/head
authorgrahambenharper <grahambenharper@gmail.com>
Tue, 27 Aug 2019 02:40:12 +0000 (20:40 -0600)
committergrahambenharper <grahambenharper@gmail.com>
Tue, 27 Aug 2019 15:25:56 +0000 (09:25 -0600)
include/deal.II/base/scalar_polynomials_base.h [new file with mode: 0644]
include/deal.II/base/tensor_polynomials_base.h
source/base/CMakeLists.txt
source/base/scalar_polynomials_base.cc [new file with mode: 0644]

diff --git a/include/deal.II/base/scalar_polynomials_base.h b/include/deal.II/base/scalar_polynomials_base.h
new file mode 100644 (file)
index 0000000..78e263c
--- /dev/null
@@ -0,0 +1,175 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2005 - 2019 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_scalar_polynomials_base_h
+#define dealii_scalar_polynomials_base_h
+
+
+#include <deal.II/base/config.h>
+
+#include <deal.II/base/exceptions.h>
+#include <deal.II/base/point.h>
+#include <deal.II/base/tensor.h>
+
+#include <vector>
+
+DEAL_II_NAMESPACE_OPEN
+
+/**
+ * This class provides a framework for the finite element polynomial
+ * classes for use with finite element classes that are derived from
+ * FE_Poly. An object of this type (or rather of a type derived from
+ * this class) is stored as a member variable in each object of type
+ * FE_Poly.
+ *
+ * <h3>Deriving classes</h3>
+ *
+ * Any derived class must provide the most basic properties for shape
+ * functions evaluated on the reference cell. This includes, but is not
+ * limited to, implementing the evaluate(), name(), and
+ * clone() member functions. These functions are necessary to store the
+ * most basic information of how the polynomials in the derived class evaluate
+ * at a given point on the reference cell. More information on each function can
+ * be found in the corresponding function's documentation.
+ *
+ * Some classes that derive from this class include
+ * <ul>
+ *   <li> <tt>PolynomialsRannacherTurek</tt>
+ *   <li> <tt>PolynomialsP</tt>
+ *   <li> <tt>PolynomialSpace</tt>
+ *   <li> <tt>TensorProductPolynomials</tt>
+ *   <li> <tt>TensorProductPolynomialsConst</tt>
+ *   <li> <tt>TensorProductPolynomialsBubbles</tt>
+ * </ul>
+ *
+ * @ingroup Polynomials
+ * @author Graham Harper
+ * @date 2019
+ */
+template <int dim>
+class ScalarPolynomialsBase
+{
+public:
+  /**
+   * Constructor. This takes the degree of the space, @p deg from the finite element
+   * class, and @p n, the number of polynomials for the space.
+   */
+  ScalarPolynomialsBase(const unsigned int deg,
+                        const unsigned int n_polynomials);
+
+  /**
+   * Move constructor.
+   */
+  ScalarPolynomialsBase(ScalarPolynomialsBase<dim> &&) = default; // NOLINT
+
+  /**
+   * Copy constructor.
+   */
+  ScalarPolynomialsBase(const ScalarPolynomialsBase<dim> &) = default;
+
+  /**
+   * Virtual destructor. Makes sure that pointers to this class are deleted
+   * properly.
+   */
+  virtual ~ScalarPolynomialsBase() = default;
+
+  /**
+   * Compute the value and the derivatives of the polynomials at
+   * @p unit_point.
+   *
+   * The size of the vectors must either be zero or equal <tt>n()</tt>.  In
+   * the first case, the function will not compute these values.
+   *
+   * If you need values or derivatives of all polynomials then use this
+   * function, rather than using any of the <tt>compute_value</tt>,
+   * <tt>compute_grad</tt> or <tt>compute_grad_grad</tt> functions, see below,
+   * in a loop over all tensor product polynomials.
+   */
+  virtual void
+  evaluate(const Point<dim> &           unit_point,
+           std::vector<double> &        values,
+           std::vector<Tensor<1, dim>> &grads,
+           std::vector<Tensor<2, dim>> &grad_grads,
+           std::vector<Tensor<3, dim>> &third_derivatives,
+           std::vector<Tensor<4, dim>> &fourth_derivatives) const = 0;
+
+  /**
+   * Return the number of polynomials.
+   */
+  unsigned int
+  n() const;
+
+  /**
+   * Return the highest polynomial degree of polynomials represented by this
+   * class. A derived class may override this if its value is different from
+   * @p my_degree.
+   */
+  unsigned int
+  degree() const;
+
+  /**
+   * A sort of virtual copy constructor, this function returns a copy of
+   * the polynomial space object. Derived classes need to override the function
+   * here in this base class and return an object of the same type as the
+   * derived class.
+   *
+   * Some places in the library, for example the constructors of FE_Poly,
+   * need to make copies of polynomial spaces without knowing their exact type.
+   * They do so through this function.
+   */
+  virtual std::unique_ptr<ScalarPolynomialsBase<dim>>
+  clone() const = 0;
+
+  /**
+   * Return the name of the space.
+   */
+  virtual std::string
+  name() const = 0;
+
+private:
+  /**
+   * The highest polynomial degree of this functions represented by this object.
+   */
+  const unsigned int polynomial_degree;
+
+  /**
+   * The number of polynomials represented by this object.
+   */
+  const unsigned int n_pols;
+};
+
+
+
+template <int dim>
+inline unsigned int
+ScalarPolynomialsBase<dim>::n() const
+{
+  return n_pols;
+}
+
+
+
+template <int dim>
+inline unsigned int
+ScalarPolynomialsBase<dim>::degree() const
+{
+  return polynomial_degree;
+}
+
+
+
+DEAL_II_NAMESPACE_CLOSE
+
+#endif
index 49cfe2c4de63dbccc2c7a62797b7b05bfef44255..d7810d08793e9debc8a88fdc2375b1a800342d02 100644 (file)
@@ -31,7 +31,7 @@ DEAL_II_NAMESPACE_OPEN
  * This class provides a framework for the finite element polynomial
  * classes for use with finite element classes that are derived from
  * FE_PolyTensor. An object of this type (or rather of a type derived
- * from this class) is stored as a member variable each object of
+ * from this class) is stored as a member variable in each object of
  * type FE_PolyTensor.
  *
  * <h3>Deriving classes</h3>
index fa6554fb48de32b3a2dd5134cceb90da91ac6ec4..6920744d5fff9ea21b236aa0e1a1c6b84c632d39 100644 (file)
@@ -75,6 +75,7 @@ SET(_unity_include_src
   quadrature.cc
   quadrature_lib.cc
   quadrature_selector.cc
+  scalar_polynomials_base.cc
   subscriptor.cc
   table_handler.cc
   tensor_function.cc
diff --git a/source/base/scalar_polynomials_base.cc b/source/base/scalar_polynomials_base.cc
new file mode 100644 (file)
index 0000000..c067980
--- /dev/null
@@ -0,0 +1,42 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2005 - 2019 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 <deal.II/base/scalar_polynomials_base.h>
+#include <deal.II/base/thread_management.h>
+
+#include <iomanip>
+#include <iostream>
+
+DEAL_II_NAMESPACE_OPEN
+
+template <int dim>
+ScalarPolynomialsBase<dim>::ScalarPolynomialsBase(
+  const unsigned int deg,
+  const unsigned int n_polynomials)
+  : polynomial_degree(deg)
+  , n_pols(n_polynomials)
+{
+  // nothing to do here for now
+}
+
+
+
+template class ScalarPolynomialsBase<1>;
+template class ScalarPolynomialsBase<2>;
+template class ScalarPolynomialsBase<3>;
+
+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.