numerics
Loading...
Searching...
No Matches
math.hpp
Go to the documentation of this file.
1/// @file math.hpp
2/// @brief Thin wrappers around C++17 `<cmath>` and `<numeric>` with readable names.
3///
4/// The standard library names for special functions (cyl_bessel_j, comp_ellint_1,
5/// cyl_neumann, etc.) are accurate but opaque. This header re-exports them under
6/// the names used in textbooks and mathematical literature.
7///
8/// Random number generation wraps the mt19937 boilerplate into a plain struct.
9///
10/// Include this header instead of pulling `<cmath>` special functions directly.
11#pragma once
12
13#include "core/types.hpp"
14#include <cmath>
15#include <numeric>
16#include <random>
17#include <vector>
18#include <cassert>
19
20// Apple libc++ does not implement C++17 special math functions.
21// Use Boost.Math as a drop-in replacement when building on macOS.
22#ifdef NUMERICS_USE_BOOST_MATH
23# include <boost/math/special_functions/bessel.hpp>
24# include <boost/math/special_functions/legendre.hpp>
25# include <boost/math/special_functions/hermite.hpp>
26# include <boost/math/special_functions/laguerre.hpp>
27# include <boost/math/special_functions/ellint_1.hpp>
28# include <boost/math/special_functions/ellint_2.hpp>
29# include <boost/math/special_functions/ellint_3.hpp>
30# include <boost/math/special_functions/expint.hpp>
31# include <boost/math/special_functions/zeta.hpp>
32# include <boost/math/special_functions/beta.hpp>
33# include <boost/math/special_functions/spherical_harmonic.hpp>
34#endif
35
36namespace num {
37
38// Mathematical constants (C++20 <numbers> not available in C++17)
39
40constexpr real pi = 3.14159265358979323846;
41constexpr real e = 2.71828182845904523536;
42constexpr real phi = 1.61803398874989484820; ///< Golden ratio
43constexpr real sqrt2 = 1.41421356237309504880;
44constexpr real sqrt3 = 1.73205080756887729353;
45constexpr real ln2 = 0.69314718055994530942;
46constexpr real inv_pi = 0.31830988618379067154; ///< 1/pi
47constexpr real two_pi = 6.28318530717958647692; ///< 2pi
48constexpr real half_pi = 1.57079632679489661923; ///< pi/2
49
50// Cylindrical Bessel functions (C++17, <cmath>)
51
52/// @brief J_nu(x) -- Bessel function of the first kind
53inline real bessel_j(real nu, real x) {
54#ifdef NUMERICS_USE_BOOST_MATH
55 return boost::math::cyl_bessel_j(nu, x);
56#else
57 return std::cyl_bessel_j(nu, x);
58#endif
59}
60
61/// @brief Y_nu(x) -- Bessel function of the second kind (Neumann function)
62inline real bessel_y(real nu, real x) {
63#ifdef NUMERICS_USE_BOOST_MATH
64 return boost::math::cyl_neumann(nu, x);
65#else
66 return std::cyl_neumann(nu, x);
67#endif
68}
69
70/// @brief I_nu(x) -- modified Bessel function of the first kind
71inline real bessel_i(real nu, real x) {
72#ifdef NUMERICS_USE_BOOST_MATH
73 return boost::math::cyl_bessel_i(nu, x);
74#else
75 return std::cyl_bessel_i(nu, x);
76#endif
77}
78
79/// @brief K_nu(x) -- modified Bessel function of the second kind
80inline real bessel_k(real nu, real x) {
81#ifdef NUMERICS_USE_BOOST_MATH
82 return boost::math::cyl_bessel_k(nu, x);
83#else
84 return std::cyl_bessel_k(nu, x);
85#endif
86}
87
88// Spherical Bessel functions (C++17, <cmath>)
89
90/// @brief j_n(x) -- spherical Bessel function of the first kind
91inline real sph_bessel_j(unsigned int n, real x) {
92#ifdef NUMERICS_USE_BOOST_MATH
93 return boost::math::sph_bessel(n, x);
94#else
95 return std::sph_bessel(n, x);
96#endif
97}
98
99/// @brief y_n(x) -- spherical Neumann function (spherical Bessel of the second kind)
100inline real sph_bessel_y(unsigned int n, real x) {
101#ifdef NUMERICS_USE_BOOST_MATH
102 return boost::math::sph_neumann(n, x);
103#else
104 return std::sph_neumann(n, x);
105#endif
106}
107
108// Orthogonal polynomials (C++17, <cmath>)
109
110/// @brief P_n(x) -- Legendre polynomial of degree n
111inline real legendre(unsigned int n, real x) {
112#ifdef NUMERICS_USE_BOOST_MATH
113 return boost::math::legendre_p(static_cast<int>(n), x);
114#else
115 return std::legendre(n, x);
116#endif
117}
118
119/// @brief P_n^m(x) -- associated Legendre polynomial
120inline real assoc_legendre(unsigned int n, unsigned int m, real x) {
121#ifdef NUMERICS_USE_BOOST_MATH
122 return boost::math::legendre_p(static_cast<int>(n), static_cast<int>(m), x);
123#else
124 return std::assoc_legendre(n, m, x);
125#endif
126}
127
128/// @brief Y_l^m(theta) -- spherical harmonic (real part, theta in radians)
129inline real sph_legendre(unsigned int l, unsigned int m, real theta) {
130#ifdef NUMERICS_USE_BOOST_MATH
131 return boost::math::spherical_harmonic_r(l, static_cast<int>(m), theta, real(0));
132#else
133 return std::sph_legendre(l, m, theta);
134#endif
135}
136
137/// @brief H_n(x) -- (physicists') Hermite polynomial
138inline real hermite(unsigned int n, real x) {
139#ifdef NUMERICS_USE_BOOST_MATH
140 return boost::math::hermite(n, x);
141#else
142 return std::hermite(n, x);
143#endif
144}
145
146/// @brief L_n(x) -- Laguerre polynomial
147inline real laguerre(unsigned int n, real x) {
148#ifdef NUMERICS_USE_BOOST_MATH
149 return boost::math::laguerre(n, x);
150#else
151 return std::laguerre(n, x);
152#endif
153}
154
155/// @brief L_n^m(x) -- associated Laguerre polynomial
156inline real assoc_laguerre(unsigned int n, unsigned int m, real x) {
157#ifdef NUMERICS_USE_BOOST_MATH
158 return boost::math::laguerre(n, m, x);
159#else
160 return std::assoc_laguerre(n, m, x);
161#endif
162}
163
164// Elliptic integrals (C++17, <cmath>)
165// Names: K(k), E(k), Pi(n,k) for complete; F(k,phi), E(k,phi), Pi(n,k,phi) for incomplete
166
167/// @brief K(k) -- complete elliptic integral of the first kind
169#ifdef NUMERICS_USE_BOOST_MATH
170 return boost::math::ellint_1(k);
171#else
172 return std::comp_ellint_1(k);
173#endif
174}
175
176/// @brief E(k) -- complete elliptic integral of the second kind
178#ifdef NUMERICS_USE_BOOST_MATH
179 return boost::math::ellint_2(k);
180#else
181 return std::comp_ellint_2(k);
182#endif
183}
184
185/// @brief Pi(n, k) -- complete elliptic integral of the third kind
187#ifdef NUMERICS_USE_BOOST_MATH
188 return boost::math::ellint_3(k, n); // boost: (k, n); std: (n, k)
189#else
190 return std::comp_ellint_3(n, k);
191#endif
192}
193
194/// @brief F(k, phi) -- incomplete elliptic integral of the first kind
196#ifdef NUMERICS_USE_BOOST_MATH
197 return boost::math::ellint_1(k, phi);
198#else
199 return std::ellint_1(k, phi);
200#endif
201}
202
203/// @brief E(k, phi) -- incomplete elliptic integral of the second kind
205#ifdef NUMERICS_USE_BOOST_MATH
206 return boost::math::ellint_2(k, phi);
207#else
208 return std::ellint_2(k, phi);
209#endif
210}
211
212/// @brief Pi(n, k, phi) -- incomplete elliptic integral of the third kind
214#ifdef NUMERICS_USE_BOOST_MATH
215 return boost::math::ellint_3(k, n, phi); // boost: (k, n, phi); std: (n, k, phi)
216#else
217 return std::ellint_3(n, k, phi);
218#endif
219}
220
221// Other special functions (C++17, <cmath>)
222
223/// @brief Ei(x) -- exponential integral
224inline real expint(real x) {
225#ifdef NUMERICS_USE_BOOST_MATH
226 return boost::math::expint(x);
227#else
228 return std::expint(x);
229#endif
230}
231
232/// @brief zeta(x) -- Riemann zeta function
233inline real zeta(real x) {
234#ifdef NUMERICS_USE_BOOST_MATH
235 return boost::math::zeta(x);
236#else
237 return std::riemann_zeta(x);
238#endif
239}
240
241/// @brief B(a, b) -- beta function
242inline real beta(real a, real b) {
243#ifdef NUMERICS_USE_BOOST_MATH
244 return boost::math::beta(a, b);
245#else
246 return std::beta(a, b);
247#endif
248}
249
250// Sequence utilities (wrapping <numeric>)
251
252/// @brief Evenly spaced values from start to stop, inclusive. MATLAB/NumPy linspace.
253inline std::vector<real> linspace(real start, real stop, idx n) {
254 assert(n >= 2);
255 std::vector<real> out(n);
256 real step = (stop - start) / static_cast<real>(n - 1);
257 for (idx i = 0; i < n; ++i)
258 out[i] = start + static_cast<real>(i) * step;
259 return out;
260}
261
262/// @brief Integer sequence [start, start+1, ..., start+n-1]. Wraps std::iota.
263inline std::vector<int> int_range(int start, int n) {
264 assert(n >= 0);
265 std::vector<int> out(static_cast<idx>(n));
266 std::iota(out.begin(), out.end(), start);
267 return out;
268}
269
270// Random number generation (wrapping the mt19937 boilerplate)
271
272/// @brief Seeded pseudo-random number generator (Mersenne Twister).
273/// Pass a pointer to rng_* functions to draw samples.
274struct Rng {
275 std::mt19937 engine;
276
277 explicit Rng(uint32_t seed) : engine(seed) {}
278
279 /// Seed from hardware entropy.
281};
282
283/// @brief Uniform real in [lo, hi).
284inline real rng_uniform(Rng* r, real lo, real hi) {
285 return std::uniform_real_distribution<real>{lo, hi}(r->engine);
286}
287
288/// @brief Normal (Gaussian) sample with given mean and standard deviation.
289inline real rng_normal(Rng* r, real mean, real stddev) {
290 return std::normal_distribution<real>{mean, stddev}(r->engine);
291}
292
293/// @brief Uniform integer in [lo, hi] (inclusive on both ends).
294inline int rng_int(Rng* r, int lo, int hi) {
295 return std::uniform_int_distribution<int>{lo, hi}(r->engine);
296}
297
298} // namespace num
Core type definitions.
double real
Definition types.hpp:10
constexpr real phi
Golden ratio.
Definition math.hpp:42
constexpr real inv_pi
1/pi
Definition math.hpp:46
constexpr real sqrt3
Definition math.hpp:44
real ellint_F(real k, real phi)
F(k, phi) – incomplete elliptic integral of the first kind.
Definition math.hpp:195
std::vector< real > linspace(real start, real stop, idx n)
Evenly spaced values from start to stop, inclusive. MATLAB/NumPy linspace.
Definition math.hpp:253
real bessel_j(real nu, real x)
J_nu(x) – Bessel function of the first kind.
Definition math.hpp:53
real expint(real x)
Ei(x) – exponential integral.
Definition math.hpp:224
int rng_int(Rng *r, int lo, int hi)
Uniform integer in [lo, hi] (inclusive on both ends).
Definition math.hpp:294
real sph_legendre(unsigned int l, unsigned int m, real theta)
Y_l^m(theta) – spherical harmonic (real part, theta in radians)
Definition math.hpp:129
real ellint_E(real k)
E(k) – complete elliptic integral of the second kind.
Definition math.hpp:177
real rng_uniform(Rng *r, real lo, real hi)
Uniform real in [lo, hi).
Definition math.hpp:284
real beta(real a, real b)
B(a, b) – beta function.
Definition math.hpp:242
constexpr T ipow(T x) noexcept
Compute x^N at compile time via repeated squaring.
std::size_t idx
Definition types.hpp:11
real legendre(unsigned int n, real x)
P_n(x) – Legendre polynomial of degree n.
Definition math.hpp:111
real rng_normal(Rng *r, real mean, real stddev)
Normal (Gaussian) sample with given mean and standard deviation.
Definition math.hpp:289
constexpr real pi
Definition math.hpp:40
real assoc_legendre(unsigned int n, unsigned int m, real x)
P_n^m(x) – associated Legendre polynomial.
Definition math.hpp:120
real bessel_k(real nu, real x)
K_nu(x) – modified Bessel function of the second kind.
Definition math.hpp:80
constexpr real ln2
Definition math.hpp:45
constexpr real e
Definition math.hpp:41
constexpr real two_pi
2pi
Definition math.hpp:47
real ellint_K(real k)
K(k) – complete elliptic integral of the first kind.
Definition math.hpp:168
std::vector< int > int_range(int start, int n)
Integer sequence [start, start+1, ..., start+n-1]. Wraps std::iota.
Definition math.hpp:263
real zeta(real x)
zeta(x) – Riemann zeta function
Definition math.hpp:233
real laguerre(unsigned int n, real x)
L_n(x) – Laguerre polynomial.
Definition math.hpp:147
constexpr real half_pi
pi/2
Definition math.hpp:48
real hermite(unsigned int n, real x)
H_n(x) – (physicists') Hermite polynomial.
Definition math.hpp:138
real bessel_i(real nu, real x)
I_nu(x) – modified Bessel function of the first kind.
Definition math.hpp:71
real bessel_y(real nu, real x)
Y_nu(x) – Bessel function of the second kind (Neumann function)
Definition math.hpp:62
real assoc_laguerre(unsigned int n, unsigned int m, real x)
L_n^m(x) – associated Laguerre polynomial.
Definition math.hpp:156
real ellint_Pi(real n, real k)
Pi(n, k) – complete elliptic integral of the third kind.
Definition math.hpp:186
real sph_bessel_j(unsigned int n, real x)
j_n(x) – spherical Bessel function of the first kind
Definition math.hpp:91
real sph_bessel_y(unsigned int n, real x)
y_n(x) – spherical Neumann function (spherical Bessel of the second kind)
Definition math.hpp:100
real ellint_Pi_inc(real n, real k, real phi)
Pi(n, k, phi) – incomplete elliptic integral of the third kind.
Definition math.hpp:213
real ellint_Ei(real k, real phi)
E(k, phi) – incomplete elliptic integral of the second kind.
Definition math.hpp:204
constexpr real sqrt2
Definition math.hpp:43
Seeded pseudo-random number generator (Mersenne Twister). Pass a pointer to rng_* functions to draw s...
Definition math.hpp:274
std::mt19937 engine
Definition math.hpp:275
Rng()
Seed from hardware entropy.
Definition math.hpp:280
Rng(uint32_t seed)
Definition math.hpp:277