1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204
|
/* specfunc/mathieu_workspace.c
*
* Copyright (C) 2003 Lowell Johnson
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 3 of the License, or (at
* your option) any later version.
*
* This program is distributed in the hope that it will be useful, but
* WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 675 Mass Ave, Cambridge, MA 02139, USA.
*/
/* Author: L. Johnson */
#include <config.h>
#include <stdlib.h>
#include <gsl/gsl_math.h>
#include <gsl/gsl_errno.h>
#include <gsl/gsl_sf_mathieu.h>
gsl_sf_mathieu_workspace *gsl_sf_mathieu_alloc(const size_t nn,
const double qq)
{
gsl_sf_mathieu_workspace *workspace;
unsigned int even_order = nn/2 + 1, odd_order = (nn + 1)/2,
extra_values;
/* Compute the maximum number of extra terms required for 10^-18 root
accuracy for a given value of q (contributed by Brian Gladman). */
extra_values = (int)(2.1*pow(fabs(qq), 0.37)) + 9;
extra_values += 20; /* additional fudge */
if (nn + 1 == 0)
{
GSL_ERROR_NULL("matrix dimension must be positive integer", GSL_EINVAL);
}
workspace =
(gsl_sf_mathieu_workspace *)malloc(sizeof(gsl_sf_mathieu_workspace));
if (workspace == NULL)
{
GSL_ERROR_NULL("failed to allocate space for workspace", GSL_ENOMEM);
}
/* Extend matrices to ensure accuracy. */
even_order += extra_values;
odd_order += extra_values;
workspace->size = nn;
workspace->even_order = even_order;
workspace->odd_order = odd_order;
workspace->extra_values = extra_values;
/* Allocate space for the characteristic values. */
workspace->aa = (double *)malloc((nn+1)*sizeof(double));
if (workspace->aa == NULL)
{
free(workspace);
GSL_ERROR_NULL("Error allocating memory for characteristic a values",
GSL_ENOMEM);
}
workspace->bb = (double *)malloc((nn+1)*sizeof(double));
if (workspace->bb == NULL)
{
free(workspace->aa);
free(workspace);
GSL_ERROR_NULL("Error allocating memory for characteristic b values",
GSL_ENOMEM);
}
/* Since even_order is always >= odd_order, dimension the arrays for
even_order. */
workspace->dd = (double *)malloc(even_order*sizeof(double));
if (workspace->dd == NULL)
{
free(workspace->aa);
free(workspace->bb);
free(workspace);
GSL_ERROR_NULL("failed to allocate space for diagonal", GSL_ENOMEM);
}
workspace->ee = (double *)malloc(even_order*sizeof(double));
if (workspace->ee == NULL)
{
free(workspace->dd);
free(workspace->aa);
free(workspace->bb);
free(workspace);
GSL_ERROR_NULL("failed to allocate space for diagonal", GSL_ENOMEM);
}
workspace->tt = (double *)malloc(3*even_order*sizeof(double));
if (workspace->tt == NULL)
{
free(workspace->ee);
free(workspace->dd);
free(workspace->aa);
free(workspace->bb);
free(workspace);
GSL_ERROR_NULL("failed to allocate space for diagonal", GSL_ENOMEM);
}
workspace->e2 = (double *)malloc(even_order*sizeof(double));
if (workspace->e2 == NULL)
{
free(workspace->tt);
free(workspace->ee);
free(workspace->dd);
free(workspace->aa);
free(workspace->bb);
free(workspace);
GSL_ERROR_NULL("failed to allocate space for diagonal", GSL_ENOMEM);
}
workspace->zz = (double *)malloc(even_order*even_order*sizeof(double));
if (workspace->zz == NULL)
{
free(workspace->e2);
free(workspace->tt);
free(workspace->ee);
free(workspace->dd);
free(workspace->aa);
free(workspace->bb);
free(workspace);
GSL_ERROR_NULL("failed to allocate space for diagonal", GSL_ENOMEM);
}
workspace->eval = gsl_vector_alloc(even_order);
if (workspace->eval == NULL)
{
free(workspace->zz);
free(workspace->e2);
free(workspace->tt);
free(workspace->ee);
free(workspace->dd);
free(workspace->aa);
free(workspace->bb);
free(workspace);
GSL_ERROR_NULL("failed to allocate space for eval", GSL_ENOMEM);
}
workspace->evec = gsl_matrix_alloc(even_order, even_order);
if (workspace->evec == NULL)
{
gsl_vector_free (workspace->eval);
free(workspace->zz);
free(workspace->e2);
free(workspace->tt);
free(workspace->ee);
free(workspace->dd);
free(workspace->aa);
free(workspace->bb);
free(workspace);
GSL_ERROR_NULL("failed to allocate space for evec", GSL_ENOMEM);
}
workspace->wmat = gsl_eigen_symmv_alloc(even_order);
if (workspace->wmat == NULL)
{
gsl_matrix_free (workspace->evec);
gsl_vector_free (workspace->eval);
free(workspace->zz);
free(workspace->e2);
free(workspace->tt);
free(workspace->ee);
free(workspace->dd);
free(workspace->aa);
free(workspace->bb);
free(workspace);
GSL_ERROR_NULL("failed to allocate space for wmat", GSL_ENOMEM);
}
return workspace;
}
void gsl_sf_mathieu_free(gsl_sf_mathieu_workspace *workspace)
{
RETURN_IF_NULL (workspace);
gsl_vector_free(workspace->eval);
gsl_matrix_free(workspace->evec);
gsl_eigen_symmv_free(workspace->wmat);
free(workspace->aa);
free(workspace->bb);
free(workspace->dd);
free(workspace->ee);
free(workspace->tt);
free(workspace->e2);
free(workspace->zz);
free(workspace);
}
|