Reference documentation for deal.II version 9.0.0
function_spherical.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/base/point.h>
17 #include <deal.II/base/function_spherical.h>
18 #include <deal.II/base/geometric_utilities.h>
19 
20 #include <cmath>
21 #include <algorithm>
22 
23 DEAL_II_NAMESPACE_OPEN
24 namespace Functions
25 {
26 
27  // other implementations/notes:
28  // https://github.com/apache/commons-math/blob/master/src/main/java/org/apache/commons/math4/geometry/euclidean/threed/SphericalCoordinates.java
29  // http://mathworld.wolfram.com/SphericalCoordinates.html
30 
31  /*derivation of Hessian in Maxima as function of tensor products of unit vectors:
32 
33  depends(ur,[theta,phi]);
34  depends(utheta,theta);
35  depends(uphi,[theta,phi]);
36  depends(f,[r,theta,phi]);
37  declare([f,r,theta,phi], scalar)@f$
38  dotscrules: true;
39  grads(a):=ur.diff(a,r)+(1/r)*uphi.diff(a,phi)+(1/(r*sin(phi)))*utheta.diff(a,theta);
40 
41 
42  H : factor(grads(grads(f)));
43  H2: subst([diff(ur,theta)=sin(phi)*utheta,
44  diff(utheta,theta)=-cos(phi)*uphi-sin(phi)*ur,
45  diff(uphi,theta)=cos(phi)*utheta,
46  diff(ur,phi)=uphi,
47  diff(uphi,phi)=-ur],
48  H);
49  H3: trigsimp(fullratsimp(H2));
50 
51 
52  srules : [diff(f,r)=sg0,
53  diff(f,theta)=sg1,
54  diff(f,phi)=sg2,
55  diff(f,r,2)=sh0,
56  diff(f,theta,2)=sh1,
57  diff(f,phi,2)=sh2,
58  diff(f,r,1,theta,1)=sh3,
59  diff(f,r,1,phi,1)=sh4,
60  diff(f,theta,1,phi,1)=sh5,
61  cos(phi)=cos_phi,
62  cos(theta)=cos_theta,
63  sin(phi)=sin_phi,
64  sin(theta)=sin_theta
65  ]@f$
66 
67  c_utheta2 : distrib(subst(srules, ratcoeff(expand(H3), utheta.utheta)));
68  c_utheta_ur : (subst(srules, ratcoeff(expand(H3), utheta.ur)));
69  (subst(srules, ratcoeff(expand(H3), ur.utheta))) - c_utheta_ur;
70  c_utheta_uphi : (subst(srules, ratcoeff(expand(H3), utheta.uphi)));
71  (subst(srules, ratcoeff(expand(H3), uphi.utheta))) - c_utheta_uphi;
72  c_ur2 : (subst(srules, ratcoeff(expand(H3), ur.ur)));
73  c_ur_uphi : (subst(srules, ratcoeff(expand(H3), ur.uphi)));
74  (subst(srules, ratcoeff(expand(H3), uphi.ur))) - c_ur_uphi;
75  c_uphi2 : (subst(srules, ratcoeff(expand(H3), uphi.uphi)));
76 
77 
78  where (used later to do tensor products):
79 
80  ur : [cos(theta)*sin(phi), sin(theta)*sin(phi), cos(phi)];
81  utheta : [-sin(theta), cos(theta), 0];
82  uphi : [cos(theta)*cos(phi), sin(theta)*cos(phi), -sin(phi)];
83 
84  with the following proof of substitution rules above:
85 
86  -diff(ur,theta)+sin(phi)*utheta;
87  trigsimp(-diff(utheta,theta)-cos(phi)*uphi-sin(phi)*ur);
88  -diff(uphi,theta)+cos(phi)*utheta;
89  -diff(ur,phi)+uphi;
90  -diff(uphi,phi)-ur;
91 
92  */
93 
94  namespace
95  {
96 
100  template <int dim>
101  void set_unit_vectors(const double &cos_theta,
102  const double &sin_theta,
103  const double &cos_phi,
104  const double &sin_phi,
105  Tensor<1,dim> &unit_r,
106  Tensor<1,dim> &unit_theta,
107  Tensor<1,dim> &unit_phi)
108  {
109  unit_r[0] = cos_theta * sin_phi;
110  unit_r[1] = sin_theta * sin_phi;
111  unit_r[2] = cos_phi;
112 
113  unit_theta[0] = -sin_theta;
114  unit_theta[1] = cos_theta;
115  unit_theta[2] = 0.;
116 
117  unit_phi[0] = cos_theta * cos_phi;
118  unit_phi[1] = sin_theta * cos_phi;
119  unit_phi[2] = -sin_phi;
120  }
121 
122 
126  template <int dim>
127  void add_outer_product(SymmetricTensor<2,dim> &out,
128  const double &val,
129  const Tensor<1,dim> &in1,
130  const Tensor<1,dim> &in2)
131  {
132  if (val != 0.)
133  for (unsigned int i = 0; i < dim; i++)
134  for (unsigned int j = i; j < dim; j++)
135  out[i][j] += (in1[i]*in2[j]+in1[j]*in2[i])*val;
136  }
137 
141  template <int dim>
142  void add_outer_product(SymmetricTensor<2,dim> &out,
143  const double &val,
144  const Tensor<1,dim> &in)
145  {
146  if (val != 0.)
147  for (unsigned int i = 0; i < dim; i++)
148  for (unsigned int j = i; j < dim; j++)
149  out[i][j] += val*in[i]*in[j];
150  }
151  }
152 
153 
154 
155  template <int dim>
157  const unsigned int n_components)
158  :
159  Function<dim>(n_components),
160  coordinate_system_offset(p)
161  {
162  AssertThrow(dim==3,
164  }
165 
166 
167 
168  template <int dim>
169  double
171  const unsigned int component) const
172  {
173  const Point<dim> p = p_ - coordinate_system_offset;
174  const std::array<double, dim> sp = GeometricUtilities::Coordinates::to_spherical(p);
175  return svalue(sp, component);
176  }
177 
178 
179 
180  template <int dim>
183  const unsigned int /*component*/) const
184 
185  {
186  Assert(false, ExcNotImplemented());
187  return {};
188  }
189 
190 
191 
192  template <>
194  Spherical<3>::gradient (const Point<3> &p_,
195  const unsigned int component) const
196  {
197  constexpr int dim = 3;
198  const Point<dim> p = p_ - coordinate_system_offset;
199  const std::array<double, dim> sp = GeometricUtilities::Coordinates::to_spherical(p);
200  const std::array<double, dim> sg = sgradient(sp, component);
201 
202  // somewhat backwards, but we need cos/sin's for unit vectors
203  const double cos_theta = std::cos(sp[1]);
204  const double sin_theta = std::sin(sp[1]);
205  const double cos_phi = std::cos(sp[2]);
206  const double sin_phi = std::sin(sp[2]);
207 
208  Tensor<1,dim> unit_r, unit_theta, unit_phi;
209  set_unit_vectors(cos_theta, sin_theta,
210  cos_phi, sin_phi,
211  unit_r, unit_theta, unit_phi);
212 
213  Tensor<1,dim> res;
214 
215  if (sg[0] != 0.)
216  {
217  res += unit_r * sg[0];
218  }
219 
220  if (sg[1] * sin_phi != 0.)
221  {
222  Assert (sp[0] != 0.,
223  ExcDivideByZero());
224  res += unit_theta * sg[1] / (sp[0] * sin_phi);
225  }
226 
227  if (sg[2] != 0.)
228  {
229  Assert (sp[0] != 0.,
230  ExcDivideByZero());
231  res += unit_phi * sg[2] / sp[0];
232  }
233 
234  return res;
235  }
236 
237 
238 
239  template <int dim>
242  const unsigned int /*component*/) const
243  {
244  Assert(false, ExcNotImplemented());
245  return {};
246  }
247 
248 
249 
250  template <>
252  Spherical<3>::hessian (const Point<3> &p_,
253  const unsigned int component) const
254 
255  {
256  constexpr int dim = 3;
257  const Point<dim> p = p_ - coordinate_system_offset;
258  const std::array<double, dim> sp = GeometricUtilities::Coordinates::to_spherical(p);
259  const std::array<double, dim> sg = sgradient(sp, component);
260  const std::array<double, 6> sh = shessian(sp, component);
261 
262  // somewhat backwards, but we need cos/sin's for unit vectors
263  const double cos_theta = std::cos(sp[1]);
264  const double sin_theta = std::sin(sp[1]);
265  const double cos_phi = std::cos(sp[2]);
266  const double sin_phi = std::sin(sp[2]);
267  const double r = sp[0];
268 
269  Tensor<1,dim> unit_r, unit_theta, unit_phi;
270  set_unit_vectors(cos_theta, sin_theta,
271  cos_phi, sin_phi,
272  unit_r, unit_theta, unit_phi);
273 
274  const double sin_phi2 = sin_phi*sin_phi;
275  const double r2 = r*r;
276  Assert (r != 0.,
277  ExcDivideByZero());
278 
279  const double c_utheta2 = sg[0]/r +
280  ((sin_phi!= 0.) ?
281  (cos_phi*sg[2])/(r2*sin_phi)+sh[1]/(r2*sin_phi2) :
282  0.);
283  const double c_utheta_ur = ((sin_phi != 0.) ?
284  (r*sh[3]-sg[1])/(r2*sin_phi) :
285  0.);
286  const double c_utheta_uphi = ((sin_phi != 0.) ?
287  (sh[5]*sin_phi-cos_phi*sg[1])/(r2*sin_phi2) :
288  0.);
289  const double c_ur2 = sh[0];
290  const double c_ur_uphi = (r*sh[4]-sg[2])/r2;
291  const double c_uphi2 = (sh[2]+r*sg[0])/r2;
292 
293  // go through each tensor product
295 
296  add_outer_product(res,
297  c_utheta2,
298  unit_theta);
299 
300  add_outer_product(res,
301  c_utheta_ur,
302  unit_theta,
303  unit_r);
304 
305  add_outer_product(res,
306  c_utheta_uphi,
307  unit_theta,
308  unit_phi);
309 
310  add_outer_product(res,
311  c_ur2,
312  unit_r);
313 
314  add_outer_product(res,
315  c_ur_uphi,
316  unit_r,
317  unit_phi);
318 
319  add_outer_product(res,
320  c_uphi2,
321  unit_phi);
322 
323  return res;
324  }
325 
326 
327 
328  template <int dim>
329  std::size_t
330  Spherical<dim>::memory_consumption () const
331  {
332  return sizeof(Spherical<dim>);
333  }
334 
335 
336 
337  template <int dim>
338  double
339  Spherical<dim>::svalue(const std::array<double, dim> & /* sp */,
340  const unsigned int /*component*/) const
341  {
342  AssertThrow(false,
344  return 0.;
345  }
346 
347 
348 
349  template <int dim>
350  std::array<double, dim>
351  Spherical<dim>::sgradient(const std::array<double, dim> & /* sp */,
352  const unsigned int /*component*/) const
353  {
354  AssertThrow(false,
356  return std::array<double,dim>();
357  }
358 
359 
360 
361  template <int dim>
362  std::array<double, 6>
363  Spherical<dim>::shessian (const std::array<double, dim> & /* sp */,
364  const unsigned int /*component*/) const
365  {
366  AssertThrow(false,
368  return std::array<double, 6>();
369  }
370 
371 
372 
373 // explicit instantiations
374  template class Spherical<1>;
375  template class Spherical<2>;
376  template class Spherical<3>;
377 
378 }
379 
380 DEAL_II_NAMESPACE_CLOSE
virtual SymmetricTensor< 2, dim > hessian(const Point< dim > &p, const unsigned int component=0) const
virtual std::array< double, 6 > shessian(const std::array< double, dim > &sp, const unsigned int component) const
virtual double value(const Point< dim > &point, const unsigned int component=0) const
Spherical(const Point< dim > &center=Point< dim >(), const unsigned int n_components=1)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
static ::ExceptionBase & ExcDivideByZero()
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const
virtual std::array< double, dim > sgradient(const std::array< double, dim > &sp, const unsigned int component) const
virtual double svalue(const std::array< double, dim > &sp, const unsigned int component) const
std::array< double, dim > to_spherical(const Point< dim > &point)
#define Assert(cond, exc)
Definition: exceptions.h:1142
static ::ExceptionBase & ExcNotImplemented()