mat_utils.txx 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518
  1. /**
  2. * \file mat_utils.txx
  3. * \author Dhairya Malhotra, dhairya.malhotra@gmail.com
  4. * \date 2-11-2011
  5. * \brief This file contains BLAS and LAPACK wrapper functions.
  6. */
  7. #include <omp.h>
  8. #include <cmath>
  9. #include <cassert>
  10. #include <cstdlib>
  11. #include <algorithm>
  12. #include <iostream>
  13. #include <vector>
  14. #include <matrix.hpp>
  15. #include <device_wrapper.hpp>
  16. #if defined(PVFMM_HAVE_CUDA)
  17. #include <cuda_runtime_api.h>
  18. #include <cublas_v2.h>
  19. #endif
  20. namespace pvfmm{
  21. namespace mat{
  22. template <class T>
  23. inline void gemm(char TransA, char TransB, int M, int N, int K, T alpha, T *A, int lda, T *B, int ldb, T beta, T *C, int ldc){
  24. if((TransA=='N' || TransA=='n') && (TransB=='N' || TransB=='n')){
  25. for(size_t n=0;n<N;n++){ // Columns of C
  26. for(size_t m=0;m<M;m++){ // Rows of C
  27. T AxB=0;
  28. for(size_t k=0;k<K;k++){
  29. AxB+=A[m+lda*k]*B[k+ldb*n];
  30. }
  31. C[m+ldc*n]=alpha*AxB+(beta==0?0:beta*C[m+ldc*n]);
  32. }
  33. }
  34. }else if(TransA=='N' || TransA=='n'){
  35. for(size_t n=0;n<N;n++){ // Columns of C
  36. for(size_t m=0;m<M;m++){ // Rows of C
  37. T AxB=0;
  38. for(size_t k=0;k<K;k++){
  39. AxB+=A[m+lda*k]*B[n+ldb*k];
  40. }
  41. C[m+ldc*n]=alpha*AxB+(beta==0?0:beta*C[m+ldc*n]);
  42. }
  43. }
  44. }else if(TransB=='N' || TransB=='n'){
  45. for(size_t n=0;n<N;n++){ // Columns of C
  46. for(size_t m=0;m<M;m++){ // Rows of C
  47. T AxB=0;
  48. for(size_t k=0;k<K;k++){
  49. AxB+=A[k+lda*m]*B[k+ldb*n];
  50. }
  51. C[m+ldc*n]=alpha*AxB+(beta==0?0:beta*C[m+ldc*n]);
  52. }
  53. }
  54. }else{
  55. for(size_t n=0;n<N;n++){ // Columns of C
  56. for(size_t m=0;m<M;m++){ // Rows of C
  57. T AxB=0;
  58. for(size_t k=0;k<K;k++){
  59. AxB+=A[k+lda*m]*B[n+ldb*k];
  60. }
  61. C[m+ldc*n]=alpha*AxB+(beta==0?0:beta*C[m+ldc*n]);
  62. }
  63. }
  64. }
  65. }
  66. template<>
  67. void gemm<float>(char TransA, char TransB, int M, int N, int K, float alpha, float *A, int lda, float *B, int ldb, float beta, float *C, int ldc);
  68. template<>
  69. void gemm<double>(char TransA, char TransB, int M, int N, int K, double alpha, double *A, int lda, double *B, int ldb, double beta, double *C, int ldc);
  70. #if defined(PVFMM_HAVE_CUDA)
  71. template<>
  72. inline void cublasgemm<float>(char TransA, char TransB, int M, int N, int K, float alpha, float *A, int lda, float *B, int ldb, float beta, float *C, int ldc) {
  73. cublasOperation_t cublasTransA, cublasTransB;
  74. cublasHandle_t *handle = CUDA_Lock::acquire_handle();
  75. if (TransA == 'T' || TransA == 't') cublasTransA = CUBLAS_OP_T;
  76. else if (TransA == 'N' || TransA == 'n') cublasTransA = CUBLAS_OP_N;
  77. if (TransB == 'T' || TransB == 't') cublasTransB = CUBLAS_OP_T;
  78. else if (TransB == 'N' || TransB == 'n') cublasTransB = CUBLAS_OP_N;
  79. cublasStatus_t status = cublasSgemm(*handle, cublasTransA, cublasTransB, M, N, K, &alpha, A, lda, B, ldb, &beta, C, ldc);
  80. }
  81. template<>
  82. inline void cublasgemm<double>(char TransA, char TransB, int M, int N, int K, double alpha, double *A, int lda, double *B, int ldb, double beta, double *C, int ldc){
  83. cublasOperation_t cublasTransA, cublasTransB;
  84. cublasHandle_t *handle = CUDA_Lock::acquire_handle();
  85. if (TransA == 'T' || TransA == 't') cublasTransA = CUBLAS_OP_T;
  86. else if (TransA == 'N' || TransA == 'n') cublasTransA = CUBLAS_OP_N;
  87. if (TransB == 'T' || TransB == 't') cublasTransB = CUBLAS_OP_T;
  88. else if (TransB == 'N' || TransB == 'n') cublasTransB = CUBLAS_OP_N;
  89. cublasStatus_t status = cublasDgemm(*handle, cublasTransA, cublasTransB, M, N, K, &alpha, A, lda, B, ldb, &beta, C, ldc);
  90. }
  91. #endif
  92. #define U(i,j) U_[(i)*dim[0]+(j)]
  93. #define S(i,j) S_[(i)*dim[1]+(j)]
  94. #define V(i,j) V_[(i)*dim[1]+(j)]
  95. //#define SVD_DEBUG
  96. template <class T>
  97. static void GivensL(T* S_, const size_t dim[2], size_t m, T a, T b){
  98. T r=pvfmm::sqrt<T>(a*a+b*b);
  99. T c=a/r;
  100. T s=-b/r;
  101. #pragma omp parallel for
  102. for(size_t i=0;i<dim[1];i++){
  103. T S0=S(m+0,i);
  104. T S1=S(m+1,i);
  105. S(m ,i)+=S0*(c-1);
  106. S(m ,i)+=S1*(-s );
  107. S(m+1,i)+=S0*( s );
  108. S(m+1,i)+=S1*(c-1);
  109. }
  110. }
  111. template <class T>
  112. static void GivensR(T* S_, const size_t dim[2], size_t m, T a, T b){
  113. T r=pvfmm::sqrt<T>(a*a+b*b);
  114. T c=a/r;
  115. T s=-b/r;
  116. #pragma omp parallel for
  117. for(size_t i=0;i<dim[0];i++){
  118. T S0=S(i,m+0);
  119. T S1=S(i,m+1);
  120. S(i,m )+=S0*(c-1);
  121. S(i,m )+=S1*(-s );
  122. S(i,m+1)+=S0*( s );
  123. S(i,m+1)+=S1*(c-1);
  124. }
  125. }
  126. template <class T>
  127. static void SVD(const size_t dim[2], T* U_, T* S_, T* V_, T eps=-1){
  128. assert(dim[0]>=dim[1]);
  129. #ifdef SVD_DEBUG
  130. Matrix<T> M0(dim[0],dim[1],S_);
  131. #endif
  132. { // Bi-diagonalization
  133. size_t n=std::min(dim[0],dim[1]);
  134. std::vector<T> house_vec(std::max(dim[0],dim[1]));
  135. for(size_t i=0;i<n;i++){
  136. // Column Householder
  137. {
  138. T x1=S(i,i);
  139. if(x1<0) x1=-x1;
  140. T x_inv_norm=0;
  141. for(size_t j=i;j<dim[0];j++){
  142. x_inv_norm+=S(j,i)*S(j,i);
  143. }
  144. if(x_inv_norm>0) x_inv_norm=1/pvfmm::sqrt<T>(x_inv_norm);
  145. T alpha=pvfmm::sqrt<T>(1+x1*x_inv_norm);
  146. T beta=x_inv_norm/alpha;
  147. house_vec[i]=-alpha;
  148. for(size_t j=i+1;j<dim[0];j++){
  149. house_vec[j]=-beta*S(j,i);
  150. }
  151. if(S(i,i)<0) for(size_t j=i+1;j<dim[0];j++){
  152. house_vec[j]=-house_vec[j];
  153. }
  154. }
  155. #pragma omp parallel for
  156. for(size_t k=i;k<dim[1];k++){
  157. T dot_prod=0;
  158. for(size_t j=i;j<dim[0];j++){
  159. dot_prod+=S(j,k)*house_vec[j];
  160. }
  161. for(size_t j=i;j<dim[0];j++){
  162. S(j,k)-=dot_prod*house_vec[j];
  163. }
  164. }
  165. #pragma omp parallel for
  166. for(size_t k=0;k<dim[0];k++){
  167. T dot_prod=0;
  168. for(size_t j=i;j<dim[0];j++){
  169. dot_prod+=U(k,j)*house_vec[j];
  170. }
  171. for(size_t j=i;j<dim[0];j++){
  172. U(k,j)-=dot_prod*house_vec[j];
  173. }
  174. }
  175. // Row Householder
  176. if(i>=n-1) continue;
  177. {
  178. T x1=S(i,i+1);
  179. if(x1<0) x1=-x1;
  180. T x_inv_norm=0;
  181. for(size_t j=i+1;j<dim[1];j++){
  182. x_inv_norm+=S(i,j)*S(i,j);
  183. }
  184. if(x_inv_norm>0) x_inv_norm=1/pvfmm::sqrt<T>(x_inv_norm);
  185. T alpha=pvfmm::sqrt<T>(1+x1*x_inv_norm);
  186. T beta=x_inv_norm/alpha;
  187. house_vec[i+1]=-alpha;
  188. for(size_t j=i+2;j<dim[1];j++){
  189. house_vec[j]=-beta*S(i,j);
  190. }
  191. if(S(i,i+1)<0) for(size_t j=i+2;j<dim[1];j++){
  192. house_vec[j]=-house_vec[j];
  193. }
  194. }
  195. #pragma omp parallel for
  196. for(size_t k=i;k<dim[0];k++){
  197. T dot_prod=0;
  198. for(size_t j=i+1;j<dim[1];j++){
  199. dot_prod+=S(k,j)*house_vec[j];
  200. }
  201. for(size_t j=i+1;j<dim[1];j++){
  202. S(k,j)-=dot_prod*house_vec[j];
  203. }
  204. }
  205. #pragma omp parallel for
  206. for(size_t k=0;k<dim[1];k++){
  207. T dot_prod=0;
  208. for(size_t j=i+1;j<dim[1];j++){
  209. dot_prod+=V(j,k)*house_vec[j];
  210. }
  211. for(size_t j=i+1;j<dim[1];j++){
  212. V(j,k)-=dot_prod*house_vec[j];
  213. }
  214. }
  215. }
  216. }
  217. size_t k0=0;
  218. size_t iter=0;
  219. if(eps<0){
  220. eps=1.0;
  221. while(eps+(T)1.0>1.0) eps*=0.5;
  222. eps*=64.0;
  223. }
  224. while(k0<dim[1]-1){ // Diagonalization
  225. iter++;
  226. T S_max=0.0;
  227. for(size_t i=0;i<dim[1];i++) S_max=(S_max>S(i,i)?S_max:S(i,i));
  228. //while(k0<dim[1]-1 && pvfmm::fabs<T>(S(k0,k0+1))<=eps*(pvfmm::fabs<T>(S(k0,k0))+pvfmm::fabs<T>(S(k0+1,k0+1)))) k0++;
  229. while(k0<dim[1]-1 && pvfmm::fabs<T>(S(k0,k0+1))<=eps*S_max) k0++;
  230. if(k0==dim[1]-1) continue;
  231. size_t n=k0+2;
  232. //while(n<dim[1] && pvfmm::fabs<T>(S(n-1,n))>eps*(pvfmm::fabs<T>(S(n-1,n-1))+pvfmm::fabs<T>(S(n,n)))) n++;
  233. while(n<dim[1] && pvfmm::fabs<T>(S(n-1,n))>eps*S_max) n++;
  234. T alpha=0;
  235. T beta=0;
  236. { // Compute mu
  237. T C[2][2];
  238. C[0][0]=S(n-2,n-2)*S(n-2,n-2);
  239. if(n-k0>2) C[0][0]+=S(n-3,n-2)*S(n-3,n-2);
  240. C[0][1]=S(n-2,n-2)*S(n-2,n-1);
  241. C[1][0]=S(n-2,n-2)*S(n-2,n-1);
  242. C[1][1]=S(n-1,n-1)*S(n-1,n-1)+S(n-2,n-1)*S(n-2,n-1);
  243. T b=-(C[0][0]+C[1][1])/2;
  244. T c= C[0][0]*C[1][1] - C[0][1]*C[1][0];
  245. T d=0;
  246. if(b*b-c>0) d=pvfmm::sqrt<T>(b*b-c);
  247. else{
  248. T b=(C[0][0]-C[1][1])/2;
  249. T c=-C[0][1]*C[1][0];
  250. if(b*b-c>0) d=pvfmm::sqrt<T>(b*b-c);
  251. }
  252. T lambda1=-b+d;
  253. T lambda2=-b-d;
  254. T d1=lambda1-C[1][1]; d1=(d1<0?-d1:d1);
  255. T d2=lambda2-C[1][1]; d2=(d2<0?-d2:d2);
  256. T mu=(d1<d2?lambda1:lambda2);
  257. alpha=S(k0,k0)*S(k0,k0)-mu;
  258. beta=S(k0,k0)*S(k0,k0+1);
  259. }
  260. for(size_t k=k0;k<n-1;k++){
  261. size_t dimU[2]={dim[0],dim[0]};
  262. size_t dimV[2]={dim[1],dim[1]};
  263. GivensR(S_,dim ,k,alpha,beta);
  264. GivensL(V_,dimV,k,alpha,beta);
  265. alpha=S(k,k);
  266. beta=S(k+1,k);
  267. GivensL(S_,dim ,k,alpha,beta);
  268. GivensR(U_,dimU,k,alpha,beta);
  269. alpha=S(k,k+1);
  270. beta=S(k,k+2);
  271. }
  272. { // Make S bi-diagonal again
  273. for(size_t i0=k0;i0<n-1;i0++){
  274. for(size_t i1=0;i1<dim[1];i1++){
  275. if(i0>i1 || i0+1<i1) S(i0,i1)=0;
  276. }
  277. }
  278. for(size_t i0=0;i0<dim[0];i0++){
  279. for(size_t i1=k0;i1<n-1;i1++){
  280. if(i0>i1 || i0+1<i1) S(i0,i1)=0;
  281. }
  282. }
  283. for(size_t i=0;i<dim[1]-1;i++){
  284. if(pvfmm::fabs<T>(S(i,i+1))<=eps*S_max){
  285. S(i,i+1)=0;
  286. }
  287. }
  288. }
  289. //std::cout<<iter<<' '<<k0<<' '<<n<<'\n';
  290. }
  291. { // Check Error
  292. #ifdef SVD_DEBUG
  293. Matrix<T> U0(dim[0],dim[0],U_);
  294. Matrix<T> S0(dim[0],dim[1],S_);
  295. Matrix<T> V0(dim[1],dim[1],V_);
  296. Matrix<T> E=M0-U0*S0*V0;
  297. T max_err=0;
  298. T max_nondiag0=0;
  299. T max_nondiag1=0;
  300. for(size_t i=0;i<E.Dim(0);i++)
  301. for(size_t j=0;j<E.Dim(1);j++){
  302. if(max_err<pvfmm::fabs<T>(E[i][j])) max_err=pvfmm::fabs<T>(E[i][j]);
  303. if((i>j+0 || i+0<j) && max_nondiag0<pvfmm::fabs<T>(S0[i][j])) max_nondiag0=pvfmm::fabs<T>(S0[i][j]);
  304. if((i>j+1 || i+1<j) && max_nondiag1<pvfmm::fabs<T>(S0[i][j])) max_nondiag1=pvfmm::fabs<T>(S0[i][j]);
  305. }
  306. std::cout<<max_err<<'\n';
  307. std::cout<<max_nondiag0<<'\n';
  308. std::cout<<max_nondiag1<<'\n';
  309. #endif
  310. }
  311. }
  312. #undef U
  313. #undef S
  314. #undef V
  315. #undef SVD_DEBUG
  316. template<class T>
  317. inline void svd(char *JOBU, char *JOBVT, int *M, int *N, T *A, int *LDA,
  318. T *S, T *U, int *LDU, T *VT, int *LDVT, T *WORK, int *LWORK,
  319. int *INFO){
  320. const size_t dim[2]={std::max(*N,*M), std::min(*N,*M)};
  321. T* U_=mem::aligned_new<T>(dim[0]*dim[0]); memset(U_, 0, dim[0]*dim[0]*sizeof(T));
  322. T* V_=mem::aligned_new<T>(dim[1]*dim[1]); memset(V_, 0, dim[1]*dim[1]*sizeof(T));
  323. T* S_=mem::aligned_new<T>(dim[0]*dim[1]);
  324. const size_t lda=*LDA;
  325. const size_t ldu=*LDU;
  326. const size_t ldv=*LDVT;
  327. if(dim[1]==*M){
  328. for(size_t i=0;i<dim[0];i++)
  329. for(size_t j=0;j<dim[1];j++){
  330. S_[i*dim[1]+j]=A[i*lda+j];
  331. }
  332. }else{
  333. for(size_t i=0;i<dim[0];i++)
  334. for(size_t j=0;j<dim[1];j++){
  335. S_[i*dim[1]+j]=A[j*lda+i];
  336. }
  337. }
  338. for(size_t i=0;i<dim[0];i++){
  339. U_[i*dim[0]+i]=1;
  340. }
  341. for(size_t i=0;i<dim[1];i++){
  342. V_[i*dim[1]+i]=1;
  343. }
  344. SVD<T>(dim, U_, S_, V_, (T)-1);
  345. for(size_t i=0;i<dim[1];i++){ // Set S
  346. S[i]=S_[i*dim[1]+i];
  347. }
  348. if(dim[1]==*M){ // Set U
  349. for(size_t i=0;i<dim[1];i++)
  350. for(size_t j=0;j<*M;j++){
  351. U[j+ldu*i]=V_[j+i*dim[1]]*(S[i]<0.0?-1.0:1.0);
  352. }
  353. }else{
  354. for(size_t i=0;i<dim[1];i++)
  355. for(size_t j=0;j<*M;j++){
  356. U[j+ldu*i]=U_[i+j*dim[0]]*(S[i]<0.0?-1.0:1.0);
  357. }
  358. }
  359. if(dim[0]==*N){ // Set V
  360. for(size_t i=0;i<*N;i++)
  361. for(size_t j=0;j<dim[1];j++){
  362. VT[j+ldv*i]=U_[j+i*dim[0]];
  363. }
  364. }else{
  365. for(size_t i=0;i<*N;i++)
  366. for(size_t j=0;j<dim[1];j++){
  367. VT[j+ldv*i]=V_[i+j*dim[1]];
  368. }
  369. }
  370. for(size_t i=0;i<dim[1];i++){
  371. S[i]=S[i]*(S[i]<0.0?-1.0:1.0);
  372. }
  373. mem::aligned_delete<T>(U_);
  374. mem::aligned_delete<T>(S_);
  375. mem::aligned_delete<T>(V_);
  376. if(0){ // Verify
  377. const size_t dim[2]={std::max(*N,*M), std::min(*N,*M)};
  378. const size_t lda=*LDA;
  379. const size_t ldu=*LDU;
  380. const size_t ldv=*LDVT;
  381. Matrix<T> A1(*M,*N);
  382. Matrix<T> S1(dim[1],dim[1]);
  383. Matrix<T> U1(*M,dim[1]);
  384. Matrix<T> V1(dim[1],*N);
  385. for(size_t i=0;i<*N;i++)
  386. for(size_t j=0;j<*M;j++){
  387. A1[j][i]=A[j+i*lda];
  388. }
  389. S1.SetZero();
  390. for(size_t i=0;i<dim[1];i++){ // Set S
  391. S1[i][i]=S[i];
  392. }
  393. for(size_t i=0;i<dim[1];i++)
  394. for(size_t j=0;j<*M;j++){
  395. U1[j][i]=U[j+ldu*i];
  396. }
  397. for(size_t i=0;i<*N;i++)
  398. for(size_t j=0;j<dim[1];j++){
  399. V1[j][i]=VT[j+ldv*i];
  400. }
  401. std::cout<<U1*S1*V1-A1<<'\n';
  402. }
  403. }
  404. template<>
  405. void svd<float>(char *JOBU, char *JOBVT, int *M, int *N, float *A, int *LDA,
  406. float *S, float *U, int *LDU, float *VT, int *LDVT, float *WORK, int *LWORK,
  407. int *INFO);
  408. template<>
  409. void svd<double>(char *JOBU, char *JOBVT, int *M, int *N, double *A, int *LDA,
  410. double *S, double *U, int *LDU, double *VT, int *LDVT, double *WORK, int *LWORK,
  411. int *INFO);
  412. /**
  413. * \brief Computes the pseudo inverse of matrix M(n1xn2) (in row major form)
  414. * and returns the output M_(n2xn1). Original contents of M are destroyed.
  415. */
  416. template <class T>
  417. void pinv(T* M, int n1, int n2, T eps, T* M_){
  418. if(n1*n2==0) return;
  419. int m = n2;
  420. int n = n1;
  421. int k = (m<n?m:n);
  422. T* tU =mem::aligned_new<T>(m*k);
  423. T* tS =mem::aligned_new<T>(k);
  424. T* tVT=mem::aligned_new<T>(k*n);
  425. //SVD
  426. int INFO=0;
  427. char JOBU = 'S';
  428. char JOBVT = 'S';
  429. //int wssize = max(3*min(m,n)+max(m,n), 5*min(m,n));
  430. int wssize = 3*(m<n?m:n)+(m>n?m:n);
  431. int wssize1 = 5*(m<n?m:n);
  432. wssize = (wssize>wssize1?wssize:wssize1);
  433. T* wsbuf = mem::aligned_new<T>(wssize);
  434. svd(&JOBU, &JOBVT, &m, &n, &M[0], &m, &tS[0], &tU[0], &m, &tVT[0], &k,
  435. wsbuf, &wssize, &INFO);
  436. if(INFO!=0)
  437. std::cout<<INFO<<'\n';
  438. assert(INFO==0);
  439. mem::aligned_delete<T>(wsbuf);
  440. T eps_=tS[0]*eps;
  441. for(int i=0;i<k;i++)
  442. if(tS[i]<eps_)
  443. tS[i]=0;
  444. else
  445. tS[i]=1.0/tS[i];
  446. for(int i=0;i<m;i++){
  447. for(int j=0;j<k;j++){
  448. tU[i+j*m]*=tS[j];
  449. }
  450. }
  451. gemm<T>('T','T',n,m,k,1.0,&tVT[0],k,&tU[0],m,0.0,M_,n);
  452. mem::aligned_delete<T>(tU);
  453. mem::aligned_delete<T>(tS);
  454. mem::aligned_delete<T>(tVT);
  455. }
  456. }//end namespace
  457. }//end namespace