1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130
|
// Copyright (c) 2017-2023, University of Tennessee. All rights reserved.
// SPDX-License-Identifier: BSD-3-Clause
// This program is free software: you can redistribute it and/or modify it under
// the terms of the BSD 3-Clause license. See the accompanying LICENSE file.
#include "test.hh"
#include "lapack.hh"
#include "lapack/flops.hh"
#include "print_matrix.hh"
#include "error.hh"
#include "lapacke_wrappers.hh"
#include <vector>
// -----------------------------------------------------------------------------
template< typename scalar_t >
void test_pbcon_work( Params& params, bool run )
{
using real_t = blas::real_type< scalar_t >;
// get & mark input values
lapack::Uplo uplo = params.uplo();
int64_t n = params.dim.n();
int64_t kd = params.kd();
int64_t align = params.align();
real_t eps = std::numeric_limits< real_t >::epsilon();
// mark non-standard output values
params.ref_time();
//params.ref_gflops();
////params.gflops();
if (! run)
return;
// ---------- setup
int64_t ldab = roundup( kd+1, align );
real_t anorm;
real_t rcond_tst;
real_t rcond_ref;
size_t size_AB = (size_t) ldab * n;
std::vector< scalar_t > AB( size_AB );
int64_t idist = 1;
int64_t iseed[4] = { 0, 1, 2, 3 };
lapack::larnv( idist, iseed, AB.size(), &AB[0] );
// diagonally dominant -> positive definite
if (uplo == lapack::Uplo::Upper) {
for (int64_t j = 0; j < n; ++j) {
AB[ kd + j*ldab ] += n;
}
}
else { // lower
for (int64_t j = 0; j < n; ++j) {
AB[ j*ldab ] += n;
}
}
anorm = lapack::lanhb( lapack::Norm::One, uplo, n, kd, &AB[0], ldab );
int64_t info = lapack::pbtrf( uplo, n, kd, &AB[0], ldab );
if (info != 0) {
fprintf( stderr, "lapack::pbtrf returned error %lld\n", llong( info ) );
}
// ---------- run test
testsweeper::flush_cache( params.cache() );
double time = testsweeper::get_wtime();
int64_t info_tst = lapack::pbcon( uplo, n, kd, &AB[0], ldab, anorm, &rcond_tst );
time = testsweeper::get_wtime() - time;
if (info_tst != 0) {
fprintf( stderr, "lapack::pbcon returned error %lld\n", llong( info_tst ) );
}
params.time() = time;
//double gflop = lapack::Gflop< scalar_t >::pbcon( n, kd );
//params.gflops() = gflop / time;
if (params.ref() == 'y' || params.check() == 'y') {
// ---------- run reference
testsweeper::flush_cache( params.cache() );
time = testsweeper::get_wtime();
int64_t info_ref = LAPACKE_pbcon( to_char( uplo ), n, kd, &AB[0], ldab, anorm, &rcond_ref );
time = testsweeper::get_wtime() - time;
if (info_ref != 0) {
fprintf( stderr, "LAPACKE_pbcon returned error %lld\n", llong( info_ref ) );
}
params.ref_time() = time;
//params.ref_gflops() = gflop / time;
// ---------- check error compared to reference
real_t error = 0;
if (info_tst != info_ref) {
error = 1;
}
error += std::abs( rcond_tst - rcond_ref );
params.error() = error;
params.okay() = (error < 3*eps);
}
}
// -----------------------------------------------------------------------------
void test_pbcon( Params& params, bool run )
{
switch (params.datatype()) {
case testsweeper::DataType::Single:
test_pbcon_work< float >( params, run );
break;
case testsweeper::DataType::Double:
test_pbcon_work< double >( params, run );
break;
case testsweeper::DataType::SingleComplex:
test_pbcon_work< std::complex<float> >( params, run );
break;
case testsweeper::DataType::DoubleComplex:
test_pbcon_work< std::complex<double> >( params, run );
break;
default:
throw std::runtime_error( "unknown datatype" );
break;
}
}
|