matrix.txx 13 KB

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