blob: 921c51569708f5a2ca5dede14162c14ce3caf2fd [file] [log] [blame]
Brian Silverman72890c22015-09-19 14:37:37 -04001// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 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
Austin Schuh189376f2018-12-20 22:11:15 +110010#include "lapack_common.h"
Brian Silverman72890c22015-09-19 14:37:37 -040011#include <Eigen/Eigenvalues>
12
Austin Schuh189376f2018-12-20 22:11:15 +110013// computes eigen values and vectors of a general N-by-N matrix A
Brian Silverman72890c22015-09-19 14:37:37 -040014EIGEN_LAPACK_FUNC(syev,(char *jobz, char *uplo, int* n, Scalar* a, int *lda, Scalar* w, Scalar* /*work*/, int* lwork, int *info))
15{
16 // TODO exploit the work buffer
17 bool query_size = *lwork==-1;
18
19 *info = 0;
20 if(*jobz!='N' && *jobz!='V') *info = -1;
21 else if(UPLO(*uplo)==INVALID) *info = -2;
22 else if(*n<0) *info = -3;
23 else if(*lda<std::max(1,*n)) *info = -5;
24 else if((!query_size) && *lwork<std::max(1,3**n-1)) *info = -8;
Austin Schuh189376f2018-12-20 22:11:15 +110025
Brian Silverman72890c22015-09-19 14:37:37 -040026 if(*info!=0)
27 {
28 int e = -*info;
29 return xerbla_(SCALAR_SUFFIX_UP"SYEV ", &e, 6);
30 }
31
32 if(query_size)
33 {
34 *lwork = 0;
35 return 0;
36 }
37
38 if(*n==0)
39 return 0;
40
41 PlainMatrixType mat(*n,*n);
42 if(UPLO(*uplo)==UP) mat = matrix(a,*n,*n,*lda).adjoint();
43 else mat = matrix(a,*n,*n,*lda);
44
45 bool computeVectors = *jobz=='V' || *jobz=='v';
46 SelfAdjointEigenSolver<PlainMatrixType> eig(mat,computeVectors?ComputeEigenvectors:EigenvaluesOnly);
47
48 if(eig.info()==NoConvergence)
49 {
Austin Schuh189376f2018-12-20 22:11:15 +110050 make_vector(w,*n).setZero();
Brian Silverman72890c22015-09-19 14:37:37 -040051 if(computeVectors)
52 matrix(a,*n,*n,*lda).setIdentity();
53 //*info = 1;
54 return 0;
55 }
56
Austin Schuh189376f2018-12-20 22:11:15 +110057 make_vector(w,*n) = eig.eigenvalues();
Brian Silverman72890c22015-09-19 14:37:37 -040058 if(computeVectors)
59 matrix(a,*n,*n,*lda) = eig.eigenvectors();
60
61 return 0;
62}