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