mirror of
https://github.com/amd/blis.git
synced 2026-04-19 15:18:52 +00:00
Description:
Added bench_trsm.c to read inputs from AOCL DTL logs to benchmark
Added sample input file
Change-Id: I6806e42244bf775cbed457553ca07fb0222ef597
445 lines
14 KiB
C
445 lines
14 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) 2020, Advanced Micro Devices, Inc. All rights reserved.
|
|
|
|
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"
|
|
|
|
/*
|
|
Benchmark application to process aocl logs generated
|
|
by BLIS library for trsm.
|
|
*/
|
|
|
|
#ifndef N_REPEAT
|
|
#define N_REPEAT 30
|
|
#endif
|
|
|
|
#define AOCL_MATRIX_INITIALISATION
|
|
|
|
//#define BLIS_ENABLE_CBLAS
|
|
|
|
/* For BLIS since logs are collected at BLAS interfaces
|
|
* we disable cblas interfaces for this benchmark application
|
|
*/
|
|
|
|
/*
|
|
#ifdef BLIS_ENABLE_CBLAS
|
|
#define CBLAS
|
|
#endif
|
|
*/
|
|
|
|
int main( int argc, char** argv )
|
|
{
|
|
obj_t a, c;
|
|
obj_t c_save;
|
|
obj_t alpha;
|
|
dim_t m, n;
|
|
num_t dt = BLIS_DOUBLE;
|
|
int r, n_repeats;
|
|
char side;
|
|
uplo_t uploa;
|
|
trans_t transa;
|
|
diag_t diaga;
|
|
|
|
f77_char f77_side;
|
|
f77_char f77_uploa;
|
|
f77_char f77_transa;
|
|
f77_char f77_diaga;
|
|
|
|
double dtime;
|
|
double dtime_save;
|
|
double gflops;
|
|
double alphaR;
|
|
double alphaI;
|
|
|
|
FILE* fin = NULL;
|
|
FILE* fout = NULL;
|
|
|
|
n_repeats = N_REPEAT;
|
|
|
|
if(argc < 3)
|
|
{
|
|
printf("Usage: ./test_trsm_XX.x input.csv output.csv\n");
|
|
exit(1);
|
|
}
|
|
fin = fopen(argv[1], "r");
|
|
if(fin == NULL)
|
|
{
|
|
printf("Error opening the file %s\n", argv[1]);
|
|
exit(1);
|
|
}
|
|
|
|
fout = fopen(argv[2], "w");
|
|
if(fout == NULL)
|
|
{
|
|
printf("Error opening the file %s\n", argv[2]);
|
|
exit(1);
|
|
}
|
|
|
|
fprintf(fout,"dt\t side\t uploa\t transa\t diaga\t m\t n\t lda\t ldb\t\
|
|
lphaR\t alphaI\t gflops\n");
|
|
|
|
inc_t cs_a,rs_a,cs_c,rs_c;
|
|
char dt_type, side_c, uploa_c, transa_c, diaga_c;
|
|
char logline[255];
|
|
|
|
while(fscanf(fin, "%s %c %c %c %ld %ld %lu %lu %lu %lu %c %c %lf %lf\n",
|
|
logline,&dt_type, &side_c,&uploa_c, &m, &n, &cs_a, &cs_c,&rs_a,&rs_c,
|
|
&transa_c,&diaga_c,&alphaR,&alphaI) == 14)
|
|
{
|
|
if((dt_type=='S')||(dt_type=='s')) dt = BLIS_FLOAT;
|
|
if((dt_type=='D')||(dt_type=='d')) dt = BLIS_DOUBLE;
|
|
if((dt_type=='C')||(dt_type=='c')) dt = BLIS_SCOMPLEX;
|
|
if((dt_type=='Z')||(dt_type=='z')) dt = BLIS_DCOMPLEX;
|
|
|
|
if( 'l' == side_c|| 'L' == side_c)
|
|
side = BLIS_LEFT;
|
|
else if('r' == side_c || 'R' == side_c)
|
|
side = BLIS_RIGHT;
|
|
else
|
|
{
|
|
printf("Invalid entry for the argument 'side':%c\n",side_c);
|
|
continue;
|
|
}
|
|
|
|
if('l' == uploa_c || 'L' == uploa_c)
|
|
uploa = BLIS_LOWER;
|
|
else if('u' == uploa_c || 'U' == uploa_c)
|
|
uploa = BLIS_UPPER;
|
|
else
|
|
{
|
|
printf("Invalid entry for the argument 'uplo':%c\n",uploa_c);
|
|
continue;
|
|
}
|
|
|
|
if('t' == transa_c || 'T' == transa_c)
|
|
transa = BLIS_TRANSPOSE;
|
|
else if('n' == transa_c || 'N' == transa_c)
|
|
transa = BLIS_NO_TRANSPOSE;
|
|
else
|
|
{
|
|
printf("Invalid entry for the argument 'transa':%c\n",transa_c);
|
|
continue;
|
|
}
|
|
|
|
if('u' == diaga_c || 'U' == diaga_c)
|
|
diaga = BLIS_UNIT_DIAG;
|
|
else if('n' == diaga_c || 'N' == diaga_c)
|
|
diaga = BLIS_NONUNIT_DIAG;
|
|
else
|
|
{
|
|
printf("Invalid entry for the argument 'diaga':%c\n", diaga_c);
|
|
continue;
|
|
}
|
|
|
|
bli_param_map_blis_to_netlib_side( side, &f77_side );
|
|
bli_param_map_blis_to_netlib_uplo( uploa, &f77_uploa );
|
|
bli_param_map_blis_to_netlib_trans( transa, &f77_transa );
|
|
bli_param_map_blis_to_netlib_diag( diaga, &f77_diaga );
|
|
|
|
if ( bli_is_left( side ) )
|
|
bli_obj_create( dt, m, m, rs_a, cs_a, &a );
|
|
else
|
|
bli_obj_create( dt, n, n, rs_a, cs_a, &a );
|
|
|
|
bli_obj_create( dt, m, n, rs_c, cs_c, &c );
|
|
bli_obj_create( dt, m, n, rs_c, cs_c, &c_save );
|
|
|
|
#ifdef AOCL_MATRIX_INITIALISATION
|
|
bli_randm( &a );
|
|
bli_randm( &c );
|
|
#endif
|
|
|
|
bli_obj_set_struc( BLIS_TRIANGULAR, &a );
|
|
bli_obj_set_uplo( uploa, &a );
|
|
bli_obj_set_conjtrans( transa, &a );
|
|
bli_obj_set_diag( diaga, &a );
|
|
|
|
// Randomize A and zero the unstored triangle to ensure the
|
|
// implementation reads only from the stored region.
|
|
bli_randm( &a );
|
|
bli_mktrim( &a );
|
|
|
|
// Load the diagonal of A to make it more likely to be invertible.
|
|
bli_shiftd( &BLIS_TWO, &a );
|
|
|
|
bli_obj_create( dt, 1, 1, 0, 0, &alpha );
|
|
bli_setsc( alphaR, alphaI, &alpha );
|
|
|
|
bli_copym( &c, &c_save );
|
|
|
|
dtime_save = DBL_MAX;
|
|
|
|
for ( r = 0; r < n_repeats; ++r )
|
|
{
|
|
bli_copym( &c_save, &c );
|
|
#ifdef PRINT
|
|
bli_printm( "a", &a, "%4.1f", "" );
|
|
bli_printm( "c", &c, "%4.1f", "" );
|
|
#endif
|
|
dtime = bli_clock();
|
|
|
|
#ifdef BLIS
|
|
|
|
bli_trsm( &side,
|
|
&alpha,
|
|
&a,
|
|
&c );
|
|
|
|
#else
|
|
|
|
#ifdef CBLAS
|
|
enum CBLAS_ORDER cblas_order;
|
|
enum CBLAS_TRANSPOSE cblas_transa;
|
|
enum CBLAS_UPLO cblas_uplo;
|
|
enum CBLAS_SIDE cblas_side;
|
|
enum CBLAS_DIAG cblas_diag;
|
|
|
|
if ( bli_obj_row_stride( &c ) == 1 )
|
|
cblas_order = CblasColMajor;
|
|
else
|
|
cblas_order = CblasRowMajor;
|
|
|
|
if( bli_is_trans( transa ) )
|
|
cblas_transa = CblasTrans;
|
|
else if( bli_is_conjtrans( transa ) )
|
|
cblas_transa = CblasConjTrans;
|
|
else
|
|
cblas_transa = CblasNoTrans;
|
|
|
|
if(bli_is_upper(uploa))
|
|
cblas_uplo = CblasUpper;
|
|
else
|
|
cblas_uplo = CblasLower;
|
|
|
|
if(bli_is_left(side))
|
|
cblas_side = CblasLeft;
|
|
else
|
|
cblas_side = CblasRight;
|
|
|
|
if(bli_is_unit_diag(diaga))
|
|
cblas_diag = CblasUnit;
|
|
else
|
|
cblas_diag = CblasNonUnit;
|
|
|
|
#else
|
|
f77_char f77_transa;
|
|
bli_param_map_blis_to_netlib_trans( transa, &f77_transa );
|
|
#endif
|
|
if ( bli_is_float( dt ) )
|
|
{
|
|
f77_int mm = bli_obj_length( &c );
|
|
f77_int nn = bli_obj_width( &c );
|
|
f77_int lda = bli_obj_col_stride( &a );
|
|
f77_int ldc = bli_obj_col_stride( &c );
|
|
|
|
float* alphap = bli_obj_buffer( &alpha );
|
|
float* ap = bli_obj_buffer( &a );
|
|
float* cp = bli_obj_buffer( &c );
|
|
|
|
#ifdef CBLAS
|
|
cblas_strsm( cblas_order,
|
|
cblas_side,
|
|
cblas_uplo,
|
|
cblas_transa,
|
|
cblas_diag,
|
|
mm,
|
|
nn,
|
|
*alphap,
|
|
ap, lda,
|
|
cp, ldc
|
|
);
|
|
#else
|
|
strsm_( &f77_side,
|
|
&f77_uploa,
|
|
&f77_transa,
|
|
&f77_diaga,
|
|
&mm,
|
|
&nn,
|
|
alphap,
|
|
ap, &lda,
|
|
cp, &ldc );
|
|
#endif
|
|
}
|
|
else if ( bli_is_double( dt ) )
|
|
{
|
|
f77_int mm = bli_obj_length( &c );
|
|
f77_int nn = bli_obj_width( &c );
|
|
f77_int lda = bli_obj_col_stride( &a );
|
|
f77_int ldc = bli_obj_col_stride( &c );
|
|
double* alphap = bli_obj_buffer( &alpha );
|
|
double* ap = bli_obj_buffer( &a );
|
|
double* cp = bli_obj_buffer( &c );
|
|
|
|
#ifdef CBLAS
|
|
cblas_dtrsm( cblas_order,
|
|
cblas_side,
|
|
cblas_uplo,
|
|
cblas_transa,
|
|
cblas_diag,
|
|
mm,
|
|
nn,
|
|
*alphap,
|
|
ap, lda,
|
|
cp, ldc
|
|
);
|
|
#else
|
|
dtrsm_( &f77_side,
|
|
&f77_uploa,
|
|
&f77_transa,
|
|
&f77_diaga,
|
|
&mm,
|
|
&nn,
|
|
alphap,
|
|
ap, &lda,
|
|
cp, &ldc );
|
|
#endif
|
|
|
|
}
|
|
else if ( bli_is_scomplex( dt ) )
|
|
{
|
|
f77_int mm = bli_obj_length( &c );
|
|
f77_int nn = bli_obj_width( &c );
|
|
f77_int lda = bli_obj_col_stride( &a );
|
|
f77_int ldc = bli_obj_col_stride( &c );
|
|
scomplex* alphap = bli_obj_buffer( &alpha );
|
|
scomplex* ap = bli_obj_buffer( &a );
|
|
scomplex* cp = bli_obj_buffer( &c );
|
|
|
|
#ifdef CBLAS
|
|
cblas_ctrsm( cblas_order,
|
|
cblas_side,
|
|
cblas_uplo,
|
|
cblas_transa,
|
|
cblas_diag,
|
|
mm,
|
|
nn,
|
|
alphap,
|
|
ap, lda,
|
|
cp, ldc
|
|
);
|
|
#else
|
|
ctrsm_( &f77_side,
|
|
&f77_uploa,
|
|
&f77_transa,
|
|
&f77_diaga,
|
|
&mm,
|
|
&nn,
|
|
alphap,
|
|
ap, &lda,
|
|
cp, &ldc );
|
|
#endif
|
|
}
|
|
else if ( bli_is_dcomplex( dt ) )
|
|
{
|
|
f77_int mm = bli_obj_length( &c );
|
|
f77_int nn = bli_obj_width( &c );
|
|
f77_int lda = bli_obj_col_stride( &a );
|
|
f77_int ldc = bli_obj_col_stride( &c );
|
|
dcomplex* alphap = bli_obj_buffer( &alpha );
|
|
dcomplex* ap = bli_obj_buffer( &a );
|
|
dcomplex* cp = bli_obj_buffer( &c );
|
|
#ifdef CBLAS
|
|
cblas_ztrsm( cblas_order,
|
|
cblas_side,
|
|
cblas_uplo,
|
|
cblas_transa,
|
|
cblas_diag,
|
|
mm,
|
|
nn,
|
|
alphap,
|
|
ap, lda,
|
|
cp, ldc
|
|
);
|
|
#else
|
|
ztrsm_( &f77_side,
|
|
&f77_uploa,
|
|
&f77_transa,
|
|
&f77_diaga,
|
|
&mm,
|
|
&nn,
|
|
alphap,
|
|
ap, &lda,
|
|
cp, &ldc );
|
|
#endif
|
|
}else{
|
|
printf("Invalid data type! Exiting!\n");
|
|
exit(1);
|
|
}
|
|
#endif
|
|
dtime_save = bli_clock_min_diff( dtime_save, dtime );
|
|
}
|
|
|
|
if ( bli_is_left( side ) )
|
|
gflops = ( 1.0 * m * m * n ) / ( dtime_save * 1.0e9 );
|
|
else
|
|
gflops = ( 1.0 * m * n * n ) / ( dtime_save * 1.0e9 );
|
|
|
|
if ( bli_is_complex( dt ) ) gflops *= 4.0;
|
|
|
|
#ifdef BLIS
|
|
printf( "data_trsm_blis\t\t");
|
|
#else
|
|
printf( "data_trsm_%s\t\t",BLAS );
|
|
#endif
|
|
|
|
printf("%c\t %c\t %c\t %c\t %c\t %4lu\t %4lu\t %4lu\t %4lu\t %6.3f\t %6.3f\t %6.3f\n",
|
|
dt_type,side_c, uploa_c, transa_c,
|
|
diaga_c, (unsigned long )m, (unsigned long ) n, (unsigned long )cs_a,
|
|
(unsigned long )cs_c, alphaR, alphaI, gflops);
|
|
|
|
fprintf(fout,"%c\t %c\t %c\t %c\t %c\t %4lu\t %4lu\t %4lu\t %4lu\t %6.3f\t %6.3f\t %6.3f\n",
|
|
dt_type,side_c, uploa_c, transa_c,
|
|
diaga_c, (unsigned long )m, (unsigned long ) n, (unsigned long )cs_a,
|
|
(unsigned long )cs_c, alphaR, alphaI, gflops);
|
|
|
|
fflush(fout);
|
|
|
|
bli_obj_free( &alpha );
|
|
bli_obj_free( &a );
|
|
bli_obj_free( &c );
|
|
bli_obj_free( &c_save );
|
|
}
|
|
|
|
fclose(fin);
|
|
fclose(fout);
|
|
|
|
//bli_finalize();
|
|
|
|
return 0;
|
|
} |