Skip to content

BLIS Usage in FORTRAN and C through BLAS and CBLAS APIs with Code Examples

Nisanth M P edited this page Feb 23, 2018 · 3 revisions

BLIS has a new BLAS-like API that can be used for BLAS-equivalent operations. This API also provides some additional operations compared to what is available with the traditional BLAS API. At the same time, BLIS also includes a BLAS compatibility layer which gives application developers access to BLIS implementations via traditional BLAS API calls, that can be used in FORTRAN as well as C code. BLIS also provides a CBLAS API (navigate LAPACK > Files > File List > CBLAS), which is a C-style interface for BLAS, that can be called from C code.

Using BLIS with FORTRAN code through standard BLAS API

BLIS can be used with FORTRAN applications through the standard BLAS API.

For example, see below, FORTRAN code that does double precision general matrix-matrix multiplication. It calls the 'DGEMM' BLAS API function to accomplish this. An example command to compile it and link with the BLIS library is also shown below the code.

! File: BLAS_DGEMM_usage.f
! Example code to demonstrate BLAS DGEMM usage

program dgemm_usage

implicit none

EXTERNAL DGEMM

DOUBLE PRECISION, ALLOCATABLE :: a(:,:)
DOUBLE PRECISION, ALLOCATABLE :: b(:,:)
DOUBLE PRECISION, ALLOCATABLE :: c(:,:)
INTEGER I, J, M, N, K, lda, ldb, ldc
DOUBLE PRECISION alpha, beta

M=2
N=M
K=M
lda=M
ldb=K
ldc=M
alpha=1.0
beta=0.0

ALLOCATE(a(lda,K), b(ldb,N), c(ldc,N))

a=RESHAPE((/ 1.0, 3.0, &
             2.0, 4.0  /), &
             (/lda,K/))
b=RESHAPE((/ 5.0, 7.0, &
             6.0, 8.0  /), &
             (/ldb,N/))

WRITE(*,*) ("a =")
DO I = LBOUND(a,1), UBOUND(a,1)
    WRITE(*,*) (a(I,J), J=LBOUND(a,2), UBOUND(a,2))
END DO
WRITE(*,*) ("b =")
DO I = LBOUND(b,1), UBOUND(b,1)
    WRITE(*,*) (b(I,J), J=LBOUND(b,2), UBOUND(b,2))
END DO

CALL DGEMM('N','N',M,N,K,alpha,a,lda,b,ldb,beta,c,ldc)

WRITE(*,*) ("c =")
DO I = LBOUND(c,1), UBOUND(c,1)
    WRITE(*,*) (c(I,J), J=LBOUND(c,2), UBOUND(c,2))
END DO

end program dgemm_usage

Example compilation command, for the above code, with gfortran compiler:

gfortran -ffree-form BLAS_DGEMM_usage.f path/to/libblis.a

Using BLIS with C code through BLAS and CBLAS APIs

There are multiple ways to use BLIS with an application written in C. While one can always use the native BLIS API for the same, BLIS also includes BLAS and CBLAS interfaces.

Using BLIS with BLAS API in C code

Shown below is the C version of the code listed above in FORTRAN. It uses the standard BLAS API. Note that (a) the matrices are transposed to account for the row-major storage of C and the column-major convention of BLAS (inherited from FORTRAN), (b) the function arguments are passed by address, again to be in line with FORTRAN conventions, (c) there is a trailing underscore in the function name ('dgemm_'), as BLIS' BLAS APIs expect that (FORTRAN compilers add a trailing underscore), and (d) "blis.h" is included as a header. An example command to compile it and link with the BLIS library is also shown below the code.

// File: BLAS_DGEMM_usage.c
// Example code to demonstrate BLAS DGEMM usage

#include<stdio.h>
#include "blis.h"

#define DIM 2

int main() {

	double a[DIM * DIM] = { 1.0, 3.0, 2.0, 4.0 };
	double b[DIM * DIM] = { 5.0, 7.0, 6.0, 8.0 };
	double c[DIM * DIM];
	int I, J, M, N, K, lda, ldb, ldc;
	double alpha, beta;

	M = DIM;
	N = M;
	K = M;
	lda = M;
	ldb = K;
	ldc = M;
	alpha = 1.0;
	beta = 0.0;

	printf("a = \n");
	for ( I = 0; I < M; I ++ ) {
		for ( J = 0; J < K; J ++ ) {
			printf("%f\t", a[J * K + I]);
		}
		printf("\n");
	}
	printf("b = \n");
	for ( I = 0; I < K; I ++ ) {
		for ( J = 0; J < N; J ++ ) {
			printf("%f\t", b[J * N + I]);
		}
		printf("\n");
	}

	dgemm_("N","N",&M,&N,&K,&alpha,a,&lda,b,&ldb,&beta,c,&ldc);

	printf("c = \n");
	for ( I = 0; I < M; I ++ ) {
		for ( J = 0; J < N; J ++ ) {
			printf("%f\t", c[J * N + I]);
		}
		printf("\n");
	}

	return 0;
}

Example compilation command, for the above code, with gcc compiler:

gcc BLAS_DGEMM_usage.c -Ipath/to/include/blis/ path/to/libblis.a

Using BLIS with CBLAS API

Shown below is the C code using CBLAS APIs for the same functionality listed above. CBLAS API can be looked up over here (navigate LAPACK > Files > File List > CBLAS). Note that (a) CBLAS Layout option allows us to choose between row-major and column-major layouts (row-major layout is used in the example, which is in line with C-style), (b) the function arguments can be passed by value also, and (c) "cblas.h" is included as a header. An example command to compile it and link with the BLIS library is also shown below the code. Also, note that, in order to get CBLAS API with BLIS, one has to supply the flag '--enable-cblas' to the 'configure' command while building the BLIS library.

// File: CBLAS_DGEMM_usage.c
// Example code to demonstrate CBLAS DGEMM usage

#include<stdio.h>
#include "cblas.h"

#define DIM 2

int main() {

	double a[DIM * DIM] = { 1.0, 2.0, 3.0, 4.0 };
	double b[DIM * DIM] = { 5.0, 6.0, 7.0, 8.0 };
	double c[DIM * DIM];
	int I, J, M, N, K, lda, ldb, ldc;
	double alpha, beta;

	M = DIM;
	N = M;
	K = M;
	lda = M;
	ldb = K;
	ldc = M;
	alpha = 1.0;
	beta = 0.0;

	printf("a = \n");
	for ( I = 0; I < M; I ++ ) {
		for ( J = 0; J < K; J ++ ) {
			printf("%f\t", a[I * K + J]);
		}
		printf("\n");
	}
	printf("b = \n");
	for ( I = 0; I < K; I ++ ) {
		for ( J = 0; J < N; J ++ ) {
			printf("%f\t", b[I * N + J]);
		}
		printf("\n");
	}

	cblas_dgemm(CblasRowMajor,  CblasNoTrans, CblasNoTrans, M, N, K, alpha, a, lda, b, ldb, beta, c, ldc);

	printf("c = \n");
	for ( I = 0; I < M; I ++ ) {
		for ( J = 0; J < N; J ++ ) {
			printf("%f\t", c[I * N + J]);
		}
		printf("\n");
	}

	return 0;
}

Example compilation command, for the above code, with gcc compiler:

gcc CBLAS_DGEMM_usage.c -Ipath/to/include/blis/ path/to/libblis.a