ROL
Functions
OrthogonalPolynomials.hpp File Reference
#include <iostream>
#include <iomanip>
#include <cmath>
#include "ROL_StdVector.hpp"
#include "Teuchos_LAPACK.hpp"
#include "LinearAlgebra.hpp"

Go to the source code of this file.

Functions

template<class Real >
void rec_jacobi (ROL::Ptr< Teuchos::LAPACK< int, Real > > lapack, const double alpha, const double beta, std::vector< Real > &a, std::vector< Real > &b)
 Generate the Jacobi polynomial recursion coeffcients \(a_k,b_k\). More...
 
template<class Real >
void vandermonde (const std::vector< Real > &a, const std::vector< Real > &b, const std::vector< Real > &x, std::vector< Real > &V)
 Construct the generalized Vandermonde matrix (in column stacked form) based upon the recurrence coefficients (a,b) on the grid x. More...
 
template<class Real >
void gauss (ROL::Ptr< Teuchos::LAPACK< int, Real > > lapack, const std::vector< Real > &a, const std::vector< Real > &b, std::vector< Real > &x, std::vector< Real > &w)
 Compute the Gauss quadrature nodes and weights for the polynomials generated by the recurrence coefficients. More...
 
template<class Real >
void rec_lobatto (ROL::Ptr< Teuchos::LAPACK< int, Real > > const lapack, const double xl1, const double xl2, std::vector< Real > &a, std::vector< Real > &b)
 Modify the given recurrence coefficients so that the set of zeros of the maximal order polynomial include the two prescribed points. More...
 

Function Documentation

template<class Real >
void rec_jacobi ( ROL::Ptr< Teuchos::LAPACK< int, Real > >  lapack,
const double  alpha,
const double  beta,
std::vector< Real > &  a,
std::vector< Real > &  b 
)

Generate the Jacobi polynomial recursion coeffcients \(a_k,b_k\).

The Jacobi polynomials satisfy the recurrence relation

\[ P^{(\alpha,\beta)}_{k+1}(x) = (x-a_k)P_{k}^{(\alpha,\beta)}(x) - b_k P_{k-1}^{(\alpha,\beta)}(x) \]

and form an orthogonal basis on \([-1,1]\) with respect to the weight function \(w(x)=(1-x)^\alpha(1+x)^\beta\).

Parameters
[in]lapackis a pointer to the Teuchos::LAPACK interface
[in]alphais a parameter that defines the weight function
[in]betais a parameter that defines the weight function
[out]apis a vector of recursion coefficients
[out]bpis a vector of recursion coefficients

Adapted from the MATLAB code by Dirk Laurie and Walter Gautschi http://www.cs.purdue.edu/archives/2002/wxg/codes/r_jacobi.m

Definition at line 29 of file OrthogonalPolynomials.hpp.

Referenced by NodalBasis< Real >::NodalBasis().

template<class Real >
void vandermonde ( const std::vector< Real > &  a,
const std::vector< Real > &  b,
const std::vector< Real > &  x,
std::vector< Real > &  V 
)

Construct the generalized Vandermonde matrix (in column stacked form) based upon the recurrence coefficients (a,b) on the grid x.

Parameters
[in]avector of recursion coefficients
[in]bvector of recursion coefficients
[in]xvector of quadrature nodes
[in]Vcolumn-stacked Vandermonde matrix

Definition at line 72 of file OrthogonalPolynomials.hpp.

template<class Real >
void gauss ( ROL::Ptr< Teuchos::LAPACK< int, Real > >  lapack,
const std::vector< Real > &  a,
const std::vector< Real > &  b,
std::vector< Real > &  x,
std::vector< Real > &  w 
)

Compute the Gauss quadrature nodes and weights for the polynomials generated by the recurrence coefficients.

Parameters
[in]lapackpointer to the Teuchos::LAPACK interface
[in]avector of recursion coefficients
[in]bvector of recursion coefficients
[out]xvector of quadrature nodes
[out]wvector of quadrature weights

Adapted from the MATLAB code by Walter Gautschi http://www.cs.purdue.edu/archives/2002/wxg/codes/gauss.m

Definition at line 104 of file OrthogonalPolynomials.hpp.

Referenced by NodalBasis< Real >::NodalBasis().

template<class Real >
void rec_lobatto ( ROL::Ptr< Teuchos::LAPACK< int, Real > > const  lapack,
const double  xl1,
const double  xl2,
std::vector< Real > &  a,
std::vector< Real > &  b 
)

Modify the given recurrence coefficients so that the set of zeros of the maximal order polynomial include the two prescribed points.

Parameters
[in]lapackpointer to the Teuchos::LAPACK interface
[in]xl1location of one pre-assigned node
[in]xl2location of another pre-assigned node
in/out]ap pointer to vector of recursion coefficients
in/out]bp pointer to vector of recursion coefficients

Based on the section 7 of the paper "Some modified matrix eigenvalue problems" by Gene Golub, SIAM Review Vol 15, No. 2, April 1973, pp.318–334

Definition at line 151 of file OrthogonalPolynomials.hpp.

References trisolve().

Referenced by NodalBasis< Real >::NodalBasis().