matrix.txx 13 KB

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