matrix.txx 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578
  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. #if defined(PVFMM_HAVE_CUDA)
  119. cudaEventCreate(&lock);
  120. cudaEventRecord(lock, 0);
  121. #endif
  122. }
  123. template <class T>
  124. void Matrix<T>::Device2HostWait(){
  125. #if defined(PVFMM_HAVE_CUDA)
  126. cudaEventSynchronize(lock);
  127. cudaEventDestroy(lock);
  128. #endif
  129. DeviceWrapper::wait(dev.lock_idx);
  130. dev.lock_idx=-1;
  131. }
  132. template <class T>
  133. void Matrix<T>::FreeDevice(bool copy){
  134. if(dev.dev_ptr==(uintptr_t)NULL) return;
  135. if(copy) DeviceWrapper::device2host((char*)data_ptr,dev.dev_ptr,(char*)data_ptr,dim[0]*dim[1]*sizeof(T));
  136. DeviceWrapper::free_device((char*)data_ptr, dev.dev_ptr);
  137. dev.dev_ptr=(uintptr_t)NULL;
  138. dev.dim[0]=0;
  139. dev.dim[1]=0;
  140. }
  141. template <class T>
  142. void Matrix<T>::Write(const char* fname){
  143. FILE* f1=fopen(fname,"wb+");
  144. if(f1==NULL){
  145. std::cout<<"Unable to open file for writing:"<<fname<<'\n';
  146. return;
  147. }
  148. size_t dim_[2]={dim[0],dim[1]};
  149. fwrite(dim_,sizeof(size_t),2,f1);
  150. fwrite(data_ptr,sizeof(T),dim[0]*dim[1],f1);
  151. fclose(f1);
  152. }
  153. template <class T>
  154. size_t Matrix<T>::Dim(size_t i) const{
  155. return dim[i];
  156. }
  157. template <class T>
  158. void Matrix<T>::Resize(size_t i, size_t j){
  159. if(dim[0]==i && dim[1]==j) return;
  160. FreeDevice(false);
  161. if(own_data){
  162. if(data_ptr!=NULL){
  163. mem::aligned_free(data_ptr);
  164. #ifndef __MIC__
  165. Profile::Add_MEM(-dim[0]*dim[1]*sizeof(T));
  166. #endif
  167. }
  168. }
  169. dim[0]=i;
  170. dim[1]=j;
  171. if(own_data){
  172. if(dim[0]*dim[1]>0){
  173. data_ptr=mem::aligned_malloc<T>(dim[0]*dim[1]);
  174. #ifndef __MIC__
  175. Profile::Add_MEM(dim[0]*dim[1]*sizeof(T));
  176. #endif
  177. }else
  178. data_ptr=NULL;
  179. }
  180. }
  181. template <class T>
  182. void Matrix<T>::SetZero(){
  183. if(dim[0]*dim[1]>0)
  184. memset(data_ptr,0,dim[0]*dim[1]*sizeof(T));
  185. }
  186. template <class T>
  187. Matrix<T>& Matrix<T>::operator=(const Matrix<T>& M){
  188. if(this!=&M){
  189. FreeDevice(false);
  190. if(own_data && dim[0]*dim[1]!=M.dim[0]*M.dim[1]){
  191. if(data_ptr!=NULL){
  192. mem::aligned_free(data_ptr); data_ptr=NULL;
  193. #ifndef __MIC__
  194. Profile::Add_MEM(-dim[0]*dim[1]*sizeof(T));
  195. #endif
  196. }
  197. if(M.dim[0]*M.dim[1]>0){
  198. data_ptr=mem::aligned_malloc<T>(M.dim[0]*M.dim[1]);
  199. #ifndef __MIC__
  200. Profile::Add_MEM(M.dim[0]*M.dim[1]*sizeof(T));
  201. #endif
  202. }
  203. }
  204. dim[0]=M.dim[0];
  205. dim[1]=M.dim[1];
  206. mem::memcopy(data_ptr,M.data_ptr,dim[0]*dim[1]*sizeof(T));
  207. }
  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>& M){
  220. assert(M.Dim(0)==Dim(0) && M.Dim(1)==Dim(1));
  221. Profile::Add_FLOP(dim[0]*dim[1]);
  222. for(size_t i=0;i<M.Dim(0)*M.Dim(1);i++)
  223. data_ptr[i]-=M.data_ptr[i];
  224. return *this;
  225. }
  226. template <class T>
  227. Matrix<T> Matrix<T>::operator+(const Matrix<T>& M2){
  228. Matrix<T>& M1=*this;
  229. assert(M2.Dim(0)==M1.Dim(0) && M2.Dim(1)==M1.Dim(1));
  230. Profile::Add_FLOP(dim[0]*dim[1]);
  231. Matrix<T> M_r(M1.Dim(0),M1.Dim(1),NULL);
  232. for(size_t i=0;i<M1.Dim(0)*M1.Dim(1);i++)
  233. M_r[0][i]=M1[0][i]+M2[0][i];
  234. return M_r;
  235. }
  236. template <class T>
  237. Matrix<T> Matrix<T>::operator-(const Matrix<T>& M2){
  238. Matrix<T>& M1=*this;
  239. assert(M2.Dim(0)==M1.Dim(0) && M2.Dim(1)==M1.Dim(1));
  240. Profile::Add_FLOP(dim[0]*dim[1]);
  241. Matrix<T> M_r(M1.Dim(0),M1.Dim(1),NULL);
  242. for(size_t i=0;i<M1.Dim(0)*M1.Dim(1);i++)
  243. M_r[0][i]=M1[0][i]-M2[0][i];
  244. return M_r;
  245. }
  246. template <class T>
  247. inline T& Matrix<T>::operator()(size_t i,size_t j) const{
  248. assert(i<dim[0] && j<dim[1]);
  249. return data_ptr[i*dim[1]+j];
  250. }
  251. template <class T>
  252. inline T* Matrix<T>::operator[](size_t i) const{
  253. assert(i<dim[0]);
  254. return &data_ptr[i*dim[1]];
  255. }
  256. template <class T>
  257. Matrix<T> Matrix<T>::operator*(const Matrix<T>& M){
  258. assert(dim[1]==M.dim[0]);
  259. Profile::Add_FLOP(2*(((long long)dim[0])*dim[1])*M.dim[1]);
  260. Matrix<T> M_r(dim[0],M.dim[1],NULL);
  261. mat::gemm('N','N',M.dim[1],dim[0],dim[1],
  262. 1.0,M.data_ptr,M.dim[1],data_ptr,dim[1],0.0,M_r.data_ptr,M_r.dim[1]);
  263. return M_r;
  264. }
  265. template <class T>
  266. void Matrix<T>::DGEMM(Matrix<T>& M_r, const Matrix<T>& A, const Matrix<T>& B, T beta){
  267. assert(A.dim[1]==B.dim[0]);
  268. assert(M_r.dim[0]==A.dim[0]);
  269. assert(M_r.dim[1]==B.dim[1]);
  270. #ifndef __MIC__
  271. Profile::Add_FLOP(2*(((long long)A.dim[0])*A.dim[1])*B.dim[1]);
  272. #endif
  273. mat::gemm('N','N',B.dim[1],A.dim[0],A.dim[1],
  274. 1.0,B.data_ptr,B.dim[1],A.data_ptr,A.dim[1],beta,M_r.data_ptr,M_r.dim[1]);
  275. }
  276. // cublasXgemm wrapper
  277. #if defined(PVFMM_HAVE_CUDA)
  278. template <class T>
  279. void Matrix<T>::CUBLASXGEMM(Matrix<T>& M_r, const Matrix<T>& A, const Matrix<T>& B, T beta){
  280. assert(A.dim[1]==B.dim[0]);
  281. assert(M_r.dim[0]==A.dim[0]);
  282. assert(M_r.dim[1]==B.dim[1]);
  283. Profile::Add_FLOP(2*(((long long)A.dim[0])*A.dim[1])*B.dim[1]);
  284. mat::cublasXgemm('N', 'N', B.dim[1], A.dim[0], A.dim[1],
  285. 1.0, B.data_ptr, B.dim[1], A.data_ptr, A.dim[1], beta, M_r.data_ptr, M_r.dim[1]);
  286. }
  287. #endif
  288. #define myswap(t,a,b) {t c=a;a=b;b=c;}
  289. template <class T>
  290. void Matrix<T>::RowPerm(const Permutation<T>& P){
  291. Matrix<T>& M=*this;
  292. if(P.Dim()==0) return;
  293. assert(M.Dim(0)==P.Dim());
  294. size_t d0=M.Dim(0);
  295. size_t d1=M.Dim(1);
  296. #pragma omp parallel for
  297. for(size_t i=0;i<d0;i++){
  298. T* M_=M[i];
  299. const T s=P.scal[i];
  300. for(size_t j=0;j<d1;j++) M_[j]*=s;
  301. }
  302. Permutation<T> P_=P;
  303. for(size_t i=0;i<d0;i++)
  304. while(P_.perm[i]!=i){
  305. size_t a=P_.perm[i];
  306. size_t b=i;
  307. T* M_a=M[a];
  308. T* M_b=M[b];
  309. myswap(size_t,P_.perm[a],P_.perm[b]);
  310. for(size_t j=0;j<d1;j++)
  311. myswap(T,M_a[j],M_b[j]);
  312. }
  313. }
  314. template <class T>
  315. void Matrix<T>::ColPerm(const Permutation<T>& P){
  316. Matrix<T>& M=*this;
  317. if(P.Dim()==0) return;
  318. assert(M.Dim(1)==P.Dim());
  319. size_t d0=M.Dim(0);
  320. size_t d1=M.Dim(1);
  321. int omp_p=omp_get_max_threads();
  322. Matrix<T> M_buff(omp_p,d1);
  323. const size_t* perm_=&(P.perm[0]);
  324. const T* scal_=&(P.scal[0]);
  325. #pragma omp parallel for
  326. for(size_t i=0;i<d0;i++){
  327. int pid=omp_get_thread_num();
  328. T* buff=&M_buff[pid][0];
  329. T* M_=M[i];
  330. for(size_t j=0;j<d1;j++)
  331. buff[j]=M_[j];
  332. for(size_t j=0;j<d1;j++){
  333. M_[j]=buff[perm_[j]]*scal_[j];
  334. }
  335. }
  336. }
  337. #undef myswap
  338. #define B1 128
  339. #define B2 32
  340. template <class T>
  341. Matrix<T> Matrix<T>::Transpose(){
  342. Matrix<T>& M=*this;
  343. size_t d0=M.dim[0];
  344. size_t d1=M.dim[1];
  345. Matrix<T> M_r(d1,d0,NULL);
  346. const size_t blk0=((d0+B1-1)/B1);
  347. const size_t blk1=((d1+B1-1)/B1);
  348. const size_t blks=blk0*blk1;
  349. // #pragma omp parallel for
  350. for(size_t k=0;k<blks;k++){
  351. size_t i=(k%blk0)*B1;
  352. size_t j=(k/blk0)*B1;
  353. // for(size_t i=0;i<d0;i+=B1)
  354. // for(size_t j=0;j<d1;j+=B1){
  355. size_t d0_=i+B1; if(d0_>=d0) d0_=d0;
  356. size_t d1_=j+B1; if(d1_>=d1) d1_=d1;
  357. for(size_t ii=i;ii<d0_;ii+=B2)
  358. for(size_t jj=j;jj<d1_;jj+=B2){
  359. size_t d0__=ii+B2; if(d0__>=d0) d0__=d0;
  360. size_t d1__=jj+B2; if(d1__>=d1) d1__=d1;
  361. for(size_t iii=ii;iii<d0__;iii++)
  362. for(size_t jjj=jj;jjj<d1__;jjj++){
  363. M_r[jjj][iii]=M[iii][jjj];
  364. }
  365. }
  366. }
  367. // for(size_t i=0;i<d0;i++)
  368. // for(size_t j=0;j<d1;j++)
  369. // M_r[j][i]=M[i][j];
  370. return M_r;
  371. }
  372. template <class T>
  373. void Matrix<T>::Transpose(Matrix<T>& M_r, const Matrix<T>& M){
  374. size_t d0=M.dim[0];
  375. size_t d1=M.dim[1];
  376. M_r.Resize(d1, d0);
  377. const size_t blk0=((d0+B1-1)/B1);
  378. const size_t blk1=((d1+B1-1)/B1);
  379. const size_t blks=blk0*blk1;
  380. #pragma omp parallel for
  381. for(size_t k=0;k<blks;k++){
  382. size_t i=(k%blk0)*B1;
  383. size_t j=(k/blk0)*B1;
  384. // for(size_t i=0;i<d0;i+=B1)
  385. // for(size_t j=0;j<d1;j+=B1){
  386. size_t d0_=i+B1; if(d0_>=d0) d0_=d0;
  387. size_t d1_=j+B1; if(d1_>=d1) d1_=d1;
  388. for(size_t ii=i;ii<d0_;ii+=B2)
  389. for(size_t jj=j;jj<d1_;jj+=B2){
  390. size_t d0__=ii+B2; if(d0__>=d0) d0__=d0;
  391. size_t d1__=jj+B2; if(d1__>=d1) d1__=d1;
  392. for(size_t iii=ii;iii<d0__;iii++)
  393. for(size_t jjj=jj;jjj<d1__;jjj++){
  394. M_r[jjj][iii]=M[iii][jjj];
  395. }
  396. }
  397. }
  398. }
  399. #undef B2
  400. #undef B1
  401. template <class T>
  402. Matrix<T> Matrix<T>::pinv(){
  403. T eps=(typeid(T)==typeid(float)?4*1e-5:4*1e-9);
  404. return pinv(eps);
  405. }
  406. template <class T>
  407. Matrix<T> Matrix<T>::pinv(T eps){
  408. Matrix<T> M_r(dim[1],dim[0]);
  409. mat::pinv(data_ptr,dim[0],dim[1],eps,M_r.data_ptr);
  410. return M_r;
  411. }
  412. template <class T>
  413. std::ostream& operator<<(std::ostream& output, const Permutation<T>& P){
  414. output<<std::setprecision(4)<<std::setiosflags(std::ios::left);
  415. size_t size=P.perm.size();
  416. for(size_t i=0;i<size;i++) output<<std::setw(10)<<P.perm[i]<<' ';
  417. output<<";\n";
  418. for(size_t i=0;i<size;i++) output<<std::setw(10)<<P.scal[i]<<' ';
  419. output<<";\n";
  420. return output;
  421. }
  422. template <class T>
  423. Permutation<T>::Permutation(size_t size){
  424. perm.Resize(size);
  425. scal.Resize(size);
  426. for(size_t i=0;i<size;i++){
  427. perm[i]=i;
  428. scal[i]=1.0;
  429. }
  430. }
  431. template <class T>
  432. Permutation<T> Permutation<T>::RandPerm(size_t size){
  433. Permutation<T> P(size);
  434. for(size_t i=0;i<size;i++){
  435. P.perm[i]=rand()%size;
  436. for(size_t j=0;j<i;j++)
  437. if(P.perm[i]==P.perm[j]){ i--; break; }
  438. P.scal[i]=((T)rand())/RAND_MAX;
  439. }
  440. return P;
  441. }
  442. template <class T>
  443. Matrix<T> Permutation<T>::GetMatrix() const{
  444. size_t size=perm.Dim();
  445. Matrix<T> M_r(size,size,NULL);
  446. for(size_t i=0;i<size;i++)
  447. for(size_t j=0;j<size;j++)
  448. M_r[i][j]=(perm[j]==i?scal[j]:0.0);
  449. return M_r;
  450. }
  451. template <class T>
  452. size_t Permutation<T>::Dim() const{
  453. return perm.Dim();
  454. }
  455. template <class T>
  456. Permutation<T> Permutation<T>::Transpose(){
  457. size_t size=perm.Dim();
  458. Permutation<T> P_r(size);
  459. Vector<PERM_INT_T>& perm_r=P_r.perm;
  460. Vector<T>& scal_r=P_r.scal;
  461. for(size_t i=0;i<size;i++){
  462. perm_r[perm[i]]=i;
  463. scal_r[perm[i]]=scal[i];
  464. }
  465. return P_r;
  466. }
  467. template <class T>
  468. Permutation<T> Permutation<T>::operator*(const Permutation<T>& P){
  469. size_t size=perm.Dim();
  470. assert(P.Dim()==size);
  471. Permutation<T> P_r(size);
  472. Vector<PERM_INT_T>& perm_r=P_r.perm;
  473. Vector<T>& scal_r=P_r.scal;
  474. for(size_t i=0;i<size;i++){
  475. perm_r[i]=perm[P.perm[i]];
  476. scal_r[i]=scal[P.perm[i]]*P.scal[i];
  477. }
  478. return P_r;
  479. }
  480. template <class T>
  481. Matrix<T> Permutation<T>::operator*(const Matrix<T>& M){
  482. if(Dim()==0) return M;
  483. assert(M.Dim(0)==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 T s=scal[i];
  489. const T* M_=M[i];
  490. T* M_r_=M_r[perm[i]];
  491. for(size_t j=0;j<d1;j++)
  492. M_r_[j]=M_[j]*s;
  493. }
  494. return M_r;
  495. }
  496. template <class T>
  497. Matrix<T> operator*(const Matrix<T>& M, const Permutation<T>& P){
  498. if(P.Dim()==0) return M;
  499. assert(M.Dim(1)==P.Dim());
  500. size_t d0=M.Dim(0);
  501. size_t d1=M.Dim(1);
  502. Matrix<T> M_r(d0,d1,NULL);
  503. for(size_t i=0;i<d0;i++){
  504. const PERM_INT_T* perm_=&(P.perm[0]);
  505. const T* scal_=&(P.scal[0]);
  506. const T* M_=M[i];
  507. T* M_r_=M_r[i];
  508. for(size_t j=0;j<d1;j++)
  509. M_r_[j]=M_[perm_[j]]*scal_[j];
  510. }
  511. return M_r;
  512. }
  513. }//end namespace