gsl_specfunc__trig.c 19 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772
  1. /* specfunc/trig.c
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000 Gerard Jungman
  4. *
  5. * This program is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU General Public License as published by
  7. * the Free Software Foundation; either version 3 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * This program is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  13. * General Public License for more details.
  14. *
  15. * You should have received a copy of the GNU General Public License
  16. * along with this program; if not, write to the Free Software
  17. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
  18. */
  19. /* Author: G. Jungman */
  20. #include "gsl__config.h"
  21. #include "gsl_math.h"
  22. #include "gsl_errno.h"
  23. #include "gsl_sf_log.h"
  24. #include "gsl_sf_trig.h"
  25. #include "gsl_specfunc__error.h"
  26. #include "gsl_specfunc__chebyshev.h"
  27. #include "gsl_specfunc__cheb_eval.c"
  28. /* sinh(x) series
  29. * double-precision for |x| < 1.0
  30. */
  31. inline
  32. static
  33. int
  34. sinh_series(const double x, double * result)
  35. {
  36. const double y = x*x;
  37. const double c0 = 1.0/6.0;
  38. const double c1 = 1.0/120.0;
  39. const double c2 = 1.0/5040.0;
  40. const double c3 = 1.0/362880.0;
  41. const double c4 = 1.0/39916800.0;
  42. const double c5 = 1.0/6227020800.0;
  43. const double c6 = 1.0/1307674368000.0;
  44. const double c7 = 1.0/355687428096000.0;
  45. *result = x*(1.0 + y*(c0+y*(c1+y*(c2+y*(c3+y*(c4+y*(c5+y*(c6+y*c7))))))));
  46. return GSL_SUCCESS;
  47. }
  48. /* cosh(x)-1 series
  49. * double-precision for |x| < 1.0
  50. */
  51. inline
  52. static
  53. int
  54. cosh_m1_series(const double x, double * result)
  55. {
  56. const double y = x*x;
  57. const double c0 = 0.5;
  58. const double c1 = 1.0/24.0;
  59. const double c2 = 1.0/720.0;
  60. const double c3 = 1.0/40320.0;
  61. const double c4 = 1.0/3628800.0;
  62. const double c5 = 1.0/479001600.0;
  63. const double c6 = 1.0/87178291200.0;
  64. const double c7 = 1.0/20922789888000.0;
  65. const double c8 = 1.0/6402373705728000.0;
  66. *result = y*(c0+y*(c1+y*(c2+y*(c3+y*(c4+y*(c5+y*(c6+y*(c7+y*c8))))))));
  67. return GSL_SUCCESS;
  68. }
  69. /* Chebyshev expansion for f(t) = sinc((t+1)/2), -1 < t < 1
  70. */
  71. static double sinc_data[17] = {
  72. 1.133648177811747875422,
  73. -0.532677564732557348781,
  74. -0.068293048346633177859,
  75. 0.033403684226353715020,
  76. 0.001485679893925747818,
  77. -0.000734421305768455295,
  78. -0.000016837282388837229,
  79. 0.000008359950146618018,
  80. 0.000000117382095601192,
  81. -0.000000058413665922724,
  82. -0.000000000554763755743,
  83. 0.000000000276434190426,
  84. 0.000000000001895374892,
  85. -0.000000000000945237101,
  86. -0.000000000000004900690,
  87. 0.000000000000002445383,
  88. 0.000000000000000009925
  89. };
  90. static cheb_series sinc_cs = {
  91. sinc_data,
  92. 16,
  93. -1, 1,
  94. 10
  95. };
  96. /* Chebyshev expansion for f(t) = g((t+1)Pi/8), -1<t<1
  97. * g(x) = (sin(x)/x - 1)/(x*x)
  98. */
  99. static double sin_data[12] = {
  100. -0.3295190160663511504173,
  101. 0.0025374284671667991990,
  102. 0.0006261928782647355874,
  103. -4.6495547521854042157541e-06,
  104. -5.6917531549379706526677e-07,
  105. 3.7283335140973803627866e-09,
  106. 3.0267376484747473727186e-10,
  107. -1.7400875016436622322022e-12,
  108. -1.0554678305790849834462e-13,
  109. 5.3701981409132410797062e-16,
  110. 2.5984137983099020336115e-17,
  111. -1.1821555255364833468288e-19
  112. };
  113. static cheb_series sin_cs = {
  114. sin_data,
  115. 11,
  116. -1, 1,
  117. 11
  118. };
  119. /* Chebyshev expansion for f(t) = g((t+1)Pi/8), -1<t<1
  120. * g(x) = (2(cos(x) - 1)/(x^2) + 1) / x^2
  121. */
  122. static double cos_data[11] = {
  123. 0.165391825637921473505668118136,
  124. -0.00084852883845000173671196530195,
  125. -0.000210086507222940730213625768083,
  126. 1.16582269619760204299639757584e-6,
  127. 1.43319375856259870334412701165e-7,
  128. -7.4770883429007141617951330184e-10,
  129. -6.0969994944584252706997438007e-11,
  130. 2.90748249201909353949854872638e-13,
  131. 1.77126739876261435667156490461e-14,
  132. -7.6896421502815579078577263149e-17,
  133. -3.7363121133079412079201377318e-18
  134. };
  135. static cheb_series cos_cs = {
  136. cos_data,
  137. 10,
  138. -1, 1,
  139. 10
  140. };
  141. /*-*-*-*-*-*-*-*-*-*-*-* Functions with Error Codes *-*-*-*-*-*-*-*-*-*-*-*/
  142. /* I would have prefered just using the library sin() function.
  143. * But after some experimentation I decided that there was
  144. * no good way to understand the error; library sin() is just a black box.
  145. * So we have to roll our own.
  146. */
  147. int
  148. gsl_sf_sin_e(double x, gsl_sf_result * result)
  149. {
  150. /* CHECK_POINTER(result) */
  151. {
  152. const double P1 = 7.85398125648498535156e-1;
  153. const double P2 = 3.77489470793079817668e-8;
  154. const double P3 = 2.69515142907905952645e-15;
  155. const double sgn_x = GSL_SIGN(x);
  156. const double abs_x = fabs(x);
  157. if(abs_x < GSL_ROOT4_DBL_EPSILON) {
  158. const double x2 = x*x;
  159. result->val = x * (1.0 - x2/6.0);
  160. result->err = fabs(x*x2*x2 / 100.0);
  161. return GSL_SUCCESS;
  162. }
  163. else {
  164. double sgn_result = sgn_x;
  165. double y = floor(abs_x/(0.25*M_PI));
  166. int octant = y - ldexp(floor(ldexp(y,-3)),3);
  167. int stat_cs;
  168. double z;
  169. if(GSL_IS_ODD(octant)) {
  170. octant += 1;
  171. octant &= 07;
  172. y += 1.0;
  173. }
  174. if(octant > 3) {
  175. octant -= 4;
  176. sgn_result = -sgn_result;
  177. }
  178. z = ((abs_x - y * P1) - y * P2) - y * P3;
  179. if(octant == 0) {
  180. gsl_sf_result sin_cs_result;
  181. const double t = 8.0*fabs(z)/M_PI - 1.0;
  182. stat_cs = cheb_eval_e(&sin_cs, t, &sin_cs_result);
  183. result->val = z * (1.0 + z*z * sin_cs_result.val);
  184. }
  185. else { /* octant == 2 */
  186. gsl_sf_result cos_cs_result;
  187. const double t = 8.0*fabs(z)/M_PI - 1.0;
  188. stat_cs = cheb_eval_e(&cos_cs, t, &cos_cs_result);
  189. result->val = 1.0 - 0.5*z*z * (1.0 - z*z * cos_cs_result.val);
  190. }
  191. result->val *= sgn_result;
  192. if(abs_x > 1.0/GSL_DBL_EPSILON) {
  193. result->err = fabs(result->val);
  194. }
  195. else if(abs_x > 100.0/GSL_SQRT_DBL_EPSILON) {
  196. result->err = 2.0 * abs_x * GSL_DBL_EPSILON * fabs(result->val);
  197. }
  198. else if(abs_x > 0.1/GSL_SQRT_DBL_EPSILON) {
  199. result->err = 2.0 * GSL_SQRT_DBL_EPSILON * fabs(result->val);
  200. }
  201. else {
  202. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  203. }
  204. return stat_cs;
  205. }
  206. }
  207. }
  208. int
  209. gsl_sf_cos_e(double x, gsl_sf_result * result)
  210. {
  211. /* CHECK_POINTER(result) */
  212. {
  213. const double P1 = 7.85398125648498535156e-1;
  214. const double P2 = 3.77489470793079817668e-8;
  215. const double P3 = 2.69515142907905952645e-15;
  216. const double abs_x = fabs(x);
  217. if(abs_x < GSL_ROOT4_DBL_EPSILON) {
  218. const double x2 = x*x;
  219. result->val = 1.0 - 0.5*x2;
  220. result->err = fabs(x2*x2/12.0);
  221. return GSL_SUCCESS;
  222. }
  223. else {
  224. double sgn_result = 1.0;
  225. double y = floor(abs_x/(0.25*M_PI));
  226. int octant = y - ldexp(floor(ldexp(y,-3)),3);
  227. int stat_cs;
  228. double z;
  229. if(GSL_IS_ODD(octant)) {
  230. octant += 1;
  231. octant &= 07;
  232. y += 1.0;
  233. }
  234. if(octant > 3) {
  235. octant -= 4;
  236. sgn_result = -sgn_result;
  237. }
  238. if(octant > 1) {
  239. sgn_result = -sgn_result;
  240. }
  241. z = ((abs_x - y * P1) - y * P2) - y * P3;
  242. if(octant == 0) {
  243. gsl_sf_result cos_cs_result;
  244. const double t = 8.0*fabs(z)/M_PI - 1.0;
  245. stat_cs = cheb_eval_e(&cos_cs, t, &cos_cs_result);
  246. result->val = 1.0 - 0.5*z*z * (1.0 - z*z * cos_cs_result.val);
  247. }
  248. else { /* octant == 2 */
  249. gsl_sf_result sin_cs_result;
  250. const double t = 8.0*fabs(z)/M_PI - 1.0;
  251. stat_cs = cheb_eval_e(&sin_cs, t, &sin_cs_result);
  252. result->val = z * (1.0 + z*z * sin_cs_result.val);
  253. }
  254. result->val *= sgn_result;
  255. if(abs_x > 1.0/GSL_DBL_EPSILON) {
  256. result->err = fabs(result->val);
  257. }
  258. else if(abs_x > 100.0/GSL_SQRT_DBL_EPSILON) {
  259. result->err = 2.0 * abs_x * GSL_DBL_EPSILON * fabs(result->val);
  260. }
  261. else if(abs_x > 0.1/GSL_SQRT_DBL_EPSILON) {
  262. result->err = 2.0 * GSL_SQRT_DBL_EPSILON * fabs(result->val);
  263. }
  264. else {
  265. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  266. }
  267. return stat_cs;
  268. }
  269. }
  270. }
  271. int
  272. gsl_sf_hypot_e(const double x, const double y, gsl_sf_result * result)
  273. {
  274. /* CHECK_POINTER(result) */
  275. if(x == 0.0 && y == 0.0) {
  276. result->val = 0.0;
  277. result->err = 0.0;
  278. return GSL_SUCCESS;
  279. }
  280. else {
  281. const double a = fabs(x);
  282. const double b = fabs(y);
  283. const double min = GSL_MIN_DBL(a,b);
  284. const double max = GSL_MAX_DBL(a,b);
  285. const double rat = min/max;
  286. const double root_term = sqrt(1.0 + rat*rat);
  287. if(max < GSL_DBL_MAX/root_term) {
  288. result->val = max * root_term;
  289. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  290. return GSL_SUCCESS;
  291. }
  292. else {
  293. OVERFLOW_ERROR(result);
  294. }
  295. }
  296. }
  297. int
  298. gsl_sf_complex_sin_e(const double zr, const double zi,
  299. gsl_sf_result * szr, gsl_sf_result * szi)
  300. {
  301. /* CHECK_POINTER(szr) */
  302. /* CHECK_POINTER(szi) */
  303. if(fabs(zi) < 1.0) {
  304. double ch_m1, sh;
  305. sinh_series(zi, &sh);
  306. cosh_m1_series(zi, &ch_m1);
  307. szr->val = sin(zr)*(ch_m1 + 1.0);
  308. szi->val = cos(zr)*sh;
  309. szr->err = 2.0 * GSL_DBL_EPSILON * fabs(szr->val);
  310. szi->err = 2.0 * GSL_DBL_EPSILON * fabs(szi->val);
  311. return GSL_SUCCESS;
  312. }
  313. else if(fabs(zi) < GSL_LOG_DBL_MAX) {
  314. double ex = exp(zi);
  315. double ch = 0.5*(ex+1.0/ex);
  316. double sh = 0.5*(ex-1.0/ex);
  317. szr->val = sin(zr)*ch;
  318. szi->val = cos(zr)*sh;
  319. szr->err = 2.0 * GSL_DBL_EPSILON * fabs(szr->val);
  320. szi->err = 2.0 * GSL_DBL_EPSILON * fabs(szi->val);
  321. return GSL_SUCCESS;
  322. }
  323. else {
  324. OVERFLOW_ERROR_2(szr, szi);
  325. }
  326. }
  327. int
  328. gsl_sf_complex_cos_e(const double zr, const double zi,
  329. gsl_sf_result * czr, gsl_sf_result * czi)
  330. {
  331. /* CHECK_POINTER(czr) */
  332. /* CHECK_POINTER(czi) */
  333. if(fabs(zi) < 1.0) {
  334. double ch_m1, sh;
  335. sinh_series(zi, &sh);
  336. cosh_m1_series(zi, &ch_m1);
  337. czr->val = cos(zr)*(ch_m1 + 1.0);
  338. czi->val = -sin(zr)*sh;
  339. czr->err = 2.0 * GSL_DBL_EPSILON * fabs(czr->val);
  340. czi->err = 2.0 * GSL_DBL_EPSILON * fabs(czi->val);
  341. return GSL_SUCCESS;
  342. }
  343. else if(fabs(zi) < GSL_LOG_DBL_MAX) {
  344. double ex = exp(zi);
  345. double ch = 0.5*(ex+1.0/ex);
  346. double sh = 0.5*(ex-1.0/ex);
  347. czr->val = cos(zr)*ch;
  348. czi->val = -sin(zr)*sh;
  349. czr->err = 2.0 * GSL_DBL_EPSILON * fabs(czr->val);
  350. czi->err = 2.0 * GSL_DBL_EPSILON * fabs(czi->val);
  351. return GSL_SUCCESS;
  352. }
  353. else {
  354. OVERFLOW_ERROR_2(czr,czi);
  355. }
  356. }
  357. int
  358. gsl_sf_complex_logsin_e(const double zr, const double zi,
  359. gsl_sf_result * lszr, gsl_sf_result * lszi)
  360. {
  361. /* CHECK_POINTER(lszr) */
  362. /* CHECK_POINTER(lszi) */
  363. if(zi > 60.0) {
  364. lszr->val = -M_LN2 + zi;
  365. lszi->val = 0.5*M_PI - zr;
  366. lszr->err = 2.0 * GSL_DBL_EPSILON * fabs(lszr->val);
  367. lszi->err = 2.0 * GSL_DBL_EPSILON * fabs(lszi->val);
  368. }
  369. else if(zi < -60.0) {
  370. lszr->val = -M_LN2 - zi;
  371. lszi->val = -0.5*M_PI + zr;
  372. lszr->err = 2.0 * GSL_DBL_EPSILON * fabs(lszr->val);
  373. lszi->err = 2.0 * GSL_DBL_EPSILON * fabs(lszi->val);
  374. }
  375. else {
  376. gsl_sf_result sin_r, sin_i;
  377. int status;
  378. gsl_sf_complex_sin_e(zr, zi, &sin_r, &sin_i); /* ok by construction */
  379. status = gsl_sf_complex_log_e(sin_r.val, sin_i.val, lszr, lszi);
  380. if(status == GSL_EDOM) {
  381. DOMAIN_ERROR_2(lszr, lszi);
  382. }
  383. }
  384. return gsl_sf_angle_restrict_symm_e(&(lszi->val));
  385. }
  386. int
  387. gsl_sf_lnsinh_e(const double x, gsl_sf_result * result)
  388. {
  389. /* CHECK_POINTER(result) */
  390. if(x <= 0.0) {
  391. DOMAIN_ERROR(result);
  392. }
  393. else if(fabs(x) < 1.0) {
  394. double eps;
  395. sinh_series(x, &eps);
  396. result->val = log(eps);
  397. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  398. return GSL_SUCCESS;
  399. }
  400. else if(x < -0.5*GSL_LOG_DBL_EPSILON) {
  401. result->val = x + log(0.5*(1.0 - exp(-2.0*x)));
  402. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  403. return GSL_SUCCESS;
  404. }
  405. else {
  406. result->val = -M_LN2 + x;
  407. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  408. return GSL_SUCCESS;
  409. }
  410. }
  411. int gsl_sf_lncosh_e(const double x, gsl_sf_result * result)
  412. {
  413. /* CHECK_POINTER(result) */
  414. if(fabs(x) < 1.0) {
  415. double eps;
  416. cosh_m1_series(x, &eps);
  417. return gsl_sf_log_1plusx_e(eps, result);
  418. }
  419. else if(x < -0.5*GSL_LOG_DBL_EPSILON) {
  420. result->val = x + log(0.5*(1.0 + exp(-2.0*x)));
  421. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  422. return GSL_SUCCESS;
  423. }
  424. else {
  425. result->val = -M_LN2 + x;
  426. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  427. return GSL_SUCCESS;
  428. }
  429. }
  430. /*
  431. inline int gsl_sf_sincos_e(const double theta, double * s, double * c)
  432. {
  433. double tan_half = tan(0.5 * theta);
  434. double den = 1. + tan_half*tan_half;
  435. double cos_theta = (1.0 - tan_half*tan_half) / den;
  436. double sin_theta = 2.0 * tan_half / den;
  437. }
  438. */
  439. int
  440. gsl_sf_polar_to_rect(const double r, const double theta,
  441. gsl_sf_result * x, gsl_sf_result * y)
  442. {
  443. double t = theta;
  444. int status = gsl_sf_angle_restrict_symm_e(&t);
  445. double c = cos(t);
  446. double s = sin(t);
  447. x->val = r * cos(t);
  448. y->val = r * sin(t);
  449. x->err = r * fabs(s * GSL_DBL_EPSILON * t);
  450. x->err += 2.0 * GSL_DBL_EPSILON * fabs(x->val);
  451. y->err = r * fabs(c * GSL_DBL_EPSILON * t);
  452. y->err += 2.0 * GSL_DBL_EPSILON * fabs(y->val);
  453. return status;
  454. }
  455. int
  456. gsl_sf_rect_to_polar(const double x, const double y,
  457. gsl_sf_result * r, gsl_sf_result * theta)
  458. {
  459. int stat_h = gsl_sf_hypot_e(x, y, r);
  460. if(r->val > 0.0) {
  461. theta->val = atan2(y, x);
  462. theta->err = 2.0 * GSL_DBL_EPSILON * fabs(theta->val);
  463. return stat_h;
  464. }
  465. else {
  466. DOMAIN_ERROR(theta);
  467. }
  468. }
  469. int gsl_sf_angle_restrict_symm_err_e(const double theta, gsl_sf_result * result)
  470. {
  471. /* synthetic extended precision constants */
  472. const double P1 = 4 * 7.8539812564849853515625e-01;
  473. const double P2 = 4 * 3.7748947079307981766760e-08;
  474. const double P3 = 4 * 2.6951514290790594840552e-15;
  475. const double TwoPi = 2*(P1 + P2 + P3);
  476. const double y = GSL_SIGN(theta) * 2 * floor(fabs(theta)/TwoPi);
  477. double r = ((theta - y*P1) - y*P2) - y*P3;
  478. if(r > M_PI) { r = (((r-2*P1)-2*P2)-2*P3); } /* r-TwoPi */
  479. else if (r < -M_PI) r = (((r+2*P1)+2*P2)+2*P3); /* r+TwoPi */
  480. result->val = r;
  481. if(fabs(theta) > 0.0625/GSL_DBL_EPSILON) {
  482. result->val = GSL_NAN;
  483. result->err = GSL_NAN;
  484. GSL_ERROR ("error", GSL_ELOSS);
  485. }
  486. else if(fabs(theta) > 0.0625/GSL_SQRT_DBL_EPSILON) {
  487. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val - theta);
  488. return GSL_SUCCESS;
  489. }
  490. else {
  491. double delta = fabs(result->val - theta);
  492. result->err = 2.0 * GSL_DBL_EPSILON * ((delta < M_PI) ? delta : M_PI);
  493. return GSL_SUCCESS;
  494. }
  495. }
  496. int gsl_sf_angle_restrict_pos_err_e(const double theta, gsl_sf_result * result)
  497. {
  498. /* synthetic extended precision constants */
  499. const double P1 = 4 * 7.85398125648498535156e-01;
  500. const double P2 = 4 * 3.77489470793079817668e-08;
  501. const double P3 = 4 * 2.69515142907905952645e-15;
  502. const double TwoPi = 2*(P1 + P2 + P3);
  503. const double y = 2*floor(theta/TwoPi);
  504. double r = ((theta - y*P1) - y*P2) - y*P3;
  505. if(r > TwoPi) {r = (((r-2*P1)-2*P2)-2*P3); } /* r-TwoPi */
  506. else if (r < 0) { /* may happen due to FP rounding */
  507. r = (((r+2*P1)+2*P2)+2*P3); /* r+TwoPi */
  508. }
  509. result->val = r;
  510. if(fabs(theta) > 0.0625/GSL_DBL_EPSILON) {
  511. result->val = GSL_NAN;
  512. result->err = fabs(result->val);
  513. GSL_ERROR ("error", GSL_ELOSS);
  514. }
  515. else if(fabs(theta) > 0.0625/GSL_SQRT_DBL_EPSILON) {
  516. result->err = GSL_DBL_EPSILON * fabs(result->val - theta);
  517. return GSL_SUCCESS;
  518. }
  519. else {
  520. double delta = fabs(result->val - theta);
  521. result->err = 2.0 * GSL_DBL_EPSILON * ((delta < M_PI) ? delta : M_PI);
  522. return GSL_SUCCESS;
  523. }
  524. }
  525. int gsl_sf_angle_restrict_symm_e(double * theta)
  526. {
  527. gsl_sf_result r;
  528. int stat = gsl_sf_angle_restrict_symm_err_e(*theta, &r);
  529. *theta = r.val;
  530. return stat;
  531. }
  532. int gsl_sf_angle_restrict_pos_e(double * theta)
  533. {
  534. gsl_sf_result r;
  535. int stat = gsl_sf_angle_restrict_pos_err_e(*theta, &r);
  536. *theta = r.val;
  537. return stat;
  538. }
  539. int gsl_sf_sin_err_e(const double x, const double dx, gsl_sf_result * result)
  540. {
  541. int stat_s = gsl_sf_sin_e(x, result);
  542. result->err += fabs(cos(x) * dx);
  543. result->err += GSL_DBL_EPSILON * fabs(result->val);
  544. return stat_s;
  545. }
  546. int gsl_sf_cos_err_e(const double x, const double dx, gsl_sf_result * result)
  547. {
  548. int stat_c = gsl_sf_cos_e(x, result);
  549. result->err += fabs(sin(x) * dx);
  550. result->err += GSL_DBL_EPSILON * fabs(result->val);
  551. return stat_c;
  552. }
  553. #if 0
  554. int
  555. gsl_sf_sin_pi_x_e(const double x, gsl_sf_result * result)
  556. {
  557. /* CHECK_POINTER(result) */
  558. if(-100.0 < x && x < 100.0) {
  559. result->val = sin(M_PI * x) / (M_PI * x);
  560. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  561. return GSL_SUCCESS;
  562. }
  563. else {
  564. const double N = floor(x + 0.5);
  565. const double f = x - N;
  566. if(N < INT_MAX && N > INT_MIN) {
  567. /* Make it an integer if we can. Saves another
  568. * call to floor().
  569. */
  570. const int intN = (int)N;
  571. const double sign = ( GSL_IS_ODD(intN) ? -1.0 : 1.0 );
  572. result->val = sign * sin(M_PI * f);
  573. result->err = GSL_DBL_EPSILON * fabs(result->val);
  574. }
  575. else if(N > 2.0/GSL_DBL_EPSILON || N < -2.0/GSL_DBL_EPSILON) {
  576. /* All integer-valued floating point numbers
  577. * bigger than 2/eps=2^53 are actually even.
  578. */
  579. result->val = 0.0;
  580. result->err = 0.0;
  581. }
  582. else {
  583. const double resN = N - 2.0*floor(0.5*N); /* 0 for even N, 1 for odd N */
  584. const double sign = ( fabs(resN) > 0.5 ? -1.0 : 1.0 );
  585. result->val = sign * sin(M_PI*f);
  586. result->err = GSL_DBL_EPSILON * fabs(result->val);
  587. }
  588. return GSL_SUCCESS;
  589. }
  590. }
  591. #endif
  592. int gsl_sf_sinc_e(double x, gsl_sf_result * result)
  593. {
  594. /* CHECK_POINTER(result) */
  595. {
  596. const double ax = fabs(x);
  597. if(ax < 0.8) {
  598. /* Do not go to the limit of the fit since
  599. * there is a zero there and the Chebyshev
  600. * accuracy will go to zero.
  601. */
  602. return cheb_eval_e(&sinc_cs, 2.0*ax-1.0, result);
  603. }
  604. else if(ax < 100.0) {
  605. /* Small arguments are no problem.
  606. * We trust the library sin() to
  607. * roughly machine precision.
  608. */
  609. result->val = sin(M_PI * ax)/(M_PI * ax);
  610. result->err = 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  611. return GSL_SUCCESS;
  612. }
  613. else {
  614. /* Large arguments must be handled separately.
  615. */
  616. const double r = M_PI*ax;
  617. gsl_sf_result s;
  618. int stat_s = gsl_sf_sin_e(r, &s);
  619. result->val = s.val/r;
  620. result->err = s.err/r + 2.0 * GSL_DBL_EPSILON * fabs(result->val);
  621. return stat_s;
  622. }
  623. }
  624. }
  625. /*-*-*-*-*-*-*-*-*-* Functions w/ Natural Prototypes *-*-*-*-*-*-*-*-*-*-*/
  626. #include "gsl_specfunc__eval.h"
  627. double gsl_sf_sin(const double x)
  628. {
  629. EVAL_RESULT(gsl_sf_sin_e(x, &result));
  630. }
  631. double gsl_sf_cos(const double x)
  632. {
  633. EVAL_RESULT(gsl_sf_cos_e(x, &result));
  634. }
  635. double gsl_sf_hypot(const double x, const double y)
  636. {
  637. EVAL_RESULT(gsl_sf_hypot_e(x, y, &result));
  638. }
  639. double gsl_sf_lnsinh(const double x)
  640. {
  641. EVAL_RESULT(gsl_sf_lnsinh_e(x, &result));
  642. }
  643. double gsl_sf_lncosh(const double x)
  644. {
  645. EVAL_RESULT(gsl_sf_lncosh_e(x, &result));
  646. }
  647. double gsl_sf_angle_restrict_symm(const double theta)
  648. {
  649. double result = theta;
  650. EVAL_DOUBLE(gsl_sf_angle_restrict_symm_e(&result));
  651. }
  652. double gsl_sf_angle_restrict_pos(const double theta)
  653. {
  654. double result = theta;
  655. EVAL_DOUBLE(gsl_sf_angle_restrict_pos_e(&result));
  656. }
  657. #if 0
  658. double gsl_sf_sin_pi_x(const double x)
  659. {
  660. EVAL_RESULT(gsl_sf_sin_pi_x_e(x, &result));
  661. }
  662. #endif
  663. double gsl_sf_sinc(const double x)
  664. {
  665. EVAL_RESULT(gsl_sf_sinc_e(x, &result));
  666. }