From: Denis Davydov Date: Fri, 3 Feb 2017 15:52:19 +0000 (+0100) Subject: Functions::Spherical add component argument to svalue/sgradient/shessian X-Git-Tag: v8.5.0-rc1~156^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=09ab92a2f9ba1655a4f144a95c1dd967f5a45a04;p=dealii.git Functions::Spherical add component argument to svalue/sgradient/shessian --- diff --git a/include/deal.II/base/function_spherical.h b/include/deal.II/base/function_spherical.h index fefac59e01..a2de2db170 100644 --- a/include/deal.II/base/function_spherical.h +++ b/include/deal.II/base/function_spherical.h @@ -48,8 +48,13 @@ namespace Functions /** * Constructor which should be provided with @p center defining the origin * of the coordinate system. + * + * Note that components of this function are treated as entirely separate + * quantities -- not as the components of a vector that will be + * re-interpreted in a different coordinate system. */ - Spherical(const Point ¢er = Point()); + Spherical(const Point ¢er = Point(), + const unsigned int n_components=1); /** * Return the value of the function at the given point. @@ -87,7 +92,8 @@ namespace Functions * Return the value at point @p sp. Here, @p sp is provided in spherical * coordinates. */ - virtual double svalue(const std_cxx11::array &sp) const; + virtual double svalue(const std_cxx11::array &sp, + const unsigned int component) const; /** * Return the gradient in spherical coordinates. @@ -95,7 +101,8 @@ namespace Functions * The returned object should contain derivatives in the following order: * $\{ f_{,r},\, f_{,\theta},\, f_{,\phi}\}$. */ - virtual std_cxx11::array sgradient(const std_cxx11::array &sp) const; + virtual std_cxx11::array sgradient(const std_cxx11::array &sp, + const unsigned int component) const; /** * Return the Hessian in spherical coordinates. @@ -103,7 +110,8 @@ namespace Functions * The returned object should contain derivatives in the following order: * $\{ f_{,rr},\, f_{,\theta\theta},\, f_{,\phi\phi},\, f_{,r\theta},\, f_{,r\phi},\, f_{,\theta\phi}\}$. */ - virtual std_cxx11::array shessian (const std_cxx11::array &sp) const; + virtual std_cxx11::array shessian (const std_cxx11::array &sp, + const unsigned int component) const; /** * A vector from the origin to the center of spherical coordinate system. diff --git a/source/base/function_spherical.cc b/source/base/function_spherical.cc index 06b31fe6d3..be08727ef2 100644 --- a/source/base/function_spherical.cc +++ b/source/base/function_spherical.cc @@ -153,9 +153,10 @@ namespace Functions template - Spherical::Spherical(const Point &p) + Spherical::Spherical(const Point &p, + const unsigned int n_components) : - Function(1), + Function(n_components), coordinate_system_offset(p) { AssertThrow(dim==3, @@ -167,11 +168,11 @@ namespace Functions template double Spherical::value (const Point &p_, - const unsigned int) const + const unsigned int component) const { const Point p = p_ - coordinate_system_offset; const std_cxx11::array sp = GeometricUtilities::Coordinates::to_spherical(p); - return svalue(sp); + return svalue(sp, component); } @@ -179,11 +180,11 @@ namespace Functions template Tensor<1,dim> Spherical::gradient (const Point &p_, - const unsigned int) const + const unsigned int component) const { const Point p = p_ - coordinate_system_offset; const std_cxx11::array sp = GeometricUtilities::Coordinates::to_spherical(p); - const std_cxx11::array sg = sgradient(sp); + const std_cxx11::array sg = sgradient(sp, component); // somewhat backwards, but we need cos/sin's for unit vectors const double cos_theta = std::cos(sp[1]); @@ -225,12 +226,12 @@ namespace Functions template SymmetricTensor<2,dim> Spherical::hessian (const Point &p_, - const unsigned int /* component */) const + const unsigned int component) const { const Point p = p_ - coordinate_system_offset; const std_cxx11::array sp = GeometricUtilities::Coordinates::to_spherical(p); - const std_cxx11::array sg = sgradient(sp); - const std_cxx11::array sh = shessian(sp); + const std_cxx11::array sg = sgradient(sp, component); + const std_cxx11::array sh = shessian(sp, component); // somewhat backwards, but we need cos/sin's for unit vectors const double cos_theta = std::cos(sp[1]); @@ -309,7 +310,8 @@ namespace Functions template double - Spherical::svalue(const std_cxx11::array & /* sp */) const + Spherical::svalue(const std_cxx11::array & /* sp */, + const unsigned int /*component*/) const { AssertThrow(false, ExcNotImplemented()); @@ -320,7 +322,8 @@ namespace Functions template std_cxx11::array - Spherical::sgradient(const std_cxx11::array & /* sp */) const + Spherical::sgradient(const std_cxx11::array & /* sp */, + const unsigned int /*component*/) const { AssertThrow(false, ExcNotImplemented()); @@ -331,7 +334,8 @@ namespace Functions template std_cxx11::array - Spherical::shessian (const std_cxx11::array & /* sp */) const + Spherical::shessian (const std_cxx11::array & /* sp */, + const unsigned int /*component*/) const { AssertThrow(false, ExcNotImplemented()); diff --git a/tests/base/functions_spherical_01.cc b/tests/base/functions_spherical_01.cc index 5196a10091..8e4cf8768c 100644 --- a/tests/base/functions_spherical_01.cc +++ b/tests/base/functions_spherical_01.cc @@ -86,12 +86,14 @@ public: {} private: - virtual double svalue(const std_cxx11::array &sp) const + virtual double svalue(const std_cxx11::array &sp, + const unsigned int) const { return std::exp(-Z*sp[0]); } - virtual std_cxx11::array sgradient(const std_cxx11::array &sp) const + virtual std_cxx11::array sgradient(const std_cxx11::array &sp, + const unsigned int) const { std_cxx11::array res; res[0] = -Z*std::exp(-Z*sp[0]); @@ -100,7 +102,8 @@ private: return res; } - virtual std_cxx11::array shessian (const std_cxx11::array &sp) const + virtual std_cxx11::array shessian (const std_cxx11::array &sp, + const unsigned int) const { std_cxx11::array res; res[0] = Z*Z*std::exp(-Z*sp[0]); diff --git a/tests/base/functions_spherical_02.cc b/tests/base/functions_spherical_02.cc index f4a0a2a859..b5f01ba1c8 100644 --- a/tests/base/functions_spherical_02.cc +++ b/tests/base/functions_spherical_02.cc @@ -128,12 +128,14 @@ public: {} private: - virtual double svalue(const std_cxx11::array &sp) const + virtual double svalue(const std_cxx11::array &sp, + const unsigned int) const { return sp[0]*sp[0]*std::cos(sp[1])*std::sin(sp[2]); } - virtual std_cxx11::array sgradient(const std_cxx11::array &sp) const + virtual std_cxx11::array sgradient(const std_cxx11::array &sp, + const unsigned int) const { std_cxx11::array res; const double r = sp[0]; @@ -145,7 +147,8 @@ private: return res; } - virtual std_cxx11::array shessian (const std_cxx11::array &sp) const + virtual std_cxx11::array shessian (const std_cxx11::array &sp, + const unsigned int) const { std_cxx11::array res; const double r = sp[0];