• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
1[/ math.qbk
2  Copyright 2006 Hubert Holin and John Maddock.
3  Distributed under the Boost Software License, Version 1.0.
4  (See accompanying file LICENSE_1_0.txt or copy at
5  http://www.boost.org/LICENSE_1_0.txt).
6]
7
8[def __form1 [^\[0;+'''∞'''\[]]
9[def __form2 [^\]-'''∞''';+1\[]]
10[def __form3 [^\]-'''∞''';-1\[]]
11[def __form4 [^\]+1;+'''∞'''\[]]
12[def __form5 [^\[-1;-1+'''ε'''\[]]
13[def __form6 '''ε''']
14[def __form7 [^\]+1-'''ε''';+1\]]]
15
16[def __effects [*Effects: ]]
17[def __formula [*Formula: ]]
18[def __exm1 '''<code>e<superscript>x</superscript> - 1</code>''']
19[def __ex '''<code>e<superscript>x</superscript></code>''']
20[def __te '''2&#x03B5;''']
21
22[section:inv_hyper Inverse Hyperbolic Functions]
23
24[section:inv_hyper_over Inverse Hyperbolic Functions Overview]
25
26The exponential function is defined, for all objects for which this makes sense,
27as the power series
28[equation special_functions_blurb1]
29with ['[^n! = 1x2x3x4x5...xn]] (and ['[^0! = 1]] by definition) being the factorial of ['[^n]].
30In particular, the exponential function is well defined for real numbers,
31complex number, quaternions, octonions, and matrices of complex numbers,
32among others.
33
34[: ['[*Graph of exp on R]] ]
35
36[: [$../graphs/exp_on_r.png] ]
37
38[: ['[*Real and Imaginary parts of exp on C]]]
39[: [$../graphs/im_exp_on_c.png]]
40
41The hyperbolic functions are defined as power series which
42can be computed (for reals, complex, quaternions and octonions) as:
43
44Hyperbolic cosine: [equation special_functions_blurb5]
45
46Hyperbolic sine: [equation special_functions_blurb6]
47
48Hyperbolic tangent: [equation special_functions_blurb7]
49
50[: ['[*Trigonometric functions on R (cos: purple; sin: red; tan: blue)]]]
51[: [$../graphs/trigonometric.png]]
52
53[: ['[*Hyperbolic functions on r (cosh: purple; sinh: red; tanh: blue)]]]
54[: [$../graphs/hyperbolic.png]]
55
56The hyperbolic sine is one to one on the set of real numbers,
57with range the full set of reals, while the hyperbolic tangent is
58also one to one on the set of real numbers but with range __form1, and
59therefore both have inverses.
60
61The hyperbolic cosine is one to one from __form2 onto __form3 (and from __form4 onto __form3).
62
63The inverse function we use here is defined on __form3 with range __form2.
64
65The inverse of the hyperbolic tangent is called the Argument hyperbolic tangent,
66and can be computed as [equation special_functions_blurb15]
67
68The inverse of the hyperbolic sine is called the Argument hyperbolic sine,
69and can be computed (for __form5) as [equation special_functions_blurb17]
70
71The inverse of the hyperbolic cosine is called the Argument hyperbolic cosine,
72and can be computed as [equation special_functions_blurb18]
73
74[endsect] [/section:inv_hyper_over Inverse Hyperbolic Functions Overview]
75
76[section:acosh acosh]
77
78``
79#include <boost/math/special_functions/acosh.hpp>
80``
81
82   template<class T>
83   ``__sf_result`` acosh(const T x);
84
85   template<class T, class ``__Policy``>
86   ``__sf_result`` acosh(const T x, const ``__Policy``&);
87
88Computes the reciprocal of (the restriction to the range of __form1)
89[link math_toolkit.inv_hyper.inv_hyper_over
90the hyperbolic cosine function], at x. Values returned are positive.
91
92If x is in the range __form2 then returns the result of __domain_error.
93
94The return type of this function is computed using the __arg_promotion_rules:
95the return type is `double` when T is an integer type, and T otherwise.
96
97[optional_policy]
98
99[graph acosh]
100
101[h4 Accuracy]
102
103Generally accuracy is to within 1 or 2 __epsilon across all supported platforms.
104
105[h4 Testing]
106
107This function is tested using a combination of random test values designed to give
108full function coverage computed at high precision using the "naive" formula:
109
110[equation acosh1]
111
112along with a selection of sanity check values
113computed using functions.wolfram.com to at least 50 decimal digits.
114
115[h4 Implementation]
116
117For sufficiently large x, we can use the
118[@http://functions.wolfram.com/ElementaryFunctions/ArcCosh/06/01/06/01/0001/
119approximation]:
120
121[equation acosh2]
122
123For x sufficiently close to 1 we can use the
124[@http://functions.wolfram.com/ElementaryFunctions/ArcCosh/06/01/04/01/0001/
125approximation]:
126
127[equation acosh4]
128
129Otherwise for x close to 1 we can use the following rearrangement of the
130primary definition to preserve accuracy:
131
132[equation acosh3]
133
134Otherwise the
135[@http://functions.wolfram.com/ElementaryFunctions/ArcCosh/02/
136primary definition] is used:
137
138[equation acosh1]
139
140[endsect] [/section:acosh acosh]
141
142[section:asinh asinh]
143
144``
145#include <boost/math/special_functions/asinh.hpp>
146``
147
148   template<class T>
149   ``__sf_result`` asinh(const T x);
150
151   template<class T, class ``__Policy``>
152   ``__sf_result`` asinh(const T x, const ``__Policy``&);
153
154Computes the reciprocal of
155[link math_toolkit.inv_hyper.inv_hyper_over
156the hyperbolic sine function].
157
158The return type of this function is computed using the __arg_promotion_rules:
159the return type is `double` when T is an integer type, and T otherwise.
160
161[graph asinh]
162
163[optional_policy]
164
165[h4 Accuracy]
166
167Generally accuracy is to within 1 or 2 __epsilon across all supported platforms.
168
169[h4 Testing]
170
171This function is tested using a combination of random test values designed to give
172full function coverage computed at high precision using the "naive" formula:
173
174[equation asinh1]
175
176along with a selection of sanity check values
177computed using functions.wolfram.com to at least 50 decimal digits.
178
179[h4 Implementation]
180
181For sufficiently large x we can use the
182[@http://functions.wolfram.com/ElementaryFunctions/ArcSinh/06/01/06/01/0001/
183approximation]:
184
185[equation asinh2]
186
187While for very small x we can use the
188[@http://functions.wolfram.com/ElementaryFunctions/ArcSinh/06/01/03/01/0001/
189approximation]:
190
191[equation asinh3]
192
193For 0.5 > x > [epsilon] the following rearrangement of the primary definition is used:
194
195[equation asinh4]
196
197Otherwise evaluation is via the
198[@http://functions.wolfram.com/ElementaryFunctions/ArcSinh/02/
199primary definition]:
200
201[equation asinh4]
202
203[endsect] [/section:asinh asinh]
204
205[section:atanh atanh]
206
207``
208#include <boost/math/special_functions/atanh.hpp>
209``
210
211   template<class T>
212   ``__sf_result`` atanh(const T x);
213
214   template<class T, class ``__Policy``>
215   ``__sf_result`` atanh(const T x, const ``__Policy``&);
216
217Computes the reciprocal of
218[link math_toolkit.inv_hyper.inv_hyper_over
219the hyperbolic tangent function], at x.
220
221[optional_policy]
222
223If x is in the range
224__form3
225or in the range
226__form4
227then returns the result of __domain_error.
228
229If x is in the range
230__form5,
231then the result of -__overflow_error is returned, with
232__form6
233denoting `std::numeric_limits<T>::epsilon()`.
234
235If x is in the range
236__form7,
237then the result of __overflow_error is returned, with
238__form6
239denoting
240`std::numeric_limits<T>::epsilon()`.
241
242The return type of this function is computed using the __arg_promotion_rules:
243the return type is `double` when T is an integer type, and T otherwise.
244
245[graph atanh]
246
247[h4 Accuracy]
248
249Generally accuracy is to within 1 or 2 __epsilon across all supported platforms.
250
251[h4 Testing]
252
253This function is tested using a combination of random test values designed to give
254full function coverage computed at high precision using the "naive" formula:
255
256[equation atanh1]
257
258along with a selection of sanity check values
259computed using functions.wolfram.com to at least 50 decimal digits.
260
261[h4 Implementation]
262
263For sufficiently small x we can use the
264[@http://functions.wolfram.com/ElementaryFunctions/ArcTanh/06/01/03/01/ approximation]:
265
266[equation atanh2]
267
268Otherwise the
269[@http://functions.wolfram.com/ElementaryFunctions/ArcTanh/02/ primary definition]:
270
271[equation atanh1]
272
273or its equivalent form:
274
275[equation atanh3]
276
277is used.
278
279[endsect] [/section:atanh atanh]
280
281[endsect] [/section:inv_hyper Inverse Hyperbolic Functions]
282
283