mat_utils.txx 12 KB

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