matrix.txx 13 KB

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