From: Giuseppe Pitton Date: Mon, 11 May 2015 19:34:07 +0000 (+0200) Subject: fixed comments and static member functions X-Git-Tag: v8.3.0-rc1~171^2~3 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=e1d9a3c1b4d952f2dfaac5345eb7a6253c5cc31f;p=dealii.git fixed comments and static member functions fixed indentation removed static qualifier from EndPoint ep removed useless comments --- diff --git a/include/deal.II/base/quadrature_lib.h b/include/deal.II/base/quadrature_lib.h index 33a633563a..ad7d5b23e6 100644 --- a/include/deal.II/base/quadrature_lib.h +++ b/include/deal.II/base/quadrature_lib.h @@ -517,9 +517,9 @@ public: * $w(x) = 1/\sqrt{1-x^2}$. The nodes and weights are known analytically, * and are exact for monomials up to the order $2n-1$, where $n$ is the number * of quadrature points. -* Here we rescale the quadrature formula so that it is defined on t -* he interval [0,1] instead of [-1,1]. So the quadrature formulas inte -* grate exactly the integral $\int_0^1 f(x) w(x) dx$ with the weight: +* Here we rescale the quadrature formula so that it is defined on +* the interval $[0,1]$ instead of $[-1,1]$. So the quadrature formulas +* integrate exactly the integral $\int_0^1 f(x) w(x) dx$ with the weight: * $w(x) = 1/sqrt{x(1-x)}$. * For details see: * M. Abramowitz & I.A. Stegun: Handbook of Mathematical Functions, par. 25.4.38 @@ -530,17 +530,17 @@ template class QGaussChebyshev : public Quadrature { public: - //Generate a formula with n quadrature points + /// Generate a formula with n quadrature points QGaussChebyshev(const unsigned int n); -protected: - // Sets the points of the quadrature formula. - std::vector - set_quadrature_points(const unsigned int n) const; +private: + /// Sets the points of the quadrature formula. + static std::vector + get_quadrature_points(const unsigned int n); - // Sets the weights of the quadrature formula. - std::vector - set_quadrature_weights(const unsigned int n) const; + /// Sets the weights of the quadrature formula. + static std::vector + get_quadrature_weights(const unsigned int n); }; @@ -552,10 +552,9 @@ protected: * lies at one of the two extrema of the interval. * The nodes and weights are known analytically, * and are exact for monomials up to the order $2n-2$, where $n$ is the number -* of quadrature points. -* Here we rescale the quadrature formula so that it is defined on t -* he interval [0,1] instead of [-1,1]. So the quadrature formulas inte -* grate exactly the integral $\int_0^1 f(x) w(x) dx$ with the weight: +* of quadrature points. Here we rescale the quadrature formula so that it is defined on +* the interval $[0,1]$ instead of $[-1,1]$. So the quadrature formulas +* integrate exactly the integral $\int_0^1 f(x) w(x) dx$ with the weight: * $w(x) = 1/sqrt{x(1-x)}$. By default the quadrature is constructed with the * left endpoint as quadrature node, but the quadrature node can be imposed at the * right endpoint through the variable ep that can assume the values left or right. @@ -566,32 +565,37 @@ template class QGaussRadauChebyshev : public Quadrature { public: + /* EndPoint is used to specify which of the two endpoints of the unit interval + * is used also as quadrature point + */ enum EndPoint { left,right }; - EndPoint ep; - //Generate a formula with n quadrature points + /// Generate a formula with n quadrature points QGaussRadauChebyshev(const unsigned int n, QGaussRadauChebyshev::EndPoint ep=QGaussRadauChebyshev::left); -protected: - // Sets the points of the quadrature formula. - std::vector - set_quadrature_points(const unsigned int n) const; +private: + const EndPoint ep; + /// Sets the points of the quadrature formula. + static std::vector + get_quadrature_points(const unsigned int n, EndPoint ep); - // Sets the weights of the quadrature formula. - std::vector - set_quadrature_weights(const unsigned int n) const; + /// Sets the weights of the quadrature formula. + static std::vector + get_quadrature_weights(const unsigned int n, EndPoint ep); }; /** * Gauss-Lobatto-Chebyshev quadrature rules integrate the weighted product * $\int_{-1}^1 f(x) w(x) dx$ with weight given by: -* $w(x) = 1/\sqrt{1-x^2}$. The nodes and weights are known analytically, -* and are exact for monomials up to the order $2n-1$, where $n$ is the number +* $w(x) = 1/\sqrt{1-x^2}$, with the additional constraint that two of the quadrature +* points are located at the endpoints of the quadrature interval. +* The nodes and weights are known analytically, +* and are exact for monomials up to the order $2n-3$, where $n$ is the number * of quadrature points. -* Here we rescale the quadrature formula so that it is defined on t -* he interval [0,1] instead of [-1,1]. So the quadrature formulas inte -* grate exactly the integral $\int_0^1 f(x) w(x) dx$ with the weight: +* Here we rescale the quadrature formula so that it is defined on +* the interval $[0,1]$ instead of $[-1,1]$. So the quadrature formulas +* integrate exactly the integral $\int_0^1 f(x) w(x) dx$ with the weight: * $w(x) = 1/sqrt{x(1-x)}$. * For details see: * M. Abramowitz & I.A. Stegun: Handbook of Mathematical Functions, par. 25.4.40 @@ -602,17 +606,17 @@ template class QGaussLobattoChebyshev : public Quadrature { public: - //Generate a formula with n quadrature points + /// Generate a formula with n quadrature points QGaussLobattoChebyshev(const unsigned int n); -protected: - // Sets the points of the quadrature formula. - std::vector - set_quadrature_points(const unsigned int n) const; +private: + /// Sets the points of the quadrature formula. + static std::vector + get_quadrature_points(const unsigned int n); - // Sets the weights of the quadrature formula. - std::vector - set_quadrature_weights(const unsigned int n) const; + /// Sets the weights of the quadrature formula. + static std::vector + get_quadrature_weights(const unsigned int n); }; diff --git a/source/base/quadrature_lib.cc b/source/base/quadrature_lib.cc index 0abe19989e..d9724b465f 100644 --- a/source/base/quadrature_lib.cc +++ b/source/base/quadrature_lib.cc @@ -1132,7 +1132,7 @@ QTelles<1>::QTelles ( template <> std::vector -QGaussChebyshev<1>::set_quadrature_points(const unsigned int n) const +QGaussChebyshev<1>::get_quadrature_points(const unsigned int n) { std::vector points(n); @@ -1149,7 +1149,7 @@ QGaussChebyshev<1>::set_quadrature_points(const unsigned int n) const template <> std::vector -QGaussChebyshev<1>::set_quadrature_weights(const unsigned int n) const +QGaussChebyshev<1>::get_quadrature_weights(const unsigned int n) { std::vector weights(n); @@ -1172,8 +1172,8 @@ QGaussChebyshev<1>::QGaussChebyshev(const unsigned int n) { Assert(n>0,ExcMessage("Need at least one point for the quadrature rule")); - std::vector p=set_quadrature_points(n); - std::vector w=set_quadrature_weights(n); + std::vector p=get_quadrature_points(n); + std::vector w=get_quadrature_weights(n); for (unsigned int i=0; isize(); ++i) { @@ -1196,7 +1196,8 @@ QGaussChebyshev::QGaussChebyshev (const unsigned int n) template <> std::vector -QGaussRadauChebyshev<1>::set_quadrature_points(const unsigned int n) const +QGaussRadauChebyshev<1>::get_quadrature_points(const unsigned int n, + QGaussRadauChebyshev::EndPoint ep) { std::vector points(n); @@ -1219,7 +1220,8 @@ QGaussRadauChebyshev<1>::set_quadrature_points(const unsigned int n) const template <> std::vector -QGaussRadauChebyshev<1>::set_quadrature_weights(const unsigned int n) const +QGaussRadauChebyshev<1>::get_quadrature_weights(const unsigned int n, + QGaussRadauChebyshev::EndPoint ep) { std::vector weights(n); @@ -1248,8 +1250,8 @@ QGaussRadauChebyshev<1>::QGaussRadauChebyshev(const unsigned int n, { Assert(n>0,ExcMessage("Need at least one point for quadrature rules")); - std::vector p=set_quadrature_points(n); - std::vector w=set_quadrature_weights(n); + std::vector p=get_quadrature_points(n,ep); + std::vector w=get_quadrature_weights(n,ep); for (unsigned int i=0; isize(); ++i) { @@ -1264,7 +1266,8 @@ QGaussRadauChebyshev<2>::QGaussRadauChebyshev (const unsigned int n, QGaussRadauChebyshev::EndPoint ep) : Quadrature<2> (QGaussRadauChebyshev<1>(n, static_cast::EndPoint>(ep)), - QGaussRadauChebyshev<1>(n, static_cast::EndPoint>(ep))) + QGaussRadauChebyshev<1>(n, static_cast::EndPoint>(ep))), + ep (ep) {} @@ -1273,13 +1276,14 @@ QGaussRadauChebyshev::QGaussRadauChebyshev (const unsigned int n, QGaussRadauChebyshev::EndPoint ep) : Quadrature (QGaussRadauChebyshev(n,static_cast::EndPoint>(ep)), - QGaussRadauChebyshev<1>(n,static_cast::EndPoint>(ep))) + QGaussRadauChebyshev<1>(n,static_cast::EndPoint>(ep))), + ep (ep) {} template <> std::vector -QGaussLobattoChebyshev<1>::set_quadrature_points(const unsigned int n) const +QGaussLobattoChebyshev<1>::get_quadrature_points(const unsigned int n) { std::vector points(n); @@ -1296,7 +1300,7 @@ QGaussLobattoChebyshev<1>::set_quadrature_points(const unsigned int n) const template <> std::vector -QGaussLobattoChebyshev<1>::set_quadrature_weights(const unsigned int n) const +QGaussLobattoChebyshev<1>::get_quadrature_weights(const unsigned int n) { std::vector weights(n); @@ -1321,8 +1325,8 @@ QGaussLobattoChebyshev<1>::QGaussLobattoChebyshev(const unsigned int n) { Assert(n>1,ExcMessage("Need at least two points for Gauss-Lobatto quadrature rule")); - std::vector p=set_quadrature_points(n); - std::vector w=set_quadrature_weights(n); + std::vector p=get_quadrature_points(n); + std::vector w=get_quadrature_weights(n); for (unsigned int i=0; isize(); ++i) {