gsl_sf_mathieu.h 3.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108
  1. /* specfunc/gsl_sf_mathieu.h
  2. *
  3. * Copyright (C) 2002 Lowell Johnson
  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., 675 Mass Ave, Cambridge, MA 02139, USA.
  18. */
  19. /* Author: L. Johnson */
  20. #ifndef __GSL_SF_MATHIEU_H__
  21. #define __GSL_SF_MATHIEU_H__
  22. #include "gsl_sf_result.h"
  23. #include "gsl_eigen.h"
  24. #undef __BEGIN_DECLS
  25. #undef __END_DECLS
  26. #ifdef __cplusplus
  27. # define __BEGIN_DECLS extern "C" {
  28. # define __END_DECLS }
  29. #else
  30. # define __BEGIN_DECLS /* empty */
  31. # define __END_DECLS /* empty */
  32. #endif
  33. __BEGIN_DECLS
  34. #define GSL_SF_MATHIEU_COEFF 100
  35. typedef struct
  36. {
  37. size_t size;
  38. size_t even_order;
  39. size_t odd_order;
  40. int extra_values;
  41. double qa; /* allow for caching of results: not implemented yet */
  42. double qb; /* allow for caching of results: not implemented yet */
  43. double *aa;
  44. double *bb;
  45. double *dd;
  46. double *ee;
  47. double *tt;
  48. double *e2;
  49. double *zz;
  50. gsl_vector *eval;
  51. gsl_matrix *evec;
  52. gsl_eigen_symmv_workspace *wmat;
  53. } gsl_sf_mathieu_workspace;
  54. /* Compute an array of characteristic (eigen) values from the recurrence
  55. matrices for the Mathieu equations. */
  56. int gsl_sf_mathieu_a_array(int order_min, int order_max, double qq, gsl_sf_mathieu_workspace *work, double result_array[]);
  57. int gsl_sf_mathieu_b_array(int order_min, int order_max, double qq, gsl_sf_mathieu_workspace *work, double result_array[]);
  58. /* Compute the characteristic value for a Mathieu function of order n and
  59. type ntype. */
  60. int gsl_sf_mathieu_a(int order, double qq, gsl_sf_result *result);
  61. int gsl_sf_mathieu_b(int order, double qq, gsl_sf_result *result);
  62. /* Compute the Fourier coefficients for a Mathieu function. */
  63. int gsl_sf_mathieu_a_coeff(int order, double qq, double aa, double coeff[]);
  64. int gsl_sf_mathieu_b_coeff(int order, double qq, double aa, double coeff[]);
  65. /* Allocate computational storage space for eigenvalue solution. */
  66. gsl_sf_mathieu_workspace *gsl_sf_mathieu_alloc(const size_t nn,
  67. const double qq);
  68. void gsl_sf_mathieu_free(gsl_sf_mathieu_workspace *workspace);
  69. /* Compute an angular Mathieu function. */
  70. int gsl_sf_mathieu_ce(int order, double qq, double zz, gsl_sf_result *result);
  71. int gsl_sf_mathieu_se(int order, double qq, double zz, gsl_sf_result *result);
  72. int gsl_sf_mathieu_ce_array(int nmin, int nmax, double qq, double zz,
  73. gsl_sf_mathieu_workspace *work,
  74. double result_array[]);
  75. int gsl_sf_mathieu_se_array(int nmin, int nmax, double qq, double zz,
  76. gsl_sf_mathieu_workspace *work,
  77. double result_array[]);
  78. /* Compute a radial Mathieu function. */
  79. int gsl_sf_mathieu_Mc(int kind, int order, double qq, double zz,
  80. gsl_sf_result *result);
  81. int gsl_sf_mathieu_Ms(int kind, int order, double qq, double zz,
  82. gsl_sf_result *result);
  83. int gsl_sf_mathieu_Mc_array(int kind, int nmin, int nmax, double qq,
  84. double zz, gsl_sf_mathieu_workspace *work,
  85. double result_array[]);
  86. int gsl_sf_mathieu_Ms_array(int kind, int nmin, int nmax, double qq,
  87. double zz, gsl_sf_mathieu_workspace *work,
  88. double result_array[]);
  89. __END_DECLS
  90. #endif /* !__GSL_SF_MATHIEU_H__ */