matrix.txx 14 KB

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