XLPack API Guide
What is XLPack API
XLPack API (C/C++ interface) is the interface for other programs than Excel to use numerical calculation module of XLPack. With this interface you will be able to do the followings:
- You can use the numerical calculation function module of XLPack from your C/C++ program
- The other programming languages which have the interface to C/C++ can also be used. The XLPack wrapper (interface) program examples for Python, Julia, C#, F#, VB.NET and Pascal (Delphi or Free Pascal) are provided
- LAPACKE/CBLAS library for XLPack API can be used for C/C++
- The program developed by using XLPack API can run on the personal computers in which XLPack (or XLPack Lite) is installed
Notes
XLPack API provides the internal interface in as-is condition. It has not been thoroughly tested under all conditions. Therefore, the function, performance, or reliability of this software are not guaranteed. And specifications are subject to change without notice. See more details in the document attached to XLPack SDK.
Downloading software
XLPack SDK can be downloaded from Download page.
SDK includes the necessary software to use API with programming languages C/C++, Python, Julia, .NET languages (C#, F#, VB.NET) and Pascal. Using interface programs included in SDK, it is possible to call XLPack from C/C++ and to call XLPack basic functions from other languages. Example programs are included, so that you can easily develop the application programs by referring those examples.
Reference manuals
Please refer to the following
Online manuals for C/C++, Python and Julia are available from the following links.
C/C++ interface
C/C++ numerical library reference manual (open in a separate tab)
C/C++ numerical library reference manual (basic functions) (open in a separate tab)
Python interface
Python numerical library reference manual (basic functions) (open in a separate tab)
Julia interface
Julia numerical library reference manual (basic functions) (open in a separate tab)
Using API
C/C++
C/C++ programs can call XLPack just by linking the following library included in SDK.
For Windows, link XLPack.lib (64 bit program) or XLPack_32.lib (32 bit program).
Header files cnumlib.h or cnumlib (for C++) can be used.
Python
The Python extension module which is used to call XLPack basic functions is included in SDK.
XLPack.pyd is used for 64 bit Python, and XLPack_32.pyd is used for 32 bit Python. These modules have to be placed in the same directory with calling script or the directory given by the variable PYTHONPATH. These modules are designed to work with Python 3.7 or later. However, depending on the software environments, these modules may not work well. In such case, please try XLPack.py (ctypes version which provides equivalent functions with extension module).
Numpy ndarray is used in the Python extension module to handle matrices. Numpy must be installed in advance.
Other languages
Other languages call DLL directly. The interface program examples are included in SDK.
The interface to call DLL directly is same with C/C++ library (see C/C++ numerical reference reference manual) except some functions described below.
The interface of functions which may set errno such as special functions and BLAS 2 and 3 functions is different from C/C++ library. x underbar (x_) is added in front of the name of function and errno is added as the last argument. For example,
func(double a1, double a2)
is changed as follows.
x_func(double a1, double a2, errno_t *xerrno)
In the case of normal exit, 0 will be returned in xerrno. When error occurred, EDOM or ERANGE (the values of errno.h of C language) will be returned. In BLAS 2 and 3 functions, errno is not set in most implementation. However, as an extension, EDOM will be set if the argument error is found.
The following functions may set errno.
List of functions which set errno
Special functions
factorial
ccbrt ccbrt_sub
laguerre alaguerre legendre legendred alegendre sharmonic sharmonic_sub
sharmonicr sharmonici hermite chebt chebtd chebu chebs
sqrt1pm1 powm1 sin_pi cos_pi cexpm1 cexpm1_sub clog1p clog1p_sub ccot ccot_sub
li ei e1 en spence
ci si chi shi
tgamma1pm1 lgammas rgamma tgammaratio tgammadratio cgamma cgamma_sub clgamma
clgamma_sub crgamma crgamma_sub poch poch1
beta lbeta cbeta cbeta_sub clbeta clbeta_sub
digamma trigamma polygamma cdigamma cdigamma_sub
gammai gammaic gammait gammap gammaq gammapi gammaqi gammapia gammaqia gammapd
betax betaxc ibeta ibetac ibetai ibetaci ibetaia ibetacia ibetaib ibetacib ibetad
zeta
erfi erfci dawson
besj0 besj1 besjn besy0 besy1 besyn besjnu besynu besjnd besjnud besynd besynud
sbesjn sbesyn sbesjnu sbesynu
besi0 besi0e besi1 besi1e besin besk0 besk0e besk1 besk1e beskn besinu besknu
besind besknd besinud besknud sbesin sbeskn sbesinu sbesknu
airyai airybi airyaid airybid
chu hyp1f1 lhyp1f1 hyp1f1r hyp2f1 hyp0f1 hyp1f0 hyp2f0
jelli jsn jcn jdn jns jnc jnd jsc jsd jdc jds jcs jcd
celli1 celli2 celli3 elli1 elli2 elli3 rc rd rg rf rj jzeta
BLAS 2
dgemv dgbmv dsymv dsbmv dspmv dtrmv dtbmv dtpmv dtrsv dtbsv dtpsv dger dsyr dspr
dsyr2 dspr2
zgemv zgbmv zhemv zhbmv zhpmv zsymv zsbmv zspmv ztrmv ztbmv ztpmv ztrsv ztbsv
ztpsv zgeru zgerc zher zhpr zsyr zspr zher2 zhpr2 zsyr2 zspr2
BLAS 3
dgemm dsymm dtrmm dtrsm dsyrk dsyr2k
zgemm zsymm zhemm ztrmm ztrsm zsyrk zherk zsyr2k zher2k
How to use LAPACKE/CBLAS
LAPACKE/CBLAS library for XLPack API for C/C++ is included in SDK. For Windows, link Lapacke.lib (for 64 bit XLPack) or Lapacke_32.lib (for 32 bit XLPack). Header files lapacke.h and cblas.h can be used.
Using these libraries and header files, you can call LAPACK/BLAS through LAPACKE/CBLAS interface (standard C interface included in LAPACK/BLAS). For the details of each function, please refer to the original site below.
Description of LAPACKE/CBLAS
LAPACK function reference manual including LAPACKE/CBLAS
The following LAPACKE and CBLAS routines corresponding to LAPACK and BLAS routines included in XLPack are supported.
List of supported LAPACKE/CBLAS functions
(*) shows the routines supported by XLPack basic functions.
cblas_drotg cblas_drotmg cblas_drot cblas_drotm cblas_dswap cblas_dscal
cblas_dcopy cblas_daxpy cblas_ddot cblas_dsdot cblas_dnrm2 cblas_dasum
cblas_idamax
cblas_dgemv cblas_dgbmv cblas_dger cblas_dsbmv cblas_dspmv cblas_dspr
cblas_dspr2 cblas_dsymv cblas_dsyr cblas_dsyr2 cblas_dtbmv cblas_dtbsv
cblas_dtpmv cblas_dtpsv cblas_dtrmv cblas_dtrsv
cblas_dgemm cblas_dsymm cblas_dsyrk cblas_dsyr2k cblas_dtrmm cblas_dtrsm
cblas_zswap cblas_zscal cblas_zdscal cblas_zcopy cblas_zaxpy cblas_zdotu_sub
cblas_zdotc_sub cblas_dznrm2 cblas_dzasum cblas_izamax
cblas_zgemv cblas_zgbmv cblas_zhemv cblas_zhbmv cblas_zhpmv cblas_ztrmv
cblas_ztbmv cblas_ztpmv cblas_ztrsv cblas_ztbsv cblas_ztpsv cblas_zgeru
cblas_zgerc cblas_zher cblas_zher2 cblas_zhpr cblas_zhpr2
cblas_zgemm cblas_zsymm cblas_zhemm cblas_zherk cblas_zher2k cblas_ztrmm
cblas_ztrsm cblas_zsyrk cblas_zsyr2k
LAPACKE_ddisna LAPACKE_dgbcon LAPACKE_dgbsv LAPACKE_dgbsvx LAPACKE_dgbtrf
LAPACKE_dgbtrs LAPACKE_dgebak LAPACKE_dgebal LAPACKE_dgecon(*) LAPACKE_dgees
LAPACKE_dgeesx LAPACKE_dgeev LAPACKE_dgeevx LAPACKE_dgejsv LAPACKE_dgelqf
LAPACKE_dgels(*) LAPACKE_dgelsd LAPACKE_dgelss LAPACKE_dgelsy LAPACKE_dgeqp3
LAPACKE_dgeqrf LAPACKE_dgesdd LAPACKE_dgesv(*) LAPACKE_dgesvd LAPACKE_dgesvdq
LAPACKE_dgesvdx LAPACKE_dgesvx LAPACKE_dgetrf LAPACKE_dgetri LAPACKE_dgetrs
LAPACKE_dgetsls LAPACKE_dgges LAPACKE_dggesx LAPACKE_dggev LAPACKE_dggevx
LAPACKE_dggglm LAPACKE_dgglse LAPACKE_dggsvd3 LAPACKE_dgtcon LAPACKE_dgtsv
LAPACKE_dgtsvx LAPACKE_dgttrf LAPACKE_dgttrs LAPACKE_dhsein LAPACKE_dhseqr
LAPACKE_dlamch(*) LAPACKE_dlange(*) LAPACKE_dlansy(*) LAPACKE_dlantr LAPACKE_dopgtr
LAPACKE_dopmtr LAPACKE_dorghr LAPACKE_dormhr LAPACKE_dorglq LAPACKE_dorgqr
LAPACKE_dorgtr LAPACKE_dormlq LAPACKE_dormqr LAPACKE_dormtr LAPACKE_dpbcon
LAPACKE_dpbsv LAPACKE_dpbsvx LAPACKE_dpbtrf LAPACKE_dpbtrs LAPACKE_dpocon(*)
LAPACKE_dposv(*) LAPACKE_dposvx LAPACKE_dpotrf LAPACKE_dpotri LAPACKE_dpotrs
LAPACKE_dppcon LAPACKE_dppsv LAPACKE_dppsvx LAPACKE_dpptrf LAPACKE_dpptri
LAPACKE_dpptrs LAPACKE_dptcon LAPACKE_dpteqr LAPACKE_dptsv LAPACKE_dptsvx
LAPACKE_dpttrf LAPACKE_dpttrs LAPACKE_dsbev LAPACKE_dsbevd LAPACKE_dsbevx
LAPACKE_dsbgv LAPACKE_dsbgvd LAPACKE_dsbgvx LAPACKE_dsbtrd LAPACKE_dsgesv
LAPACKE_dspcon LAPACKE_dspev LAPACKE_dspevd LAPACKE_dspevx LAPACKE_dspgv
LAPACKE_dspgvd LAPACKE_dspgvx LAPACKE_dsposv LAPACKE_dspsv LAPACKE_dspsvx
LAPACKE_dsptrd LAPACKE_dsptrf LAPACKE_dsptri LAPACKE_dsptrs LAPACKE_dstebz
LAPACKE_dstedc LAPACKE_dstein LAPACKE_dstemr LAPACKE_dsteqr LAPACKE_dsterf
LAPACKE_dstev LAPACKE_dstevd LAPACKE_dstevr LAPACKE_dstevx LAPACKE_dsycon
LAPACKE_dsyev(*) LAPACKE_dsyevd LAPACKE_dsyevr LAPACKE_dsyevx LAPACKE_dsygv
LAPACKE_dsygvd LAPACKE_dsygvx LAPACKE_dsysv LAPACKE_dsysvx LAPACKE_dsytrd
LAPACKE_dsytrf LAPACKE_dsytri LAPACKE_dsytrs LAPACKE_dtbcon LAPACKE_dtbtrs
LAPACKE_dtpcon LAPACKE_dtptri LAPACKE_dtptrs LAPACKE_dtrcon LAPACKE_dtrexc
LAPACKE_dtrsen LAPACKE_dtrsna LAPACKE_dtrsyl LAPACKE_dtrtri LAPACKE_dtrtrs
LAPACKE_dtrttf LAPACKE_dlatms
LAPACKE_slamch
LAPACKE_zcgesv LAPACKE_zcposv LAPACKE_zgbcon LAPACKE_zgbsv LAPACKE_zgbsvx
LAPACKE_zgbtrf LAPACKE_zgbtrs LAPACKE_zgebak LAPACKE_zgebal LAPACKE_zgecon
LAPACKE_zgees LAPACKE_zgeesx LAPACKE_zgeev LAPACKE_zgeevx LAPACKE_zgejsv
LAPACKE_zgelqf LAPACKE_zgels LAPACKE_zgelsd LAPACKE_zgelss LAPACKE_zgelsy
LAPACKE_zgeqp3 LAPACKE_zgeqrf LAPACKE_zgesdd LAPACKE_zgesv LAPACKE_zgesvd
LAPACKE_zgesvdq LAPACKE_zgesvdx LAPACKE_zgesvx LAPACKE_zgetrf LAPACKE_zgetri
LAPACKE_zgetrs LAPACKE_zgetsls LAPACKE_zgges LAPACKE_zggesx LAPACKE_zggev
LAPACKE_zggevx LAPACKE_zggglm LAPACKE_zgglse LAPACKE_zggsvd3 LAPACKE_zgtcon
LAPACKE_zgtsv LAPACKE_zgtsvx LAPACKE_zgttrf LAPACKE_zgttrs LAPACKE_zhbev
LAPACKE_zhbevd LAPACKE_zhbevx LAPACKE_zhbgv LAPACKE_zhbgvd LAPACKE_zhbgvx
LAPACKE_zhbtrd LAPACKE_zhecon LAPACKE_zheev LAPACKE_zheevd LAPACKE_zheevr
LAPACKE_zheevx LAPACKE_zhegv LAPACKE_zhegvd LAPACKE_zhegvx LAPACKE_zhesv
LAPACKE_zhesvx LAPACKE_zhetrd LAPACKE_zhetrf LAPACKE_zhetri LAPACKE_zhetrs
LAPACKE_zhpcon LAPACKE_zhpev LAPACKE_zhpevd LAPACKE_zhpevx LAPACKE_zhpgv
LAPACKE_zhpgvd LAPACKE_zhpgvx LAPACKE_zhpsv LAPACKE_zhpsvx LAPACKE_zhptrd
LAPACKE_zhptrf LAPACKE_zhptri LAPACKE_zhptrs LAPACKE_zhsein LAPACKE_zhseqr
LAPACKE_zlange LAPACKE_zlanhe LAPACKE_zlansy LAPACKE_zlantr LAPACKE_zpbcon
LAPACKE_zpbsv LAPACKE_zpbsvx LAPACKE_zpbtrf LAPACKE_zpbtrs LAPACKE_zpocon
LAPACKE_zposv LAPACKE_zposvx LAPACKE_zpotrf LAPACKE_zpotri LAPACKE_zpotrs
LAPACKE_zppcon LAPACKE_zppsv LAPACKE_zppsvx LAPACKE_zpptrf LAPACKE_zpptri
LAPACKE_zpptrs LAPACKE_zptcon LAPACKE_zpteqr LAPACKE_zptsv LAPACKE_zptsvx
LAPACKE_zpttrf LAPACKE_zpttrs LAPACKE_zspcon LAPACKE_zspsv LAPACKE_zspsvx
LAPACKE_zsptrf LAPACKE_zsptri LAPACKE_zsptrs LAPACKE_zstedc LAPACKE_zstein
LAPACKE_zstemr LAPACKE_zsteqr LAPACKE_zsycon LAPACKE_zsyr LAPACKE_zsysv
LAPACKE_zsysvx LAPACKE_zsytrf LAPACKE_zsytri LAPACKE_zsytrs LAPACKE_ztbcon
LAPACKE_ztbtrs LAPACKE_ztpcon LAPACKE_ztptri LAPACKE_ztptrs LAPACKE_ztpttf
LAPACKE_ztrcon LAPACKE_ztrevc3 LAPACKE_ztrexe LAPACKE_ztrsen LAPACKE_ztrsna
LAPACKE_ztrsyl LAPACKE_ztrtri LAPACKE_ztrtrs LAPACKE_ztrttf LAPACKE_zunghr
LAPACKE_zunglq LAPACKE_zungqr LAPACKE_zungtr LAPACKE_zunmhr LAPACKE_zunmlq
LAPACKE_zunmtr LAPACKE_zunmqr LAPACKE_zupgtr LAPACKE_zupmtr LAPACKE_zlatms
Different from the original version, error messages for some errors are output by the LAPACK/BLAS base routines instead of LAPACKE/CBLAS routines. In those cases, the displayed argument number may be different from LAPACKE/CBLAS. Please replace it as necessary.
The following functions which are not included in the original library are added as an extension.
LAPACKE additional functions
double LAPACKE_dlangb(int matrix_layout, char norm, int n, int kl, int ku, const double *ab, int ldab)
double LAPACKE_dlangt(char norm, int n, const double *dl, const double *d, const double *du)
double LAPACKE_dlansb(int matrix_layout, char norm, char uplo, int n, int k, const double *ab, int ldab)
double LAPACKE_dlansp(int matrix_layout, char norm, char uplo, int n, const double *ap)
double LAPACKE_dlanst(char norm, int n, const double *d, const double *e)
int LAPACKE_dtrevc3(int matrix_layout, char side, char howmny, int* select, int n, const double* t, int ldt, double* vl, int ldvl, double* vr, int ldvr, int mm, int* m)
double LAPACKE_zlangb(int matrix_layout, char norm, int n, int kl, int ku, const doublecomplex *ab, int ldab)
double LAPACKE_zlangt(char norm, int n, const doublecomplex *dl, const doublecomplex *d, const doublecomplex *du)
double LAPACKE_zlanhb(int matrix_layout, char norm, char uplo, int n, int k, const doublecomplex *ab, int ldab)
double LAPACKE_zlanhp(int matrix_layout, char norm, char uplo, int n, const doublecomplex *ap)
double LAPACKE_zlanht(char norm, int n, const double *d, const doublecomplex *e)
double LAPACKE_zlansp(int matrix_layout, char norm, char uplo, int n, const doublecomplex *ap)
int LAPACKE_ztrevc3(int matrix_layout, char side, char howmny, const int* select, int n, doublecomplex* t, int ldt, doublecomplex* vl, int ldvl, doublecomplex* vr, int ldvr, int mm, int* m)
Example programs using API
Example programs using API are described below.
The following two example problems (using dgesv and qk15) are described below. For other functions, please refer to sample programs included in SDK.
Example (1)
Solve the system of linear equations Ax = B. Coefficient matrix A, right hand matrix B and solution x are as follows.
( 0.2 -0.11 -0.93 ) ( -0.3727 ) ( 0.86 )
A = ( -0.32 0.81 0.37 ), B = ( 0.4319 ), x = ( 0.64 )
( -0.8 -0.92 -0.29 ) ( -1.4247 ) ( 0.51 )
Example (2)
Computes the following integral. The true analytic solution is atan(4) = 1.32581766366803.
∫ 1/(1 + x^2) dx [0, 4]
C/C++ example programs
C example program (Example (1), Visual C)
This is the example using Visual C (Visual Studio (Windows version)). Note that matrices are stored in column-major order (same with VBA and Fortran, rows and columns are opposite to the C/C++ storage order).
One dimensional array is used for a and b (In the reference manual, the variable length two dimensional array is used like Clang example below. However, please read with replacing it to one dimensional array in the case of Visual C).
#include <stdio.h>
#include "cnumlib.h"
int main(void)
{
double a[3][3] = {
{ 0.2, -0.32, -0.8 },
{ -0.11, 0.81, -0.92 },
{ -0.93, 0.37, -0.29 } };
double b[] = { -0.3727, 0.4319, -1.4247 };
int n = 3, nrhs = 1, lda = 3, ldb = 3, ipiv[3], info;
dgesv(n, nrhs, lda, (double *)a, ipiv, ldb, b, &info);
printf("x = %f %f %f\n", b[0], b[1], b[2]);
printf("info = %d\n", info);
return 0;
}
Result
x = 0.860000 0.640000 0.510000
info = 0
C example program (Example (1), Clang)
Two dimensional arrays are used for both a and b since Clang suports the variable length array which is standardized in C99.
#include <stdio.h>
#include "cnumlib.h"
int main(void)
{
double a[3][3] = {
{ 0.2, -0.32, -0.8 },
{ -0.11, 0.81, -0.92 },
{ -0.93, 0.37, -0.29 } };
double b[] = { -0.3727, 0.4319, -1.4247 };
int n = 3, nrhs = 1, lda = 3, ldb = 3, ipiv[3], info;
dgesv(n, nrhs, lda, a, ipiv, ldb, (double (*)[ldb])b, &info);
printf("x = %f %f %f\n", b[0], b[1], b[2]);
printf("info = %d\n", info);
return 0;
}
C example program (Example (1), LAPACKE Row-major)
It is possible to specify to use the row-major order array by the first parameter of LAPACKE routine. Using this feature, a matrix can be coded in the same order with mathematical notation. Internally, arrays are transposed before computation and transposed again after computation. Note that ldb = 1. This is because the horizontal vector is expected as the input b since not only a but also b is transposed.
#include <stdio.h>
#include "lapacke.h"
int main(void)
{
double a[3][3] = {
{ 0.2, -0.11, -0.93 },
{ -0.32, 0.81, 0.37 },
{ -0.8, -0.92, -0.29 } };
double b[] = { -0.3727, 0.4319, -1.4247 };
int n = 3, nrhs = 1, lda = 3, ldb = 1, ipiv[3], info;
info = LAPACKE_dgesv(LAPACK_ROW_MAJOR, n, nrhs, (double *)a, lda, ipiv, b, ldb);
printf("x = %f %f %f\n", b[0], b[1], b[2]);
printf("info = %d\n", info);
return 0;
}
C++ example program (Example (1))
cnumlib header file can be used (note that "info" argument is used instead of "&info"). One dimensional array is used for matrices even if Clang is used to compile.
#include <iostream>
#include "cnumlib"
using namespace std;
int main(void)
{
double a[3][3] = {
{ 0.2, -0.32, -0.8 },
{ -0.11, 0.81, -0.92 },
{ -0.93, 0.37, -0.29 } };
double b[] = { -0.3727, 0.4319, -1.4247 };
int n = 3, nrhs = 1, lda = 3, ldb = 3, ipiv[3], info;
dgesv(n, nrhs, lda, (double *)a, ipiv, ldb, b, info);
cout << "x = " << b[0] << ", " << b[1] << ", " << b[2] << endl;
cout << "info = " << info << endl;
return 0;
}
C example program (Example (2))
In Example (2), the integral of f(x) is computed using qk15. qk15 requires the external function defining f(x). It can be coded in C and qk15 will call it when necessary.
#include <stdio.h>
#include "cnumlib.h"
double f(double x)
{
return 1/(1 + x*x);
}
void test_qk15()
{
double a, b, result, abserr, resabs, resasc;
a = 0; b = 4;
qk15(f, a, b, &result, &abserr, &resabs, &resasc);
printf("result = %g, abserr = %g\n", result, abserr);
}
Result
result = 1.32582, abserr = 0.00148272
Python example programs
The XLPack Python extension module is included in SDK, so that the calling program can easily call XLPack like built-in functions by declaring “import XLPack” (or “import XLPack_32” for 32 bit Python).
Numpy ndarray is used to represent the matrix. Either two dimensional or one dimensional array can be used for both coefficient matrix a and right hand matrix b.
Note – For Python 3.8 or later, an error occurs if DLL’s path is not explicitly specified. Refer to the example program included in SDK.
Python example program (Example (1))
import numpy as np
from XLPack import *
def TestDgesv():
n = 3
a = np.array([
[0.2, -0.32, -0.8],
[-0.11, 0.81, -0.92],
[-0.93, 0.37, -0.29]])
b = np.array([-0.3727, 0.4319, -1.4247])
ipiv = np.empty(n, dtype=int)
info = dgesv(n, a, ipiv, b)
print(b, info)
Result
Python example program (Example (2))
In Example (2), the integral of f(x) is computed using qk15. qk15 requires the external function defining f(x). It can be coded in Python and qk15 will call it when necessary.
import numpy as np
from XLPack import *
def f(x):
return 1/(1 + x**2)
def TestQk15():
a = 0
b = 4
s, abserr = qk15(f, a, b)
print(s, abserr)
TestQk15()
Result
1.3258176613637855 0.0014827239412162237
Julia example programs
The XLPack Julia interface program example is included in SDK, so that the calling program can easily call XLPack like built-in functions by declaring “using XLPack”.
Float64 is used for floating point variables and Int32 is used for integer variables. Either two dimensional or one dimensional array can be used for both coefficient matrix a and right hand matrix b.
Julia example program (Example (1))
using XLPack
function TestDgesv()
n = 3
a = [ 0.2 -0.11 -0.93;
-0.32 0.81 0.37;
-0.8 -0.92 -0.29 ]
b = [ -0.3727, 0.4319, -1.4247 ]
ipiv = Vector{Cint}(undef, n)
info = dgesv(n, a, ipiv, b)
println("x = ", b, ", info = ", info)
end
TestDgesv()
Result
x = [0.8600000000000002, 0.64, 0.5099999999999999], info = 0
Julia example program (Example (2))
In Example (2), the integral of f(x) is computed using qk15. qk15 requires the external function defining f(x). It can be coded in Julia and qk15 will call it when necessary.
using XLPack
function TestQk15()
f(x) = 1/(1 + x^2)
a = 0
b = 4
s, abserr = qk15(f, a, b)
println(s, " ", abserr)
end
TestQk15()
Result
1.3258176613637855 0.0014827239412162237
C# example programs
To call XLPack from C#, it is necessary to use DLL directly. DLL interface program example xlpack.cs is included in SDK. By using it, the calling program can easily call XLPack like built-in functions without considering the detail DLL interface by declaring “using static XLPack;”.
Two dimensional array is used for the coefficient matrix a. One dimensional array is used for the right hand matrix b.
C# example program (Example (1))
using System;
using System.IO;
using static XLPack;
public class Test
{
static void Test_dgesv()
{
const int n = 3;
double[,] a =
{{ 0.2, -0.32, -0.8 },
{ -0.11, 0.81, -0.92 },
{ -0.93, 0.37, -0.29 }};
double[] b = { -0.3727, 0.4319, -1.4247 };
int[] ipiv = new int[n];
int info;
dgesv(n, a, ipiv, b, out info);
Console.Write("x = {0} {1} {2}\n", b[0], b[1], b[2]);
Console.Write("info = {0}\n", info);
}
public static int Main(string[] args)
{
Test_dgesv();
return 0;
}
}
Result
x = 0.86 0.64 0.51
info = 0
C# example program (Example (2))
In Example (2), the integral of f(x) is computed using qk15. qk15 requires the external function defining f(x). It can be coded in C# and qk15 will call it when necessary.
using System;
using System.IO;
using static XLPack;
public class Test
{
static double f_qk15(double x)
{
return 1/(1 + x*x);
}
static void Test_qk15()
{
double a, b, result, abserr;
a = 0; b = 4;
qk15(f_qk15, a, b, out result, out abserr);
Console.Write("result = {0}, abserr = {1}\n", result, abserr);
}
public static int Main(string[] args)
{
Test_qk15();
return 0;
}
}
Result
result = 1.32581766136379, abserr = 0.00148272394121622
F# example programs
To call XLPack from F#, it is necessary to use DLL directly. DLL interface program example xlpack_fs.fs is included in SDK. By using it, the calling program can easily call XLPack like built-in functions without considering the detail DLL interface by declaring “open XLPack”.
Two dimensional array is used for the coefficient matrix a. One dimensional array is used for the right hand matrix b.
F# example program (Example (1))
open XLPack
let TestDgesv() =
let n = 3
let a = array2D [
[ 0.2; -0.32; -0.8 ];
[ -0.11; 0.81; -0.92 ];
[ -0.93; 0.37; -0.29 ] ]
let b = [| -0.3727; 0.4319; -1.4247 |]
let ipiv = Array.create n 0
let info = XLPack.Dgesv(n, a, ipiv, b)
printfn "x = %A, info = %d" b info
TestDgesv()
Result
x = [|0.86; 0.64; 0.51|], info = 0
F# example program (Example (2))
In Example (2), the integral of f(x) is computed using qk15. qk15 requires the external function defining f(x). It can be coded in F# and qk15 will call it when necessary.
open XLPack
let TestQk15() =
let f(x: double) = 1.0/(1.0 + x*x)
let a, b = (0.0, 4.0)
let result, abserr = XLPack.Qk15(f, a, b)
printfn "result = %f, abserr = %g" result abserr
TestQk15()
Result
result = 1.325818, abserr = 0.00148272
VB.NET example programs
To call XLPack from VB.NET, it is necessary to use DLL directly. DLL interface program example xlpack.vb is included in SDK. By using it, the calling program can easily call XLPack like built-in functions without considering the detail DLL interface by declaring “Imports XLPack”.
Two dimensional array is used for the coefficient matrix a. One dimensional array is used for the right hand matrix b.
Different from Excel VBA, two dimensional array of VB.NET is stored in row-major order like C/C++, and the 32 bit integer data type of VB.NET is “Integer” but not “Long”.
VB.NET example program (Example (1))
Imports XLPack
Class Test
Shared Sub TestDgesv
Const N = 3
Dim A(,) As Double = {
{ 0.2, -0.32, -0.8 },
{ -0.11, 0.81, -0.92 },
{ -0.93, 0.37, -0.29 } }
Dim B() As Double = { -0.3727, 0.4319, -1.4247 }
Dim IPiv(N - 1) As Integer
Dim Info As Integer
Call Dgesv(N, A, IPiv, B, Info)
Console.WriteLine("X = {0} {1} {2}", B(0), B(1), B(2))
Console.WriteLine("Info = {0}", Info)
End Sub
Public Shared Sub Main(ByVal args() As String)
Call TestDgesv
End Sub
End Class
Result
x = 0.86 0.64 0.51
info = 0
VB.NET example program (Example (2))
In Example (2), the integral of f(x) is computed using qk15. qk15 requires the external function defining f(x). It can be coded in VB.NET and qk15 will call it when necessary.
Imports XLPack
Class Test
Shared Function F(X As Double) As Double
F = 1/(1 + X ^ 2)
End Function
Shared Sub TestQk15
Dim A As Double, B As Double, Result As Double, AbsErr As Double
A = 0: B = 4
Call Qk15(AddressOf F, A, B, Result, AbsErr)
Console.WriteLine("Result = {0}, AbsErr = {1}", Result, AbsErr)
End Sub
Public Shared Sub Main(ByVal args() As String)
Call TestQk15
End Sub
End Class
Result
Result = 1.32581766136379, AbsErr = 0.00148272394121622
Pascal example programs
To call XLPack from Pascal, it is necessary to use DLL directly. DLL interface program example xlpack.pas is included in SDK. By using it, the calling program can easily call XLPack like built-in functions without considering the detail DLL interface by declaring “uses XLPack;”.
One dimensional array is used for both coefficient matrix a and right hand matrix b.
In this example, Free Pascal is used in Delphi mode.
Pascal example program (Example (1))
{$Mode Delphi}
{$AppType CONSOLE}
program test_xlpack_1;
uses XLPack;
procedure test_dgesv;
const
n = 3; lda = n; ldb = n;
a: array[1..n*n] of Double =
(0.2, -0.32, -0.8,
-0.11, 0.81, -0.92,
-0.93, 0.37, -0.29);
b: array[1..n] of Double =
(-0.3727, 0.4319, -1.4247);
var
ipiv: array[1..n] of Integer;
info: Integer;
begin
Dgesv(n, lda, a, ipiv, ldb, b, info);
Writeln('x =');
Writeln(b[1], ' ', b[2], ' ', b[3]);
Writeln('info = ', info);
end;
begin
test_dgesv;
end.
Result
x =
8.60000000000000E-0001 6.40000000000000E-0001 5.10000000000000E-0001
info = 0
Pascal example program (Example (2))
In Example (2), the integral of f(x) is computed using qk15. qk15 requires the external function defining f(x). It can be coded in Pascal and qk15 will call it when necessary.
{$Mode Delphi}
{$AppType CONSOLE}
program test_xlpack_2;
uses XLPack;
function f(x: Double): Double;
begin
f := 1/(1 + x*x);
end;
procedure test_qk15;
var
a, b, result, abserr: Double;
begin
a := 0; b := 4;
Qk15(f, a, b, result, abserr);
Writeln('result = ', result, ', abserr = ', abserr);
end;
begin
test_qk15;
end.
Result
result = 1.32581766136379E+0000, abserr = 1.48272394121622E-0003
Source libraries and references
Refer to "Main Page" tab of each reference manual.
Published: Jul 26, 2022 (Last update: Feb 27, 2023)