Reference documentation for deal.II version 9.4.1
\(\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\}}\)
Loading...
Searching...
No Matches
geometric_utilities.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
16
19
20
22
23
24namespace GeometricUtilities
25{
26 namespace Coordinates
27 {
29 double,
30 << "The radius <" << arg1 << "> can not be negative.");
31
33 double,
34 << "The azimuth angle <" << arg1 << "> is not in [0,2Pi).");
35
37 double,
38 << "The polar angle <" << arg1 << "> is not in [0,Pi].");
39
40
41 template <int dim>
42 std::array<double, dim>
43 to_spherical(const Point<dim> &position)
44 {
45 std::array<double, dim> scoord;
46
47 // radius
48 scoord[0] = position.norm();
49 // azimuth angle \theta:
50 scoord[1] = std::atan2(position(1), position(0));
51 // correct to [0,2*pi)
52 if (scoord[1] < 0.0)
53 scoord[1] += 2.0 * numbers::PI;
54
55 // polar angle \phi:
56 if (dim == 3)
57 {
58 // acos returns the angle in the range [0,\pi]
59 if (scoord[0] > std::numeric_limits<double>::min())
60 scoord[2] = std::acos(position(2) / scoord[0]);
61 else
62 scoord[2] = 0.0;
63 }
64 return scoord;
65 }
66
67 template <std::size_t dim>
69 from_spherical(const std::array<double, dim> &scoord)
70 {
71 Point<dim> ccoord;
72
73 Assert(scoord[0] >= 0., NegativeRadius(scoord[0]));
74
75 Assert(scoord[1] >= 0. && scoord[1] < 2. * numbers::PI,
76 SphericalAzimuth(scoord[1]));
77
78 switch (dim)
79 {
80 case 2:
81 {
82 ccoord[0] = scoord[0] * std::cos(scoord[1]);
83 ccoord[1] = scoord[0] * std::sin(scoord[1]);
84 break;
85 }
86 case 3:
87 {
88 Assert(scoord[2] >= 0. && scoord[2] <= numbers::PI,
89 SphericalPolar(scoord[2]));
90
91 ccoord[0] = scoord[0] * std::sin(scoord[2]) * std::cos(scoord[1]);
92 ccoord[1] = scoord[0] * std::sin(scoord[2]) * std::sin(scoord[1]);
93 ccoord[2] = scoord[0] * std::cos(scoord[2]);
94 break;
95 }
96 default:
97 Assert(false, ExcNotImplemented());
98 break;
99 }
100
101 return ccoord;
102 }
103
104 // explicit instantiations
105#include "geometric_utilities.inst"
106
107 } // namespace Coordinates
108} // namespace GeometricUtilities
109
Definition: point.h:111
numbers::NumberTraits< Number >::real_type norm() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
static ::ExceptionBase & SphericalAzimuth(double arg1)
static ::ExceptionBase & NegativeRadius(double arg1)
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
Definition: exceptions.h:1473
static ::ExceptionBase & SphericalPolar(double arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:509
Point< dim > from_spherical(const std::array< double, dim > &scoord)
std::array< double, dim > to_spherical(const Point< dim > &point)
static constexpr double PI
Definition: numbers.h:233
::VectorizedArray< Number, width > cos(const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > sin(const ::VectorizedArray< Number, width > &)