2 Copyright (c) 2012 John Maddock
3 Use, modification and distribution are subject to the
4 Boost Software License, Version 1.0. (See accompanying file
5 LICENSE_1_0.txt or copy at http://www.boost.org/LICENSE_1_0.txt)
8 [section:jacobi Jacobi Elliptic Functions]
10 [section:jac_over Overvew of the Jacobi Elliptic Functions]
12 There are twelve Jacobi Elliptic functions, of which the three copolar functions ['sn], ['cn] and ['dn] are the most important
13 as the other nine can be computed from these three
14 [footnote [@http://en.wikipedia.org/wiki/Jacobi_elliptic_functions Wikipedia: Jacobi elliptic functions]]
15 [footnote [@http://mathworld.wolfram.com/JacobiEllipticFunctions.html Weisstein, Eric W. "Jacobi Elliptic Functions." From MathWorld - A Wolfram Web Resource.]]
16 [footnote [@http://dlmf.nist.gov/22 Digital Library of Mathematical Functions: Jacobian Elliptic Functions]].
18 These functions each take two arguments: a parameter, and a variable as described below.
20 Like all elliptic functions these can be parameterised in a number of ways:
22 * In terms of a parameter ['m].
23 * In terms of the elliptic modulus ['k] where ['m = k[super 2]].
24 * In terms of the modular angle [alpha], where ['m = sin[super 2][alpha]].
26 In our implementation, these functions all take the elliptic modulus /k/ as the parameter.
28 In addition the variable /u/ is sometimes expressed as an amplitude [phi], in our implementation we always use /u/.
30 Finally note that our functions all take the elliptic modulus as the first argument - this is for alignment with
31 the Elliptic Integrals.
33 There are twenve functions for computing the twelve individual Jacobi elliptic functions: __jacobi_cd, __jacobi_cn, __jacobi_cs,
34 __jacobi_dc, __jacobi_dn, __jacobi_ds, __jacobi_nc, __jacobi_nd, __jacobi_ns, __jacobi_sc, __jacobi_sd and __jacobi_sn.
36 They are all called as for example:
40 Note however that these individual functions are all really thin wrappers around the function __jacobi_elliptic which calculates
41 the three copolar functions ['sn], ['cn] and ['dn] in a single function call. Thus if you need more than one of these functions
42 for a given set of arguments, it's most efficient to use __jacobi_elliptic.
46 [section:jacobi_elliptic Jacobi Elliptic SN, CN and DN]
51 #include <boost/math/special_functions/jacobi_elliptic.hpp>
54 namespace boost { namespace math {
56 template <class T, class U, class V>
57 ``__sf_result`` jacobi_elliptic(T k, U u, V* pcn, V* pdn);
59 template <class T, class U, class V, class Policy>
60 ``__sf_result`` jacobi_elliptic(T k, U u, V* pcn, V* pdn, const Policy&);
66 The function __jacobi_elliptic calculates the three copolar Jacobi elliptic functions
67 ['sn(u, k)], ['cn(u, k)] and ['dn(u, k)]. The returned value is
68 ['sn(u, k)], and if provided, ['*pcn] is
69 set to ['cn(u, k)], and ['*pdn] is set to ['dn(u, k)].
71 The functions are defined as follows, given:
75 The the angle [phi] is called the ['amplitude] and:
79 [note ['[phi]] is called the amplitude.
81 ['k] is called the modulus.
85 [caution Rather like other elliptic functions, the Jacobi functions
86 are expressed in a variety of different ways. In particular,
87 the parameter /k/ (the modulus) may also be expressed using a modular
88 angle [alpha], or a parameter /m/. These are related by:
92 m = k[super 2] = sin[super 2][alpha]
94 So that the function ['sn] (for example) may be expressed as
103 To further complicate matters, some texts refer to the ['complement
104 of the parameter m], or 1 - m, where:
106 1 - m = 1 - k[super 2] = cos[super 2][alpha]
108 This implementation uses /k/ throughout, and makes this the first argument
109 to the functions: this is for alignment with the elliptic integrals which match the requirements
110 of the [@http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2005/n1836.pdf
111 Technical Report on C++ Library Extensions]. However, you should
112 be extra careful when using these functions!]
116 The following graphs illustrate how these functions change as /k/ changes: for small /k/
117 these are sine waves, while as /k/ tends to 1 they become hyperbolic functions:
127 These functions are computed using only basic arithmetic operations and trigomometric functions, so
128 there isn't much variation in accuracy over differing platforms.
129 Typically errors are trivially small for small angles, and as is typical for cyclic
130 functions, grow as the angle increases.
131 Note that only results for the widest floating point type on the
132 system are given as narrower types have __zero_error. All values
133 are relative errors in units of epsilon.
143 The tests use a mixture of spot test values calculated using the online
144 calculator at [@http://functions.wolfram.com/ functions.wolfram.com],
145 and random test data generated using
146 MPFR at 1000-bit precision and this implementation.
148 [heading Implementation]
150 For ['k > 1] we apply the relations:
154 Then filter off the special cases:
156 ['sn(0, k) = 0] and ['cn(0, k) = dn(0, k) = 1].
158 ['sn(u, 0) = sin(u), cn(u, 0) = cos(u) and dn(u, 0) = 1].
160 ['sn(u, 1) = tanh(u), cn(u, 1) = dn(u, 1) = 1 / cosh(u)].
162 And for ['k[super 4] < [epsilon]] we have:
166 Otherwise the values are calculated using the method of [@http://dlmf.nist.gov/22.20#SS2 arithmetic geometric means].
170 [section:jacobi_cd Jacobi Elliptic Function cd]
175 #include <boost/math/special_functions/jacobi_elliptic.hpp>
178 namespace boost { namespace math {
180 template <class T, class U>
181 ``__sf_result`` jacobi_cd(T k, U u);
183 template <class T, class U, class Policy>
184 ``__sf_result`` jacobi_cd(T k, U u, const Policy& pol);
188 [heading Description]
190 This function returns the Jacobi elliptic function ['cd].
194 This function is a trivial wrapper around __jacobi_elliptic, with:
196 ['cd(u, k) = cn(u, k) / dn(u, k)]
202 [section:jacobi_cn Jacobi Elliptic Function cn]
207 #include <boost/math/special_functions/jacobi_elliptic.hpp>
210 namespace boost { namespace math {
212 template <class T, class U>
213 ``__sf_result`` jacobi_cn(T k, U u);
215 template <class T, class U, class Policy>
216 ``__sf_result`` jacobi_cn(T k, U u, const Policy& pol);
220 [heading Description]
222 This function returns the Jacobi elliptic function ['cn].
226 This function is a trivial wrapper around __jacobi_elliptic.
232 [section:jacobi_cs Jacobi Elliptic Function cs]
237 #include <boost/math/special_functions/jacobi_elliptic.hpp>
240 namespace boost { namespace math {
242 template <class T, class U>
243 ``__sf_result`` jacobi_cs(T k, U u);
245 template <class T, class U, class Policy>
246 ``__sf_result`` jacobi_cs(T k, U u, const Policy& pol);
250 [heading Description]
252 This function returns the Jacobi elliptic function ['cs].
256 This function is a trivial wrapper around __jacobi_elliptic, with:
258 ['cs(u, k) = cn(u, k) / sn(u, k)]
264 [section:jacobi_dc Jacobi Elliptic Function dc]
269 #include <boost/math/special_functions/jacobi_elliptic.hpp>
272 namespace boost { namespace math {
274 template <class T, class U>
275 ``__sf_result`` jacobi_dc(T k, U u);
277 template <class T, class U, class Policy>
278 ``__sf_result`` jacobi_dc(T k, U u, const Policy& pol);
282 [heading Description]
284 This function returns the Jacobi elliptic function ['dc].
288 This function is a trivial wrapper around __jacobi_elliptic, with:
290 ['dc(u, k) = dn(u, k) / cn(u, k)]
296 [section:jacobi_dn Jacobi Elliptic Function dn]
301 #include <boost/math/special_functions/jacobi_elliptic.hpp>
304 namespace boost { namespace math {
306 template <class T, class U>
307 ``__sf_result`` jacobi_dn(T k, U u);
309 template <class T, class U, class Policy>
310 ``__sf_result`` jacobi_dn(T k, U u, const Policy& pol);
314 [heading Description]
316 This function returns the Jacobi elliptic function ['dn].
320 This function is a trivial wrapper around __jacobi_elliptic.
326 [section:jacobi_ds Jacobi Elliptic Function ds]
331 #include <boost/math/special_functions/jacobi_elliptic.hpp>
334 namespace boost { namespace math {
336 template <class T, class U>
337 ``__sf_result`` jacobi_ds(T k, U u);
339 template <class T, class U, class Policy>
340 ``__sf_result`` jacobi_ds(T k, U u, const Policy& pol);
344 [heading Description]
346 This function returns the Jacobi elliptic function ['ds].
350 This function is a trivial wrapper around __jacobi_elliptic, with:
352 ['ds(u, k) = dn(u, k) / sn(u, k)]
358 [section:jacobi_nc Jacobi Elliptic Function nc]
363 #include <boost/math/special_functions/jacobi_elliptic.hpp>
366 namespace boost { namespace math {
368 template <class T, class U>
369 ``__sf_result`` jacobi_nc(T k, U u);
371 template <class T, class U, class Policy>
372 ``__sf_result`` jacobi_nc(T k, U u, const Policy& pol);
376 [heading Description]
378 This function returns the Jacobi elliptic function ['nc].
382 This function is a trivial wrapper around __jacobi_elliptic, with:
384 ['nc(u, k) = 1 / cn(u, k)]
390 [section:jacobi_nd Jacobi Elliptic Function nd]
395 #include <boost/math/special_functions/jacobi_elliptic.hpp>
398 namespace boost { namespace math {
400 template <class T, class U>
401 ``__sf_result`` jacobi_nd(T k, U u);
403 template <class T, class U, class Policy>
404 ``__sf_result`` jacobi_nd(T k, U u, const Policy& pol);
408 [heading Description]
410 This function returns the Jacobi elliptic function ['nd].
414 This function is a trivial wrapper around __jacobi_elliptic, with:
416 ['nd(u, k) = 1 / dn(u, k)]
422 [section:jacobi_ns Jacobi Elliptic Function ns]
427 #include <boost/math/special_functions/jacobi_elliptic.hpp>
430 namespace boost { namespace math {
432 template <class T, class U>
433 ``__sf_result`` jacobi_ns(T k, U u);
435 template <class T, class U, class Policy>
436 ``__sf_result`` jacobi_ns(T k, U u, const Policy& pol);
440 [heading Description]
442 This function returns the Jacobi elliptic function ['ns].
446 This function is a trivial wrapper around __jacobi_elliptic, with:
448 ['ns(u, k) = 1 / sn(u, k)]
454 [section:jacobi_sc Jacobi Elliptic Function sc]
459 #include <boost/math/special_functions/jacobi_elliptic.hpp>
462 namespace boost { namespace math {
464 template <class T, class U>
465 ``__sf_result`` jacobi_sc(T k, U u);
467 template <class T, class U, class Policy>
468 ``__sf_result`` jacobi_sc(T k, U u, const Policy& pol);
472 [heading Description]
474 This function returns the Jacobi elliptic function ['sc].
478 This function is a trivial wrapper around __jacobi_elliptic, with:
480 ['sc(u, k) = sn(u, k) / cn(u, k)]
486 [section:jacobi_sd Jacobi Elliptic Function sd]
491 #include <boost/math/special_functions/jacobi_elliptic.hpp>
494 namespace boost { namespace math {
496 template <class T, class U>
497 ``__sf_result`` jacobi_sd(T k, U u);
499 template <class T, class U, class Policy>
500 ``__sf_result`` jacobi_sd(T k, U u, const Policy& pol);
504 [heading Description]
506 This function returns the Jacobi elliptic function ['sd].
510 This function is a trivial wrapper around __jacobi_elliptic, with:
512 ['sd(u, k) = sn(u, k) / dn(u, k)]
518 [section:jacobi_sn Jacobi Elliptic Function sn]
523 #include <boost/math/special_functions/jacobi_elliptic.hpp>
526 namespace boost { namespace math {
528 template <class T, class U>
529 ``__sf_result`` jacobi_sn(T k, U u);
531 template <class T, class U, class Policy>
532 ``__sf_result`` jacobi_sn(T k, U u, const Policy& pol);
536 [heading Description]
538 This function returns the Jacobi elliptic function ['sn].
542 This function is a trivial wrapper around __jacobi_elliptic.