]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Update PolynomialsBDM degree and pass the degree to the base class constructor
authorgrahambenharper <grahambenharper@gmail.com>
Sun, 1 Sep 2019 21:06:02 +0000 (15:06 -0600)
committergrahambenharper <grahambenharper@gmail.com>
Tue, 3 Sep 2019 00:31:27 +0000 (18:31 -0600)
include/deal.II/base/polynomials_bdm.h
source/base/polynomials_bdm.cc
tests/base/polynomials_bdm_01.cc
tests/base/polynomials_bdm_02.cc

index 4e8219ecfd051613d0022245ae028afb24ecc50b..7e4cd085a79764f25327f468257fdacb1e529782 100644 (file)
@@ -131,13 +131,6 @@ public:
            std::vector<Tensor<4, dim>> &third_derivatives,
            std::vector<Tensor<5, dim>> &fourth_derivatives) const override;
 
-  /**
-   * Return the degree of the BDM space, which is one less than the highest
-   * polynomial degree.
-   */
-  unsigned int
-  degree() const;
-
   /**
    * Return the name of the space, which is <tt>BDM</tt>.
    */
@@ -203,14 +196,6 @@ private:
 };
 
 
-template <int dim>
-inline unsigned int
-PolynomialsBDM<dim>::degree() const
-{
-  return polynomial_space.degree() - 1;
-}
-
-
 template <int dim>
 inline std::string
 PolynomialsBDM<dim>::name() const
index 80a2abe0c1d67925662eb6e047b545de201a5622..953201b086aff6f730fbbe4a6dbcef6adabffa71 100644 (file)
@@ -28,7 +28,7 @@ DEAL_II_NAMESPACE_OPEN
 
 template <int dim>
 PolynomialsBDM<dim>::PolynomialsBDM(const unsigned int k)
-  : TensorPolynomialsBase<dim>(k, n_polynomials(k))
+  : TensorPolynomialsBase<dim>(k + 1, n_polynomials(k))
   , polynomial_space(Polynomials::Legendre::generate_complete_basis(k))
   , monomials((dim == 2) ? (1) : (k + 2))
   , p_values(polynomial_space.n())
@@ -193,7 +193,8 @@ PolynomialsBDM<dim>::evaluate(
               // p(t) = t^(i+1)
               monomials[i + 1].value(unit_point(d), monovali[d]);
               // q(t) = t^(k-i)
-              monomials[degree() - i].value(unit_point(d), monovalk[d]);
+              monomials[this->degree() - 1 - i].value(unit_point(d),
+                                                      monovalk[d]);
             }
           if (values.size() != 0)
             {
index 093185cf833246aaa464582dc9b0549a420acff9..bc4043b4593532200f59974445baa8538c04640d 100644 (file)
@@ -32,7 +32,7 @@ void
 plot(const PolynomialsBDM<dim> &poly)
 {
   QTrapez<1>                  base_quadrature;
-  QIterated<dim>              quadrature(base_quadrature, poly.degree() + 4);
+  QIterated<dim>              quadrature(base_quadrature, poly.degree() + 3);
   std::vector<Tensor<1, dim>> values(poly.n());
   std::vector<Tensor<2, dim>> grads;
   std::vector<Tensor<3, dim>> grads2;
@@ -41,10 +41,10 @@ plot(const PolynomialsBDM<dim> &poly)
 
   for (unsigned int k = 0; k < quadrature.size(); ++k)
     {
-      if (k % (poly.degree() + 5) == 0)
-        deallog << "BDM" << poly.degree() << '<' << dim << '>' << std::endl;
+      if (k % (poly.degree() + 4) == 0)
+        deallog << "BDM" << poly.degree() - 1 << '<' << dim << '>' << std::endl;
 
-      deallog << "BDM" << poly.degree() << '<' << dim << '>' << '\t'
+      deallog << "BDM" << poly.degree() - 1 << '<' << dim << '>' << '\t'
               << quadrature.point(k);
       poly.evaluate(
         quadrature.point(k), values, grads, grads2, thirds, fourths);
index 68834f1cdb7139dcd9c2283e63414cc998aa8f8a..6c3a0dcbb617b22fc098f153427b2fabadbc184e 100644 (file)
@@ -32,7 +32,7 @@ void
 plot(const PolynomialsBDM<dim> &poly)
 {
   const PolynomialSpace<dim> legendre_poly_space =
-    Polynomials::Legendre::generate_complete_basis(poly.degree());
+    Polynomials::Legendre::generate_complete_basis(poly.degree() - 1);
 
   const Point<3> p0(0, 0, 0);
   const Point<3> p1(0.25, 0.5, 0.75);
@@ -57,16 +57,16 @@ plot(const PolynomialsBDM<dim> &poly)
 
   for (unsigned int k = 0; k < points.size(); ++k)
     {
-      if (k % (poly.degree() + 4) == 0)
-        deallog << "BDM" << poly.degree() << '<' << dim << '>' << std::endl;
+      if (k % (poly.degree() + 3) == 0)
+        deallog << "BDM" << poly.degree() - 1 << '<' << dim << '>' << std::endl;
 
       unsigned int start = dim * n_sub;
 
-      deallog << "BDM" << poly.degree() << '<' << dim << '>' << points[k]
+      deallog << "BDM" << poly.degree() - 1 << '<' << dim << '>' << points[k]
               << std::endl;
       poly.evaluate(points[k], values, grads, grads2, thirds, fourths);
 
-      for (unsigned int i = 0; i < poly.degree() + 1; ++i, start += dim)
+      for (unsigned int i = 0; i < poly.degree(); ++i, start += dim)
         for (unsigned int j = 0; j < dim; ++j)
           {
             for (unsigned int d1 = 0; d1 < dim; ++d1)

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.