From: Wolfgang Bangerth Date: Tue, 28 May 2002 06:57:24 +0000 (+0000) Subject: Fix a number of copy-paste documentation errors where someone seems to have left... X-Git-Tag: v8.0.0~17998 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=6afdd6e27adbfdf222193c95aa2aec84337abc4a;p=dealii.git Fix a number of copy-paste documentation errors where someone seems to have left in comments that related to other classes and where plain wrong in the present context. Remove default constructors of two polynomial classes since they left the object in a generally unusable state. Add static functions to Legendre and LagrangeEquidistant classes to generate a full basis of polynomials. git-svn-id: https://svn.dealii.org/trunk@5888 0785d39b-7218-0410-832d-ea1e28bc413d --- diff --git a/deal.II/base/include/base/polynomial.h b/deal.II/base/include/base/polynomial.h index ef5fda2220..f15e4a2c26 100644 --- a/deal.II/base/include/base/polynomial.h +++ b/deal.II/base/include/base/polynomial.h @@ -51,11 +51,6 @@ class Polynomial : public Subscriptor * minus one. */ Polynomial (const typename std::vector &coefficients); - - /** - * Default-Constructor. - */ - Polynomial (); /** * Return the value of this @@ -116,6 +111,10 @@ class Polynomial : public Subscriptor * of this class and may be * passed down by derived * classes. + * + * This vector cannot be constant + * since we want to allow copying + * of polynomials. */ typename std::vector coefficients; }; @@ -123,16 +122,18 @@ class Polynomial : public Subscriptor /** - * Lagrange polynomials with equidistant interpolation - * points in [0,1]. The polynomial of degree @p{n} has got @p{n+1} interpolation + * Lagrange polynomials with equidistant interpolation points in + * [0,1]. The polynomial of degree @p{n} has got @p{n+1} interpolation * points. The interpolation points are sorted in ascending * order. This order gives an index to each interpolation point. A - * Lagrangian polynomial equals to 1 at its `support point', - * and 0 at all other interpolation - * points. For example, if the degree is 3, and the support point is 1, - * then the polynomial represented by this object is cubic and its - * value is 1 at the point @p{x=1/3}, and zero at the point @p{x=0}, - * @p{x=2/3}, and @p{x=1}. + * Lagrangian polynomial equals to 1 at its `support point', and 0 at + * all other interpolation points. For example, if the degree is 3, + * and the support point is 1, then the polynomial represented by this + * object is cubic and its value is 1 at the point @p{x=1/3}, and zero + * at the point @p{x=0}, @p{x=2/3}, and @p{x=1}. All the polynomials + * have polynomial order equal to @p{degree}, but together they span + * the entire space of polynomials of degree less than or equal + * @p{degree}. * * The Lagrange polynomials are implemented up to degree 10. * @@ -153,11 +154,26 @@ class LagrangeEquidistant: public Polynomial LagrangeEquidistant (const unsigned int n, const unsigned int support_point); - /** - * Default-constructor. + /** + * Return a vector of polynomial + * objects of order @p{degree}, + * which then spans the full + * space of polynomials up to the + * given degree. The polynomials + * are generated by calling the + * destructor of this class with + * the same degree but support + * point running from zero to + * @p{degree}. This function may + * be used to initialize the + * @ref{TensorProductPolynomials} + * and @ref{PolynomialSpace} + * classes. */ - LagrangeEquidistant (); - + static + std::vector > + generate_complete_basis (const unsigned int degree); + private: /** @@ -195,6 +211,22 @@ class Legendre : public Polynomial */ Legendre (const unsigned int k); + /** + * Return a vector of Legendre + * polynomial objects of orders + * zero through @p{degree}, which + * then spans the full space of + * polynomials up to the given + * degree. This function may be + * used to initialize the + * @ref{TensorProductPolynomials} + * and @ref{PolynomialSpace} + * classes. + */ + static + std::vector > + generate_complete_basis (const unsigned int degree); + private: /** * Vector with already computed @@ -225,6 +257,8 @@ class Legendre : public Polynomial }; +/* -------------------------- inline functions --------------------- */ + template inline unsigned int diff --git a/deal.II/base/include/base/polynomial_space.h b/deal.II/base/include/base/polynomial_space.h index 6fd8423abe..2cd7a37a1f 100644 --- a/deal.II/base/include/base/polynomial_space.h +++ b/deal.II/base/include/base/polynomial_space.h @@ -25,12 +25,13 @@ /** - * Polynomial space of degree at most n in higher dimensions. + * Representation of the space of polynomials of degree at most n in + * higher dimensions. * * Given a vector of @{n} one-dimensional polynomials @{P0} to @{Pn}, * where @{Pi} has degree @p{i}, this class generates all polynomials - * the form @p{ Pijk(x,y,z) = Pi(x)Pj(y)Pk(z)}, where the sum of - * @p{i}, @p{j} and @p{k} is below/equal @p{n}. + * of the form @p{ Pijk(x,y,z) = Pi(x)Pj(y)Pk(z)}, where the sum of + * @p{i}, @p{j} and @p{k} is less than or equal @p{n}. * * @author Guido Kanschat, 2002 */ @@ -43,7 +44,15 @@ class PolynomialSpace * vector of pointers to * one-dimensional polynomials * and will be copied into the - * member variable @p{polynomials}. + * member variable + * @p{polynomials}. The static + * type of the template argument + * @p{pols} needs to be + * convertible to + * @p{Polynomial}, + * i.e. should usually be a + * derived class of + * @p{Polynomial}. */ template PolynomialSpace(const typename std::vector &pols); @@ -70,10 +79,10 @@ class PolynomialSpace * functions, see below, in a * loop over all polynomials. */ - void compute(const Point &unit_point, - std::vector &values, - typename std::vector > &grads, - typename std::vector > &grad_grads) const; + void compute (const Point &unit_point, + std::vector &values, + typename std::vector > &grads, + typename std::vector > &grad_grads) const; /** * Computes the value of the @@ -107,9 +116,16 @@ class PolynomialSpace const Point &p) const; /** - * Returns the number of tensor - * product polynomials. For $n$ - * 1d polynomials this is $n^dim$. + * Return the number of + * polynomials spanning the space + * represented by this + * class. Here, if @p{N} is the + * number of one-dimensional + * polynomials given, then the + * result of this function is + * @p{N} in 1d, @p{N(N+1)/2} in + * 2d, and @p{N(N+1)(N+2)/6 in + * 3d. */ unsigned int n() const; @@ -127,42 +143,33 @@ class PolynomialSpace * polynomials given to the * constructor. */ - std::vector > polynomials; + const std::vector > polynomials; /** - * Number of tensor product - * polynomials. For $n$ 1d - * polynomials this is $n^dim$. + * Store the precomputed value + * which the @p{n()} function + * returns. */ - unsigned int n_pols; - + const unsigned int n_pols; + /** - * Computes @p{x} to the power of - * @p{y} for unsigned int @p{x} - * and @p{y}. It is a private - * function as it is only used in - * this class. + * Static function used in the + * constructor to compute the + * number of polynomials. */ - static unsigned int power(const unsigned int x, const unsigned int y); + static unsigned int compute_n_pols (const unsigned int n); }; template template -PolynomialSpace::PolynomialSpace( - const typename std::vector &pols): - polynomials (pols.begin(), pols.end()) -{ - const unsigned int n=polynomials.size(); - - n_pols = n; - for (unsigned int i=1;i:: +PolynomialSpace (const typename std::vector &pols) + : + polynomials (pols.begin(), pols.end()), + n_pols (compute_n_pols(polynomials.size())) +{} diff --git a/deal.II/base/include/base/tensor_product_polynomials.h b/deal.II/base/include/base/tensor_product_polynomials.h index e56e1a2c1c..32e140eea2 100644 --- a/deal.II/base/include/base/tensor_product_polynomials.h +++ b/deal.II/base/include/base/tensor_product_polynomials.h @@ -182,7 +182,7 @@ class TensorProductPolynomials * polynomials given to the * constructor. */ - std::vector > polynomials; + const std::vector > polynomials; /** * Number of tensor product @@ -213,11 +213,12 @@ class TensorProductPolynomials template template -TensorProductPolynomials::TensorProductPolynomials( - const typename std::vector &pols): +TensorProductPolynomials:: +TensorProductPolynomials(const typename std::vector &pols) + : polynomials (pols.begin(), pols.end()), n_tensor_pols(power(pols.size(), dim)), - n_pols_to(dim+1) + n_pols_to(dim+1, 0) { const unsigned int n_pols=polynomials.size(); diff --git a/deal.II/base/source/legendre.cc b/deal.II/base/source/legendre.cc index 230dcc9a5c..2344a45f0f 100644 --- a/deal.II/base/source/legendre.cc +++ b/deal.II/base/source/legendre.cc @@ -165,5 +165,18 @@ Legendre::Legendre (const unsigned int k) +template +std::vector > +Legendre::generate_complete_basis (const unsigned int degree) +{ + std::vector > v; + v.reserve(degree+1); + for (unsigned int i=0; i<=degree; ++i) + v.push_back (Legendre(i)); + return v; +}; + + + // explicit instantiations template class Legendre; diff --git a/deal.II/base/source/polynomial.cc b/deal.II/base/source/polynomial.cc index 8952b6ecfb..71b4269030 100644 --- a/deal.II/base/source/polynomial.cc +++ b/deal.II/base/source/polynomial.cc @@ -21,14 +21,6 @@ Polynomial::Polynomial (const typename std::vector &a): -template -Polynomial::Polynomial () - : - coefficients(0) -{} - - - template number Polynomial::value (const number x) const @@ -110,11 +102,6 @@ LagrangeEquidistant::LagrangeEquidistant (const unsigned int n, -LagrangeEquidistant::LagrangeEquidistant () -{} - - - std::vector LagrangeEquidistant::compute_coefficients (const unsigned int n, const unsigned int support_point) @@ -336,6 +323,28 @@ LagrangeEquidistant::compute_coefficients (const unsigned int n, return a; } + +std::vector > +LagrangeEquidistant:: +generate_complete_basis (const unsigned int degree) +{ + if (degree==0) + // create constant polynomial + return std::vector > + (1, Polynomial (std::vector (1,1.))); + else + { + // create array of Lagrange + // polynomials + std::vector > v; + for (unsigned int i=0; i<=degree; ++i) + v.push_back(LagrangeEquidistant(degree,i)); + return v; + }; +}; + + + template class Polynomial; template class Polynomial; template class Polynomial; diff --git a/deal.II/base/source/polynomial_space.cc b/deal.II/base/source/polynomial_space.cc index 15b1000abe..460a9f9d9d 100644 --- a/deal.II/base/source/polynomial_space.cc +++ b/deal.II/base/source/polynomial_space.cc @@ -16,16 +16,18 @@ #include - template -unsigned int PolynomialSpace::power(const unsigned int x, - const unsigned int y) +unsigned int +PolynomialSpace::compute_n_pols (const unsigned int n) { - unsigned int value=1; - for (unsigned int i=0; ibase
    +
  1. + New: The Legendre and + LagrangeEquidistant classes now have + static member functions generate_complete_basis which returns an + array of polynomial objects spanning the complete space up to a + specified order in 1d. This may be used to generate the + respective polynomial spaces in higher space dimensions. +
    + (WB 2002/05/27) +

    + +
  2. + Changed: The Polynomial and + LagrangeEquidistant classes have lost + their default constructor, as that did not make much sense + anyway. +
    + (WB 2002/05/27) +

    +
  3. Fixed: When forward declaring the Tensor class, we now also forward declare diff --git a/tests/base/polynomial1d.cc b/tests/base/polynomial1d.cc index 1e9636d5f3..6e00b697c2 100644 --- a/tests/base/polynomial1d.cc +++ b/tests/base/polynomial1d.cc @@ -2,7 +2,7 @@ // $Id$ // Version: $Name$ // -// Copyright (C) 2000, 2001 by the deal.II authors +// Copyright (C) 2000, 2001, 2002 by the deal.II authors // // This file is subject to QPL and may not be distributed // without copyright and license information. Please refer @@ -44,12 +44,12 @@ int main () deallog.attach(logfile); deallog.depth_console(0); - std::vector > p (15); + std::vector > p; deallog << "Legendre" << endl; - for (unsigned int i=0;i(i); + for (unsigned int i=0;i<15;++i) + p.push_back (Legendre(i)); for (unsigned int i=0;i > p(3); - for (unsigned int i=0;i > p; + for (unsigned int i=0;i<3;++i) + p.push_back (LagrangeEquidistant(3, i)); check_dimensions(p); - for (unsigned int i=0;i(i); + p.clear (); + for (unsigned int i=0;i<3;++i) + p.push_back (Legendre(i)); check_dimensions(p); }