]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Functions::Spherical add component argument to svalue/sgradient/shessian 3921/head
authorDenis Davydov <davydden@gmail.com>
Fri, 3 Feb 2017 15:52:19 +0000 (16:52 +0100)
committerDenis Davydov <davydden@gmail.com>
Fri, 3 Feb 2017 15:52:19 +0000 (16:52 +0100)
include/deal.II/base/function_spherical.h
source/base/function_spherical.cc
tests/base/functions_spherical_01.cc
tests/base/functions_spherical_02.cc

index fefac59e01e8d6164610cf0f501c82326b193ed0..a2de2db170296220f7f9cf5e7826e3bb5f4a6be6 100644 (file)
@@ -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<dim> &center = Point<dim>());
+    Spherical(const Point<dim> &center = Point<dim>(),
+              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<double, dim> &sp) const;
+    virtual double svalue(const std_cxx11::array<double, dim> &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<double, dim> sgradient(const std_cxx11::array<double, dim> &sp) const;
+    virtual std_cxx11::array<double, dim> sgradient(const std_cxx11::array<double, dim> &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<double, 6> shessian (const std_cxx11::array<double, dim> &sp) const;
+    virtual std_cxx11::array<double, 6> shessian (const std_cxx11::array<double, dim> &sp,
+                                                  const unsigned int component) const;
 
     /**
      * A vector from the origin to the center of spherical coordinate system.
index 06b31fe6d3330593eddeaf01c35467ae269931a6..be08727ef22935f0310584e0835713870aa7bfbe 100644 (file)
@@ -153,9 +153,10 @@ namespace Functions
 
 
   template <int dim>
-  Spherical<dim>::Spherical(const Point<dim> &p)
+  Spherical<dim>::Spherical(const Point<dim> &p,
+                            const unsigned int n_components)
     :
-    Function<dim>(1),
+    Function<dim>(n_components),
     coordinate_system_offset(p)
   {
     AssertThrow(dim==3,
@@ -167,11 +168,11 @@ namespace Functions
   template <int dim>
   double
   Spherical<dim>::value (const Point<dim>   &p_,
-                         const unsigned int) const
+                         const unsigned int component) const
   {
     const Point<dim> p = p_ - coordinate_system_offset;
     const std_cxx11::array<double, dim> sp = GeometricUtilities::Coordinates::to_spherical(p);
-    return svalue(sp);
+    return svalue(sp, component);
   }
 
 
@@ -179,11 +180,11 @@ namespace Functions
   template <int dim>
   Tensor<1,dim>
   Spherical<dim>::gradient (const Point<dim>   &p_,
-                            const unsigned int) const
+                            const unsigned int component) const
   {
     const Point<dim> p = p_ - coordinate_system_offset;
     const std_cxx11::array<double, dim> sp = GeometricUtilities::Coordinates::to_spherical(p);
-    const std_cxx11::array<double, dim> sg = sgradient(sp);
+    const std_cxx11::array<double, dim> 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 <int dim>
   SymmetricTensor<2,dim>
   Spherical<dim>::hessian (const Point<dim> &p_,
-                           const unsigned int /* component */) const
+                           const unsigned int component) const
   {
     const Point<dim> p = p_ - coordinate_system_offset;
     const std_cxx11::array<double, dim> sp = GeometricUtilities::Coordinates::to_spherical(p);
-    const std_cxx11::array<double, dim> sg = sgradient(sp);
-    const std_cxx11::array<double, 6>   sh = shessian(sp);
+    const std_cxx11::array<double, dim> sg = sgradient(sp, component);
+    const std_cxx11::array<double, 6>   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 <int dim>
   double
-  Spherical<dim>::svalue(const std_cxx11::array<double, dim> & /* sp */) const
+  Spherical<dim>::svalue(const std_cxx11::array<double, dim> & /* sp */,
+                         const unsigned int /*component*/) const
   {
     AssertThrow(false,
                 ExcNotImplemented());
@@ -320,7 +322,8 @@ namespace Functions
 
   template <int dim>
   std_cxx11::array<double, dim>
-  Spherical<dim>::sgradient(const std_cxx11::array<double, dim> & /* sp */) const
+  Spherical<dim>::sgradient(const std_cxx11::array<double, dim> & /* sp */,
+                            const unsigned int /*component*/) const
   {
     AssertThrow(false,
                 ExcNotImplemented());
@@ -331,7 +334,8 @@ namespace Functions
 
   template <int dim>
   std_cxx11::array<double, 6>
-  Spherical<dim>::shessian (const std_cxx11::array<double, dim> & /* sp */) const
+  Spherical<dim>::shessian (const std_cxx11::array<double, dim> & /* sp */,
+                            const unsigned int /*component*/) const
   {
     AssertThrow(false,
                 ExcNotImplemented());
index 5196a10091b026ebaca5aca9d45b686ade536803..8e4cf8768cc63fee43be8e8406cb772ddccb7714 100644 (file)
@@ -86,12 +86,14 @@ public:
   {}
 
 private:
-  virtual double svalue(const std_cxx11::array<double, dim> &sp) const
+  virtual double svalue(const std_cxx11::array<double, dim> &sp,
+                        const unsigned int) const
   {
     return std::exp(-Z*sp[0]);
   }
 
-  virtual std_cxx11::array<double, dim> sgradient(const std_cxx11::array<double, dim> &sp) const
+  virtual std_cxx11::array<double, dim> sgradient(const std_cxx11::array<double, dim> &sp,
+                                                  const unsigned int) const
   {
     std_cxx11::array<double, dim> res;
     res[0] = -Z*std::exp(-Z*sp[0]);
@@ -100,7 +102,8 @@ private:
     return res;
   }
 
-  virtual std_cxx11::array<double, 6> shessian (const std_cxx11::array<double, dim> &sp) const
+  virtual std_cxx11::array<double, 6> shessian (const std_cxx11::array<double, dim> &sp,
+                                                const unsigned int) const
   {
     std_cxx11::array<double, 6> res;
     res[0] = Z*Z*std::exp(-Z*sp[0]);
index f4a0a2a8597f578dabad893058d2db54dd2cf61f..b5f01ba1c8e526931033af943160f61724226811 100644 (file)
@@ -128,12 +128,14 @@ public:
   {}
 
 private:
-  virtual double svalue(const std_cxx11::array<double, dim> &sp) const
+  virtual double svalue(const std_cxx11::array<double, dim> &sp,
+                        const unsigned int) const
   {
     return sp[0]*sp[0]*std::cos(sp[1])*std::sin(sp[2]);
   }
 
-  virtual std_cxx11::array<double, dim> sgradient(const std_cxx11::array<double, dim> &sp) const
+  virtual std_cxx11::array<double, dim> sgradient(const std_cxx11::array<double, dim> &sp,
+                                                  const unsigned int) const
   {
     std_cxx11::array<double, dim> res;
     const double r     = sp[0];
@@ -145,7 +147,8 @@ private:
     return res;
   }
 
-  virtual std_cxx11::array<double, 6> shessian (const std_cxx11::array<double, dim> &sp) const
+  virtual std_cxx11::array<double, 6> shessian (const std_cxx11::array<double, dim> &sp,
+                                                const unsigned int) const
   {
     std_cxx11::array<double, 6> res;
     const double r = sp[0];

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.