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.

72 lines
2.1 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2010-2011 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/Cholesky>
  11. // POTRF computes the Cholesky factorization of a real symmetric positive definite matrix A.
  12. EIGEN_LAPACK_FUNC(potrf,(char* uplo, int *n, RealScalar *pa, int *lda, int *info))
  13. {
  14. *info = 0;
  15. if(UPLO(*uplo)==INVALID) *info = -1;
  16. else if(*n<0) *info = -2;
  17. else if(*lda<std::max(1,*n)) *info = -4;
  18. if(*info!=0)
  19. {
  20. int e = -*info;
  21. return xerbla_(SCALAR_SUFFIX_UP"POTRF", &e, 6);
  22. }
  23. Scalar* a = reinterpret_cast<Scalar*>(pa);
  24. MatrixType A(a,*n,*n,*lda);
  25. int ret;
  26. if(UPLO(*uplo)==UP) ret = internal::llt_inplace<Scalar, Upper>::blocked(A);
  27. else ret = internal::llt_inplace<Scalar, Lower>::blocked(A);
  28. if(ret>=0)
  29. *info = ret+1;
  30. return 0;
  31. }
  32. // POTRS solves a system of linear equations A*X = B with a symmetric
  33. // positive definite matrix A using the Cholesky factorization
  34. // A = U**T*U or A = L*L**T computed by DPOTRF.
  35. EIGEN_LAPACK_FUNC(potrs,(char* uplo, int *n, int *nrhs, RealScalar *pa, int *lda, RealScalar *pb, int *ldb, int *info))
  36. {
  37. *info = 0;
  38. if(UPLO(*uplo)==INVALID) *info = -1;
  39. else if(*n<0) *info = -2;
  40. else if(*nrhs<0) *info = -3;
  41. else if(*lda<std::max(1,*n)) *info = -5;
  42. else if(*ldb<std::max(1,*n)) *info = -7;
  43. if(*info!=0)
  44. {
  45. int e = -*info;
  46. return xerbla_(SCALAR_SUFFIX_UP"POTRS", &e, 6);
  47. }
  48. Scalar* a = reinterpret_cast<Scalar*>(pa);
  49. Scalar* b = reinterpret_cast<Scalar*>(pb);
  50. MatrixType A(a,*n,*n,*lda);
  51. MatrixType B(b,*n,*nrhs,*ldb);
  52. if(UPLO(*uplo)==UP)
  53. {
  54. A.triangularView<Upper>().adjoint().solveInPlace(B);
  55. A.triangularView<Upper>().solveInPlace(B);
  56. }
  57. else
  58. {
  59. A.triangularView<Lower>().solveInPlace(B);
  60. A.triangularView<Lower>().adjoint().solveInPlace(B);
  61. }
  62. return 0;
  63. }