Files
blis/testsuite/src/test_trsm.c
Field G. Van Zee 29b0e1ef4e Code review + tweaks to AMD's AOCL 2.0 PR (#349).
Details:
- NOTE: This is a merge commit of 'master' of git://github.com/amd/blis
  into 'amd-master' of flame/blis.
- Fixed a bug in the downstream value of BLIS_NUM_ARCHS, which was
  inadvertantly not incremented when the Zen2 subconfiguration was
  added.
- In bli_gemm_front(), added a missing conditional constraint around the
  call to bli_gemm_small() that ensures that the computation precision
  of C matches the storage precision of C.
- In bli_syrk_front(), reorganized and relocated the notrans/trans logic
  that existed around the call to bli_syrk_small() into bli_syrk_small()
  to minimize the calling code footprint and also to bring that code
  into stylistic harmony with similar code in bli_gemm_front() and
  bli_trsm_front(). Also, replaced direct accessing of obj_t fields with
  proper accessor static functions (e.g. 'a->dim[0]' becomes
  'bli_obj_length( a )').
- Added #ifdef BLIS_ENABLE_SMALL_MATRIX guard around prototypes for
  bli_gemm_small(), bli_syrk_small(), and bli_trsm_small(). This is
  strictly speaking unnecessary, but it serves as a useful visual cue to
  those who may be reading the files.
- Removed cpp macro-protected small matrix debugging code from
  bli_trsm_front.c.
- Added a GCC_OT_9_1_0 variable to build/config.mk.in to facilitate gcc
  version check for availability of -march=znver2, and added appropriate
  support to configure script.
- Cleanups to compiler flags common to recent AMD microarchitectures in
  config/zen/amd_config.mk, including: removal of -march=znver1 et al.
  from CKVECFLAGS (since the -march flag is added within make_defs.mk);
  setting CRVECFLAGS similarly to CKVECFLAGS.
- Cleanups to config/zen/bli_cntx_init_zen.c.
- Cleanups, added comments to config/zen/make_defs.mk.
- Cleanups to config/zen2/make_defs.mk, including making use of newly-
  added GCC_OT_9_1_0 and existing GCC_OT_6_1_0 to choose the correct
  set of compiler flags based on the version of gcc being used.
- Reverted downstream changes to test/test_gemm.c.
- Various whitespace/comment changes.
2019-10-11 10:24:24 -05:00

388 lines
10 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) 2018 - 2019, 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(s) of the copyright holder(s) 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.
*/
#include "blis.h"
#include "test_libblis.h"
// Static variables.
static char* op_str = "trsm";
static char* o_types = "mm"; // a b
static char* p_types = "suhd"; // side uploa transa diaga
static thresh_t thresh[BLIS_NUM_FP_TYPES] = { { 1e-04, 1e-05 }, // warn, pass for s
{ 1e-04, 1e-05 }, // warn, pass for c
{ 1e-13, 1e-14 }, // warn, pass for d
{ 1e-13, 1e-14 } }; // warn, pass for z
// Local prototypes.
void libblis_test_trsm_deps
(
thread_data_t* tdata,
test_params_t* params,
test_op_t* op
);
void libblis_test_trsm_experiment
(
test_params_t* params,
test_op_t* op,
iface_t iface,
char* dc_str,
char* pc_str,
char* sc_str,
unsigned int p_cur,
double* perf,
double* resid
);
void libblis_test_trsm_impl
(
iface_t iface,
side_t side,
obj_t* alpha,
obj_t* a,
obj_t* b
);
void libblis_test_trsm_check
(
test_params_t* params,
side_t side,
obj_t* alpha,
obj_t* a,
obj_t* b,
obj_t* b_orig,
double* resid
);
void libblis_test_trsm_deps
(
thread_data_t* tdata,
test_params_t* params,
test_op_t* op
)
{
libblis_test_randv( tdata, params, &(op->ops->randv) );
libblis_test_randm( tdata, params, &(op->ops->randm) );
libblis_test_setv( tdata, params, &(op->ops->setv) );
libblis_test_normfv( tdata, params, &(op->ops->normfv) );
libblis_test_subv( tdata, params, &(op->ops->subv) );
libblis_test_scalv( tdata, params, &(op->ops->scalv) );
libblis_test_copym( tdata, params, &(op->ops->copym) );
libblis_test_scalm( tdata, params, &(op->ops->scalm) );
libblis_test_gemv( tdata, params, &(op->ops->gemv) );
libblis_test_trsv( tdata, params, &(op->ops->trsv) );
}
void libblis_test_trsm
(
thread_data_t* tdata,
test_params_t* params,
test_op_t* op
)
{
// Return early if this test has already been done.
if ( libblis_test_op_is_done( op ) ) return;
// Return early if operation is disabled.
if ( libblis_test_op_is_disabled( op ) ||
libblis_test_l3_is_disabled( op ) ) return;
// Call dependencies first.
if ( TRUE ) libblis_test_trsm_deps( tdata, params, op );
// Execute the test driver for each implementation requested.
//if ( op->front_seq == ENABLE )
{
libblis_test_op_driver( tdata,
params,
op,
BLIS_TEST_SEQ_FRONT_END,
op_str,
p_types,
o_types,
thresh,
libblis_test_trsm_experiment );
}
}
void libblis_test_trsm_experiment
(
test_params_t* params,
test_op_t* op,
iface_t iface,
char* dc_str,
char* pc_str,
char* sc_str,
unsigned int p_cur,
double* perf,
double* resid
)
{
unsigned int n_repeats = params->n_repeats;
unsigned int i;
double time_min = DBL_MAX;
double time;
num_t datatype;
dim_t m, n;
dim_t mn_side;
side_t side;
uplo_t uploa;
trans_t transa;
diag_t diaga;
obj_t alpha, a, b;
obj_t b_save;
// Use the datatype of the first char in the datatype combination string.
bli_param_map_char_to_blis_dt( dc_str[0], &datatype );
// Map the dimension specifier to actual dimensions.
m = libblis_test_get_dim_from_prob_size( op->dim_spec[0], p_cur );
n = libblis_test_get_dim_from_prob_size( op->dim_spec[1], p_cur );
// Map parameter characters to BLIS constants.
bli_param_map_char_to_blis_side( pc_str[0], &side );
bli_param_map_char_to_blis_uplo( pc_str[1], &uploa );
bli_param_map_char_to_blis_trans( pc_str[2], &transa );
bli_param_map_char_to_blis_diag( pc_str[3], &diaga );
// Create test scalars.
bli_obj_scalar_init_detached( datatype, &alpha );
// Create test operands (vectors and/or matrices).
bli_set_dim_with_side( side, m, n, &mn_side );
libblis_test_mobj_create( params, datatype, transa,
sc_str[1], mn_side, mn_side, &a );
libblis_test_mobj_create( params, datatype, BLIS_NO_TRANSPOSE,
sc_str[0], m, n, &b );
libblis_test_mobj_create( params, datatype, BLIS_NO_TRANSPOSE,
sc_str[0], m, n, &b_save );
// Set alpha.
if ( bli_obj_is_real( &b ) )
{
bli_setsc( 2.0, 0.0, &alpha );
}
else
{
bli_setsc( 2.0, 0.0, &alpha );
}
// Set the structure and uplo properties of A.
bli_obj_set_struc( BLIS_TRIANGULAR, &a );
bli_obj_set_uplo( uploa, &a );
// Randomize A, load the diagonal, make it densely triangular.
libblis_test_mobj_randomize( params, TRUE, &a );
libblis_test_mobj_load_diag( params, &a );
bli_mktrim( &a );
// Randomize B and save B.
libblis_test_mobj_randomize( params, TRUE, &b );
bli_copym( &b, &b_save );
// Apply the remaining parameters.
bli_obj_set_conjtrans( transa, &a );
bli_obj_set_diag( diaga, &a );
// Repeat the experiment n_repeats times and record results.
for ( i = 0; i < n_repeats; ++i )
{
bli_copym( &b_save, &b );
time = bli_clock();
libblis_test_trsm_impl( iface, side, &alpha, &a, &b );
time_min = bli_clock_min_diff( time_min, time );
}
// Estimate the performance of the best experiment repeat.
*perf = ( 1.0 * mn_side * m * n ) / time_min / FLOPS_PER_UNIT_PERF;
if ( bli_obj_is_complex( &b ) ) *perf *= 4.0;
// Perform checks.
libblis_test_trsm_check( params, side, &alpha, &a, &b, &b_save, resid );
// Zero out performance and residual if output matrix is empty.
libblis_test_check_empty_problem( &b, perf, resid );
// Free the test objects.
bli_obj_free( &a );
bli_obj_free( &b );
bli_obj_free( &b_save );
}
void libblis_test_trsm_impl
(
iface_t iface,
side_t side,
obj_t* alpha,
obj_t* a,
obj_t* b
)
{
switch ( iface )
{
case BLIS_TEST_SEQ_FRONT_END:
#if 0
bli_printm( "a", a, "%5.2f", "" );
bli_printm( "b", b, "%5.2f", "" );
//bli_printm( "alpha", alpha, "%5.2f", "" );
#endif
bli_trsm( side, alpha, a, b );
#if 0
bli_printm( "b after", b, "%5.2f", "" );
#endif
break;
default:
libblis_test_printf_error( "Invalid interface type.\n" );
}
}
void libblis_test_trsm_check
(
test_params_t* params,
side_t side,
obj_t* alpha,
obj_t* a,
obj_t* b,
obj_t* b_orig,
double* resid
)
{
num_t dt = bli_obj_dt( b );
num_t dt_real = bli_obj_dt_proj_to_real( b );
dim_t m = bli_obj_length( b );
dim_t n = bli_obj_width( b );
obj_t norm;
obj_t t, v, w, z;
double junk;
//
// Pre-conditions:
// - a is randomized and triangular.
// - b_orig is randomized.
// Note:
// - alpha should have a non-zero imaginary component in the
// complex cases in order to more fully exercise the implementation.
//
// Under these conditions, we assume that the implementation for
//
// B := alpha * inv(transa(A)) * B_orig (side = left)
// B := alpha * B_orig * inv(transa(A)) (side = right)
//
// is functioning correctly if
//
// normfv( v - z )
//
// is negligible, where
//
// v = B * t
//
// z = ( alpha * inv(transa(A)) * B ) * t (side = left)
// = alpha * inv(transa(A)) * B * t
// = alpha * inv(transa(A)) * w
//
// z = ( alpha * B * inv(transa(A)) ) * t (side = right)
// = alpha * B * tinv(ransa(A)) * t
// = alpha * B * w
bli_obj_scalar_init_detached( dt_real, &norm );
if ( bli_is_left( side ) )
{
bli_obj_create( dt, n, 1, 0, 0, &t );
bli_obj_create( dt, m, 1, 0, 0, &v );
bli_obj_create( dt, m, 1, 0, 0, &w );
bli_obj_create( dt, m, 1, 0, 0, &z );
}
else // else if ( bli_is_left( side ) )
{
bli_obj_create( dt, n, 1, 0, 0, &t );
bli_obj_create( dt, m, 1, 0, 0, &v );
bli_obj_create( dt, n, 1, 0, 0, &w );
bli_obj_create( dt, m, 1, 0, 0, &z );
}
libblis_test_vobj_randomize( params, TRUE, &t );
bli_gemv( &BLIS_ONE, b, &t, &BLIS_ZERO, &v );
if ( bli_is_left( side ) )
{
bli_gemv( alpha, b_orig, &t, &BLIS_ZERO, &w );
bli_trsv( &BLIS_ONE, a, &w );
bli_copyv( &w, &z );
}
else
{
bli_copyv( &t, &w );
bli_trsv( &BLIS_ONE, a, &w );
bli_gemv( alpha, b_orig, &w, &BLIS_ZERO, &z );
}
bli_subv( &z, &v );
bli_normfv( &v, &norm );
bli_getsc( &norm, resid, &junk );
bli_obj_free( &t );
bli_obj_free( &v );
bli_obj_free( &w );
bli_obj_free( &z );
}