Main Page   File List   File Members  

polylib.c File Reference

#include <stdlib.h>
#include <sys/types.h>
#include <stdio.h>
#include <float.h>
#include <math.h>
#include "polylib.h"

Go to the source code of this file.

Defines

#define STOP   30
#define EPS   100*DBL_EPSILON
#define sign(a, b)   ((b)<0 ? -fabs(a) : fabs(a))
#define jacobz(n, z, alpha, beta)   JacZeros(n,z,alpha,beta)
 zero determination using eigenvalues of tridiagaonl matrix


Functions

void Jacobz (int n, double *z, double alpha, double beta)
 Calculate the n zeros, z, of the Jacobi polynomial, i.e. .

void JacZeros (int n, double *a, double alpha, double beta)
 Zero determination through the eigenvalues of a tridiagonal matrix from teh three term recursion relationship.

void TriQL (int n, double *d, double *e)
 QL algorithm for symmetric tridiagonal matrix.

double gammaF (double)
 Calculate the Gamma function , , for integer values and halves.

void zwgj (double *z, double *w, int np, double alpha, double beta)
 Gauss-Jacobi zeros and weights.

void zwgrjm (double *z, double *w, int np, double alpha, double beta)
 Gauss-Radau-Jacobi zeros and weights with end point at z=-1.

void zwgrjp (double *z, double *w, int np, double alpha, double beta)
 Gauss-Radau-Jacobi zeros and weights with end point at z=1.

void zwglj (double *z, double *w, int np, double alpha, double beta)
 Gauss-Lobatto-Jacobi zeros and weights with end point at z=-1,1.

void Dgj (double *D, double *Dt, double *z, int np, double alpha, double beta)
 Compute the Derivative Matrix and its transpose associated with the Gauss-Jacobi zeros.

void Dgrjm (double *D, double *Dt, double *z, int np, double alpha, double beta)
 Compute the Derivative Matrix and its transpose associated with the Gauss-Radau-Jacobi zeros with a zero at z=-1.

void Dgrjp (double *D, double *Dt, double *z, int np, double alpha, double beta)
 Compute the Derivative Matrix and its transpose associated with the Gauss-Radau-Jacobi zeros with a zero at z=1.

void Dglj (double *D, double *Dt, double *z, int np, double alpha, double beta)
 Compute the Derivative Matrix and its transpose associated with the Gauss-Lobatto-Jacobi zeros.

double hgj (int i, double z, double *zgj, int np, double alpha, double beta)
 Compute the value of the i th Lagrangian interpolant through the np Gauss-Jacobi points zgj at the arbitrary location z.

double hgrjm (int i, double z, double *zgrj, int np, double alpha, double beta)
 Compute the value of the i th Lagrangian interpolant through the np Gauss-Radau-Jacobi points zgrj at the arbitrary location z. This routine assumes zgrj includes the point -1.

double hgrjp (int i, double z, double *zgrj, int np, double alpha, double beta)
 Compute the value of the i th Lagrangian interpolant through the np Gauss-Radau-Jacobi points zgrj at the arbitrary location z. This routine assumes zgrj includes the point +1.

double hglj (int i, double z, double *zglj, int np, double alpha, double beta)
 Compute the value of the i th Lagrangian interpolant through the np Gauss-Lobatto-Jacobi points zgrj at the arbitrary location z.

void Imgj (double *im, double *zgj, double *zm, int nz, int mz, double alpha, double beta)
 Interpolation Operator from Gauss-Jacobi points to an arbitrary distrubtion at points zm.

void Imgrjm (double *im, double *zgrj, double *zm, int nz, int mz, double alpha, double beta)
 Interpolation Operator from Gauss-Radau-Jacobi points (including z=-1) to an arbitrary distrubtion at points zm.

void Imgrjp (double *im, double *zgrj, double *zm, int nz, int mz, double alpha, double beta)
 Interpolation Operator from Gauss-Radau-Jacobi points (including z=1) to an arbitrary distrubtion at points zm.

void Imglj (double *im, double *zglj, double *zm, int nz, int mz, double alpha, double beta)
 Interpolation Operator from Gauss-Lobatto-Jacobi points to an arbitrary distrubtion at points zm.

void jacobfd (int np, double *z, double *poly_in, double *polyd, int n, double alpha, double beta)
 Routine to calculate Jacobi polynomials, , and their first derivative, .

void jacobd (int np, double *z, double *polyd, int n, double alpha, double beta)
 Calculate the derivative of Jacobi polynomials.


Define Documentation

#define EPS   100*DBL_EPSILON
 

Definition at line 9 of file polylib.c.

Referenced by hgj(), hglj(), hgrjm(), hgrjp(), and Jacobz().

#define jacobz n,
z,
alpha,
beta   )     JacZeros(n,z,alpha,beta)
 

zero determination using eigenvalues of tridiagaonl matrix

Definition at line 164 of file polylib.c.

Referenced by zwgj(), zwglj(), zwgrjm(), and zwgrjp().

#define sign a,
 )     ((b)<0 ? -fabs(a) : fabs(a))
 

Definition at line 10 of file polylib.c.

Referenced by TriQL().

#define STOP   30
 

Definition at line 8 of file polylib.c.

Referenced by Jacobz(), and TriQL().


Function Documentation

void Dgj double *  D,
double *  Dt,
double *  z,
int  np,
double  alpha,
double  beta
 

Compute the Derivative Matrix and its transpose associated with the Gauss-Jacobi zeros.

  • Compute the derivative matrix, d, and its transpose, dt, associated with the n_th order Lagrangian interpolants through the np Gauss-Jacobi points z such that
  • d and dt are both square matrices.

Definition at line 321 of file polylib.c.

References jacobd().

void Dglj double *  D,
double *  Dt,
double *  z,
int  np,
double  alpha,
double  beta
 

Compute the Derivative Matrix and its transpose associated with the Gauss-Lobatto-Jacobi zeros.

  • Compute the derivative matrix, d, and its transpose, dt, associated with the n_th order Lagrangian interpolants through the np Gauss-Lobatto-Jacobi points z such that
  • d and dt are both square matrices.

Definition at line 469 of file polylib.c.

References gammaF(), and jacobd().

void Dgrjm double *  D,
double *  Dt,
double *  z,
int  np,
double  alpha,
double  beta
 

Compute the Derivative Matrix and its transpose associated with the Gauss-Radau-Jacobi zeros with a zero at z=-1.

  • Compute the derivative matrix, d, and its transpose, dt, associated with the n_th order Lagrangian interpolants through the np Gauss-Radau-Jacobi points z such that
  • d and dt are both square matrices.

Definition at line 365 of file polylib.c.

References gammaF(), and jacobd().

void Dgrjp double *  D,
double *  Dt,
double *  z,
int  np,
double  alpha,
double  beta
 

Compute the Derivative Matrix and its transpose associated with the Gauss-Radau-Jacobi zeros with a zero at z=1.

  • Compute the derivative matrix, d, and its transpose, dt, associated with the n_th order Lagrangian interpolants through the np Gauss-Radau-Jacobi points z such that
  • d and dt are both square matrices.

Definition at line 417 of file polylib.c.

References gammaF(), and jacobd().

double gammaF double  x  )  [static]
 

Calculate the Gamma function , , for integer values and halves.

Determine the value of using:

where

Definition at line 947 of file polylib.c.

Referenced by Dglj(), Dgrjm(), Dgrjp(), JacZeros(), zwgj(), zwglj(), zwgrjm(), and zwgrjp().

double hgj int  i,
double  z,
double *  zgj,
int  np,
double  alpha,
double  beta
 

Compute the value of the i th Lagrangian interpolant through the np Gauss-Jacobi points zgj at the arbitrary location z.

  • Uses the defintion of the Lagrangian interpolant:

Definition at line 532 of file polylib.c.

References EPS, jacobd(), and jacobfd().

Referenced by Imgj().

double hglj int  i,
double  z,
double *  zglj,
int  np,
double  alpha,
double  beta
 

Compute the value of the i th Lagrangian interpolant through the np Gauss-Lobatto-Jacobi points zgrj at the arbitrary location z.

  • Uses the defintion of the Lagrangian interpolant:

Definition at line 651 of file polylib.c.

References EPS, jacobd(), and jacobfd().

Referenced by Imglj().

double hgrjm int  i,
double  z,
double *  zgrj,
int  np,
double  alpha,
double  beta
 

Compute the value of the i th Lagrangian interpolant through the np Gauss-Radau-Jacobi points zgrj at the arbitrary location z. This routine assumes zgrj includes the point -1.

  • Uses the defintion of the Lagrangian interpolant:

Definition at line 569 of file polylib.c.

References EPS, jacobd(), and jacobfd().

Referenced by Imgrjm().

double hgrjp int  i,
double  z,
double *  zgrj,
int  np,
double  alpha,
double  beta
 

Compute the value of the i th Lagrangian interpolant through the np Gauss-Radau-Jacobi points zgrj at the arbitrary location z. This routine assumes zgrj includes the point +1.

  • Uses the defintion of the Lagrangian interpolant:

Definition at line 610 of file polylib.c.

References EPS, jacobd(), and jacobfd().

Referenced by Imgrjp().

void Imgj double *  im,
double *  zgj,
double *  zm,
int  nz,
int  mz,
double  alpha,
double  beta
 

Interpolation Operator from Gauss-Jacobi points to an arbitrary distrubtion at points zm.

  • Computes the one-dimensional interpolation matrix, im, to interpolate a function from at Gauss-Jacobi distribution of nz zeros zgrj to an arbitrary distribution of mz points zm, i.e.

Definition at line 684 of file polylib.c.

References hgj().

void Imglj double *  im,
double *  zglj,
double *  zm,
int  nz,
int  mz,
double  alpha,
double  beta
 

Interpolation Operator from Gauss-Lobatto-Jacobi points to an arbitrary distrubtion at points zm.

  • Computes the one-dimensional interpolation matrix, im, to interpolate a function from at Gauss-Lobatto-Jacobi distribution of nz zeros zgrj (where zgrj[0]=-1) to an arbitrary distribution of mz points zm, i.e.

Definition at line 766 of file polylib.c.

References hglj().

void Imgrjm double *  im,
double *  zgrj,
double *  zm,
int  nz,
int  mz,
double  alpha,
double  beta
 

Interpolation Operator from Gauss-Radau-Jacobi points (including z=-1) to an arbitrary distrubtion at points zm.

  • Computes the one-dimensional interpolation matrix, im, to interpolate a function from at Gauss-Radau-Jacobi distribution of nz zeros zgrj (where zgrj[0]=-1) to an arbitrary distribution of mz points zm, i.e.

Definition at line 711 of file polylib.c.

References hgrjm().

void Imgrjp double *  im,
double *  zgrj,
double *  zm,
int  nz,
int  mz,
double  alpha,
double  beta
 

Interpolation Operator from Gauss-Radau-Jacobi points (including z=1) to an arbitrary distrubtion at points zm.

  • Computes the one-dimensional interpolation matrix, im, to interpolate a function from at Gauss-Radau-Jacobi distribution of nz zeros zgrj (where zgrj[nz-1]=1) to an arbitrary distribution of mz points zm, i.e.

Definition at line 738 of file polylib.c.

References hgrjp().

void jacobd int  np,
double *  z,
double *  polyd,
int  n,
double  alpha,
double  beta
 

Calculate the derivative of Jacobi polynomials.

  • Generates a vector poly of values of the derivative of the n th order Jacobi polynomial at the np points z.
  • To do this we have used the relation
  • This formulation is valid for

Definition at line 921 of file polylib.c.

References jacobfd().

Referenced by Dgj(), Dglj(), Dgrjm(), Dgrjp(), hgj(), hglj(), hgrjm(), hgrjp(), and zwgj().

void jacobfd int  np,
double *  z,
double *  poly_in,
double *  polyd,
int  n,
double  alpha,
double  beta
 

Routine to calculate Jacobi polynomials, , and their first derivative, .

  • This function returns the vectors poly_in and poly_d containing the value of the order Jacobi polynomial and its derivative at the np points in z[i]
  • If poly_in = NULL then only calculate derivatice

  • If polyd = NULL then only calculate polynomial

  • To calculate the polynomial this routine uses the recursion relationship (see appendix A ref [4]) :

  • To calculate the derivative of the polynomial this routine uses the relationship (see appendix A ref [4]) :

  • Note the derivative from this routine is only valid for -1 < z < 1.

Definition at line 821 of file polylib.c.

Referenced by hgj(), hglj(), hgrjm(), hgrjp(), jacobd(), Jacobz(), zwglj(), zwgrjm(), and zwgrjp().

void Jacobz int  n,
double *  z,
double  alpha,
double  beta
[static]
 

Calculate the n zeros, z, of the Jacobi polynomial, i.e. .

This routine is only value for and uses polynomial deflation in a Newton iteration

Definition at line 984 of file polylib.c.

References EPS, jacobfd(), and STOP.

void JacZeros int  n,
double *  a,
double  alpha,
double  beta
[static]
 

Zero determination through the eigenvalues of a tridiagonal matrix from teh three term recursion relationship.

Set up a symmetric tridiagonal matrix

Where the coefficients a[n], b[n] come from the recurrence relation

where and are the Jacobi (normalized) orthogonal polynomials ( integer values and halves). Since the polynomials are orthonormalized, the tridiagonal matrix is guaranteed to be symmetric. The eigenvalues of this matrix are the zeros of the Jacobi polynomial.

Definition at line 1038 of file polylib.c.

References gammaF(), and TriQL().

void TriQL int  n,
double *  d,
double *  e
[static]
 

QL algorithm for symmetric tridiagonal matrix.

This subroutine is a translation of an algol procedure, num. math. 12, 377-383(1968) by martin and wilkinson, as modified in num. math. 15, 450(1970) by dubrulle. Handbook for auto. comp., vol.ii-linear algebra, 241-248(1971). This is a modified version from numerical recipes.

This subroutine finds the eigenvalues and first components of the eigenvectors of a symmetric tridiagonal matrix by the implicit QL method.

on input:

  • n is the order of the matrix;
  • d contains the diagonal elements of the input matrix;
  • e contains the subdiagonal elements of the input matrix in its first n-1 positions. e(n) is arbitrary;

on output:

  • d contains the eigenvalues in ascending order.
  • e has been destroyed;

Definition at line 1096 of file polylib.c.

References sign, and STOP.

Referenced by JacZeros().

void zwgj double *  z,
double *  w,
int  np,
double  alpha,
double  beta
 

Gauss-Jacobi zeros and weights.

  • Generate np Gauss Jacobi zeros, z, and weights,w, associated with the Jacobi polynomial ,
  • Exact for polynomials of order 2np-1 or less

Definition at line 186 of file polylib.c.

References gammaF(), jacobd(), and jacobz.

void zwglj double *  z,
double *  w,
int  np,
double  alpha,
double  beta
 

Gauss-Lobatto-Jacobi zeros and weights with end point at z=-1,1.

  • Generate np Gauss-Lobatto-Jacobi points, z, and weights, w, associated with polynomial
  • Exact for polynomials of order 2np-3 or less

Definition at line 282 of file polylib.c.

References gammaF(), jacobfd(), and jacobz.

void zwgrjm double *  z,
double *  w,
int  np,
double  alpha,
double  beta
 

Gauss-Radau-Jacobi zeros and weights with end point at z=-1.

  • Generate np Gauss-Radau-Jacobi zeros, z, and weights,w, associated with the polynomial .
  • Exact for polynomials of order 2np-2 or less

Definition at line 212 of file polylib.c.

References gammaF(), jacobfd(), and jacobz.

void zwgrjp double *  z,
double *  w,
int  np,
double  alpha,
double  beta
 

Gauss-Radau-Jacobi zeros and weights with end point at z=1.

  • Generate np Gauss-Radau-Jacobi zeros, z, and weights,w, associated with the polynomial .
  • Exact for polynomials of order 2np-2 or less

Definition at line 248 of file polylib.c.

References gammaF(), jacobfd(), and jacobz.


Generated on Sun Jul 13 22:18:36 2003 for Polylib by doxygen1.3