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 Overview 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, Reinhardt, W. P., Walker, O. L.]].
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][thin][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 /k/ as the *first* argument - this is for alignment with
31 the Elliptic Integrals (but is different from other implementations, for example Mathworks).
33 A simple example comparing use of __WolframAlpha with Boost.Math (including much higher precision using Boost.Multiprecision)
34 is [@../../example/jacobi_zeta_example.cpp jacobi_zeta_example.cpp].
36 There are twelve functions for computing the twelve individual Jacobi elliptic functions: __jacobi_cd, __jacobi_cn, __jacobi_cs,
37 __jacobi_dc, __jacobi_dn, __jacobi_ds, __jacobi_nc, __jacobi_nd, __jacobi_ns, __jacobi_sc, __jacobi_sd and __jacobi_sn.
39 They are all called as for example:
43 Note however that these individual functions are all really thin wrappers around the function __jacobi_elliptic which calculates
44 the three copolar functions ['sn], ['cn] and ['dn] in a single function call.
46 [tip If you need more than one of these functions
47 for a given set of arguments, it's most efficient to use __jacobi_elliptic.]
49 [endsect] [/section:jac_over Overvew of the Jacobi Elliptic Functions]
51 [section:jacobi_elliptic Jacobi Elliptic SN, CN and DN]
56 #include <boost/math/special_functions/jacobi_elliptic.hpp>
59 namespace boost { namespace math {
61 template <class T, class U, class V>
62 ``__sf_result`` jacobi_elliptic(T k, U u, V* pcn, V* pdn);
64 template <class T, class U, class V, class Policy>
65 ``__sf_result`` jacobi_elliptic(T k, U u, V* pcn, V* pdn, const Policy&);
71 The function __jacobi_elliptic calculates the three copolar Jacobi elliptic functions
72 ['sn(u, k)], ['cn(u, k)] and ['dn(u, k)]. The returned value is
73 ['sn(u, k)], and if provided, `*pcn` is
74 set to ['cn(u, k)], and `*pdn` is set to ['dn(u, k)].
76 The functions are defined as follows, given:
80 The the angle ['[phi]] is called the ['amplitude] and:
85 ['[phi]] is called the amplitude.
86 ['k] is called the elliptic modulus.
89 [caution Rather like other elliptic functions, the Jacobi functions
90 are expressed in a variety of different ways. In particular,
91 the parameter /k/ (the modulus) may also be expressed using a modular
92 angle [alpha], or a parameter /m/. These are related by:
94 [expression k = sin [alpha]]
96 [expression m = k[super 2] = sin[super 2][alpha]]
98 So that the function ['sn] (for example) may be expressed as
101 [expression sn(u, k)]
103 [expression sn(u \\ [alpha])]
105 [expression sn(u | m)]
107 To further complicate matters, some texts refer to the ['complement
108 of the parameter m], or 1 - m, where:
110 [expression 1 - m = 1 - k[super 2] = cos[super 2][alpha]]
112 This implementation uses /k/ throughout, and makes this the first argument
113 to the functions: this is for alignment with the elliptic integrals which match the requirements
114 of the [@http://www.open-std.org/jtc1/sc22/wg21/docs/papers/2005/n1836.pdf
115 Technical Report on C++ Library Extensions]. However, you should
116 be extra careful when using these functions!]
120 The following graphs illustrate how these functions change as /k/ changes: for small /k/
121 these are sine waves, while as /k/ tends to 1 they become hyperbolic functions:
131 These functions are computed using only basic arithmetic operations and trigomometric functions, so
132 there isn't much variation in accuracy over differing platforms.
133 Typically errors are trivially small for small angles, and as is typical for cyclic
134 functions, grow as the angle increases.
135 Note that only results for the widest floating-point type on the
136 system are given as narrower types have __zero_error. All values
137 are relative errors in units of epsilon.
147 The tests use a mixture of spot test values calculated using the online
148 calculator at [@http://functions.wolfram.com/ functions.wolfram.com],
149 and random test data generated using
150 MPFR at 1000-bit precision and this implementation.
152 [heading Implementation]
154 For ['k > 1] we apply the relations:
158 Then filter off the special cases:
160 [expression ['sn(0, k) = 0] and ['cn(0, k) = dn(0, k) = 1]]
162 [expression ['sn(u, 0) = sin(u), cn(u, 0) = cos(u) and dn(u, 0) = 1]]
164 [expression ['sn(u, 1) = tanh(u), cn(u, 1) = dn(u, 1) = 1 / cosh(u)]]
166 And for ['k[super 4] < [epsilon]] we have:
170 Otherwise the values are calculated using the method of [@http://dlmf.nist.gov/22.20#SS2 arithmetic geometric means].
172 [endsect] [/section:jacobi_elliptic Jacobi Elliptic SN, CN and DN]
175 [section:jacobi_cd Jacobi Elliptic Function cd]
180 #include <boost/math/special_functions/jacobi_elliptic.hpp>
183 namespace boost { namespace math {
185 template <class T, class U>
186 ``__sf_result`` jacobi_cd(T k, U u);
188 template <class T, class U, class Policy>
189 ``__sf_result`` jacobi_cd(T k, U u, const Policy& pol);
193 [heading Description]
195 This function returns the Jacobi elliptic function ['cd].
199 This function is a trivial wrapper around __jacobi_elliptic, with:
201 [expression ['cd(u, k) = cn(u, k) / dn(u, k)]]
205 [endsect] [/section:jacobi_cd Jacobi Elliptic Function cd]
208 [section:jacobi_cn Jacobi Elliptic Function cn]
213 #include <boost/math/special_functions/jacobi_elliptic.hpp>
216 namespace boost { namespace math {
218 template <class T, class U>
219 ``__sf_result`` jacobi_cn(T k, U u);
221 template <class T, class U, class Policy>
222 ``__sf_result`` jacobi_cn(T k, U u, const Policy& pol);
226 [heading Description]
228 This function returns the Jacobi elliptic function ['cn].
232 This function is a trivial wrapper around __jacobi_elliptic.
236 [endsect] [/section:jacobi_cn Jacobi Elliptic Function cn]
239 [section:jacobi_cs Jacobi Elliptic Function cs]
244 #include <boost/math/special_functions/jacobi_elliptic.hpp>
247 namespace boost { namespace math {
249 template <class T, class U>
250 ``__sf_result`` jacobi_cs(T k, U u);
252 template <class T, class U, class Policy>
253 ``__sf_result`` jacobi_cs(T k, U u, const Policy& pol);
257 [heading Description]
259 This function returns the Jacobi elliptic function ['cs].
263 This function is a trivial wrapper around __jacobi_elliptic, with:
265 [expression ['cs(u, k) = cn(u, k) / sn(u, k)]]
269 [endsect] [/section:jacobi_cs Jacobi Elliptic Function cs]
272 [section:jacobi_dc Jacobi Elliptic Function dc]
277 #include <boost/math/special_functions/jacobi_elliptic.hpp>
280 namespace boost { namespace math {
282 template <class T, class U>
283 ``__sf_result`` jacobi_dc(T k, U u);
285 template <class T, class U, class Policy>
286 ``__sf_result`` jacobi_dc(T k, U u, const Policy& pol);
290 [heading Description]
292 This function returns the Jacobi elliptic function ['dc].
296 This function is a trivial wrapper around __jacobi_elliptic, with:
298 [expression ['dc(u, k) = dn(u, k) / cn(u, k)]]
302 [endsect] [/section:jacobi_dc Jacobi Elliptic Function dc]
304 [section:jacobi_dn Jacobi Elliptic Function dn]
309 #include <boost/math/special_functions/jacobi_elliptic.hpp>
312 namespace boost { namespace math {
314 template <class T, class U>
315 ``__sf_result`` jacobi_dn(T k, U u);
317 template <class T, class U, class Policy>
318 ``__sf_result`` jacobi_dn(T k, U u, const Policy& pol);
322 [heading Description]
324 This function returns the Jacobi elliptic function ['dn].
328 This function is a trivial wrapper around __jacobi_elliptic.
334 [section:jacobi_ds Jacobi Elliptic Function ds]
339 #include <boost/math/special_functions/jacobi_elliptic.hpp>
342 namespace boost { namespace math {
344 template <class T, class U>
345 ``__sf_result`` jacobi_ds(T k, U u);
347 template <class T, class U, class Policy>
348 ``__sf_result`` jacobi_ds(T k, U u, const Policy& pol);
352 [heading Description]
354 This function returns the Jacobi elliptic function ['ds].
358 This function is a trivial wrapper around __jacobi_elliptic, with:
360 [expression ['ds(u, k) = dn(u, k) / sn(u, k)]]
364 [endsect] [/section:jacobi_dn Jacobi Elliptic Function dn]
366 [section:jacobi_nc Jacobi Elliptic Function nc]
371 #include <boost/math/special_functions/jacobi_elliptic.hpp>
374 namespace boost { namespace math {
376 template <class T, class U>
377 ``__sf_result`` jacobi_nc(T k, U u);
379 template <class T, class U, class Policy>
380 ``__sf_result`` jacobi_nc(T k, U u, const Policy& pol);
384 [heading Description]
386 This function returns the Jacobi elliptic function ['nc].
390 This function is a trivial wrapper around __jacobi_elliptic, with:
392 [expression ['nc(u, k) = 1 / cn(u, k)]]
396 [endsect] [/section:jacobi_nc Jacobi Elliptic Function nc]
398 [section:jacobi_nd Jacobi Elliptic Function nd]
403 #include <boost/math/special_functions/jacobi_elliptic.hpp>
406 namespace boost { namespace math {
408 template <class T, class U>
409 ``__sf_result`` jacobi_nd(T k, U u);
411 template <class T, class U, class Policy>
412 ``__sf_result`` jacobi_nd(T k, U u, const Policy& pol);
416 [heading Description]
418 This function returns the Jacobi elliptic function ['nd].
422 This function is a trivial wrapper around __jacobi_elliptic, with:
424 [expression ['nd(u, k) = 1 / dn(u, k)]]
428 [endsect] [/section:jacobi_nd Jacobi Elliptic Function nd]
430 [section:jacobi_ns Jacobi Elliptic Function ns]
435 #include <boost/math/special_functions/jacobi_elliptic.hpp>
438 namespace boost { namespace math {
440 template <class T, class U>
441 ``__sf_result`` jacobi_ns(T k, U u);
443 template <class T, class U, class Policy>
444 ``__sf_result`` jacobi_ns(T k, U u, const Policy& pol);
448 [heading Description]
450 This function returns the Jacobi elliptic function ['ns].
454 This function is a trivial wrapper around __jacobi_elliptic, with:
456 [expression ['ns(u, k) = 1 / sn(u, k)]]
460 [endsect] [/section:jacobi_ns Jacobi Elliptic Function ns]
462 [section:jacobi_sc Jacobi Elliptic Function sc]
467 #include <boost/math/special_functions/jacobi_elliptic.hpp>
470 namespace boost { namespace math {
472 template <class T, class U>
473 ``__sf_result`` jacobi_sc(T k, U u);
475 template <class T, class U, class Policy>
476 ``__sf_result`` jacobi_sc(T k, U u, const Policy& pol);
480 [heading Description]
482 This function returns the Jacobi elliptic function ['sc].
486 This function is a trivial wrapper around __jacobi_elliptic, with:
488 [expression ['sc(u, k) = sn(u, k) / cn(u, k)]]
492 [endsect] [/section:jacobi_sc Jacobi Elliptic Function sc]
494 [section:jacobi_sd Jacobi Elliptic Function sd]
499 #include <boost/math/special_functions/jacobi_elliptic.hpp>
502 namespace boost { namespace math {
504 template <class T, class U>
505 ``__sf_result`` jacobi_sd(T k, U u);
507 template <class T, class U, class Policy>
508 ``__sf_result`` jacobi_sd(T k, U u, const Policy& pol);
512 [heading Description]
514 This function returns the Jacobi elliptic function ['sd].
518 This function is a trivial wrapper around __jacobi_elliptic, with:
520 [expression ['sd(u, k) = sn(u, k) / dn(u, k)]]
524 [endsect] [/section:jacobi_sd Jacobi Elliptic Function sd]
526 [section:jacobi_sn Jacobi Elliptic Function sn]
531 #include <boost/math/special_functions/jacobi_elliptic.hpp>
534 namespace boost { namespace math {
536 template <class T, class U>
537 ``__sf_result`` jacobi_sn(T k, U u);
539 template <class T, class U, class Policy>
540 ``__sf_result`` jacobi_sn(T k, U u, const Policy& pol);
544 [heading Description]
546 This function returns the Jacobi elliptic function ['sn].
550 This function is a trivial wrapper around __jacobi_elliptic.
554 [endsect] [/section:jacobi_sn Jacobi Elliptic Function sn]
556 [endsect] [/section:jacobi Jacobi Elliptic Functions]