mat_utils.txx 14 KB

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