MODgvec_Basis1D Module

Module ** Basis 1D **

Routines to provide and evaluate 1D polynomial Lagrange basis functions, interpolation and integration points


Uses

  • module~~modgvec_basis1d~~UsesGraph module~modgvec_basis1d MODgvec_Basis1D module~modgvec_globals MODgvec_Globals module~modgvec_basis1d->module~modgvec_globals iso_fortran_env iso_fortran_env module~modgvec_globals->iso_fortran_env

Used by

  • module~~modgvec_basis1d~~UsedByGraph module~modgvec_basis1d MODgvec_Basis1D proc~sbase_eval t_sBase%sBase_eval proc~sbase_eval->module~modgvec_basis1d proc~sbase_init t_sBase%sBase_init proc~sbase_init->module~modgvec_basis1d

Variables

Type Visibility Attributes Name Initial
real(kind=wp), private, parameter :: PP_RealTolerance = EPSILON(1.0_wp)

machine precision

real(kind=wp), private, parameter :: PP_Pi = ACOS(-1.0_wp)

Pi up to machine accuracy


Interfaces

private interface BuildLegendreVdm

private interface InitializeVandermonde

private interface ChebyshevGaussNodesAndWeights

private interface ChebyGaussLobNodesAndWeights

private interface ClenshawCurtisNodesAndWeights

private interface LegendreGaussNodesAndWeights

private interface LegGaussLobNodesAndWeights

private interface PolynomialDerivativeMatrix

private interface MthPolynomialDerivativeMatrix

private interface BarycentricWeights

private interface LagrangeInterpolationPolys

private interface EQUALTOTOLERANCE


Functions

private function ALMOSTEQUAL(x, y) result(AlmostEqual)

Determines if two REAL(wp) numbers are equal up to a specified tolerance (=PP_RealTolerance, normaly set to machine precision) Takes into account that x,y are located in-between [-1;1] Based on Algorithm 139, Kopriva

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: x

(IN) first scalar to be compared

real(kind=wp), intent(in) :: y

(IN) second scalar to be compared

Return Value logical

(OUT) TRUE if |x-y| < 2*PP_RealTolerance

public function EQUALTOTOLERANCE(x, y, tolerance) result(EqualToTolerance)

Determines if two REAL(wp) numbers are equal up to a given tolerance. Routine requires: x,y > tolerance

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: x

(IN) first scalar to be compared

real(kind=wp), intent(in) :: y

(IN) second scalar to be compared

real(kind=wp), intent(in) :: tolerance

(IN) Tolerance to be checked against

Return Value logical

(OUT) TRUE if x and y are closer than tolerance


Subroutines

public subroutine buildLegendreVdm(N_In, xi_In, Vdm_Leg, sVdm_Leg)

Build a 1D Vandermonde matrix from an orthonormal Legendre basis to a nodal basis and reverse

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_In

input polynomial degree

real(kind=wp), intent(in) :: xi_In(0:N_In)

nodal positions [-1,1]

real(kind=wp), intent(out) :: Vdm_Leg(0:N_In,0:N_In)

Vandermonde from Legendre to nodal basis

real(kind=wp), intent(out) :: sVdm_Leg(0:N_In,0:N_In)

Vandermonde from nodal basis to Legendre

public subroutine InitializeVandermonde(N_In, N_Out, wBary_In, xi_In, xi_Out, Vdm)

Build a 1D Vandermonde matrix using the Lagrange basis functions of degree N_In, evaluated at the interpolation points xi_Out

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_In

(IN) input polynomial degree

integer, intent(in) :: N_Out

(IN) output polynomial degree

real(kind=wp), intent(in) :: wBary_In(0:N_In)

(IN) input interpolation weights

real(kind=wp), intent(in) :: xi_In(0:N_In)

(IN) input nodal positions [-1,1]

real(kind=wp), intent(in) :: xi_Out(0:N_Out)

(IN) outout nodal positions [-1,1]

real(kind=wp), intent(out) :: Vdm(0:N_Out,0:N_In)

(OUT) nodal Vandermonde from N_In to N_out

public subroutine LegendrePolynomialAndDerivative(N_in, x, L, Lder)

Evaluate the Legendre polynomial L_N and its derivative at position x[-1,1] recursive algorithm using the N_in-1 N_in-2 Legendre polynomials algorithm 22, Kopriva book

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_in

(IN) polynomial degree, (N+1) CLpoints

real(kind=wp), intent(in) :: x

(IN) coordinate value in the interval [-1,1]

real(kind=wp), intent(out) :: L

(OUT) Legedre polynomial evaluated at \f$ \xi: L_N(\xi), \partial/\partial\xi L_N(\xi) \f$

real(kind=wp), intent(out) :: Lder

(OUT) Legedre polynomial deriv. evaluated at \f$ \xi: L_N(\xi), \partial/\partial\xi L_N(\xi) \f$

public subroutine ChebyshevGaussNodesAndWeights(N_in, xGP, wGP)

Compute Chebychev-Gauss nodes and integration weights (algorithm 27, Kopriva book)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_in

polynomial degree, (N_in+1) CLpoints

real(kind=wp), intent(out) :: xGP(0:N_in)

Gauss point positions for the reference interval [-1,1]

real(kind=wp), intent(out), optional :: wGP(0:N_in)

Gauss point integration weights

public subroutine ChebyGaussLobNodesAndWeights(N_in, xGP, wGP)

Compute Chebychev-Gauss-Lobatto nodes and integration weights (algorithm 27, Kopriva book)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_in

polynomial degree, (N_in+1) CLpoints

real(kind=wp), intent(out) :: xGP(0:N_in)

Gauss point positions for the reference interval [-1,1]

real(kind=wp), intent(out), optional :: wGP(0:N_in)

Gauss point weights

public subroutine ClenshawCurtisNodesAndWeights(N_in, xGP, wGP)

Compute Clenshaw-Curtis nodes and integration weights

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_in

polynomial degree, (N_in+1) CLpoints

real(kind=wp), intent(out) :: xGP(0:N_in)

Gauss point positions for the reference interval [-1,1]

real(kind=wp), intent(out), optional :: wGP(0:N_in)

Gauss point weights

public subroutine LegendreGaussNodesAndWeights(N_in, xGP, wGP)

Compute Legendre-Gauss nodes and integration weights (algorithm 23, Kopriva book)

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_in

polynomial degree, (N_in+1) Gausspoints

real(kind=wp), intent(out) :: xGP(0:N_in)

Gauss point positions for the reference interval [-1,1]

real(kind=wp), intent(out), optional :: wGP(0:N_in)

Gauss point weights

private subroutine qAndLEvaluation(N_in, x, q, qder, L)

Evaluate the polynomial q=L_{N_in+1}-L_{N_in-1} and its derivative at position x in [-1,1] Recursive algorithm using the N_in-1 N_in-2 Legendre polynomials. (Algorithm 24, Kopriva book)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_in

polynomial degree

real(kind=wp), intent(in) :: x

coordinate value in the interval [-1,1]

real(kind=wp), intent(out) :: q

\f$ q_N(\xi) \f$

real(kind=wp), intent(out) :: qder

\f$ \partial/\partial\xi \; L_N(\xi) \f$

real(kind=wp), intent(out) :: L

\f$ L_N(\xi) \f$

public subroutine LegGaussLobNodesAndWeights(N_in, xGP, wGP)

Starting with initial guess by Parter Relation, a Newton method is used to find the roots of the Legendre Polynomial Lder_(N_in), which are the positions of Gauss-Lobatto points. Uses qAndLEvaluation subroutine. algorithm 25, Kopriva

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_in

polynomial degree (N_in+1) Gausspoints

real(kind=wp), intent(out) :: xGP(0:N_in)

Gauss point positions for the reference interval [-1,1]

real(kind=wp), intent(out), optional :: wGP(0:N_in)

Gauss point weights

public subroutine BarycentricWeights(N_in, xGP, wBary)

Computes barycentric (interpolation) weights for interpolation polynomial given by set of nodes. (Algorithm 30, Kopriva book)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_in

polynomial degree

real(kind=wp), intent(in) :: xGP(0:N_in)

Gauss point positions for the reference interval [-1,1]

real(kind=wp), intent(out) :: wBary(0:N_in)

barycentric weights

public subroutine PolynomialDerivativeMatrix(N_in, xGP, D)

Computes polynomial differentiation matrix for interpolation polynomial given by set of nodes. (Algorithm 37, Kopriva book)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_in

polynomial degree

real(kind=wp), intent(in) :: xGP(0:N_in)

Gauss point positions for the reference interval [-1,1]

real(kind=wp), intent(out) :: D(0:N_in,0:N_in)

differentiation Matrix

public subroutine MthPolynomialDerivativeMatrix(N_in, xGP, deriv, D)

Computes mth polynomial differentiation matrix for interpolation polynomial given by set of nodes. (Algorithm 38, Kopriva book)

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_in

polynomial degree

real(kind=wp), intent(in) :: xGP(0:N_in)

Gauss point positions for the reference interval [-1,1]

integer, intent(in) :: deriv

derivative (starting at 1)

real(kind=wp), intent(out) :: D(0:N_in,0:N_in)

differentiation Matrix

public subroutine LagrangeInterpolationPolys(x, N_in, xGP, wBary, L)

Computes all Lagrange functions evaluated at position x in [-1,1] For details see paper Barycentric Lagrange Interpolation by Berrut and Trefethen (SIAM 2004) Uses function ALMOSTEQUAL Algorithm 34, Kopriva book

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in) :: x

Coordinate

integer, intent(in) :: N_in

polynomial degree

real(kind=wp), intent(in) :: xGP(0:N_in)

Gauss point positions for the reference interval [-1,1]

real(kind=wp), intent(in) :: wBary(0:N_in)

Barycentric weights

real(kind=wp), intent(out) :: L(0:N_in)

Lagrange basis functions evaluated at x

private subroutine GaussRadauNodesAndWeights(N_in, xGR, wGR)

This routine was taken fom QUADRULE (http://people.sc.fsu.edu/~jburkardt/f_src/quadrule/quadrule.html)

Read more…

Arguments

Type IntentOptional Attributes Name
integer, intent(in) :: N_in

polynomial degree (N_in+1) Gausspoints

real(kind=wp), intent(out) :: xGR(0:N_in)

Gauss point positions for the reference interval [-1,1]

real(kind=wp), intent(out), optional :: wGR(0:N_in)

Gauss point weights