]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Implement compute_*_derivative
authorDaniel Arndt <arndtd@ornl.gov>
Thu, 21 May 2020 16:04:10 +0000 (12:04 -0400)
committerDaniel Arndt <arndtd@ornl.gov>
Thu, 21 May 2020 16:04:10 +0000 (12:04 -0400)
include/deal.II/base/polynomial_space.h
include/deal.II/base/polynomials_adini.h
include/deal.II/base/polynomials_rannacher_turek.h
include/deal.II/base/scalar_polynomials_base.h
include/deal.II/base/tensor_product_polynomials.h
include/deal.II/base/tensor_product_polynomials_bubbles.h
include/deal.II/base/tensor_product_polynomials_const.h

index 4e050b6225c3bdfa265ecfd036a757fda14655c4..e1260dd53ab11ce5fed944a957e16fddb3086e17 100644 (file)
@@ -171,6 +171,34 @@ public:
   Tensor<order, dim>
   compute_derivative(const unsigned int i, const Point<dim> &p) const;
 
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_1st_derivative()
+   */
+  virtual Tensor<1, dim>
+  compute_1st_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_2nd_derivative()
+   */
+  virtual Tensor<2, dim>
+  compute_2nd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_3rd_derivative()
+   */
+  virtual Tensor<3, dim>
+  compute_3rd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_4th_derivative()
+   */
+  virtual Tensor<4, dim>
+  compute_4th_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
   /**
    * Compute the gradient of the <tt>i</tt>th polynomial at unit point
    * <tt>p</tt>.
@@ -425,6 +453,45 @@ PolynomialSpace<dim>::compute_derivative(const unsigned int i,
 }
 
 
+
+template <int dim>
+inline Tensor<1, dim>
+PolynomialSpace<dim>::compute_1st_derivative(const unsigned int i,
+                                             const Point<dim> & p) const
+{
+  return compute_derivative<1>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<2, dim>
+PolynomialSpace<dim>::compute_2nd_derivative(const unsigned int i,
+                                             const Point<dim> & p) const
+{
+  return compute_derivative<2>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<3, dim>
+PolynomialSpace<dim>::compute_3rd_derivative(const unsigned int i,
+                                             const Point<dim> & p) const
+{
+  return compute_derivative<3>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<4, dim>
+PolynomialSpace<dim>::compute_4th_derivative(const unsigned int i,
+                                             const Point<dim> & p) const
+{
+  return compute_derivative<4>(i, p);
+}
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif
index 64daffad89a0ee75628ffcb223e6cef9150361e0..49eb1033f0d57c58a9ad590f95b258e75071ad1a 100644 (file)
@@ -80,6 +80,34 @@ public:
   double
   compute_value(const unsigned int i, const Point<dim> &p) const override;
 
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_1st_derivative()
+   */
+  virtual Tensor<1, dim>
+  compute_1st_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_2nd_derivative()
+   */
+  virtual Tensor<2, dim>
+  compute_2nd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_3rd_derivative()
+   */
+  virtual Tensor<3, dim>
+  compute_3rd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_4th_derivative()
+   */
+  virtual Tensor<4, dim>
+  compute_4th_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
   /**
    * Compute the gradient of the <tt>i</tt>th polynomial at
    * <tt>unit_point</tt>.
@@ -150,6 +178,50 @@ private:
 
 
 
+template <int dim>
+inline Tensor<1, dim>
+PolynomialsAdini<dim>::compute_1st_derivative(const unsigned int /*i*/,
+                                              const Point<dim> & /*p*/) const
+{
+  Assert(false, ExcNotImplemented());
+  return {};
+}
+
+
+
+template <int dim>
+inline Tensor<2, dim>
+PolynomialsAdini<dim>::compute_2nd_derivative(const unsigned int /*i*/,
+                                              const Point<dim> & /*p*/) const
+{
+  Assert(false, ExcNotImplemented());
+  return {};
+}
+
+
+
+template <int dim>
+inline Tensor<3, dim>
+PolynomialsAdini<dim>::compute_3rd_derivative(const unsigned int /*i*/,
+                                              const Point<dim> & /*p*/) const
+{
+  Assert(false, ExcNotImplemented());
+  return {};
+}
+
+
+
+template <int dim>
+inline Tensor<4, dim>
+PolynomialsAdini<dim>::compute_4th_derivative(const unsigned int /*i*/,
+                                              const Point<dim> & /*p*/) const
+{
+  Assert(false, ExcNotImplemented());
+  return {};
+}
+
+
+
 template <int dim>
 inline std::string
 PolynomialsAdini<dim>::name() const
index 4fb378458b4ea53b55e9027cf6901f9125a541ee..b411c6e74ca7012956e0d2a5b62a4cbb9ae68065 100644 (file)
@@ -69,6 +69,34 @@ public:
   Tensor<order, dim>
   compute_derivative(const unsigned int i, const Point<dim> &p) const;
 
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_1st_derivative()
+   */
+  virtual Tensor<1, dim>
+  compute_1st_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_2nd_derivative()
+   */
+  virtual Tensor<2, dim>
+  compute_2nd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_3rd_derivative()
+   */
+  virtual Tensor<3, dim>
+  compute_3rd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_4th_derivative()
+   */
+  virtual Tensor<4, dim>
+  compute_4th_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
   /**
    * Gradient of basis function @p i at @p p.
    */
@@ -220,6 +248,50 @@ PolynomialsRannacherTurek<dim>::compute_derivative(const unsigned int i,
 
 
 
+template <int dim>
+inline Tensor<1, dim>
+PolynomialsRannacherTurek<dim>::compute_1st_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<1>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<2, dim>
+PolynomialsRannacherTurek<dim>::compute_2nd_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<2>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<3, dim>
+PolynomialsRannacherTurek<dim>::compute_3rd_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<3>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<4, dim>
+PolynomialsRannacherTurek<dim>::compute_4th_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<4>(i, p);
+}
+
+
+
 template <int dim>
 inline std::string
 PolynomialsRannacherTurek<dim>::name() const
index dc702bbc06d2741d8c8396f0a8ee2f29787f4a96..bd56e819364037b55f76618401a34c98ea8d9bd4 100644 (file)
@@ -113,11 +113,7 @@ public:
    * Consider using compute() instead.
    */
   virtual double
-  compute_value(const unsigned int /*i*/, const Point<dim> & /*p*/) const
-  {
-    Assert(false, ExcNotImplemented());
-    return 0;
-  }
+  compute_value(const unsigned int /*i*/, const Point<dim> & /*p*/) const = 0;
 
   /**
    * Compute the <tt>order</tt>th derivative of the <tt>i</tt>th polynomial
@@ -129,12 +125,43 @@ public:
    */
   template <int order>
   Tensor<order, dim>
-  compute_derivative(const unsigned int /*i*/, const Point<dim> & /*p*/) const
-  {
-    Assert(false, ExcNotImplemented());
-    Tensor<order, dim> empty;
-    return empty;
-  }
+  compute_derivative(const unsigned int i, const Point<dim> &p) const;
+
+  /**
+   * Compute the first derivative of the <tt>i</tt>th polynomial
+   * at unit point <tt>p</tt>.
+   *
+   * Consider using evaluate() instead.
+   */
+  virtual Tensor<1, dim>
+  compute_1st_derivative(const unsigned int i, const Point<dim> &p) const = 0;
+
+  /**
+   * Compute the second derivative of the <tt>i</tt>th polynomial
+   * at unit point <tt>p</tt>.
+   *
+   * Consider using evaluate() instead.
+   */
+  virtual Tensor<2, dim>
+  compute_2nd_derivative(const unsigned int i, const Point<dim> &p) const = 0;
+
+  /**
+   * Compute the third derivative of the <tt>i</tt>th polynomial
+   * at unit point <tt>p</tt>.
+   *
+   * Consider using evaluate() instead.
+   */
+  virtual Tensor<3, dim>
+  compute_3rd_derivative(const unsigned int i, const Point<dim> &p) const = 0;
+
+  /**
+   * Compute the fourth derivative of the <tt>i</tt>th polynomial
+   * at unit point <tt>p</tt>.
+   *
+   * Consider using evaluate() instead.
+   */
+  virtual Tensor<4, dim>
+  compute_4th_derivative(const unsigned int i, const Point<dim> &p) const = 0;
 
   /**
    * Compute the gradient of the <tt>i</tt>th polynomial at unit point
@@ -143,12 +170,7 @@ public:
    * Consider using compute() instead.
    */
   virtual Tensor<1, dim>
-  compute_grad(const unsigned int /*i*/, const Point<dim> & /*p*/) const
-  {
-    Assert(false, ExcNotImplemented());
-    Tensor<1, dim> empty;
-    return empty;
-  }
+  compute_grad(const unsigned int /*i*/, const Point<dim> & /*p*/) const = 0;
 
   /**
    * Compute the second derivative (grad_grad) of the <tt>i</tt>th polynomial
@@ -157,12 +179,8 @@ public:
    * Consider using compute() instead.
    */
   virtual Tensor<2, dim>
-  compute_grad_grad(const unsigned int /*i*/, const Point<dim> & /*p*/) const
-  {
-    Assert(false, ExcNotImplemented());
-    Tensor<2, dim> empty;
-    return empty;
-  }
+  compute_grad_grad(const unsigned int /*i*/,
+                    const Point<dim> & /*p*/) const = 0;
 
   /**
    * Return the number of polynomials.
@@ -235,6 +253,37 @@ ScalarPolynomialsBase<dim>::degree() const
 
 
 
+template <int dim>
+template <int order>
+inline Tensor<order, dim>
+ScalarPolynomialsBase<dim>::compute_derivative(const unsigned int i,
+                                               const Point<dim> & p) const
+{
+  if (order == 1)
+    {
+      auto derivative = compute_1st_derivative(i, p);
+      return *reinterpret_cast<Tensor<order, dim> *>(&derivative);
+    }
+  if (order == 2)
+    {
+      auto derivative = compute_2nd_derivative(i, p);
+      return *reinterpret_cast<Tensor<order, dim> *>(&derivative);
+    }
+  if (order == 3)
+    {
+      auto derivative = compute_3rd_derivative(i, p);
+      return *reinterpret_cast<Tensor<order, dim> *>(&derivative);
+    }
+  if (order == 4)
+    {
+      auto derivative = compute_4th_derivative(i, p);
+      return *reinterpret_cast<Tensor<order, dim> *>(&derivative);
+    }
+  Assert(false, ExcNotImplemented());
+  Tensor<order, dim> empty;
+  return empty;
+}
+
 DEAL_II_NAMESPACE_CLOSE
 
 #endif
index e024553d6934c584718328299ba9b93332533b69..653fbcbd8bfe7fb0baa2da5070e66bbca0c82360 100644 (file)
@@ -167,6 +167,34 @@ public:
   Tensor<order, dim>
   compute_derivative(const unsigned int i, const Point<dim> &p) const;
 
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_1st_derivative()
+   */
+  virtual Tensor<1, dim>
+  compute_1st_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_2nd_derivative()
+   */
+  virtual Tensor<2, dim>
+  compute_2nd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_3rd_derivative()
+   */
+  virtual Tensor<3, dim>
+  compute_3rd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_4th_derivative()
+   */
+  virtual Tensor<4, dim>
+  compute_4th_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
   /**
    * Compute the grad of the <tt>i</tt>th tensor product polynomial at
    * <tt>unit_point</tt>. Here <tt>i</tt> is given in tensor product
@@ -360,6 +388,34 @@ public:
   Tensor<order, dim>
   compute_derivative(const unsigned int i, const Point<dim> &p) const;
 
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_1st_derivative()
+   */
+  virtual Tensor<1, dim>
+  compute_1st_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_2nd_derivative()
+   */
+  virtual Tensor<2, dim>
+  compute_2nd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_3rd_derivative()
+   */
+  virtual Tensor<3, dim>
+  compute_3rd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_4th_derivative()
+   */
+  virtual Tensor<4, dim>
+  compute_4th_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
   /**
    * Compute the grad of the <tt>i</tt>th tensor product polynomial at
    * <tt>unit_point</tt>. Here <tt>i</tt> is given in tensor product
@@ -606,6 +662,52 @@ TensorProductPolynomials<dim, PolynomialType>::compute_derivative(
     }
 }
 
+
+
+template <int dim, typename PolynomialType>
+inline Tensor<1, dim>
+TensorProductPolynomials<dim, PolynomialType>::compute_1st_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<1>(i, p);
+}
+
+
+
+template <int dim, typename PolynomialType>
+inline Tensor<2, dim>
+TensorProductPolynomials<dim, PolynomialType>::compute_2nd_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<2>(i, p);
+}
+
+
+
+template <int dim, typename PolynomialType>
+inline Tensor<3, dim>
+TensorProductPolynomials<dim, PolynomialType>::compute_3rd_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<3>(i, p);
+}
+
+
+
+template <int dim, typename PolynomialType>
+inline Tensor<4, dim>
+TensorProductPolynomials<dim, PolynomialType>::compute_4th_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<4>(i, p);
+}
+
+
+
 template <int dim>
 template <int order>
 Tensor<order, dim>
@@ -725,6 +827,47 @@ AnisotropicPolynomials<dim>::compute_derivative(const unsigned int i,
 }
 
 
+
+template <int dim>
+inline Tensor<1, dim>
+AnisotropicPolynomials<dim>::compute_1st_derivative(const unsigned int i,
+                                                    const Point<dim> & p) const
+{
+  return compute_derivative<1>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<2, dim>
+AnisotropicPolynomials<dim>::compute_2nd_derivative(const unsigned int i,
+                                                    const Point<dim> & p) const
+{
+  return compute_derivative<2>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<3, dim>
+AnisotropicPolynomials<dim>::compute_3rd_derivative(const unsigned int i,
+                                                    const Point<dim> & p) const
+{
+  return compute_derivative<3>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<4, dim>
+AnisotropicPolynomials<dim>::compute_4th_derivative(const unsigned int i,
+                                                    const Point<dim> & p) const
+{
+  return compute_derivative<4>(i, p);
+}
+
+
+
 template <int dim>
 inline std::string
 AnisotropicPolynomials<dim>::name() const
index 0e8d28823f03a32182ef1cca99ecd0d27b2e5c50..59baa0366ecdb85d5e339f2c9a9912182e081b3a 100644 (file)
@@ -138,6 +138,34 @@ public:
   Tensor<order, dim>
   compute_derivative(const unsigned int i, const Point<dim> &p) const;
 
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_1st_derivative()
+   */
+  virtual Tensor<1, dim>
+  compute_1st_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_2nd_derivative()
+   */
+  virtual Tensor<2, dim>
+  compute_2nd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_3rd_derivative()
+   */
+  virtual Tensor<3, dim>
+  compute_3rd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_4th_derivative()
+   */
+  virtual Tensor<4, dim>
+  compute_4th_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
   /**
    * Compute the grad of the <tt>i</tt>th tensor product polynomial at
    * <tt>unit_point</tt>. Here <tt>i</tt> is given in tensor product
@@ -427,6 +455,49 @@ TensorProductPolynomialsBubbles<dim>::compute_derivative(
 }
 
 
+
+template <int dim>
+inline Tensor<1, dim>
+TensorProductPolynomialsBubbles<dim>::compute_1st_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<1>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<2, dim>
+TensorProductPolynomialsBubbles<dim>::compute_2nd_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<2>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<3, dim>
+TensorProductPolynomialsBubbles<dim>::compute_3rd_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<3>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<4, dim>
+TensorProductPolynomialsBubbles<dim>::compute_4th_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<4>(i, p);
+}
+
 #endif // DOXYGEN
 DEAL_II_NAMESPACE_CLOSE
 
index 0c75e2f59c16fd5839c00be0cec893a81933ff85..274f247492758af21eea44481b22ba074571198a 100644 (file)
@@ -141,6 +141,34 @@ public:
   Tensor<order, dim>
   compute_derivative(const unsigned int i, const Point<dim> &p) const;
 
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_1st_derivative()
+   */
+  virtual Tensor<1, dim>
+  compute_1st_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_2nd_derivative()
+   */
+  virtual Tensor<2, dim>
+  compute_2nd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_3rd_derivative()
+   */
+  virtual Tensor<3, dim>
+  compute_3rd_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
+  /**
+   * @copydoc ScalarPolynomialsBase::compute_4th_derivative()
+   */
+  virtual Tensor<4, dim>
+  compute_4th_derivative(const unsigned int i,
+                         const Point<dim> & p) const override;
+
   /**
    * Compute the grad of the <tt>i</tt>th tensor product polynomial at
    * <tt>unit_point</tt>. Here <tt>i</tt> is given in tensor product
@@ -287,6 +315,49 @@ TensorProductPolynomialsConst<dim>::compute_derivative(
 }
 
 
+
+template <int dim>
+inline Tensor<1, dim>
+TensorProductPolynomialsConst<dim>::compute_1st_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<1>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<2, dim>
+TensorProductPolynomialsConst<dim>::compute_2nd_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<2>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<3, dim>
+TensorProductPolynomialsConst<dim>::compute_3rd_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<3>(i, p);
+}
+
+
+
+template <int dim>
+inline Tensor<4, dim>
+TensorProductPolynomialsConst<dim>::compute_4th_derivative(
+  const unsigned int i,
+  const Point<dim> & p) const
+{
+  return compute_derivative<4>(i, p);
+}
+
 #endif // DOXYGEN
 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.