Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
function_spherical.cc
Go to the documentation of this file.
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.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
18#include <deal.II/base/point.h>
19
20#include <algorithm>
21#include <cmath>
22
24namespace Functions
25{
26 // other implementations/notes:
27 // https://github.com/apache/commons-math/blob/master/src/main/java/org/apache/commons/math4/geometry/euclidean/threed/SphericalCoordinates.java
28 // http://mathworld.wolfram.com/SphericalCoordinates.html
29
30 /*derivation of Hessian in Maxima as function of tensor products of unit
31 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 {
99 template <int dim>
100 void
101 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 } // namespace
152
153
154
155 template <int dim>
157 const unsigned int n_components)
158 : Function<dim>(n_components)
159 , coordinate_system_offset(p)
160 {
161 AssertThrow(dim == 3, ExcNotImplemented());
162 }
163
164
165
166 template <int dim>
167 double
169 const unsigned int component) const
170 {
171 const Point<dim> p = p_ - coordinate_system_offset;
172 const std::array<double, dim> sp =
174 return svalue(sp, component);
175 }
176
177
178
179 template <int dim>
182 const unsigned int /*component*/) const
183
184 {
185 Assert(false, ExcNotImplemented());
186 return {};
187 }
188
189
190
191 template <>
193 Spherical<3>::gradient(const Point<3> &p_, const unsigned int component) const
194 {
195 constexpr int dim = 3;
196 const Point<dim> p = p_ - coordinate_system_offset;
197 const std::array<double, dim> sp =
199 const std::array<double, dim> sg = sgradient(sp, component);
200
201 // somewhat backwards, but we need cos/sin's for unit vectors
202 const double cos_theta = std::cos(sp[1]);
203 const double sin_theta = std::sin(sp[1]);
204 const double cos_phi = std::cos(sp[2]);
205 const double sin_phi = std::sin(sp[2]);
206
207 Tensor<1, dim> unit_r, unit_theta, unit_phi;
208 set_unit_vectors(
209 cos_theta, sin_theta, cos_phi, sin_phi, unit_r, unit_theta, unit_phi);
210
211 Tensor<1, dim> res;
212
213 if (sg[0] != 0.)
214 {
215 res += unit_r * sg[0];
216 }
217
218 if (sg[1] * sin_phi != 0.)
219 {
220 Assert(sp[0] != 0., ExcDivideByZero());
221 res += unit_theta * sg[1] / (sp[0] * sin_phi);
222 }
223
224 if (sg[2] != 0.)
225 {
226 Assert(sp[0] != 0., ExcDivideByZero());
227 res += unit_phi * sg[2] / sp[0];
228 }
229
230 return res;
231 }
232
233
234
235 template <int dim>
238 const unsigned int /*component*/) const
239 {
240 Assert(false, ExcNotImplemented());
241 return {};
242 }
243
244
245
246 template <>
248 Spherical<3>::hessian(const Point<3> &p_, const unsigned int component) const
249
250 {
251 constexpr int dim = 3;
252 const Point<dim> p = p_ - coordinate_system_offset;
253 const std::array<double, dim> sp =
255 const std::array<double, dim> sg = sgradient(sp, component);
256 const std::array<double, 6> sh = shessian(sp, component);
257
258 // somewhat backwards, but we need cos/sin's for unit vectors
259 const double cos_theta = std::cos(sp[1]);
260 const double sin_theta = std::sin(sp[1]);
261 const double cos_phi = std::cos(sp[2]);
262 const double sin_phi = std::sin(sp[2]);
263 const double r = sp[0];
264
265 Tensor<1, dim> unit_r, unit_theta, unit_phi;
266 set_unit_vectors(
267 cos_theta, sin_theta, cos_phi, sin_phi, unit_r, unit_theta, unit_phi);
268
269 const double sin_phi2 = sin_phi * sin_phi;
270 const double r2 = r * r;
271 Assert(r != 0., ExcDivideByZero());
272
273 const double c_utheta2 =
274 sg[0] / r + ((sin_phi != 0.) ? (cos_phi * sg[2]) / (r2 * sin_phi) +
275 sh[1] / (r2 * sin_phi2) :
276 0.);
277 const double c_utheta_ur =
278 ((sin_phi != 0.) ? (r * sh[3] - sg[1]) / (r2 * sin_phi) : 0.);
279 const double c_utheta_uphi =
280 ((sin_phi != 0.) ? (sh[5] * sin_phi - cos_phi * sg[1]) / (r2 * sin_phi2) :
281 0.);
282 const double c_ur2 = sh[0];
283 const double c_ur_uphi = (r * sh[4] - sg[2]) / r2;
284 const double c_uphi2 = (sh[2] + r * sg[0]) / r2;
285
286 // go through each tensor product
288
289 add_outer_product(res, c_utheta2, unit_theta);
290
291 add_outer_product(res, c_utheta_ur, unit_theta, unit_r);
292
293 add_outer_product(res, c_utheta_uphi, unit_theta, unit_phi);
294
295 add_outer_product(res, c_ur2, unit_r);
296
297 add_outer_product(res, c_ur_uphi, unit_r, unit_phi);
298
299 add_outer_product(res, c_uphi2, unit_phi);
300
301 return res;
302 }
303
304
305
306 template <int dim>
307 std::size_t
309 {
310 return sizeof(Spherical<dim>);
311 }
312
313
314
315 template <int dim>
316 double
317 Spherical<dim>::svalue(const std::array<double, dim> & /* sp */,
318 const unsigned int /*component*/) const
319 {
321 return 0.;
322 }
323
324
325
326 template <int dim>
327 std::array<double, dim>
328 Spherical<dim>::sgradient(const std::array<double, dim> & /* sp */,
329 const unsigned int /*component*/) const
330 {
332 return std::array<double, dim>();
333 }
334
335
336
337 template <int dim>
338 std::array<double, 6>
339 Spherical<dim>::shessian(const std::array<double, dim> & /* sp */,
340 const unsigned int /*component*/) const
341 {
343 return std::array<double, 6>();
344 }
345
346
347
348 // explicit instantiations
349 template class Spherical<1>;
350 template class Spherical<2>;
351 template class Spherical<3>;
352
353} // namespace Functions
354
virtual SymmetricTensor< 2, dim > hessian(const Point< dim > &p, const unsigned int component=0) const override
Spherical(const Point< dim > &center=Point< dim >(), const unsigned int n_components=1)
virtual Tensor< 1, dim > gradient(const Point< dim > &p, const unsigned int component=0) const override
virtual std::array< double, dim > sgradient(const std::array< double, dim > &sp, const unsigned int component) const
virtual double value(const Point< dim > &point, const unsigned int component=0) const override
virtual std::size_t memory_consumption() const override
virtual std::array< double, 6 > shessian(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
Definition: point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcDivideByZero()
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
std::array< double, dim > to_spherical(const Point< dim > &point)
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)