summaryrefslogtreecommitdiff
diff options
context:
space:
mode:
-rw-r--r--include/internal/routine.h4
-rw-r--r--include/internal/utilities.h4
-rw-r--r--test/correctness/testabc.cc2
-rw-r--r--test/correctness/testabc.h1
-rw-r--r--test/correctness/testaxy.cc2
-rw-r--r--test/correctness/testaxy.h1
-rw-r--r--test/correctness/tester.cc25
-rw-r--r--test/correctness/tester.h3
-rw-r--r--test/correctness/testxy.cc2
-rw-r--r--test/correctness/testxy.h3
10 files changed, 41 insertions, 6 deletions
diff --git a/include/internal/routine.h b/include/internal/routine.h
index 1b2e0dbb..a65ced20 100644
--- a/include/internal/routine.h
+++ b/include/internal/routine.h
@@ -29,10 +29,6 @@ namespace clblast {
class Routine {
public:
- // Khronos OpenCL extensions
- const std::string kKhronosHalfPrecision = "cl_khr_fp16";
- const std::string kKhronosDoublePrecision = "cl_khr_fp64";
-
// The cache of compiled OpenCL programs, along with some meta-data
struct ProgramCache {
Program program;
diff --git a/include/internal/utilities.h b/include/internal/utilities.h
index 4cefdc4a..6ad17a6a 100644
--- a/include/internal/utilities.h
+++ b/include/internal/utilities.h
@@ -31,6 +31,10 @@ namespace clblast {
using float2 = std::complex<float>;
using double2 = std::complex<double>;
+// Khronos OpenCL extensions
+const std::string kKhronosHalfPrecision = "cl_khr_fp16";
+const std::string kKhronosDoublePrecision = "cl_khr_fp64";
+
// =================================================================================================
// The routine-specific arguments in string form
diff --git a/test/correctness/testabc.cc b/test/correctness/testabc.cc
index 107ffb70..eed17560 100644
--- a/test/correctness/testabc.cc
+++ b/test/correctness/testabc.cc
@@ -46,6 +46,7 @@ TestABC<T>::TestABC(int argc, char *argv[], const bool silent,
// Tests the routine for a wide variety of parameters
template <typename T>
void TestABC<T>::TestRegular(Arguments<T> &args, const std::string &name) {
+ if (!PrecisionSupported()) { return; }
TestStart("regular behaviour", name);
// Computes whether or not the matrices are transposed. Note that we assume a default of
@@ -161,6 +162,7 @@ void TestABC<T>::TestRegular(Arguments<T> &args, const std::string &name) {
// does not test for results (if any).
template <typename T>
void TestABC<T>::TestInvalidBufferSizes(Arguments<T> &args, const std::string &name) {
+ if (!PrecisionSupported()) { return; }
TestStart("invalid buffer sizes", name);
// Sets example test parameters
diff --git a/test/correctness/testabc.h b/test/correctness/testabc.h
index b7601024..2c44d532 100644
--- a/test/correctness/testabc.h
+++ b/test/correctness/testabc.h
@@ -42,6 +42,7 @@ class TestABC: public Tester<T> {
using Tester<T>::TestErrorCodes;
using Tester<T>::GetExampleScalars;
using Tester<T>::GetOffsets;
+ using Tester<T>::PrecisionSupported;
// Test settings for the regular test. Append to this list in case more tests are required.
const std::vector<size_t> kMatrixDims = { 7, 64 };
diff --git a/test/correctness/testaxy.cc b/test/correctness/testaxy.cc
index 3588573a..cb5e9923 100644
--- a/test/correctness/testaxy.cc
+++ b/test/correctness/testaxy.cc
@@ -47,6 +47,7 @@ TestAXY<T>::TestAXY(int argc, char *argv[], const bool silent,
// Tests the routine for a wide variety of parameters
template <typename T>
void TestAXY<T>::TestRegular(Arguments<T> &args, const std::string &name) {
+ if (!PrecisionSupported()) { return; }
TestStart("regular behaviour", name);
// Iterates over the dimension for the matrix and vectors
@@ -151,6 +152,7 @@ void TestAXY<T>::TestRegular(Arguments<T> &args, const std::string &name) {
// does not test for results (if any).
template <typename T>
void TestAXY<T>::TestInvalidBufferSizes(Arguments<T> &args, const std::string &name) {
+ if (!PrecisionSupported()) { return; }
TestStart("invalid buffer sizes", name);
// Sets example test parameters
diff --git a/test/correctness/testaxy.h b/test/correctness/testaxy.h
index 00d4f6ec..fa2c4a98 100644
--- a/test/correctness/testaxy.h
+++ b/test/correctness/testaxy.h
@@ -42,6 +42,7 @@ class TestAXY: public Tester<T> {
using Tester<T>::TestErrorCodes;
using Tester<T>::GetExampleScalars;
using Tester<T>::GetOffsets;
+ using Tester<T>::PrecisionSupported;
// Test settings for the regular test. Append to this list in case more tests are required.
const std::vector<size_t> kMatrixVectorDims = { 61, 512 };
diff --git a/test/correctness/tester.cc b/test/correctness/tester.cc
index 09b0b5bd..74b6679d 100644
--- a/test/correctness/tester.cc
+++ b/test/correctness/tester.cc
@@ -56,8 +56,18 @@ Tester<T>::Tester(int argc, char *argv[], const bool silent,
// Prints the header
fprintf(stdout, "* Running on OpenCL device '%s'.\n", device_.Name().c_str());
- fprintf(stdout, "* Starting tests for the %s'%s'%s routine. Legend:\n",
+ fprintf(stdout, "* Starting tests for the %s'%s'%s routine.",
kPrintMessage.c_str(), name.c_str(), kPrintEnd.c_str());
+
+ // Checks whether the precision is supported
+ if (!PrecisionSupported()) {
+ fprintf(stdout, "\n* Tests skipped: %sUnsupported precision%s\n",
+ kPrintWarning.c_str(), kPrintEnd.c_str());
+ return;
+ }
+
+ // Prints the legend
+ fprintf(stdout, " Legend:\n");
fprintf(stdout, " %s -> Test produced correct results\n", kSuccessData.c_str());
fprintf(stdout, " %s -> Test returned the correct error code\n", kSuccessStatus.c_str());
fprintf(stdout, " %s -> Test produced incorrect results\n", kErrorData.c_str());
@@ -290,6 +300,19 @@ const std::vector<size_t> Tester<T>::GetOffsets() {
// =================================================================================================
+template <> bool Tester<float>::PrecisionSupported() const { return true; }
+template <> bool Tester<float2>::PrecisionSupported() const { return true; }
+template <> bool Tester<double>::PrecisionSupported() const {
+ auto extensions = device_.Extensions();
+ return (extensions.find(kKhronosDoublePrecision) == std::string::npos) ? false : true;
+}
+template <> bool Tester<double2>::PrecisionSupported() const {
+ auto extensions = device_.Extensions();
+ return (extensions.find(kKhronosDoublePrecision) == std::string::npos) ? false : true;
+}
+
+// =================================================================================================
+
// A test can either pass, be skipped, or fail
template <typename T>
void Tester<T>::ReportPass() {
diff --git a/test/correctness/tester.h b/test/correctness/tester.h
index 1085c472..b627d8f1 100644
--- a/test/correctness/tester.h
+++ b/test/correctness/tester.h
@@ -100,6 +100,9 @@ class Tester {
// Retrieves a list of offset values to test
const std::vector<size_t> GetOffsets();
+ // Returns false is this precision is not supported by the device
+ bool PrecisionSupported() const;
+
// The help-message
std::string help_;
diff --git a/test/correctness/testxy.cc b/test/correctness/testxy.cc
index b2d6c147..5c86182c 100644
--- a/test/correctness/testxy.cc
+++ b/test/correctness/testxy.cc
@@ -44,6 +44,7 @@ TestXY<T>::TestXY(int argc, char *argv[], const bool silent,
// Tests the routine for a wide variety of parameters
template <typename T>
void TestXY<T>::TestRegular(Arguments<T> &args, const std::string &name) {
+ if (!PrecisionSupported()) { return; }
TestStart("regular behaviour", name);
// Iterates over the vector dimension
@@ -122,6 +123,7 @@ void TestXY<T>::TestRegular(Arguments<T> &args, const std::string &name) {
// does not test for results (if any).
template <typename T>
void TestXY<T>::TestInvalidBufferSizes(Arguments<T> &args, const std::string &name) {
+ if (!PrecisionSupported()) { return; }
TestStart("invalid buffer sizes", name);
// Sets example test parameters
diff --git a/test/correctness/testxy.h b/test/correctness/testxy.h
index 8f2b6876..ec2cbcc7 100644
--- a/test/correctness/testxy.h
+++ b/test/correctness/testxy.h
@@ -40,11 +40,12 @@ class TestXY: public Tester<T> {
using Tester<T>::TestErrorCodes;
using Tester<T>::GetExampleScalars;
using Tester<T>::GetOffsets;
+ using Tester<T>::PrecisionSupported;
// Test settings for the regular test. Append to this list in case more tests are required.
const std::vector<size_t> kVectorDims = { 7, 93, 4096 };
const std::vector<size_t> kOffsets = GetOffsets();
- const std::vector<size_t> kIncrements = { 1, 2 };
+ const std::vector<size_t> kIncrements = { 1, 2, 7 };
const std::vector<T> kAlphaValues = GetExampleScalars();
// Test settings for the invalid test