Reference documentation for deal.II version GIT 9042b9283b 2023-12-02 14:50:02+00:00
\(\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\}}\)
geometric_utilities.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2016 - 2023 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 #include <array>
21 #include <cmath>
22 #include <limits>
23 
24 
26 
27 
28 namespace GeometricUtilities
29 {
30  namespace Coordinates
31  {
33  double,
34  << "The radius <" << arg1 << "> can not be negative.");
35 
37  double,
38  << "The azimuth angle <" << arg1 << "> is not in [0,2Pi).");
39 
41  double,
42  << "The polar angle <" << arg1 << "> is not in [0,Pi].");
43 
44 
45  template <int dim>
46  std::array<double, dim>
47  to_spherical(const Point<dim> &position)
48  {
49  std::array<double, dim> scoord{};
50  Assert(dim > 1, ExcNotImplemented());
51 
52  // radius
53  if constexpr (dim > 1)
54  {
55  scoord[0] = position.norm();
56  // azimuth angle \theta:
57  scoord[1] = std::atan2(position(1), position(0));
58  // correct to [0,2*pi)
59  if (scoord[1] < 0.0)
60  scoord[1] += 2.0 * numbers::PI;
61  }
62 
63  // polar angle \phi:
64  if constexpr (dim == 3)
65  {
66  // acos returns the angle in the range [0,\pi]
67  if (scoord[0] > std::numeric_limits<double>::min())
68  scoord[2] = std::acos(position(2) / scoord[0]);
69  else
70  scoord[2] = 0.0;
71  }
72  return scoord;
73  }
74 
75  template <std::size_t dim>
77  from_spherical(const std::array<double, dim> &scoord)
78  {
79  Point<dim> ccoord;
80 
81  Assert(scoord[0] >= 0., NegativeRadius(scoord[0]));
82 
83  Assert(scoord[1] >= 0. && scoord[1] < 2. * numbers::PI,
84  SphericalAzimuth(scoord[1]));
85 
86  switch (dim)
87  {
88  case 2:
89  {
90  ccoord[0] = scoord[0] * std::cos(scoord[1]);
91  ccoord[1] = scoord[0] * std::sin(scoord[1]);
92  break;
93  }
94  case 3:
95  {
96  Assert(scoord[2] >= 0. && scoord[2] <= numbers::PI,
97  SphericalPolar(scoord[2]));
98 
99  ccoord[0] = scoord[0] * std::sin(scoord[2]) * std::cos(scoord[1]);
100  ccoord[1] = scoord[0] * std::sin(scoord[2]) * std::sin(scoord[1]);
101  ccoord[2] = scoord[0] * std::cos(scoord[2]);
102  break;
103  }
104  default:
105  Assert(false, ExcNotImplemented());
106  break;
107  }
108 
109  return ccoord;
110  }
111 
112  // explicit instantiations
113 #include "geometric_utilities.inst"
114 
115  } // namespace Coordinates
116 } // namespace GeometricUtilities
117 
Definition: point.h:112
numbers::NumberTraits< Number >::real_type norm() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:477
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:478
static ::ExceptionBase & NegativeRadius(double arg1)
static ::ExceptionBase & SphericalAzimuth(double arg1)
static ::ExceptionBase & SphericalPolar(double arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1631
static ::ExceptionBase & ExcNotImplemented()
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:517
Expression atan2(const Expression &y, const Expression &x)
Expression acos(const Expression &x)
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:260