You can not select more than 25 topics Topics must start with a chinese character,a letter or number, can include dashes ('-') and can be up to 35 characters long.

dgemm_thread_safety.cpp 4.2 kB

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