gaussjacobiquadccompat Module Reference

Compatibility layer for Gauss-Jacobi quadrature with C bindings. More...

Functions/Subroutines

subroutine gauss_jacobi_rec_c (npts, alpha, beta, x, wts)
 C-compatible wrapper for gauss_jacobi_rec subroutine. More...
 
subroutine gauss_jacobi_gw_c (npts, alpha, beta, x, wts)
 C-compatible wrapper for gauss_jacobi_gw subroutine. More...
 
subroutine gauss_jacobi_algo665_c (npts, alpha, beta, x, wts)
 C-compatible wrapper for gauss_jacobi_algo665 subroutine. More...
 

Detailed Description

Compatibility layer for Gauss-Jacobi quadrature with C bindings.

Function/Subroutine Documentation

◆ gauss_jacobi_algo665_c()

subroutine gaussjacobiquadccompat::gauss_jacobi_algo665_c ( integer(c_int), intent(in)  npts,
real(c_double), intent(in)  alpha,
real(c_double), intent(in)  beta,
real(c_double), dimension(npts), intent(out)  x,
real(c_double), dimension(npts), intent(out)  wts 
)

C-compatible wrapper for gauss_jacobi_algo665 subroutine.

This subroutine is a C-compatible wrapper that calls the original gauss_jacobi_algo665 subroutine for calculating Gauss-Jacobi quadrature nodes and weights using the "algo665" method.

Parameters
[in]nptsNumber of quadrature points.
[in]alphaParameter alpha in the Jacobi polynomial. Must be greater than -1.
[in]betaParameter beta in the Jacobi polynomial. Must be greater than -1.
[out]xQuadrature nodes.
[out]wtsQuadrature weights.

Definition at line 75 of file GaussJacobiQuadCCompat.f90.

76  integer(c_int), intent(in) :: npts
77  real(c_double), intent(in) :: alpha, beta
78  real(c_double), intent(out) :: x(npts), wts(npts)
79  call gauss_jacobi_algo665(npts, alpha, beta, x, wts)

◆ gauss_jacobi_gw_c()

subroutine gaussjacobiquadccompat::gauss_jacobi_gw_c ( integer(c_int), intent(in)  npts,
real(c_double), intent(in)  alpha,
real(c_double), intent(in)  beta,
real(c_double), dimension(npts), intent(out)  x,
real(c_double), dimension(npts), intent(out)  wts 
)

C-compatible wrapper for gauss_jacobi_gw subroutine.

This subroutine is a C-compatible wrapper that calls the original gauss_jacobi_gw subroutine for calculating Gauss-Jacobi quadrature nodes and weights using the "gw" method.

Parameters
[in]nptsNumber of quadrature points.
[in]alphaParameter alpha in the Jacobi polynomial. Must be greater than -1.
[in]betaParameter beta in the Jacobi polynomial. Must be greater than -1.
[out]xQuadrature nodes.
[out]wtsQuadrature weights.

Definition at line 58 of file GaussJacobiQuadCCompat.f90.

59  integer(c_int), intent(in) :: npts
60  real(c_double), intent(in) :: alpha, beta
61  real(c_double), intent(out) :: x(npts), wts(npts)
62  call gauss_jacobi_gw(npts, alpha, beta, x, wts)

◆ gauss_jacobi_rec_c()

subroutine gaussjacobiquadccompat::gauss_jacobi_rec_c ( integer(c_int), intent(in)  npts,
real(c_double), intent(in)  alpha,
real(c_double), intent(in)  beta,
real(c_double), dimension(npts), intent(out)  x,
real(c_double), dimension(npts), intent(out)  wts 
)

C-compatible wrapper for gauss_jacobi_rec subroutine.

This subroutine is a C-compatible wrapper that calls the original gauss_jacobi_rec subroutine for calculating Gauss-Jacobi quadrature nodes and weights using the recurrence method.

Parameters
[in]nptsNumber of quadrature points.
[in]alphaParameter alpha in the Jacobi polynomial. Must be greater than -1.
[in]betaParameter beta in the Jacobi polynomial. Must be greater than -1.
[out]xQuadrature nodes.
[out]wtsQuadrature weights.

Definition at line 41 of file GaussJacobiQuadCCompat.f90.

42  integer(c_int), intent(in) :: npts
43  real(c_double), intent(in) :: alpha, beta
44  real(c_double), intent(out) :: x(npts), wts(npts)
45  call gauss_jacobi_rec(npts, alpha, beta, x, wts)