mirror of
https://github.com/amd/blis.git
synced 2026-05-11 17:50:00 +00:00
Merged contributions from AMD's AOCL BLIS (#448). Details: - Added support for level-3 operation gemmt, which performs a gemm on only the lower or upper triangle of a square matrix C. For now, only the conventional/large code path will be supported (in vanilla BLIS). This was accomplished by leveraging the existing variant logic for herk. However, some of the infrastructure to support a gemmtsup is included in this commit, including - A bli_gemmtsup() front-end, similar to bli_gemmsup(). - A bli_gemmtsup_ref() reference handler function. - A bli_gemmtsup_int() variant chooser function (with variant calls commented out). - Added support for inducing complex domain gemmt via the 1m method. - Added gemmt APIs to the BLAS and CBLAS compatiblity layers. - Added gemmt test module to testsuite. - Added standalone gemmt test driver to 'test' directory. - Documented gemmt APIs in BLISObjectAPI.md and BLISTypedAPI.md. - Added a C++ template header (blis.hh) containing a BLAS-inspired wrapper to a set of polymorphic CBLAS-like function wrappers defined in another header (cblas.hh). These two headers are installed if running the 'install' target with INSTALL_HH is set to 'yes'. (Also added a set of unit tests that exercise blis.hh, although they are disabled for now because they aren't compatible with out-of-tree builds.) These files now live in the 'vendor' top-level directory. - Various updates to 'zen' and 'zen2' subconfigurations, particularly within the context initialization functions. - Added s and d copyv, setv, and swapv kernels to kernels/zen/1, and various minor updates to dotv and scalv kernels. Also added various sup kernels contributed by AMD to kernels/zen/3. However, these kernels are (for now) not yet used, in part because they caused AppVeyor clang failures, and also because I have not found time to review and vet them. - Output the python found during configure into the definition of PYTHON in build/config.mk (via build/config.mk.in). - Added early-return checks (A, B, or C with zero dimension; alpha = 0) to bli_gemm_front.c. - Implemented explicit beta = 0 handling in for the sgemm ukernel in bli_gemm_armv7a_int_d4x4.c, which was previously missing. This latent bug surfaced because the gemmt module verifies its computation using gemm with its beta parameter set to zero, which, on a cortexa15 system caused the gemm kernel code to unconditionally multiply the uninitialized C data by beta. The C matrix likely contained non-numeric values such as NaN, which then would have resulted in a false failure. - Fixed a bug whereby the implementation for bli_herk_determine_kc(), in bli_l3_blocksize.c, was inadvertantly being defined in terms of helper functions meant for trmm. This bug was probably harmless since the trmm code should have also done the right thing for herk. - Used cpp macros to neutralize the various AOCL_DTL_TRACE_ macros in kernels/zen/3/bli_gemm_small.c since those macros are not used in vanilla BLIS. - Added cpp guard to definition of bli_mem_clear() in bli_mem.h to accommodate C++'s stricter type checking. - Added cpp guard to test/*.c drivers that facilitate compilation on Windows systems. - Various whitespace changes.
484 lines
12 KiB
C
484 lines
12 KiB
C
/*
|
|
|
|
BLIS
|
|
An object-based framework for developing high-performance BLAS-like
|
|
libraries.
|
|
|
|
Copyright (C) 2014, The University of Texas at Austin
|
|
Copyright (C) 2019 - 2020, Advanced Micro Devices, Inc.
|
|
|
|
Redistribution and use in source and binary forms, with or without
|
|
modification, are permitted provided that the following conditions are
|
|
met:
|
|
- Redistributions of source code must retain the above copyright
|
|
notice, this list of conditions and the following disclaimer.
|
|
- Redistributions in binary form must reproduce the above copyright
|
|
notice, this list of conditions and the following disclaimer in the
|
|
documentation and/or other materials provided with the distribution.
|
|
- Neither the name of The University of Texas nor the names of its
|
|
contributors may be used to endorse or promote products derived
|
|
from this software without specific prior written permission.
|
|
|
|
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
|
|
"AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
|
|
LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
|
|
A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
|
|
HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
|
|
SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
|
|
LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
|
|
DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
|
|
THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
|
|
(INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
|
|
OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
|
|
|
|
*/
|
|
|
|
#ifdef WIN32
|
|
#include <io.h>
|
|
#else
|
|
#include <unistd.h>
|
|
#endif
|
|
#include "blis.h"
|
|
|
|
//#define CBLAS
|
|
//#define C_STOR_R
|
|
|
|
//#define PRINT
|
|
|
|
int main( int argc, char** argv )
|
|
{
|
|
obj_t a, b, c;
|
|
obj_t c_save;
|
|
obj_t alpha, beta;
|
|
dim_t m, k;
|
|
dim_t p;
|
|
dim_t p_begin, p_end, p_inc;
|
|
int m_input, k_input;
|
|
num_t dt;
|
|
int r, n_repeats;
|
|
uplo_t uploc;
|
|
trans_t transa;
|
|
trans_t transb;
|
|
f77_char f77_uploc;
|
|
f77_char f77_transa;
|
|
f77_char f77_transb;
|
|
|
|
double dtime;
|
|
double dtime_save;
|
|
double gflops;
|
|
|
|
//bli_init();
|
|
|
|
//bli_error_checking_level_set( BLIS_NO_ERROR_CHECKING );
|
|
|
|
n_repeats = 3;
|
|
|
|
#ifndef PRINT
|
|
p_begin = 200;
|
|
p_end = 2000;
|
|
p_inc = 200;
|
|
|
|
m_input = -1;
|
|
k_input = -1;
|
|
#else
|
|
p_begin = 16;
|
|
p_end = 16;
|
|
p_inc = 1;
|
|
|
|
m_input = 5;
|
|
k_input = 4;
|
|
#endif
|
|
|
|
#if 1
|
|
//dt = BLIS_FLOAT;
|
|
dt = BLIS_DOUBLE;
|
|
#else
|
|
//dt = BLIS_SCOMPLEX;
|
|
dt = BLIS_DCOMPLEX;
|
|
#endif
|
|
|
|
uploc = BLIS_LOWER;
|
|
//uploc = BLIS_UPPER;
|
|
|
|
transa = BLIS_NO_TRANSPOSE;
|
|
transb = BLIS_NO_TRANSPOSE;
|
|
|
|
bli_param_map_blis_to_netlib_uplo( uploc, &f77_uploc );
|
|
bli_param_map_blis_to_netlib_trans( transa, &f77_transa );
|
|
bli_param_map_blis_to_netlib_trans( transb, &f77_transb );
|
|
|
|
char uplocl = tolower( f77_uploc );
|
|
char transal = tolower( f77_transa );
|
|
char transbl = tolower( f77_transb );
|
|
|
|
f77_int cbla_uploc = ( uplocl == 'l' ? CblasLower : CblasUpper );
|
|
f77_int cbla_transa = ( transal == 'n' ? CblasNoTrans : CblasTrans );
|
|
f77_int cbla_transb = ( transbl == 'n' ? CblasNoTrans : CblasTrans );
|
|
|
|
( void )cbla_uploc;
|
|
( void )cbla_transa;
|
|
( void )cbla_transb;
|
|
|
|
// Begin with initializing the last entry to zero so that
|
|
// matlab allocates space for the entire array once up-front.
|
|
for ( p = p_begin; p + p_inc <= p_end; p += p_inc ) ;
|
|
#ifdef BLIS
|
|
printf( "data_gemmt_blis" );
|
|
#else
|
|
printf( "data_gemmt_%s", BLAS );
|
|
#endif
|
|
printf( "( %2lu, 1:3 ) = [ %4lu %4lu %7.2f ];\n",
|
|
( unsigned long )(p - p_begin)/p_inc + 1,
|
|
( unsigned long )0,
|
|
( unsigned long )0, 0.0 );
|
|
|
|
//for ( p = p_begin; p <= p_end; p += p_inc )
|
|
for ( p = p_end; p_begin <= p; p -= p_inc )
|
|
{
|
|
if ( m_input < 0 ) m = p * ( dim_t )abs(m_input);
|
|
else m = ( dim_t ) m_input;
|
|
if ( k_input < 0 ) k = p * ( dim_t )abs(k_input);
|
|
else k = ( dim_t ) k_input;
|
|
|
|
bli_obj_create( dt, 1, 1, 0, 0, &alpha );
|
|
bli_obj_create( dt, 1, 1, 0, 0, &beta );
|
|
|
|
#ifndef C_STOR_R
|
|
if ( bli_does_trans( transa ) )
|
|
bli_obj_create( dt, k, m, 0, 0, &a );
|
|
else
|
|
bli_obj_create( dt, m, k, 0, 0, &a );
|
|
|
|
if ( bli_does_trans( transb ) )
|
|
bli_obj_create( dt, m, k, 0, 0, &b );
|
|
else
|
|
bli_obj_create( dt, k, m, 0, 0, &b );
|
|
|
|
bli_obj_create( dt, m, m, 0, 0, &c );
|
|
bli_obj_create( dt, m, m, 0, 0, &c_save );
|
|
#else
|
|
if ( bli_does_trans( transa ) )
|
|
bli_obj_create( dt, k, m, -1, -1, &a );
|
|
else
|
|
bli_obj_create( dt, m, k, -1, -1, &a );
|
|
|
|
if ( bli_does_trans( transb ) )
|
|
bli_obj_create( dt, m, k, -1, -1, &b );
|
|
else
|
|
bli_obj_create( dt, k, m, -1, -1, &b );
|
|
|
|
bli_obj_create( dt, m, m, -1, -1, &c );
|
|
bli_obj_create( dt, m, m, -1, -1, &c_save );
|
|
#endif
|
|
|
|
bli_randm( &a );
|
|
bli_randm( &b );
|
|
bli_randm( &c );
|
|
|
|
bli_obj_set_uplo( uploc, &c );
|
|
|
|
bli_obj_set_conjtrans( transa, &a );
|
|
bli_obj_set_conjtrans( transb, &b );
|
|
|
|
bli_setsc( (0.9/1.0), 0.2, &alpha );
|
|
bli_setsc( -(1.1/1.0), 0.3, &beta );
|
|
|
|
|
|
bli_copym( &c, &c_save );
|
|
|
|
dtime_save = DBL_MAX;
|
|
|
|
for ( r = 0; r < n_repeats; ++r )
|
|
{
|
|
bli_copym( &c_save, &c );
|
|
|
|
|
|
dtime = bli_clock();
|
|
|
|
|
|
#ifdef PRINT
|
|
bli_printm( "a", &a, "%4.1f", "" );
|
|
bli_printm( "b", &b, "%4.1f", "" );
|
|
bli_printm( "c", &c, "%4.1f", "" );
|
|
#endif
|
|
|
|
#ifdef BLIS
|
|
|
|
bli_gemmt( &alpha,
|
|
&a,
|
|
&b,
|
|
&beta,
|
|
&c );
|
|
|
|
#else
|
|
|
|
#ifndef CBLAS
|
|
|
|
if ( bli_is_float( dt ) )
|
|
{
|
|
f77_int mm = bli_obj_length( &c );
|
|
f77_int kk = bli_obj_width_after_trans( &a );
|
|
f77_int lda = bli_obj_col_stride( &a );
|
|
f77_int ldb = bli_obj_col_stride( &b );
|
|
f77_int ldc = bli_obj_col_stride( &c );
|
|
float* alphap = bli_obj_buffer( &alpha );
|
|
float* ap = bli_obj_buffer( &a );
|
|
float* bp = bli_obj_buffer( &b );
|
|
float* betap = bli_obj_buffer( &beta );
|
|
float* cp = bli_obj_buffer( &c );
|
|
|
|
sgemmt_( &f77_uploc,
|
|
&f77_transa,
|
|
&f77_transb,
|
|
&mm,
|
|
&kk,
|
|
alphap,
|
|
ap, &lda,
|
|
bp, &ldb,
|
|
betap,
|
|
cp, &ldc );
|
|
}
|
|
else if ( bli_is_double( dt ) )
|
|
{
|
|
f77_int mm = bli_obj_length( &c );
|
|
f77_int kk = bli_obj_width_after_trans( &a );
|
|
f77_int lda = bli_obj_col_stride( &a );
|
|
f77_int ldb = bli_obj_col_stride( &b );
|
|
f77_int ldc = bli_obj_col_stride( &c );
|
|
double* alphap = bli_obj_buffer( &alpha );
|
|
double* ap = bli_obj_buffer( &a );
|
|
double* bp = bli_obj_buffer( &b );
|
|
double* betap = bli_obj_buffer( &beta );
|
|
double* cp = bli_obj_buffer( &c );
|
|
|
|
dgemmt_( &f77_uploc,
|
|
&f77_transa,
|
|
&f77_transb,
|
|
&mm,
|
|
&kk,
|
|
alphap,
|
|
ap, &lda,
|
|
bp, &ldb,
|
|
betap,
|
|
cp, &ldc );
|
|
}
|
|
else if ( bli_is_scomplex( dt ) )
|
|
{
|
|
f77_int mm = bli_obj_length( &c );
|
|
f77_int kk = bli_obj_width_after_trans( &a );
|
|
f77_int lda = bli_obj_col_stride( &a );
|
|
f77_int ldb = bli_obj_col_stride( &b );
|
|
f77_int ldc = bli_obj_col_stride( &c );
|
|
scomplex* alphap = bli_obj_buffer( &alpha );
|
|
scomplex* ap = bli_obj_buffer( &a );
|
|
scomplex* bp = bli_obj_buffer( &b );
|
|
scomplex* betap = bli_obj_buffer( &beta );
|
|
scomplex* cp = bli_obj_buffer( &c );
|
|
|
|
cgemmt_( &f77_uploc,
|
|
&f77_transa,
|
|
&f77_transb,
|
|
&mm,
|
|
&kk,
|
|
alphap,
|
|
ap, &lda,
|
|
bp, &ldb,
|
|
betap,
|
|
cp, &ldc );
|
|
}
|
|
else if ( bli_is_dcomplex( dt ) )
|
|
{
|
|
f77_int mm = bli_obj_length( &c );
|
|
f77_int kk = bli_obj_width_after_trans( &a );
|
|
f77_int lda = bli_obj_col_stride( &a );
|
|
f77_int ldb = bli_obj_col_stride( &b );
|
|
f77_int ldc = bli_obj_col_stride( &c );
|
|
dcomplex* alphap = bli_obj_buffer( &alpha );
|
|
dcomplex* ap = bli_obj_buffer( &a );
|
|
dcomplex* bp = bli_obj_buffer( &b );
|
|
dcomplex* betap = bli_obj_buffer( &beta );
|
|
dcomplex* cp = bli_obj_buffer( &c );
|
|
|
|
zgemmt_( &f77_uploc,
|
|
&f77_transa,
|
|
&f77_transb,
|
|
&mm,
|
|
&kk,
|
|
alphap,
|
|
ap, &lda,
|
|
bp, &ldb,
|
|
betap,
|
|
cp, &ldc );
|
|
}
|
|
|
|
#else // #ifdef CBLAS
|
|
|
|
f77_int cbla_storage = ( bli_obj_is_row_stored( &c ) ? CblasRowMajor
|
|
: CblasColMajor );
|
|
|
|
if ( bli_is_float( dt ) )
|
|
{
|
|
f77_int mm = bli_obj_length( &c );
|
|
f77_int kk = bli_obj_width_after_trans( &a );
|
|
#ifdef C_STOR_R
|
|
f77_int lda = bli_obj_row_stride( &a );
|
|
f77_int ldb = bli_obj_row_stride( &b );
|
|
f77_int ldc = bli_obj_row_stride( &c );
|
|
#else
|
|
f77_int lda = bli_obj_col_stride( &a );
|
|
f77_int ldb = bli_obj_col_stride( &b );
|
|
f77_int ldc = bli_obj_col_stride( &c );
|
|
#endif
|
|
float* alphap = bli_obj_buffer( &alpha );
|
|
float* ap = bli_obj_buffer( &a );
|
|
float* bp = bli_obj_buffer( &b );
|
|
float* betap = bli_obj_buffer( &beta );
|
|
float* cp = bli_obj_buffer( &c );
|
|
|
|
cblas_sgemmt( cbla_storage,
|
|
cbla_uploc,
|
|
cbla_transa,
|
|
cbla_transb,
|
|
mm,
|
|
kk,
|
|
*alphap,
|
|
ap, lda,
|
|
bp, ldb,
|
|
*betap,
|
|
cp, ldc );
|
|
}
|
|
else if ( bli_is_double( dt ) )
|
|
{
|
|
f77_int mm = bli_obj_length( &c );
|
|
f77_int kk = bli_obj_width_after_trans( &a );
|
|
#ifdef C_STOR_R
|
|
f77_int lda = bli_obj_row_stride( &a );
|
|
f77_int ldb = bli_obj_row_stride( &b );
|
|
f77_int ldc = bli_obj_row_stride( &c );
|
|
#else
|
|
f77_int lda = bli_obj_col_stride( &a );
|
|
f77_int ldb = bli_obj_col_stride( &b );
|
|
f77_int ldc = bli_obj_col_stride( &c );
|
|
#endif
|
|
double* alphap = bli_obj_buffer( &alpha );
|
|
double* ap = bli_obj_buffer( &a );
|
|
double* bp = bli_obj_buffer( &b );
|
|
double* betap = bli_obj_buffer( &beta );
|
|
double* cp = bli_obj_buffer( &c );
|
|
|
|
cblas_dgemmt( cbla_storage,
|
|
cbla_uploc,
|
|
cbla_transa,
|
|
cbla_transb,
|
|
mm,
|
|
kk,
|
|
*alphap,
|
|
ap, lda,
|
|
bp, ldb,
|
|
*betap,
|
|
cp, ldc );
|
|
}
|
|
else if ( bli_is_scomplex( dt ) )
|
|
{
|
|
f77_int mm = bli_obj_length( &c );
|
|
f77_int kk = bli_obj_width_after_trans( &a );
|
|
#ifdef C_STOR_R
|
|
f77_int lda = bli_obj_row_stride( &a );
|
|
f77_int ldb = bli_obj_row_stride( &b );
|
|
f77_int ldc = bli_obj_row_stride( &c );
|
|
#else
|
|
f77_int lda = bli_obj_col_stride( &a );
|
|
f77_int ldb = bli_obj_col_stride( &b );
|
|
f77_int ldc = bli_obj_col_stride( &c );
|
|
#endif
|
|
scomplex* alphap = bli_obj_buffer( &alpha );
|
|
scomplex* ap = bli_obj_buffer( &a );
|
|
scomplex* bp = bli_obj_buffer( &b );
|
|
scomplex* betap = bli_obj_buffer( &beta );
|
|
scomplex* cp = bli_obj_buffer( &c );
|
|
|
|
cblas_cgemmt( cbla_storage,
|
|
cbla_uploc,
|
|
cbla_transa,
|
|
cbla_transb,
|
|
mm,
|
|
kk,
|
|
alphap,
|
|
ap, lda,
|
|
bp, ldb,
|
|
betap,
|
|
cp, ldc );
|
|
}
|
|
else if ( bli_is_dcomplex( dt ) )
|
|
{
|
|
f77_int mm = bli_obj_length( &c );
|
|
f77_int kk = bli_obj_width_after_trans( &a );
|
|
#ifdef C_STOR_R
|
|
f77_int lda = bli_obj_row_stride( &a );
|
|
f77_int ldb = bli_obj_row_stride( &b );
|
|
f77_int ldc = bli_obj_row_stride( &c );
|
|
#else
|
|
f77_int lda = bli_obj_col_stride( &a );
|
|
f77_int ldb = bli_obj_col_stride( &b );
|
|
f77_int ldc = bli_obj_col_stride( &c );
|
|
#endif
|
|
dcomplex* alphap = bli_obj_buffer( &alpha );
|
|
dcomplex* ap = bli_obj_buffer( &a );
|
|
dcomplex* bp = bli_obj_buffer( &b );
|
|
dcomplex* betap = bli_obj_buffer( &beta );
|
|
dcomplex* cp = bli_obj_buffer( &c );
|
|
|
|
cblas_zgemmt( cbla_storage,
|
|
cbla_uploc,
|
|
cbla_transa,
|
|
cbla_transb,
|
|
mm,
|
|
kk,
|
|
alphap,
|
|
ap, lda,
|
|
bp, ldb,
|
|
betap,
|
|
cp, ldc );
|
|
}
|
|
#endif
|
|
|
|
#endif
|
|
|
|
#ifdef PRINT
|
|
bli_printm( "c after", &c, "%4.1f", "" );
|
|
exit(1);
|
|
#endif
|
|
|
|
|
|
dtime_save = bli_clock_min_diff( dtime_save, dtime );
|
|
}
|
|
|
|
gflops = ( 1.0 * m * k * m ) / ( dtime_save * 1.0e9 );
|
|
|
|
if ( bli_is_complex( dt ) ) gflops *= 4.0;
|
|
|
|
#ifdef BLIS
|
|
printf( "data_gemmt_blis" );
|
|
#else
|
|
printf( "data_gemmt_%s", BLAS );
|
|
#endif
|
|
printf( "( %2lu, 1:3 ) = [ %4lu %4lu %7.2f ];\n",
|
|
( unsigned long )(p - p_begin)/p_inc + 1,
|
|
( unsigned long )m,
|
|
( unsigned long )k, gflops );
|
|
|
|
bli_obj_free( &alpha );
|
|
bli_obj_free( &beta );
|
|
|
|
bli_obj_free( &a );
|
|
bli_obj_free( &b );
|
|
bli_obj_free( &c );
|
|
bli_obj_free( &c_save );
|
|
}
|
|
|
|
//bli_finalize();
|
|
|
|
return 0;
|
|
}
|
|
|