Merge branch 'develop' into msvc
[platform/upstream/openblas.git] / cpp_thread_test / dgemv_thread_safety.cpp
1 #include <iostream>
2 #include <vector>
3 #include <random>
4 #include <future>
5 #include <omp.h>
6 #include "../cblas.h"
7 #include "cpp_thread_safety_common.h"
8
9 void launch_cblas_dgemv(double* A, double* x, double* y, const blasint randomMatSize){
10         const blasint inc = 1;
11         cblas_dgemv(CblasColMajor, CblasNoTrans, randomMatSize, randomMatSize, 1.0, A, randomMatSize, x, inc, 0.1, y, inc);
12 }
13
14 int main(int argc, char* argv[]){
15         blasint randomMatSize = 1024; //dimension of the random square matrices and vectors being used
16         uint32_t numConcurrentThreads = 52; //number of concurrent calls of the functions being tested
17         uint32_t numTestRounds = 16; //number of testing rounds before success exit
18         uint32_t maxHwThreads = omp_get_max_threads();
19         
20         if (maxHwThreads < 52)
21                 numConcurrentThreads = maxHwThreads;
22         
23         if (argc > 4){
24                 std::cout<<"ERROR: too many arguments for thread safety tester"<<std::endl;
25                 abort();
26         }
27         if(argc == 4){
28                 std::vector<std::string> cliArgs;
29                 for (int i = 1; i < argc; i++){
30                         cliArgs.push_back(argv[i]);
31                         std::cout<<argv[i]<<std::endl;
32                 }
33                 randomMatSize = std::stoul(cliArgs.at(0));
34                 numConcurrentThreads = std::stoul(cliArgs.at(1));
35                 numTestRounds = std::stoul(cliArgs.at(2));
36         }
37         
38         std::uniform_real_distribution<double> rngdist{-1.0, 1.0};
39         std::vector<std::vector<double>> matBlock(numConcurrentThreads);
40         std::vector<std::vector<double>> vecBlock(numConcurrentThreads*2);
41         std::vector<std::future<void>> futureBlock(numConcurrentThreads);
42         
43         std::cout<<"*----------------------------*\n";
44         std::cout<<"| DGEMV thread safety tester |\n";
45         std::cout<<"*----------------------------*\n";
46         std::cout<<"Size of random matrices and vectors(N=M): "<<randomMatSize<<'\n';
47         std::cout<<"Number of concurrent calls into OpenBLAS : "<<numConcurrentThreads<<'\n';
48         std::cout<<"Number of testing rounds : "<<numTestRounds<<'\n';
49         std::cout<<"This test will need "<<((static_cast<uint64_t>(randomMatSize*randomMatSize)*numConcurrentThreads*8)+(static_cast<uint64_t>(randomMatSize)*numConcurrentThreads*8*2))/static_cast<double>(1024*1024)<<" MiB of RAM\n"<<std::endl;
50
51         FailIfThreadsAreZero(numConcurrentThreads);
52         
53         std::cout<<"Initializing random number generator..."<<std::flush;
54         std::mt19937_64 PRNG = InitPRNG();
55         std::cout<<"done\n";
56         
57         std::cout<<"Preparing to test CBLAS DGEMV thread safety\n";
58         std::cout<<"Allocating matrices..."<<std::flush;
59         for(uint32_t i=0; i<numConcurrentThreads; i++){
60                 matBlock.at(i).resize(randomMatSize*randomMatSize);
61         }
62         std::cout<<"done\n";
63         std::cout<<"Allocating vectors..."<<std::flush;
64         for(uint32_t i=0; i<(numConcurrentThreads*2); i++){
65                 vecBlock.at(i).resize(randomMatSize);
66         }
67         std::cout<<"done\n";
68         //pauser();
69         
70         std::cout<<"Filling matrices with random numbers..."<<std::flush;
71         FillMatrices(matBlock, PRNG, rngdist, randomMatSize, numConcurrentThreads, 1);
72         //PrintMatrices(matBlock, randomMatSize, numConcurrentThreads);
73         std::cout<<"done\n";
74         std::cout<<"Filling vectors with random numbers..."<<std::flush;
75         FillVectors(vecBlock, PRNG, rngdist, randomMatSize, numConcurrentThreads, 2);
76         std::cout<<"done\n";
77         
78         std::cout<<"Testing CBLAS DGEMV thread safety"<<std::endl;
79         omp_set_num_threads(numConcurrentThreads);
80         for(uint32_t R=0; R<numTestRounds; R++){
81                 std::cout<<"DGEMV round #"<<R<<std::endl;
82                 std::cout<<"Launching "<<numConcurrentThreads<<" threads simultaneously using OpenMP..."<<std::flush;
83                 #pragma omp parallel for default(none) shared(futureBlock, matBlock, vecBlock, randomMatSize, numConcurrentThreads)
84                 for(uint32_t i=0; i<numConcurrentThreads; i++){
85                         futureBlock[i] = std::async(std::launch::async, launch_cblas_dgemv, &matBlock[i][0], &vecBlock[i*2][0], &vecBlock[i*2+1][0], randomMatSize);
86                 }
87                 std::cout<<"done\n";
88                 std::cout<<"Waiting for threads to finish..."<<std::flush;
89                 for(uint32_t i=0; i<numConcurrentThreads; i++){
90                         futureBlock[i].get();
91                 }
92                 std::cout<<"done\n";
93                 std::cout<<"Comparing results from different threads..."<<std::flush;
94                 for(uint32_t i=2; i<(numConcurrentThreads*2); i+=2){ //i is the index of vector x, for a given thread
95                         for(uint32_t j = 0; j < static_cast<uint32_t>(randomMatSize); j++){
96                                 if (std::abs(vecBlock[i+1][j] - vecBlock[1][j]) > 1.0E-13){ //i+1 is the index of vector y, for a given thread
97                                         std::cout<<"ERROR: one of the threads returned a different result! Index : "<<i+1<<std::endl;
98                                         std::cout<<"CBLAS DGEMV thread safety test FAILED!"<<std::endl;
99                                         return -1;
100                                 }
101                         }
102                 }
103                 std::cout<<"OK!\n"<<std::endl;
104         }
105         std::cout<<"CBLAS DGEMV thread safety test PASSED!\n"<<std::endl;
106         return 0;
107 }