You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

137 lines
4.8 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2014 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #include "lapack_common.h"
  10. #include <Eigen/SVD>
  11. // computes the singular values/vectors a general M-by-N matrix A using divide-and-conquer
  12. EIGEN_LAPACK_FUNC(gesdd,(char *jobz, int *m, int* n, Scalar* a, int *lda, RealScalar *s, Scalar *u, int *ldu, Scalar *vt, int *ldvt, Scalar* /*work*/, int* lwork,
  13. EIGEN_LAPACK_ARG_IF_COMPLEX(RealScalar */*rwork*/) int * /*iwork*/, int *info))
  14. {
  15. // TODO exploit the work buffer
  16. bool query_size = *lwork==-1;
  17. int diag_size = (std::min)(*m,*n);
  18. *info = 0;
  19. if(*jobz!='A' && *jobz!='S' && *jobz!='O' && *jobz!='N') *info = -1;
  20. else if(*m<0) *info = -2;
  21. else if(*n<0) *info = -3;
  22. else if(*lda<std::max(1,*m)) *info = -5;
  23. else if(*lda<std::max(1,*m)) *info = -8;
  24. else if(*ldu <1 || (*jobz=='A' && *ldu <*m)
  25. || (*jobz=='O' && *m<*n && *ldu<*m)) *info = -8;
  26. else if(*ldvt<1 || (*jobz=='A' && *ldvt<*n)
  27. || (*jobz=='S' && *ldvt<diag_size)
  28. || (*jobz=='O' && *m>=*n && *ldvt<*n)) *info = -10;
  29. if(*info!=0)
  30. {
  31. int e = -*info;
  32. return xerbla_(SCALAR_SUFFIX_UP"GESDD ", &e, 6);
  33. }
  34. if(query_size)
  35. {
  36. *lwork = 0;
  37. return 0;
  38. }
  39. if(*n==0 || *m==0)
  40. return 0;
  41. PlainMatrixType mat(*m,*n);
  42. mat = matrix(a,*m,*n,*lda);
  43. int option = *jobz=='A' ? ComputeFullU|ComputeFullV
  44. : *jobz=='S' ? ComputeThinU|ComputeThinV
  45. : *jobz=='O' ? ComputeThinU|ComputeThinV
  46. : 0;
  47. BDCSVD<PlainMatrixType> svd(mat,option);
  48. make_vector(s,diag_size) = svd.singularValues().head(diag_size);
  49. if(*jobz=='A')
  50. {
  51. matrix(u,*m,*m,*ldu) = svd.matrixU();
  52. matrix(vt,*n,*n,*ldvt) = svd.matrixV().adjoint();
  53. }
  54. else if(*jobz=='S')
  55. {
  56. matrix(u,*m,diag_size,*ldu) = svd.matrixU();
  57. matrix(vt,diag_size,*n,*ldvt) = svd.matrixV().adjoint();
  58. }
  59. else if(*jobz=='O' && *m>=*n)
  60. {
  61. matrix(a,*m,*n,*lda) = svd.matrixU();
  62. matrix(vt,*n,*n,*ldvt) = svd.matrixV().adjoint();
  63. }
  64. else if(*jobz=='O')
  65. {
  66. matrix(u,*m,*m,*ldu) = svd.matrixU();
  67. matrix(a,diag_size,*n,*lda) = svd.matrixV().adjoint();
  68. }
  69. return 0;
  70. }
  71. // computes the singular values/vectors a general M-by-N matrix A using two sided jacobi algorithm
  72. EIGEN_LAPACK_FUNC(gesvd,(char *jobu, char *jobv, int *m, int* n, Scalar* a, int *lda, RealScalar *s, Scalar *u, int *ldu, Scalar *vt, int *ldvt, Scalar* /*work*/, int* lwork,
  73. EIGEN_LAPACK_ARG_IF_COMPLEX(RealScalar */*rwork*/) int *info))
  74. {
  75. // TODO exploit the work buffer
  76. bool query_size = *lwork==-1;
  77. int diag_size = (std::min)(*m,*n);
  78. *info = 0;
  79. if( *jobu!='A' && *jobu!='S' && *jobu!='O' && *jobu!='N') *info = -1;
  80. else if((*jobv!='A' && *jobv!='S' && *jobv!='O' && *jobv!='N')
  81. || (*jobu=='O' && *jobv=='O')) *info = -2;
  82. else if(*m<0) *info = -3;
  83. else if(*n<0) *info = -4;
  84. else if(*lda<std::max(1,*m)) *info = -6;
  85. else if(*ldu <1 || ((*jobu=='A' || *jobu=='S') && *ldu<*m)) *info = -9;
  86. else if(*ldvt<1 || (*jobv=='A' && *ldvt<*n)
  87. || (*jobv=='S' && *ldvt<diag_size)) *info = -11;
  88. if(*info!=0)
  89. {
  90. int e = -*info;
  91. return xerbla_(SCALAR_SUFFIX_UP"GESVD ", &e, 6);
  92. }
  93. if(query_size)
  94. {
  95. *lwork = 0;
  96. return 0;
  97. }
  98. if(*n==0 || *m==0)
  99. return 0;
  100. PlainMatrixType mat(*m,*n);
  101. mat = matrix(a,*m,*n,*lda);
  102. int option = (*jobu=='A' ? ComputeFullU : *jobu=='S' || *jobu=='O' ? ComputeThinU : 0)
  103. | (*jobv=='A' ? ComputeFullV : *jobv=='S' || *jobv=='O' ? ComputeThinV : 0);
  104. JacobiSVD<PlainMatrixType> svd(mat,option);
  105. make_vector(s,diag_size) = svd.singularValues().head(diag_size);
  106. if(*jobu=='A') matrix(u,*m,*m,*ldu) = svd.matrixU();
  107. else if(*jobu=='S') matrix(u,*m,diag_size,*ldu) = svd.matrixU();
  108. else if(*jobu=='O') matrix(a,*m,diag_size,*lda) = svd.matrixU();
  109. if(*jobv=='A') matrix(vt,*n,*n,*ldvt) = svd.matrixV().adjoint();
  110. else if(*jobv=='S') matrix(vt,diag_size,*n,*ldvt) = svd.matrixV().adjoint();
  111. else if(*jobv=='O') matrix(a,diag_size,*n,*lda) = svd.matrixV().adjoint();
  112. return 0;
  113. }