matrix.txx 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570
  1. /**
  2. * \file matrix.txx
  3. * \author Dhairya Malhotra, dhairya.malhotra@gmail.com
  4. * \date 2-11-2011
  5. * \brief This file contains inplementation of the class Matrix.
  6. */
  7. #include <cstring>
  8. #include <cassert>
  9. #include <iomanip>
  10. #include <typeinfo>
  11. #include <profile.hpp>
  12. #include <mat_utils.hpp>
  13. namespace pvfmm{
  14. template <class T>
  15. std::ostream& operator<<(std::ostream& output, const Matrix<T>& M){
  16. output<<std::fixed<<std::setprecision(4)<<std::setiosflags(std::ios::left);
  17. for(size_t i=0;i<M.Dim(0);i++){
  18. for(size_t j=0;j<M.Dim(1);j++)
  19. output<<std::setw(10)<<M(i,j)<<' ';
  20. output<<";\n";
  21. }
  22. return output;
  23. }
  24. template <class T>
  25. Matrix<T>::Matrix(){
  26. dim[0]=0;
  27. dim[1]=0;
  28. own_data=true;
  29. data_ptr=NULL;
  30. dev.dev_ptr=(uintptr_t)NULL;
  31. }
  32. template <class T>
  33. Matrix<T>::Matrix(size_t dim1, size_t dim2, T* data_, bool own_data_){
  34. dim[0]=dim1;
  35. dim[1]=dim2;
  36. own_data=own_data_;
  37. if(own_data){
  38. if(dim[0]*dim[1]>0){
  39. data_ptr=mem::aligned_malloc<T>(dim[0]*dim[1]);
  40. #ifndef __MIC__
  41. Profile::Add_MEM(dim[0]*dim[1]*sizeof(T));
  42. #endif
  43. if(data_!=NULL) mem::memcopy(data_ptr,data_,dim[0]*dim[1]*sizeof(T));
  44. }else data_ptr=NULL;
  45. }else
  46. data_ptr=data_;
  47. dev.dev_ptr=(uintptr_t)NULL;
  48. }
  49. template <class T>
  50. Matrix<T>::Matrix(const Matrix<T>& M){
  51. dim[0]=M.dim[0];
  52. dim[1]=M.dim[1];
  53. own_data=true;
  54. if(dim[0]*dim[1]>0){
  55. data_ptr=mem::aligned_malloc<T>(dim[0]*dim[1]);
  56. #ifndef __MIC__
  57. Profile::Add_MEM(dim[0]*dim[1]*sizeof(T));
  58. #endif
  59. mem::memcopy(data_ptr,M.data_ptr,dim[0]*dim[1]*sizeof(T));
  60. }else
  61. data_ptr=NULL;
  62. dev.dev_ptr=(uintptr_t)NULL;
  63. }
  64. template <class T>
  65. Matrix<T>::~Matrix(){
  66. FreeDevice(false);
  67. if(own_data){
  68. if(data_ptr!=NULL){
  69. mem::aligned_free(data_ptr);
  70. #ifndef __MIC__
  71. Profile::Add_MEM(-dim[0]*dim[1]*sizeof(T));
  72. #endif
  73. }
  74. }
  75. data_ptr=NULL;
  76. dim[0]=0;
  77. dim[1]=0;
  78. }
  79. template <class T>
  80. void Matrix<T>::Swap(Matrix<T>& M){
  81. size_t dim_[2]={dim[0],dim[1]};
  82. T* data_ptr_=data_ptr;
  83. bool own_data_=own_data;
  84. Device dev_=dev;
  85. Vector<char> dev_sig_=dev_sig;
  86. dim[0]=M.dim[0];
  87. dim[1]=M.dim[1];
  88. data_ptr=M.data_ptr;
  89. own_data=M.own_data;
  90. dev=M.dev;
  91. dev_sig=M.dev_sig;
  92. M.dim[0]=dim_[0];
  93. M.dim[1]=dim_[1];
  94. M.data_ptr=data_ptr_;
  95. M.own_data=own_data_;
  96. M.dev=dev_;
  97. M.dev_sig=dev_sig_;
  98. }
  99. template <class T>
  100. void Matrix<T>::ReInit(size_t dim1, size_t dim2, T* data_, bool own_data_){
  101. Matrix<T> tmp(dim1,dim2,data_,own_data_);
  102. this->Swap(tmp);
  103. }
  104. template <class T>
  105. typename Matrix<T>::Device& Matrix<T>::AllocDevice(bool copy){
  106. size_t len=dim[0]*dim[1];
  107. if(dev.dev_ptr==(uintptr_t)NULL && len>0) // Allocate data on device.
  108. dev.dev_ptr=DeviceWrapper::alloc_device((char*)data_ptr, len*sizeof(T));
  109. if(dev.dev_ptr!=(uintptr_t)NULL && copy) // Copy data to device
  110. dev.lock_idx=DeviceWrapper::host2device((char*)data_ptr,(char*)data_ptr,dev.dev_ptr,len*sizeof(T));
  111. dev.dim[0]=dim[0];
  112. dev.dim[1]=dim[1];
  113. return dev;
  114. }
  115. template <class T>
  116. void Matrix<T>::Device2Host(T* host_ptr){
  117. dev.lock_idx=DeviceWrapper::device2host((char*)data_ptr,dev.dev_ptr,(char*)(host_ptr==NULL?data_ptr:host_ptr),dim[0]*dim[1]*sizeof(T));
  118. }
  119. template <class T>
  120. void Matrix<T>::Device2HostWait(){
  121. DeviceWrapper::wait(dev.lock_idx);
  122. dev.lock_idx=-1;
  123. }
  124. template <class T>
  125. void Matrix<T>::FreeDevice(bool copy){
  126. if(dev.dev_ptr==(uintptr_t)NULL) return;
  127. if(copy) DeviceWrapper::device2host((char*)data_ptr,dev.dev_ptr,(char*)data_ptr,dim[0]*dim[1]*sizeof(T));
  128. DeviceWrapper::free_device((char*)data_ptr, dev.dev_ptr);
  129. dev.dev_ptr=(uintptr_t)NULL;
  130. dev.dim[0]=0;
  131. dev.dim[1]=0;
  132. }
  133. template <class T>
  134. void Matrix<T>::Write(const char* fname){
  135. FILE* f1=fopen(fname,"wb+");
  136. if(f1==NULL){
  137. std::cout<<"Unable to open file for writing:"<<fname<<'\n';
  138. return;
  139. }
  140. size_t dim_[2]={dim[0],dim[1]};
  141. fwrite(dim_,sizeof(size_t),2,f1);
  142. fwrite(data_ptr,sizeof(T),dim[0]*dim[1],f1);
  143. fclose(f1);
  144. }
  145. template <class T>
  146. size_t Matrix<T>::Dim(size_t i) const{
  147. return dim[i];
  148. }
  149. template <class T>
  150. void Matrix<T>::Resize(size_t i, size_t j){
  151. if(dim[0]==i && dim[1]==j) return;
  152. FreeDevice(false);
  153. if(own_data){
  154. if(data_ptr!=NULL){
  155. mem::aligned_free(data_ptr);
  156. #ifndef __MIC__
  157. Profile::Add_MEM(-dim[0]*dim[1]*sizeof(T));
  158. #endif
  159. }
  160. }
  161. dim[0]=i;
  162. dim[1]=j;
  163. if(own_data){
  164. if(dim[0]*dim[1]>0){
  165. data_ptr=mem::aligned_malloc<T>(dim[0]*dim[1]);
  166. #ifndef __MIC__
  167. Profile::Add_MEM(dim[0]*dim[1]*sizeof(T));
  168. #endif
  169. }else
  170. data_ptr=NULL;
  171. }
  172. }
  173. template <class T>
  174. void Matrix<T>::SetZero(){
  175. if(dim[0]*dim[1]>0)
  176. memset(data_ptr,0,dim[0]*dim[1]*sizeof(T));
  177. }
  178. template <class T>
  179. Matrix<T>& Matrix<T>::operator=(const Matrix<T>& M){
  180. if(this!=&M){
  181. FreeDevice(false);
  182. if(own_data && dim[0]*dim[1]!=M.dim[0]*M.dim[1]){
  183. if(data_ptr!=NULL){
  184. mem::aligned_free(data_ptr); data_ptr=NULL;
  185. #ifndef __MIC__
  186. Profile::Add_MEM(-dim[0]*dim[1]*sizeof(T));
  187. #endif
  188. }
  189. if(M.dim[0]*M.dim[1]>0){
  190. data_ptr=mem::aligned_malloc<T>(M.dim[0]*M.dim[1]);
  191. #ifndef __MIC__
  192. Profile::Add_MEM(M.dim[0]*M.dim[1]*sizeof(T));
  193. #endif
  194. }
  195. }
  196. dim[0]=M.dim[0];
  197. dim[1]=M.dim[1];
  198. mem::memcopy(data_ptr,M.data_ptr,dim[0]*dim[1]*sizeof(T));
  199. }
  200. return *this;
  201. }
  202. template <class T>
  203. Matrix<T>& Matrix<T>::operator+=(const Matrix<T>& M){
  204. assert(M.Dim(0)==Dim(0) && M.Dim(1)==Dim(1));
  205. Profile::Add_FLOP(dim[0]*dim[1]);
  206. for(size_t i=0;i<M.Dim(0)*M.Dim(1);i++)
  207. data_ptr[i]+=M.data_ptr[i];
  208. return *this;
  209. }
  210. template <class T>
  211. Matrix<T>& Matrix<T>::operator-=(const Matrix<T>& M){
  212. assert(M.Dim(0)==Dim(0) && M.Dim(1)==Dim(1));
  213. Profile::Add_FLOP(dim[0]*dim[1]);
  214. for(size_t i=0;i<M.Dim(0)*M.Dim(1);i++)
  215. data_ptr[i]-=M.data_ptr[i];
  216. return *this;
  217. }
  218. template <class T>
  219. Matrix<T> Matrix<T>::operator+(const Matrix<T>& M2){
  220. Matrix<T>& M1=*this;
  221. assert(M2.Dim(0)==M1.Dim(0) && M2.Dim(1)==M1.Dim(1));
  222. Profile::Add_FLOP(dim[0]*dim[1]);
  223. Matrix<T> M_r(M1.Dim(0),M1.Dim(1),NULL);
  224. for(size_t i=0;i<M1.Dim(0)*M1.Dim(1);i++)
  225. M_r[0][i]=M1[0][i]+M2[0][i];
  226. return M_r;
  227. }
  228. template <class T>
  229. Matrix<T> Matrix<T>::operator-(const Matrix<T>& M2){
  230. Matrix<T>& M1=*this;
  231. assert(M2.Dim(0)==M1.Dim(0) && M2.Dim(1)==M1.Dim(1));
  232. Profile::Add_FLOP(dim[0]*dim[1]);
  233. Matrix<T> M_r(M1.Dim(0),M1.Dim(1),NULL);
  234. for(size_t i=0;i<M1.Dim(0)*M1.Dim(1);i++)
  235. M_r[0][i]=M1[0][i]-M2[0][i];
  236. return M_r;
  237. }
  238. template <class T>
  239. inline T& Matrix<T>::operator()(size_t i,size_t j) const{
  240. assert(i<dim[0] && j<dim[1]);
  241. return data_ptr[i*dim[1]+j];
  242. }
  243. template <class T>
  244. inline T* Matrix<T>::operator[](size_t i) const{
  245. assert(i<dim[0]);
  246. return &data_ptr[i*dim[1]];
  247. }
  248. template <class T>
  249. Matrix<T> Matrix<T>::operator*(const Matrix<T>& M){
  250. assert(dim[1]==M.dim[0]);
  251. Profile::Add_FLOP(2*(((long long)dim[0])*dim[1])*M.dim[1]);
  252. Matrix<T> M_r(dim[0],M.dim[1],NULL);
  253. mat::gemm('N','N',M.dim[1],dim[0],dim[1],
  254. 1.0,M.data_ptr,M.dim[1],data_ptr,dim[1],0.0,M_r.data_ptr,M_r.dim[1]);
  255. return M_r;
  256. }
  257. template <class T>
  258. void Matrix<T>::DGEMM(Matrix<T>& M_r, const Matrix<T>& A, const Matrix<T>& B, T beta){
  259. assert(A.dim[1]==B.dim[0]);
  260. assert(M_r.dim[0]==A.dim[0]);
  261. assert(M_r.dim[1]==B.dim[1]);
  262. #ifndef __MIC__
  263. Profile::Add_FLOP(2*(((long long)A.dim[0])*A.dim[1])*B.dim[1]);
  264. #endif
  265. mat::gemm('N','N',B.dim[1],A.dim[0],A.dim[1],
  266. 1.0,B.data_ptr,B.dim[1],A.data_ptr,A.dim[1],beta,M_r.data_ptr,M_r.dim[1]);
  267. }
  268. // cublasXgemm wrapper
  269. #if defined(PVFMM_HAVE_CUDA)
  270. template <class T>
  271. void Matrix<T>::CUBLASXGEMM(Matrix<T>& M_r, const Matrix<T>& A, const Matrix<T>& B, T beta){
  272. assert(A.dim[1]==B.dim[0]);
  273. assert(M_r.dim[0]==A.dim[0]);
  274. assert(M_r.dim[1]==B.dim[1]);
  275. Profile::Add_FLOP(2*(((long long)A.dim[0])*A.dim[1])*B.dim[1]);
  276. mat::cublasXgemm('N', 'N', B.dim[1], A.dim[0], A.dim[1],
  277. 1.0, B.data_ptr, B.dim[1], A.data_ptr, A.dim[1], beta, M_r.data_ptr, M_r.dim[1]);
  278. }
  279. #endif
  280. #define myswap(t,a,b) {t c=a;a=b;b=c;}
  281. template <class T>
  282. void Matrix<T>::RowPerm(const Permutation<T>& P){
  283. Matrix<T>& M=*this;
  284. if(P.Dim()==0) return;
  285. assert(M.Dim(0)==P.Dim());
  286. size_t d0=M.Dim(0);
  287. size_t d1=M.Dim(1);
  288. #pragma omp parallel for
  289. for(size_t i=0;i<d0;i++){
  290. T* M_=M[i];
  291. const T s=P.scal[i];
  292. for(size_t j=0;j<d1;j++) M_[j]*=s;
  293. }
  294. Permutation<T> P_=P;
  295. for(size_t i=0;i<d0;i++)
  296. while(P_.perm[i]!=i){
  297. size_t a=P_.perm[i];
  298. size_t b=i;
  299. T* M_a=M[a];
  300. T* M_b=M[b];
  301. myswap(size_t,P_.perm[a],P_.perm[b]);
  302. for(size_t j=0;j<d1;j++)
  303. myswap(T,M_a[j],M_b[j]);
  304. }
  305. }
  306. template <class T>
  307. void Matrix<T>::ColPerm(const Permutation<T>& P){
  308. Matrix<T>& M=*this;
  309. if(P.Dim()==0) return;
  310. assert(M.Dim(1)==P.Dim());
  311. size_t d0=M.Dim(0);
  312. size_t d1=M.Dim(1);
  313. int omp_p=omp_get_max_threads();
  314. Matrix<T> M_buff(omp_p,d1);
  315. const size_t* perm_=&(P.perm[0]);
  316. const T* scal_=&(P.scal[0]);
  317. #pragma omp parallel for
  318. for(size_t i=0;i<d0;i++){
  319. int pid=omp_get_thread_num();
  320. T* buff=&M_buff[pid][0];
  321. T* M_=M[i];
  322. for(size_t j=0;j<d1;j++)
  323. buff[j]=M_[j];
  324. for(size_t j=0;j<d1;j++){
  325. M_[j]=buff[perm_[j]]*scal_[j];
  326. }
  327. }
  328. }
  329. #undef myswap
  330. #define B1 128
  331. #define B2 32
  332. template <class T>
  333. Matrix<T> Matrix<T>::Transpose(){
  334. Matrix<T>& M=*this;
  335. size_t d0=M.dim[0];
  336. size_t d1=M.dim[1];
  337. Matrix<T> M_r(d1,d0,NULL);
  338. const size_t blk0=((d0+B1-1)/B1);
  339. const size_t blk1=((d1+B1-1)/B1);
  340. const size_t blks=blk0*blk1;
  341. // #pragma omp parallel for
  342. for(size_t k=0;k<blks;k++){
  343. size_t i=(k%blk0)*B1;
  344. size_t j=(k/blk0)*B1;
  345. // for(size_t i=0;i<d0;i+=B1)
  346. // for(size_t j=0;j<d1;j+=B1){
  347. size_t d0_=i+B1; if(d0_>=d0) d0_=d0;
  348. size_t d1_=j+B1; if(d1_>=d1) d1_=d1;
  349. for(size_t ii=i;ii<d0_;ii+=B2)
  350. for(size_t jj=j;jj<d1_;jj+=B2){
  351. size_t d0__=ii+B2; if(d0__>=d0) d0__=d0;
  352. size_t d1__=jj+B2; if(d1__>=d1) d1__=d1;
  353. for(size_t iii=ii;iii<d0__;iii++)
  354. for(size_t jjj=jj;jjj<d1__;jjj++){
  355. M_r[jjj][iii]=M[iii][jjj];
  356. }
  357. }
  358. }
  359. // for(size_t i=0;i<d0;i++)
  360. // for(size_t j=0;j<d1;j++)
  361. // M_r[j][i]=M[i][j];
  362. return M_r;
  363. }
  364. template <class T>
  365. void Matrix<T>::Transpose(Matrix<T>& M_r, const Matrix<T>& M){
  366. size_t d0=M.dim[0];
  367. size_t d1=M.dim[1];
  368. M_r.Resize(d1, d0);
  369. const size_t blk0=((d0+B1-1)/B1);
  370. const size_t blk1=((d1+B1-1)/B1);
  371. const size_t blks=blk0*blk1;
  372. #pragma omp parallel for
  373. for(size_t k=0;k<blks;k++){
  374. size_t i=(k%blk0)*B1;
  375. size_t j=(k/blk0)*B1;
  376. // for(size_t i=0;i<d0;i+=B1)
  377. // for(size_t j=0;j<d1;j+=B1){
  378. size_t d0_=i+B1; if(d0_>=d0) d0_=d0;
  379. size_t d1_=j+B1; if(d1_>=d1) d1_=d1;
  380. for(size_t ii=i;ii<d0_;ii+=B2)
  381. for(size_t jj=j;jj<d1_;jj+=B2){
  382. size_t d0__=ii+B2; if(d0__>=d0) d0__=d0;
  383. size_t d1__=jj+B2; if(d1__>=d1) d1__=d1;
  384. for(size_t iii=ii;iii<d0__;iii++)
  385. for(size_t jjj=jj;jjj<d1__;jjj++){
  386. M_r[jjj][iii]=M[iii][jjj];
  387. }
  388. }
  389. }
  390. }
  391. #undef B2
  392. #undef B1
  393. template <class T>
  394. Matrix<T> Matrix<T>::pinv(){
  395. T eps=(typeid(T)==typeid(float)?4*1e-5:4*1e-9);
  396. return pinv(eps);
  397. }
  398. template <class T>
  399. Matrix<T> Matrix<T>::pinv(T eps){
  400. Matrix<T> M_r(dim[1],dim[0]);
  401. mat::pinv(data_ptr,dim[0],dim[1],eps,M_r.data_ptr);
  402. return M_r;
  403. }
  404. template <class T>
  405. std::ostream& operator<<(std::ostream& output, const Permutation<T>& P){
  406. output<<std::setprecision(4)<<std::setiosflags(std::ios::left);
  407. size_t size=P.perm.size();
  408. for(size_t i=0;i<size;i++) output<<std::setw(10)<<P.perm[i]<<' ';
  409. output<<";\n";
  410. for(size_t i=0;i<size;i++) output<<std::setw(10)<<P.scal[i]<<' ';
  411. output<<";\n";
  412. return output;
  413. }
  414. template <class T>
  415. Permutation<T>::Permutation(size_t size){
  416. perm.Resize(size);
  417. scal.Resize(size);
  418. for(size_t i=0;i<size;i++){
  419. perm[i]=i;
  420. scal[i]=1.0;
  421. }
  422. }
  423. template <class T>
  424. Permutation<T> Permutation<T>::RandPerm(size_t size){
  425. Permutation<T> P(size);
  426. for(size_t i=0;i<size;i++){
  427. P.perm[i]=rand()%size;
  428. for(size_t j=0;j<i;j++)
  429. if(P.perm[i]==P.perm[j]){ i--; break; }
  430. P.scal[i]=((T)rand())/RAND_MAX;
  431. }
  432. return P;
  433. }
  434. template <class T>
  435. Matrix<T> Permutation<T>::GetMatrix() const{
  436. size_t size=perm.Dim();
  437. Matrix<T> M_r(size,size,NULL);
  438. for(size_t i=0;i<size;i++)
  439. for(size_t j=0;j<size;j++)
  440. M_r[i][j]=(perm[j]==i?scal[j]:0.0);
  441. return M_r;
  442. }
  443. template <class T>
  444. size_t Permutation<T>::Dim() const{
  445. return perm.Dim();
  446. }
  447. template <class T>
  448. Permutation<T> Permutation<T>::Transpose(){
  449. size_t size=perm.Dim();
  450. Permutation<T> P_r(size);
  451. Vector<PERM_INT_T>& perm_r=P_r.perm;
  452. Vector<T>& scal_r=P_r.scal;
  453. for(size_t i=0;i<size;i++){
  454. perm_r[perm[i]]=i;
  455. scal_r[perm[i]]=scal[i];
  456. }
  457. return P_r;
  458. }
  459. template <class T>
  460. Permutation<T> Permutation<T>::operator*(const Permutation<T>& P){
  461. size_t size=perm.Dim();
  462. assert(P.Dim()==size);
  463. Permutation<T> P_r(size);
  464. Vector<PERM_INT_T>& perm_r=P_r.perm;
  465. Vector<T>& scal_r=P_r.scal;
  466. for(size_t i=0;i<size;i++){
  467. perm_r[i]=perm[P.perm[i]];
  468. scal_r[i]=scal[P.perm[i]]*P.scal[i];
  469. }
  470. return P_r;
  471. }
  472. template <class T>
  473. Matrix<T> Permutation<T>::operator*(const Matrix<T>& M){
  474. if(Dim()==0) return M;
  475. assert(M.Dim(0)==Dim());
  476. size_t d0=M.Dim(0);
  477. size_t d1=M.Dim(1);
  478. Matrix<T> M_r(d0,d1,NULL);
  479. for(size_t i=0;i<d0;i++){
  480. const T s=scal[i];
  481. const T* M_=M[i];
  482. T* M_r_=M_r[perm[i]];
  483. for(size_t j=0;j<d1;j++)
  484. M_r_[j]=M_[j]*s;
  485. }
  486. return M_r;
  487. }
  488. template <class T>
  489. Matrix<T> operator*(const Matrix<T>& M, const Permutation<T>& P){
  490. if(P.Dim()==0) return M;
  491. assert(M.Dim(1)==P.Dim());
  492. size_t d0=M.Dim(0);
  493. size_t d1=M.Dim(1);
  494. Matrix<T> M_r(d0,d1,NULL);
  495. for(size_t i=0;i<d0;i++){
  496. const PERM_INT_T* perm_=&(P.perm[0]);
  497. const T* scal_=&(P.scal[0]);
  498. const T* M_=M[i];
  499. T* M_r_=M_r[i];
  500. for(size_t j=0;j<d1;j++)
  501. M_r_[j]=M_[perm_[j]]*scal_[j];
  502. }
  503. return M_r;
  504. }
  505. }//end namespace