This wiki summarizes the user-level BLIS APIs, including BLIS-specific type definitions, header files, and function prototypes. This document also includes APIs to key [micro-]kernels which are used to accelerate and optimize various level-2 and level-3 operations.
There are many functions that BLIS implements that are not listed here, either because they are lower-level functions, or they are considered for use primarily by developers and experts.
The following tables list various types used throughout the BLIS API.
BLIS type | Definition | Used to represent... |
---|---|---|
gint_t |
int32_t or int64_t |
general-purpose signed integer; used to define signed integer types. |
dim_t |
gint_t |
matrix and vector dimensions. |
inc_t |
gint_t |
matrix row/column strides and vector increments. |
doff_t |
gint_t |
matrix diagonal offset: if k < 0, diagonal begins at element (-k,0); otherwise diagonal begins at element (0,k). |
BLIS type | BLIS char | Definition | Used to represent... |
---|---|---|---|
float |
s |
N/A | single-precision real numbers |
double |
d |
N/A | double-precision real numbers |
scomplex |
c |
struct { float real; float imag; } |
single-precision complex numbers |
dcomplex |
z |
struct { double real; double imag; } |
double-precision complex numbers |
trans_t |
Semantic meaning: Corresponding matrix operand... |
---|---|
BLIS_NO_TRANSPOSE |
will be used as given. |
BLIS_TRANSPOSE |
will be implicitly transposed. |
BLIS_CONJ_NO_TRANSPOSE |
will be implicitly conjugated. |
BLIS_CONJ_TRANSPOSE |
will be implicitly transposed and conjugated. |
conj_t |
Semantic meaning: Corresponding matrix/vector operand... |
---|---|
BLIS_NO_CONJUGATE |
will be used as given. |
BLIS_CONJUGATE |
will be implicitly conjugated. |
side_t |
Semantic meaning: Corresponding matrix operand... |
---|---|
BLIS_LEFT |
appears on the left. |
BLIS_RIGHT |
appears on the right. |
uplo_t |
Semantic meaning: Corresponding matrix operand... |
---|---|
BLIS_LOWER |
is stored in (and will be accessed only from) the lower triangle. |
BLIS_UPPER |
is stored in (and will be accessed only from) the upper triangle. |
BLIS_DENSE |
is stored as a full matrix (ie: in both triangles). |
diag_t |
Semantic meaning: Corresponding matrix operand... |
---|---|
BLIS_NONUNIT_DIAG |
has a non-unit diagonal that should be explicitly read from. |
BLIS_UNIT_DIAG |
has a unit diagonal that should be implicitly assumed (and not read from). |
The functions listed in this document belong to the basic interface subset of the BLIS typed API. There is a companion "expert" interface that mirrors the basic interface, except that it also contains at least one additional parameter that is only of interest to experts and library developers. The expert interfaces use the same name as the basic function names, except for an additional "_ex" suffix. For example, the basic interface for gemm
is
void bli_?gemm
(
trans_t transa,
trans_t transb,
dim_t m,
dim_t n,
dim_t k,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb,
ctype* beta,
ctype* c, inc_t rsc, inc_t csc
);
while the expert interface is:
void bli_?gemm_ex
(
trans_t transa,
trans_t transb,
dim_t m,
dim_t n,
dim_t k,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb,
ctype* beta,
ctype* c, inc_t rsc, inc_t csc,
cntx_t* cntx
);
The expert interface contains an additional cntx_t*
parameter. Note that calling a function from the expert interface with the cntx_t*
argument set to NULL
is equivalent to calling the corresponding basic interface.
The cntx_t*
argument appears in the interfaces to the gemm
, trsm
, and gemmtrsm
micro-kernels. However, assembly implementations for micro-kernels generally do not make use of the context structure, and therefore NULL
may be passed in by casual developers.
All BLIS definitions and prototypes may be included in your C source file by including a single header file:
#include "blis.h"
As of 9804adf, BLIS no longer requires explicit initialization and finalization at runtime. In other words, users do not need to call bli_init()
before the application can make use of the library (and bli_finalize()
after the application is finished with the library). Instead, all computational operations (and some non-computational functions) in BLIS will initialize the library on behalf of the user if it has not already been initialized. This change was made to simplify the user experience.
Application developers should keep in mind, however, that this new self-initialization regime implies the following: unless the library is explicitly finalized via bli_finalize()
, it will, once initialized, remain initialized for the life of the application. This is likely not a problem in the vast majority of cases. However, a memory-constrained application that performs all of its DLA up-front, for example, may wish to explicitly finalize the library after BLIS is no longer needed in order to free up memory for other purposes.
Similarly, an expert user may call bli_init()
manually in order to control when the overhead of library initialization is incurred, even though the library would have self-initialized.
The interfaces to bli_init()
and bli_finalize()
are quite simple; they require no arguments and return no values:
void bli_init( void );
void bli_finalize( void );
Notes for interpreting the following prototypes:
- Any occurrence of
?
should be replaced withs
,d
,c
, orz
to form an actual function name. - Any occurrence of
ctype
should be replaced with the actual C type corresponding to the datatype instance in question, whilertype
should be replaced by the real projection ofctype
. For example:- If we consider the prototype for
bli_zaxpyv()
below,ctype
refers todcomplex
. - If we consider the prototype for
bli_zfnormv()
below,ctype
refers todcomplex
whilertype
refers todouble
.
- If we consider the prototype for
- Any occurrence of
itype
should be replaced with the general-purpose signed integer type,gint_t
. - All vector arguments have associated increments that proceed them, typically listed as
incX
for a given vectorx
. The semantic meaning of a vector increment is "the distance, in units of elements, between any two adjacent elements in the vector." - All matrix arguments have associated row and column strides arguments that proceed them, typically listed as
rsX
andcsX
for a given matrixX
. Row strides are always listed first, and column strides are always listed second. The semantic meaning of a row stride is "the distance, in units of elements, to the next row (within a column)," and the meaning of a column stride is "the distance, in units of elements, to the next column (within a row)." Thus, unit row stride implies column-major storage and unit column stride implies row-major storage.
Notes for interpreting function descriptions:
conjX()
andtransX()
should be interpreted as predicates that capture the operand X with any value ofconj_t
ortrans_t
applied. For example:conjx(x)
refers to a vectorx
that is either conjugated or used as given.transa(A)
refers to a matrixA
that is either transposed, conjugated and transposed, conjugated only, or used as given.
- Any operand marked with
conj()
is unconditionally conjugated. - Any operand marked with
^T
is unconditionally transposed. Similarly, any operand that is marked with^H
is unconditionally conjugate-transposed. - All occurrences of
alpha
,beta
, andrho
parameters are scalars.
- Level-1v: Operations on vectors:
- Level-1d: Element-wise operations on matrix diagonals:
- Level-1m: Element-wise operations on matrices:
- Level-1f: Fused operations on multiple vectors:
- Level-2: Operations with one matrix and (at least) one vector operand:
- Level-3: Operations with matrices that are multiplication-like:
- Utility: Miscellaneous operations on matrices and vectors:
Level-1v operations perform various level-1 BLAS-like operations on vectors (hence the v). Note: Each level-1v operation has a corresponding level-1v kernel through which it is primarily implemented.
void bli_?addv
(
conj_t conjx,
dim_t m,
ctype* x, inc_t incx,
ctype* y, inc_t incy
);
Perform
y := y + conjx(x)
where y
and x
are vectors of length m.
void bli_?amaxv
(
dim_t n,
ctype* x, inc_t incx,
dim_t* index
);
Find the element of vector x
which contains the maximum absolute value. The index of the element found is stored to index
.
Note: This function attempts to mimic the algorithm for finding the element with the maximum absolute value in the netlib BLAS routines i?amax()
.
void bli_?axpyv
(
conj_t conjx,
dim_t m,
ctype* alpha,
ctype* x, inc_t incx,
ctype* y, inc_t incy
);
Perform
y := y + alpha * conjx(x)
where y
and x
are vectors of length m.
void bli_?copyv
(
conj_t conjx,
dim_t m,
ctype* x, inc_t incx,
ctype* y, inc_t incy
);
Perform
y := conjx(x)
where y
and x
are vectors of length m.
void bli_?dotv
(
conj_t conjx,
conj_t conjy,
dim_t m,
ctype* x, inc_t incx,
ctype* y, inc_t incy,
ctype* rho
);
Perform
rho := conjx(x)^T * conjy(y)
where y
and x
are vectors of length m and rho
is a scalar.
void bli_?dotxv
(
conj_t conjx,
conj_t conjy,
dim_t m,
ctype* alpha,
ctype* x, inc_t incx,
ctype* y, inc_t incy,
ctype* beta,
ctype* rho
);
Perform
rho := beta * rho + alpha * conjx(x)^T * conjy(y)
where y
and x
are vectors of length m and rho
is a scalar.
void bli_?invertv
(
dim_t m,
ctype* x, inc_t incx
);
Invert all elements of an m-length vector x
.
void bli_?scalv
(
conj_t conjalpha,
dim_t m,
ctype* alpha,
ctype* x, inc_t incx
);
Perform
x := conjalpha(alpha) * x
where x
is a vector of length m.
void bli_?scal2v
(
conj_t conjx,
dim_t m,
ctype* alpha,
ctype* x, inc_t incx,
ctype* y, inc_t incy
);
Perform
y := alpha * conjx(x)
where y
and x
are vectors of length m.
void bli_?setv
(
conj_t conjalpha,
dim_t m,
ctype* alpha,
ctype* x, inc_t incx
);
Set all elements of an m-length vector x
to conjalpha(alpha)
.
void bli_?subv
(
conj_t conjx,
dim_t m,
ctype* x, inc_t incx,
ctype* y, inc_t incy
);
Perform
y := y - conjx(x)
where y
and x
are vectors of length m.
void bli_?swapv
(
dim_t m,
ctype* x, inc_t incx,
ctype* y, inc_t incy
);
Swap corresponding elements of two m-length vectors x
and y
.
Level-1d operations perform various level-1 BLAS-like operations on matrix diagonals (hence the d).
These operations are similar to their level-1m counterparts, except they only read and update matrix diagonals. Please see the descriptions for the corresponding level-1m operation for a description of the arguments.
void bli_?addd
(
doff_t diagoffa,
diag_t diaga,
trans_t transa,
dim_t m,
dim_t n,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb
);
void bli_?axpyd
(
doff_t diagoffa,
diag_t diaga,
trans_t transa,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb
);
void bli_?copyd
(
doff_t diagoffa,
diag_t diaga,
trans_t transa,
dim_t m,
dim_t n,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb
);
void bli_?invertd
(
doff_t diagoffa,
dim_t m,
dim_t n,
ctype* a, inc_t rsa, inc_t csa
);
void bli_?scald
(
conj_t conjalpha,
doff_t diagoffa,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa
);
void bli_?scal2d
(
doff_t diagoffa,
diag_t diaga,
trans_t transa,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb
);
void bli_?setd
(
conj_t conjalpha,
doff_t diagoffa,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa
);
void bli_?setid
(
doff_t diagoffa,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa
);
Set the imaginary components of a matrix diagonal.
void bli_?subd
(
doff_t diagoffa,
diag_t diaga,
trans_t transa,
dim_t m,
dim_t n,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb
);
Level-1m operations perform various level-1 BLAS-like operations on matrices (hence the m).
void bli_?addm
(
doff_t diagoffa,
diag_t diaga,
uplo_t uploa,
trans_t transa,
dim_t m,
dim_t n,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb
);
Perform
B := B + transa(A)
where B
is an m x n matrix, A
is stored as a dense matrix, or lower- or upper-triangular/trapezoidal matrix, as specified by uploa
, with the diagonal offset of A
specified by diagoffa
and unit/non-unit nature of the diagonal specified by diaga
. If uploa
indicates lower or upper storage, only that part of matrix A
will be referenced and used to update B
.
void bli_?axpym
(
doff_t diagoffa,
diag_t diaga,
uplo_t uploa,
trans_t transa,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb
);
Perform
B := B + alpha * transa(A)
where B
is an m x n matrix, A
is stored as a dense matrix, or lower- or upper-triangular/trapezoidal matrix, as specified by uploa
, with the diagonal offset of A
specified by diagoffa
and unit/non-unit nature of the diagonal specified by diaga
. If uploa
indicates lower or upper storage, only that part of matrix A
will be referenced and used to update B
.
void bli_?copym
(
doff_t diagoffa,
diag_t diaga,
uplo_t uploa,
trans_t transa,
dim_t m,
dim_t n,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb
);
Perform
B := transa(A)
where B
is an m x n matrix, A
is stored as a dense matrix, or lower- or upper-triangular/trapezoidal matrix, as specified by uploa
, with the diagonal offset of A
specified by diagoffa
and unit/non-unit nature of the diagonal specified by diaga
. If uploa
indicates lower or upper storage, only that part of matrix A
will be referenced and used to update B
.
void bli_?scalm
(
conj_t conjalpha,
doff_t diagoffa,
uplo_t uploa,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa
);
Perform
A := conjalpha(alpha) * A
where A
is an m x n matrix stored as a dense matrix, or lower- or upper-triangular/trapezoidal matrix, as specified by uploa
, with the diagonal offset of A
specified by diagoffa
. If uploa
indicates lower or upper storage, only that part of matrix A
will be updated.
void bli_?scal2m
(
doff_t diagoffa,
diag_t diaga,
uplo_t uploa,
trans_t transa,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb
);
Perform
B := alpha * transa(A)
where B
is an m x n matrix, A
is stored as a dense matrix, or lower- or upper-triangular/trapezoidal matrix, as specified by uploa
, with the diagonal offset of A
specified by diagoffa
and unit/non-unit nature of the diagonal specified by diaga
. If uploa
indicates lower or upper storage, only that part of matrix A
will be referenced and used to update B
.
void bli_?setm
(
conj_t conjalpha,
doff_t diagoffa,
diag_t diaga,
uplo_t uploa,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa
);
Set all elements of an m x n matrix A
to conjalpha(alpha)
, where A
is stored as a dense matrix, or lower- or upper- triangular/trapezoidal matrix, as specified by uploa
, with the diagonal offset of A
specified by diagoffa
and unit/non-unit nature of the diagonal specified by diaga
. If uploa
indicates lower or upper storage, only that part of matrix A
will be updated.
void bli_?subm
(
doff_t diagoffa,
diag_t diaga,
uplo_t uploa,
trans_t transa,
dim_t m,
dim_t n,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb
);
Perform
B := B - transa(A)
where B
is an m x n matrix, A
is stored as a dense matrix, or lower- or upper-triangular/trapezoidal matrix, as specified by uploa
, with the diagonal offset of A
specified by diagoffa
and unit/non-unit nature of the diagonal specified by diaga
. If uploa
indicates lower or upper storage, only that part of matrix A
will be referenced and used to update B
.
Level-1f operations implement various fused combinations of level-1 operations (hence the f). Note: Each level-1f operation has a corresponding level-1f kernel through which it is primarily implemented.
Level-1f kernels are employed when optimizing level-2 operations.
void bli_?axpy2v
(
conj_t conjx,
conj_t conjy,
dim_t m,
ctype* alphax,
ctype* alphay,
ctype* x, inc_t incx,
ctype* y, inc_t incy,
ctype* z, inc_t incz
);
Perform
y := y + alphax * conjx(x) + alphay * conjy(y)
where x
, y
, and z
are vectors of length m. The kernel, if optimized, is implemented as a fused pair of calls to axpyv.
void bli_?dotaxpyv
(
conj_t conjxt,
conj_t conjx,
conj_t conjy,
dim_t m,
ctype* alpha,
ctype* x, inc_t incx,
ctype* y, inc_t incy,
ctype* rho,
ctype* z, inc_t incz
);
Perform
rho := conjxt(x^T) * conjy(y)
y := y + alpha * conjx(x)
where x
, y
, and z
are vectors of length m and alpha
and rho
are scalars. The kernel, if optimized, is implemented as a fusion of calls to dotv and axpyv.
void bli_?axpyf
(
conj_t conja,
conj_t conjx,
dim_t m,
dim_t nf,
ctype* alpha,
ctype* a, inc_t inca, inc_t lda,
ctype* x, inc_t incx,
ctype* y, inc_t incy
);
Perform
y := y + alpha * conja(A) * conjx(x)
where A
is an m x nf matrix, and y
and x
are vectors. The kernel, if optimized, is implemented as a fused series of calls to axpyv where nf is less than or equal to an implementation-dependent fusing factor specific to axpyf
.
void bli_?dotxf
(
conj_t conjat,
conj_t conjx,
dim_t m,
dim_t nf,
ctype* alpha,
ctype* a, inc_t inca, inc_t lda,
ctype* x, inc_t incx,
ctype* beta,
ctype* y, inc_t incy
);
Perform
y := y + alpha * conjat(A^T) * conjx(x)
where A
is an m x nf matrix, and y
and x
are vectors. The kernel, if optimized, is implemented as a fused series of calls to dotxv where nf is less than or equal to an implementation-dependent fusing factor specific to dotxf
.
void bli_?dotxaxpyf
(
conj_t conjat,
conj_t conja,
conj_t conjw,
conj_t conjx,
dim_t m,
dim_t nf,
ctype* alpha,
ctype* a, inc_t inca, inc_t lda,
ctype* w, inc_t incw,
ctype* x, inc_t incx,
ctype* beta,
ctype* y, inc_t incy,
ctype* z, inc_t incz
);
Perform
y := beta * y + alpha * conjat(A^T) * conjw(w)
z := z + alpha * conja(A) * conjx(x)
where A
is an m x nf matrix, w
and z
are vectors of length m, x
and y
are vectors of length nf
, and alpha
and beta
are scalars. The kernel, if optimized, is implemented as a fusion of calls to dotxf and axpyf.
Level-2 operations perform various level-2 BLAS-like operations.
void bli_?gemv
(
trans_t transa,
conj_t conjx,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* x, inc_t incx,
ctype* beta,
ctype* y, inc_t incy
);
Perform
y := beta * y + alpha * transa(A) * conjx(x)
where transa(A)
is an m x n matrix, and y
and x
are vectors.
void bli_?ger
(
conj_t conjx,
conj_t conjy,
dim_t m,
dim_t n,
ctype* alpha,
ctype* x, inc_t incx,
ctype* y, inc_t incy,
ctype* a, inc_t rsa, inc_t csa
);
Perform
A := A + alpha * conjx(x) * conjy(y)^T
where A
is an m x n matrix, and x
and y
are vectors of length m and n, respectively.
void bli_?hemv
(
uplo_t uploa,
conj_t conja,
conj_t conjx,
dim_t m,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* x, inc_t incx,
ctype* beta,
ctype* y, inc_t incy
);
Perform
y := beta * y + alpha * conja(A) * conjx(x)
where A
is an m x m Hermitian matrix stored in the lower or upper triangle as specified by uploa
, and y
and x
are vectors of length m.
void bli_?her
(
uplo_t uploa,
conj_t conjx,
dim_t m,
rtype* alpha,
ctype* x, inc_t incx,
ctype* a, inc_t rsa, inc_t csa
);
Perform
A := A + alpha * conjx(x) * conjx(x)^H
where A
is an m x m Hermitian matrix stored in the lower or upper triangle as specified by uploa
, and x
is a vector of length m.
Note: The floating-point type of alpha
is always the real projection of the floating-point types of x
and A
.
void bli_?her2
(
uplo_t uploa,
conj_t conjx,
dim_t m,
ctype* alpha,
ctype* x, inc_t incx,
ctype* y, inc_t incy,
ctype* a, inc_t rsa, inc_t csa
);
Perform
A := A + alpha * conjx(x) * conjy(y)^H + conj(alpha) * conjy(y) * conjx(x)^H
where A
is an m x m Hermitian matrix stored in the lower or upper triangle as specified by uploa
, and x
and y
are vectors of length m.
void bli_?symv
(
uplo_t uploa,
conj_t conja,
conj_t conjx,
dim_t m,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* x, inc_t incx,
ctype* beta,
ctype* y, inc_t incy
);
Perform
y := beta * y + alpha * conja(A) * conjx(x)
where A
is an m x m symmetric matrix stored in the lower or upper triangle as specified by uploa
, and y
and x
are vectors of length m.
void bli_?syr
(
uplo_t uploa,
conj_t conjx,
dim_t m,
ctype* alpha,
ctype* x, inc_t incx,
ctype* a, inc_t rsa, inc_t csa
);
Perform
A := A + alpha * conjx(x) * conjx(x)^T
where A
is an m x m symmetric matrix stored in the lower or upper triangle as specified by uploa
, and x
is a vector of length m.
void bli_?syr2
(
uplo_t uploa,
conj_t conjx,
dim_t m,
ctype* alpha,
ctype* x, inc_t incx,
ctype* y, inc_t incy,
ctype* a, inc_t rsa, inc_t csa
);
Perform
A := A + alpha * conjx(x) * conjy(y)^T + conj(alpha) * conjy(y) * conjx(x)^T
where A
is an m x m symmetric matrix stored in the lower or upper triangle as specified by uploa
, and x
and y
are vectors of length m.
void bli_?trmv
(
uplo_t uploa,
trans_t transa,
diag_t diaga,
dim_t m,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* x, inc_t incx
);
Perform
x := alpha * transa(A) * x
where A
is an m x m triangular matrix stored in the lower or upper triangle as specified by uploa
with unit/non-unit nature specified by diaga
, and x
is a vector of length m.
void bli_?trsv
(
uplo_t uploa,
trans_t transa,
diag_t diaga,
dim_t m,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* y, inc_t incy
);
Solve the linear system
transa(A) * x = alpha * y
where A
is an m x m triangular matrix stored in the lower or upper triangle as specified by uploa
with unit/non-unit nature specified by diaga
, and x
and y
are vectors of length m. The right-hand side vector operand y
is overwritten with the solution vector x
.
Level-3 operations perform various level-3 BLAS-like operations. Note: Each All level-3 operations are implemented through a handful of level-3 micro-kernels. Please see the Kernels wiki for more details.
void bli_?gemm
(
trans_t transa,
trans_t transb,
dim_t m,
dim_t n,
dim_t k,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb,
ctype* beta,
ctype* c, inc_t rsc, inc_t csc
);
Perform
C := beta * C + alpha * transa(A) * transb(B)
where C is an m x n matrix, transa(A)
is an m x k matrix, and transb(B)
is a k x n matrix.
void bli_?hemm
(
side_t sidea,
uplo_t uploa,
conj_t conja,
trans_t transb,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb,
ctype* beta,
ctype* c, inc_t rsc, inc_t csc
);
Perform
C := beta * C + alpha * conja(A) * transb(B)
if sidea
is BLIS_LEFT
, or
C := beta * C + alpha * transb(B) * conja(A)
if sidea
is BLIS_RIGHT
, where C
and B
are m x n matrices and A
is a Hermitian matrix stored in the lower or upper triangle as specified by uploa
. When sidea
is BLIS_LEFT
, A
is m x m, and when sidea
is BLIS_RIGHT
, A
is n x n.
void bli_?herk
(
uplo_t uploc,
trans_t transa,
dim_t m,
dim_t k,
rtype* alpha,
ctype* a, inc_t rsa, inc_t csa,
rtype* beta,
ctype* c, inc_t rsc, inc_t csc
);
Perform
C := beta * C + alpha * transa(A) * transa(A)^H
where C is an m x m Hermitian matrix stored in the lower or upper triangle as specified by uploa
and transa(A)
is an m x k matrix.
Note: The floating-point types of alpha
and beta
are always the real projection of the floating-point types of A
and C
.
void bli_?her2k
(
uplo_t uploc,
trans_t transa,
dim_t m,
dim_t k,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb,
rtype* beta,
ctype* c, inc_t rsc, inc_t csc
);
Perform
C := beta * C + alpha * transa(A) * transb(B)^H + conj(alpha) * transb(B) * transa(A)^H
where C is an m x m Hermitian matrix stored in the lower or upper triangle as specified by uploa
and transa(A)
and transb(B)
are m x k matrices.
Note: The floating-point type of beta
is always the real projection of the floating-point types of A
and C
.
void bli_?symm
(
side_t sidea,
uplo_t uploa,
conj_t conja,
trans_t transb,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb,
ctype* beta,
ctype* c, inc_t rsc, inc_t csc
);
Perform
C := beta * C + alpha * conja(A) * transb(B)
if sidea
is BLIS_LEFT
, or
C := beta * C + alpha * transb(B) * conja(A)
if sidea
is BLIS_RIGHT
, where C
and B
are m x n matrices and A
is a symmetric matrix stored in the lower or upper triangle as specified by uploa
. When sidea
is BLIS_LEFT
, A
is m x m, and when sidea
is BLIS_RIGHT
, A
is n x n.
void bli_?syrk
(
uplo_t uploc,
trans_t transa,
dim_t m,
dim_t k,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* beta,
ctype* c, inc_t rsc, inc_t csc
);
Perform
C := beta * C + alpha * transa(A) * transa(A)^T
where C is an m x m symmetric matrix stored in the lower or upper triangle as specified by uploa
and transa(A)
is an m x k matrix.
void bli_?syr2k
(
uplo_t uploc,
trans_t transa,
dim_t m,
dim_t k,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb,
ctype* beta,
ctype* c, inc_t rsc, inc_t csc
);
Perform
C := beta * C + alpha * transa(A) * transb(B)^T + alpha * transb(B) * transa(A)^T
where C is an m x m symmetric matrix stored in the lower or upper triangle as specified by uploa
and transa(A)
and transb(B)
are m x k matrices.
void bli_?trmm
(
side_t sidea,
uplo_t uploa,
trans_t transa,
diag_t diaga,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb
);
Perform
B := alpha * transa(A) * B
if sidea
is BLIS_LEFT
, or
B := alpha * B * transa(A)
if sidea
is BLIS_RIGHT
, where B
is an m x n matrix and A
is a triangular matrix stored in the lower or upper triangle as specified by uploa
with unit/non-unit nature specified by diaga
. When sidea
is BLIS_LEFT
, A
is m x m, and when sidea
is BLIS_RIGHT
, A
is n x n.
void bli_?trmm3
(
side_t sidea,
uplo_t uploa,
trans_t transa,
diag_t diaga,
trans_t transb,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb,
ctype* beta,
ctype* c, inc_t rsc, inc_t csc
);
Perform
C := beta * C + alpha * transa(A) * transb(B)
if sidea
is BLIS_LEFT
, or
C := beta * C + alpha * transb(B) * transa(A)
if sidea
is BLIS_RIGHT
, where C
and transb(B)
are m x n matrices and A
is a triangular matrix stored in the lower or upper triangle as specified by uploa
with unit/non-unit nature specified by diaga
. When sidea
is BLIS_LEFT
, A
is m x m, and when sidea
is BLIS_RIGHT
, A
is n x n.
void bli_?trsm
(
side_t sidea,
uplo_t uploa,
trans_t transa,
diag_t diaga,
dim_t m,
dim_t n,
ctype* alpha,
ctype* a, inc_t rsa, inc_t csa,
ctype* b, inc_t rsb, inc_t csb
);
Solve the linear system with multiple right-hand sides
transa(A) * X = alpha * B
if sidea
is BLIS_LEFT
, or
X * transa(A) = alpha * B
if sidea
is BLIS_RIGHT
, where X
and B
are an m x n matrices and A
is a triangular matrix stored in the lower or upper triangle as specified by uploa
with unit/non-unit nature specified by diaga
. When sidea
is BLIS_LEFT
, A
is m x m, and when sidea
is BLIS_RIGHT
, A
is n x n. The right-hand side matrix operand B
is overwritten with the solution matrix X
.
void bli_?asumv
(
dim_t n,
ctype* x, inc_t incx,
rtype* asum
);
Compute the su of the absolute values of the fundamental elements of vector x
. The resulting sum is stored to asum
.
Note: The floating-point type of asum
is always the real projection of the floating-point type of x
.
Note: This function attempts to mimic the algorithm for computing the absolute vector sum in the netlib BLAS routines *asum()
.
void bli_?norm[1fi]m
(
doff_t diagoffa,
doff_t diaga,
uplo_t uploa,
dim_t m,
dim_t n,
ctype* a, inc_t rs_a, inc_t cs_a,
rtype* norm
);
Compute the one-norm (norm1m
), Frobenius norm (normfm
), or infinity norm (normim
) of the elements in an m x n matrix A
. If uploa
is BLIS_LOWER
or BLIS_UPPER
then A
is assumed to be lower or upper triangular, respectively, with the main diagonal located at offset diagoffa
. The resulting norm is stored to norm
.
Note: The floating-point type of norm
is always the real projection of the floating-point type of x
.
void bli_?norm[1fi]v
(
dim_t n,
ctype* x, inc_t incx,
rtype* norm
);
Compute the one-norm (norm1v
), Frobenius norm (normfv
), or infinity norm (normiv
) of the elements in a vector x
of length n. The resulting norm is stored to norm
.
Note: The floating-point type of norm
is always the real projection of the floating-point type of x
.
void bli_?mkherm
(
uplo_t uploa,
dim_t m,
ctype* a, inc_t rs_a, inc_t cs_a
);
Make an m x m matrix A
explicitly Hermitian by copying the conjugate of the triangle specified by uploa
to the opposite triangle. Imaginary components of diagonal elements are explicitly set to zero. It is assumed that the diagonal offset of A
is zero.
void bli_?mksymm
(
uplo_t uploa,
dim_t m,
ctype* a, inc_t rs_a, inc_t cs_a
);
Make an m x m matrix A
explicitly symmetric by copying the triangle specified by uploa
to the opposite triangle. It is assumed that the diagonal offset of A
is zero.
void bli_?mktrim
(
uplo_t uploa,
dim_t m,
ctype* a, inc_t rs_a, inc_t cs_a
);
Make an m x m matrix A
explicitly triangular by preserving the triangle specified by uploa
and zeroing the elements in the opposite triangle. It is assumed that the diagonal offset of A
is zero.
void bli_?fprintv
(
FILE* file,
char* s1,
dim_t m,
ctype* x, inc_t incx,
char* format,
char* s2
);
Print a vector x
of length m to file stream file
, where file
is a file pointer returned by the standard C library function fopen()
. The caller may also pass in a global file pointer such as stdout
or stderr
. The strings s1
and s2
are printed immediately before and after the output (respectively), and the format specifier format
is used to format the individual elements. For valid format specifiers, please see documentation for the standard C library function printf()
.
Note: For complex datatypes, the format specifier is applied to both the real and imaginary components individually. Therefore, you should use format specifiers such as "%5.2f"
, but not "%5.2f + %5.2f"
.
void bli_?fprintm
(
FILE* file,
char* s1,
dim_t m,
dim_t n,
ctype* a, inc_t rs_a, inc_t cs_a,
char* format,
char* s2
);
Print an m x n matrix A
to file stream file
, where file
is a file pointer returned by the standard C library function fopen()
. The caller may also pass in a global file pointer such as stdout
or stderr
. The strings s1
and s2
are printed immediately before and after the output (respectively), and the format specifier format
is used to format the individual elements. For valid format specifiers, please see documentation for the standard C library function printf()
.
Note: For complex datatypes, the format specifier is applied to both the real and imaginary components individually. Therefore, you should use format specifiers such as "%5.2f"
, but not "%5.2f + %5.2f"
.
void bli_?printv
(
char* s1,
dim_t m,
ctype* x, inc_t incx,
char* format,
char* s2
);
Print a vector x
of length m to standard output. This function call is equivalent to calling bli_?fprintv()
with stdout
as the file pointer.
void bli_?printm
(
char* s1,
dim_t m,
dim_t n,
ctype* a, inc_t rs_a, inc_t cs_a,
char* format,
char* s2
);
Print an m x n matrix a
to standard output. This function call is equivalent to calling bli_?fprintm()
with stdout
as the file pointer.
void bli_?randv
(
dim_t m,
ctype* x, inc_t incx
);
Set the elements of a vector x
of length m to random values on the interval [-1,1)
.
Note: For complex datatypes, the real and imaginary components of each element are randomized individually and independently of one another.
void bli_?randm
(
doff_t diagoffa,
uplo_t uploa,
dim_t m,
dim_t n,
ctype* a, inc_t rs_a, inc_t cs_a
);
Set the elements of an m x n matrix A
to random values on the interval [-1,1)
. If uploa
is BLIS_LOWER
or BLIS_UPPER
, then additional scaling occurs so that the resulting matrix is diagonally dominant. Specifically, the diagonal elements (identified by diagonal offset diagoffa
) are shifted so that they lie on the interval [1,2)
and the off-diagonal elements (in the triangle specified by uploa
) are scaled by 1.0/max(m,n)
.
Note: For complex datatypes, the real and imaginary components of each off-diagonal element are randomized individually and independently of one another.
void bli_?sumsqv
(
dim_t m,
ctype* x, inc_t incx,
rtype* scale,
rtype* sumsq
);
Compute the sum of the squares of the elements in a vector x
of length m. The result is computed in scaled form, and in such a way that it may be used repeatedly to accumulate the sum of the squares of several vectors.
The function computes scale_new and sumsq_new such that
scale_new^2 * sumsq_new = x[0]^2 + x[1]^2 + ... x[m-1]^2 + scale_old^2 * sumsq_old
where, on entry, scale
and sumsq
contain scale_old
and sumsq_old
, respectively, and on exit, scale
and sumsq
contain scale_new
and sumsq_new
, respectively.
Note: This function attempts to mimic the algorithm for computing the Frobenius norm in the netlib LAPACK routine ?lassq()
.
Note: The *
in level-3 micro-kernel function names shown below reflect that there is no exact naming convention required for the micro-kernels, except that they must begin with bli_?
. We strongly recommend, however, that the micro-kernel function names include the name of the micro-kernel itself. For example, the gemm
micro-kernel should be named with the prefix bli_?gemm_
and the trsm
micro-kernels should be named with the prefixes bli_?trsm_l_
(lower triangular) and bli_?trsm_u_
(upper triangular).
void bli_?gemm_*
(
dim_t k,
ctype* restrict alpha,
ctype* restrict a1,
ctype* restrict b1,
ctype* restrict beta,
ctype* restrict c11, inc_t rsc, inc_t csc,
auxinfo_t* restrict data,
cntx_t* restrict cntx
);
Perform
C11 := beta * C11 + alpha * A1 * B1
where C11
is an MR x NR matrix, A1
is an MR x k "micro-panel" matrix stored in packed (column-stored) format, B1
is a k x NR "micro-panel" matrix in packed (row-stored) format, and alpha and beta are scalars. The storage of C11
is specified by its row and column strides, rsc
and csc
.
Please see the BLIS Kernel wiki for more information on the gemm
micro-kernel.
void bli_?trsm_l_*
(
ctype* restrict a11,
ctype* restrict b11,
ctype* restrict c11, inc_t rsc, inc_t csc
auxinfo_t* restrict data,
cntx_t* restrict cntx
);
void bli_?trsm_u_*
(
ctype* restrict a11,
ctype* restrict b11,
ctype* restrict c11, inc_t rsc, inc_t csc
auxinfo_t* restrict data,
cntx_t* restrict cntx
);
Perform
B11 := inv(A11) * B11
C11 := B11
where A11
is an MR x MR lower or upper triangular matrix stored in packed (column-stored) format, B11
is an MR x NR matrix stored in packed (row-stored) format, and C11
is an MR x NR matrix stored according to row and column strides rsc
and csc
.
Please see the BLIS Kernel wiki for more information on the trsm
micro-kernel.
void bli_?gemmtrsm_l_*
(
dim_t k,
ctype* restrict alpha,
ctype* restrict a10,
ctype* restrict a11,
ctype* restrict b01,
ctype* restrict b11,
ctype* restrict c11, inc_t rs_c, inc_t cs_c,
auxinfo_t* restrict data,
cntx_t* restrict cntx
);
void bli_?gemmtrsm_u_*
(
dim_t k,
ctype* restrict alpha,
ctype* restrict a12,
ctype* restrict a11,
ctype* restrict b21,
ctype* restrict b11,
ctype* restrict c11, inc_t rs_c, inc_t cs_c,
auxinfo_t* restrict data,
cntx_t* restrict cntx
);
Perform
B11 := alpha * B11 - A10 * B01
B11 := inv(A11) * B11
C11 := B11
if A11
is lower triangular, or
B11 := alpha * B11 - A12 * B21
B11 := inv(A11) * B11
C11 := B11
if A11
is upper triangular.
Please see the BLIS Kernel wiki for more information on the gemmtrsm
micro-kernel.
BLIS allows applications to query information about how BLIS was configured. The bli_info_
API provides several categories of query routines. Most values are returned as a gint_t
, which is a signed integer. The size of this integer can be queried through a special routine that returns the size in a character string:
char* bli_info_get_int_type_size_str( void );
Note: All of the bli_info_
functions are always thread-safe, no matter how BLIS was configured.
The following routine returns the address the full BLIS version string:
char* bli_info_get_version_str( void );
The following routine returns a unique ID of type arch_t
that identifies the current current active configuration:
arch_t bli_arch_query_id( void );
This is most useful when BLIS is configured with multiple configurations. (When linking to multi-configuration builds of BLIS, you don't know for sure which configuration will be used until runtime since the configuration-specific parameters are not loaded until after calling a hueristic to detect the hardware--usually based the CPUID
instruction.)
Once the configuration's ID is known, it can be used to query a string that contains the name of the configuration:
char* bli_arch_string( arch_t id );
The following routines return various general-purpose constants that affect the entire framework. All of these settings default to sane values, which can then be overridden by the configuration in bli_config.h. If they are absent from a particular configuration's bli_config.h
header file, then the default value is used, as specified in frame/include/bli_config_macro_defs.h.
gint_t bli_info_get_int_type_size( void );
gint_t bli_info_get_num_fp_types( void );
gint_t bli_info_get_max_type_size( void );
gint_t bli_info_get_page_size( void );
gint_t bli_info_get_simd_num_registers( void );
gint_t bli_info_get_simd_size( void );
gint_t bli_info_get_simd_align_size( void );
gint_t bli_info_get_stack_buf_max_size( void );
gint_t bli_info_get_stack_buf_align_size( void );
gint_t bli_info_get_heap_addr_align_size( void );
gint_t bli_info_get_heap_stride_align_size( void );
gint_t bli_info_get_pool_addr_align_size( void );
gint_t bli_info_get_enable_stay_auto_init( void );
gint_t bli_info_get_enable_blas( void );
gint_t bli_info_get_blas_int_type_size( void );
The following routines allow the caller to obtain a string that identifies the implementation type of each micro-kernel that is currently active (ie: part of the current active configuration, as identified bi bli_arch_query_id()
).
char* bli_info_get_gemm_ukr_impl_string( ind_t method, num_t dt )
char* bli_info_get_gemmtrsm_l_ukr_impl_string( ind_t method, num_t dt )
char* bli_info_get_gemmtrsm_u_ukr_impl_string( ind_t method, num_t dt )
char* bli_info_get_trsm_l_ukr_impl_string( ind_t method, num_t dt )
char* bli_info_get_trsm_u_ukr_impl_string( ind_t method, num_t dt )
Possible implementation (ie: the ind_t method
argument) types are:
BLIS_3MH
: Implementation based on the 3m method applied at the highest level, outside the 5th loop around the micro-kernel.BLIS_3M1
: Implementation based on the 3m method applied within the 1st loop around the micro-kernel.BLIS_4MH
: Implementation based on the 4m method applied at the highest level, outside the 5th loop around the micro-kernel.BLIS_4M1B
: Implementation based on the 4m method applied within the 1st loop around the micro-kernel. Computation is ordered such that the 1st loop is fissured into two loops, the first of which multiplies the real part of the current micro-panel of packed matrix B (against all real and imaginary parts of packed matrix A), and the second of which multiplies the imaginary part of the current micro-panel of packed matrix B.BLIS_4M1A
: Implementation based on the 4m method applied within the 1st loop around the micro-kernel. Computation is ordered such that real and imaginary components of the current micro-panels are completely used before proceeding to the next virtual micro-kernel invocation.BLIS_1M
: Implementation based on the 1m method. (This is the default induced method when real domain kernels are present but complex kernels are missing.)BLIS_NAT
: Implementation based on "native" execution (ie: NOT an induced method).
NOTE: BLIS_3M3
and BLIS_3M2
have been deprecated from the typedef enum
of ind_t
, and BLIS_4M1B
is also effectively no longer available, though the typedef enum
value still exists.
Possible micro-kernel types (ie: the return values for bli_info_get_*_ukr_impl_string()
) are:
BLIS_REFERENCE_UKERNEL
("refrnce"
): This value is returned when the queried micro-kernel is provided by the reference implementation.BLIS_VIRTUAL_UKERNEL
("virtual"
): This value is returned when the queried micro-kernel is driven by a the "virtual" micro-kernel provided by an induced method. This happens for anymethod
value that is notBLIS_NAT
(ie: native), but only applies to the complex domain.BLIS_OPTIMIZED_UKERNEL
("optimzd"
): This value is returned when the queried micro-kernel is provided by an implementation that is neither reference nor virtual, and thus we assume the kernel author would deem it to be "optimized". Such a micro-kernel may not be optimal in the literal sense of the word, but nonetheless is intended to be optimized, at least relative to the reference micro-kernels.BLIS_NOTAPPLIC_UKERNEL
("notappl"
): This value is returned usually when performing agemmtrsm
ortrsm
micro-kernel type query for anymethod
value that is notBLIS_NAT
(ie: native). That is, induced methods cannot be (purely) used ontrsm
-based micro-kernels because these micro-kernels perform more a triangular inversion, which is not matrix multiplication.
The following routines allow the caller to obtain a string that identifies the implementation (ind_t
) that is currently active (ie: implemented and enabled) for each level-3 operation. Possible implementation types are listed in the section above covering micro-kernel implemenation query.
char* bli_info_get_gemm_impl_string( num_t dt );
char* bli_info_get_hemm_impl_string( num_t dt );
char* bli_info_get_herk_impl_string( num_t dt );
char* bli_info_get_her2k_impl_string( num_t dt );
char* bli_info_get_symm_impl_string( num_t dt );
char* bli_info_get_syrk_impl_string( num_t dt );
char* bli_info_get_syr2k_impl_string( num_t dt );
char* bli_info_get_trmm_impl_string( num_t dt );
char* bli_info_get_trmm3_impl_string( num_t dt );
char* bli_info_get_trsm_impl_string( num_t dt );