GaussJacobiQuadCCompat.f90
Go to the documentation of this file.
1 ! BEGIN_HEADER
2 ! -----------------------------------------------------------------------------
3 ! Gauss-Jacobi Quadrature Implementation
4 ! Authors: Rohit Goswami <rgoswami[at]ieee.org>
5 ! Source: GaussJacobiQuad Library
6 ! License: MIT
7 ! GitHub Repository: https://github.com/HaoZeke/GaussJacobiQuad
8 ! Date: 2023-08-28
9 ! Commit: c442f77
10 ! -----------------------------------------------------------------------------
11 ! This code is part of the GaussJacobiQuad library, providing an efficient
12 ! implementation for Gauss-Jacobi quadrature nodes and weights computation.
13 ! -----------------------------------------------------------------------------
14 ! To cite this software:
15 ! Rohit Goswami (2023). HaoZeke/GaussJacobiQuad: v0.1.0.
16 ! Zenodo: https://doi.org/10.5281/ZENODO.8285112
17 ! ---------------------------------------------------------------------
18 ! END_HEADER
20 !
21 !! @details This module provides a C-compatible interface to the Gauss-Jacobi quadrature methods
22 !! implemented in the GaussJacobiQuad module. It's designed to make the Fortran implementations
23 !! accessible from C/C++ code. It wraps the original Fortran subroutines with ones that are C-compatible,
24 !! using the `iso_c_binding` module for the necessary type conversions.
26 use iso_c_binding, only: c_double, c_int
27 use gaussjacobiquad, only: gauss_jacobi_rec, gauss_jacobi_gw, gauss_jacobi_algo665
28 implicit none
29 contains
30 
41 subroutine gauss_jacobi_rec_c(npts, alpha, beta, x, wts) bind(C, name="gauss_jacobi_rec_c")
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)
46 end subroutine gauss_jacobi_rec_c
47 
58 subroutine gauss_jacobi_gw_c(npts, alpha, beta, x, wts) bind(C, name="gauss_jacobi_gw_c")
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)
63 end subroutine gauss_jacobi_gw_c
64 
75 subroutine gauss_jacobi_algo665_c(npts, alpha, beta, x, wts) bind(C, name="gauss_jacobi_algo665_c")
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)
80 end subroutine gauss_jacobi_algo665_c
81 
82 end module gaussjacobiquadccompat
void gauss_jacobi_algo665_c(int *npts, double *alpha, double *beta, double x[], double wts[])
void gauss_jacobi_gw_c(int *npts, double *alpha, double *beta, double x[], double wts[])
void gauss_jacobi_rec_c(int *npts, double *alpha, double *beta, double x[], double wts[])
Overall driver for Gauss-Jacobi quadrature.
Compatibility layer for Gauss-Jacobi quadrature with C bindings.