Expose gauss_legendre_tbl publicly as gauss_legendre_roots

This commit is contained in:
Dag Sverre Seljebotn 2015-04-24 09:06:20 +02:00
parent 765831ea2b
commit 7ecd1ddc93
7 changed files with 152 additions and 65 deletions

View file

@ -32,6 +32,7 @@
#include <math.h>
#include "sharp_geomhelpers.h"
#include "sharp_legendre_roots.h"
#include "c_utils.h"
#include "ls_fft.h"
#include <stdio.h>
@ -106,69 +107,6 @@ void sharp_make_weighted_healpix_geom_info (int nside, int stride,
sharp_make_subset_healpix_geom_info(nside, stride, 4 * nside - 1, NULL, weight, geom_info);
}
static inline double one_minus_x2 (double x)
{ return (fabs(x)>0.1) ? (1.+x)*(1.-x) : 1.-x*x; }
/* Function adapted from GNU GSL file glfixed.c
Original author: Pavel Holoborodko (http://www.holoborodko.com)
Adjustments by M. Reinecke
- adjusted interface (keep epsilon internal, return full number of points)
- removed precomputed tables
- tweaked Newton iteration to obtain higher accuracy */
static void gauss_legendre_tbl(int n, double *x, double *w)
{
const double pi = 3.141592653589793238462643383279502884197;
const double eps = 3e-14;
int m = (n+1)>>1;
double t0 = 1 - (1-1./n) / (8.*n*n);
double t1 = 1./(4.*n+2.);
#pragma omp parallel
{
int i;
#pragma omp for schedule(dynamic,100)
for (i=1; i<=m; ++i)
{
double x0 = cos(pi * ((i<<2)-1) * t1) * t0;
int dobreak=0;
int j=0;
double dpdx;
while(1)
{
double P_1 = 1.0;
double P0 = x0;
double dx, x1;
for (int k=2; k<=n; k++)
{
double P_2 = P_1;
P_1 = P0;
// P0 = ((2*k-1)*x0*P_1-(k-1)*P_2)/k;
P0 = x0*P_1 + (k-1.)/k * (x0*P_1-P_2);
}
dpdx = (P_1 - x0*P0) * n / one_minus_x2(x0);
/* Newton step */
x1 = x0 - P0/dpdx;
dx = x0-x1;
x0 = x1;
if (dobreak) break;
if (fabs(dx)<=eps) dobreak=1;
UTIL_ASSERT(++j<100,"convergence problem");
}
x[i-1] = -x0;
x[n-i] = x0;
w[i-1] = w[n-i] = 2. / (one_minus_x2(x0) * dpdx * dpdx);
}
} // end of parallel region
}
void sharp_make_gauss_geom_info (int nrings, int nphi, double phi0,
int stride_lon, int stride_lat, sharp_geom_info **geom_info)
{
@ -181,7 +119,7 @@ void sharp_make_gauss_geom_info (int nrings, int nphi, double phi0,
ptrdiff_t *ofs=RALLOC(ptrdiff_t,nrings);
int *stride_=RALLOC(int,nrings);
gauss_legendre_tbl(nrings,theta,weight);
sharp_legendre_roots(nrings,theta,weight);
for (int m=0; m<nrings; ++m)
{
theta[m] = acos(-theta[m]);