hypre/lapack/dsyev.c

292 lines
9.1 KiB
C

/*BHEADER**********************************************************************
* Copyright (c) 2006 The Regents of the University of California.
* Produced at the Lawrence Livermore National Laboratory.
* Written by the HYPRE team. UCRL-CODE-222953.
* All rights reserved.
*
* This file is part of HYPRE (see http://www.llnl.gov/CASC/hypre/).
* Please see the COPYRIGHT_and_LICENSE file for the copyright notice,
* disclaimer, contact information and the GNU Lesser General Public License.
*
* HYPRE is free software; you can redistribute it and/or modify it under the
* terms of the GNU General Public License (as published by the Free Software
* Foundation) version 2.1 dated February 1999.
*
* HYPRE is distributed in the hope that it will be useful, but WITHOUT ANY
* WARRANTY; without even the IMPLIED WARRANTY OF MERCHANTABILITY or FITNESS
* FOR A PARTICULAR PURPOSE. See the terms and conditions of the GNU General
* Public License for more details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with this program; if not, write to the Free Software Foundation,
* Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
*
* $Revision$
***********************************************************************EHEADER*/
#include "../blas/hypre_blas.h"
#include "hypre_lapack.h"
#include "f2c.h"
/* Subroutine */ int dsyev_(char *jobz, char *uplo, integer *n, doublereal *a,
integer *lda, doublereal *w, doublereal *work, integer *lwork,
integer *info)
{
/* -- LAPACK driver routine (version 3.0) --
Univ. of Tennessee, Univ. of California Berkeley, NAG Ltd.,
Courant Institute, Argonne National Lab, and Rice University
June 30, 1999
Purpose
=======
DSYEV computes all eigenvalues and, optionally, eigenvectors of a
real symmetric matrix A.
Arguments
=========
JOBZ (input) CHARACTER*1
= 'N': Compute eigenvalues only;
= 'V': Compute eigenvalues and eigenvectors.
UPLO (input) CHARACTER*1
= 'U': Upper triangle of A is stored;
= 'L': Lower triangle of A is stored.
N (input) INTEGER
The order of the matrix A. N >= 0.
A (input/output) DOUBLE PRECISION array, dimension (LDA, N)
On entry, the symmetric matrix A. If UPLO = 'U', the
leading N-by-N upper triangular part of A contains the
upper triangular part of the matrix A. If UPLO = 'L',
the leading N-by-N lower triangular part of A contains
the lower triangular part of the matrix A.
On exit, if JOBZ = 'V', then if INFO = 0, A contains the
orthonormal eigenvectors of the matrix A.
If JOBZ = 'N', then on exit the lower triangle (if UPLO='L')
or the upper triangle (if UPLO='U') of A, including the
diagonal, is destroyed.
LDA (input) INTEGER
The leading dimension of the array A. LDA >= max(1,N).
W (output) DOUBLE PRECISION array, dimension (N)
If INFO = 0, the eigenvalues in ascending order.
WORK (workspace/output) DOUBLE PRECISION array, dimension (LWORK)
On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
LWORK (input) INTEGER
The length of the array WORK. LWORK >= max(1,3*N-1).
For optimal efficiency, LWORK >= (NB+2)*N,
where NB is the blocksize for DSYTRD returned by ILAENV.
If LWORK = -1, then a workspace query is assumed; the routine
only calculates the optimal size of the WORK array, returns
this value as the first entry of the WORK array, and no error
message related to LWORK is issued by XERBLA.
INFO (output) INTEGER
= 0: successful exit
< 0: if INFO = -i, the i-th argument had an illegal value
> 0: if INFO = i, the algorithm failed to converge; i
off-diagonal elements of an intermediate tridiagonal
form did not converge to zero.
=====================================================================
Test the input parameters.
Parameter adjustments */
/* Table of constant values */
static integer c__1 = 1;
static integer c_n1 = -1;
static integer c__0 = 0;
static doublereal c_b17 = 1.;
/* System generated locals */
integer a_dim1, a_offset, i__1, i__2;
doublereal d__1;
/* Builtin functions */
double sqrt(doublereal);
/* Local variables */
static integer inde;
static doublereal anrm;
static integer imax;
static doublereal rmin, rmax;
/***static integer lopt;***/
extern /* Subroutine */ int dscal_(integer *, doublereal *, doublereal *,
integer *);
static doublereal sigma;
extern logical lsame_(char *, char *);
static integer iinfo;
static logical lower, wantz;
static integer nb;
extern doublereal dlamch_(char *);
static integer iscale;
extern /* Subroutine */ int dlascl_(char *, integer *, integer *,
doublereal *, doublereal *, integer *, integer *, doublereal *,
integer *, integer *);
static doublereal safmin;
extern integer ilaenv_(integer *, char *, char *, integer *, integer *,
integer *, integer *, ftnlen, ftnlen);
extern /* Subroutine */ int xerbla_(char *, integer *);
static doublereal bignum;
static integer indtau;
extern /* Subroutine */ int dsterf_(integer *, doublereal *, doublereal *,
integer *);
extern doublereal dlansy_(char *, char *, integer *, doublereal *,
integer *, doublereal *);
static integer indwrk;
extern /* Subroutine */ int dorgtr_(char *, integer *, doublereal *,
integer *, doublereal *, doublereal *, integer *, integer *), dsteqr_(char *, integer *, doublereal *, doublereal *,
doublereal *, integer *, doublereal *, integer *),
dsytrd_(char *, integer *, doublereal *, integer *, doublereal *,
doublereal *, doublereal *, doublereal *, integer *, integer *);
static integer llwork;
static doublereal smlnum;
static integer lwkopt;
static logical lquery;
static doublereal eps;
#define a_ref(a_1,a_2) a[(a_2)*a_dim1 + a_1]
a_dim1 = *lda;
a_offset = 1 + a_dim1 * 1;
a -= a_offset;
--w;
--work;
/* Function Body */
wantz = lsame_(jobz, "V");
lower = lsame_(uplo, "L");
lquery = *lwork == -1;
*info = 0;
if (! (wantz || lsame_(jobz, "N"))) {
*info = -1;
} else if (! (lower || lsame_(uplo, "U"))) {
*info = -2;
} else if (*n < 0) {
*info = -3;
} else if (*lda < max(1,*n)) {
*info = -5;
} else /* if(complicated condition) */ {
/* Computing MAX */
i__1 = 1, i__2 = *n * 3 - 1;
if (*lwork < max(i__1,i__2) && ! lquery) {
*info = -8;
}
}
if (*info == 0) {
nb = ilaenv_(&c__1, "DSYTRD", uplo, n, &c_n1, &c_n1, &c_n1, (ftnlen)6,
(ftnlen)1);
/* Computing MAX */
i__1 = 1, i__2 = (nb + 2) * *n;
lwkopt = max(i__1,i__2);
work[1] = (doublereal) lwkopt;
}
if (*info != 0) {
i__1 = -(*info);
xerbla_("DSYEV ", &i__1);
return 0;
} else if (lquery) {
return 0;
}
/* Quick return if possible */
if (*n == 0) {
work[1] = 1.;
return 0;
}
if (*n == 1) {
w[1] = a_ref(1, 1);
work[1] = 3.;
if (wantz) {
a_ref(1, 1) = 1.;
}
return 0;
}
/* Get machine constants. */
safmin = dlamch_("Safe minimum");
eps = dlamch_("Precision");
smlnum = safmin / eps;
bignum = 1. / smlnum;
rmin = sqrt(smlnum);
rmax = sqrt(bignum);
/* Scale matrix to allowable range, if necessary. */
anrm = dlansy_("M", uplo, n, &a[a_offset], lda, &work[1]);
iscale = 0;
if (anrm > 0. && anrm < rmin) {
iscale = 1;
sigma = rmin / anrm;
} else if (anrm > rmax) {
iscale = 1;
sigma = rmax / anrm;
}
if (iscale == 1) {
dlascl_(uplo, &c__0, &c__0, &c_b17, &sigma, n, n, &a[a_offset], lda,
info);
}
/* Call DSYTRD to reduce symmetric matrix to tridiagonal form. */
inde = 1;
indtau = inde + *n;
indwrk = indtau + *n;
llwork = *lwork - indwrk + 1;
dsytrd_(uplo, n, &a[a_offset], lda, &w[1], &work[inde], &work[indtau], &
work[indwrk], &llwork, &iinfo);
/***lopt = (integer) ((*n << 1) + work[indwrk]);***/
/* For eigenvalues only, call DSTERF. For eigenvectors, first call
DORGTR to generate the orthogonal matrix, then call DSTEQR. */
if (! wantz) {
dsterf_(n, &w[1], &work[inde], info);
} else {
dorgtr_(uplo, n, &a[a_offset], lda, &work[indtau], &work[indwrk], &
llwork, &iinfo);
dsteqr_(jobz, n, &w[1], &work[inde], &a[a_offset], lda, &work[indtau],
info);
}
/* If matrix was scaled, then rescale eigenvalues appropriately. */
if (iscale == 1) {
if (*info == 0) {
imax = *n;
} else {
imax = *info - 1;
}
d__1 = 1. / sigma;
dscal_(&imax, &d__1, &w[1], &c__1);
}
/* Set WORK(1) to optimal workspace size. */
work[1] = (doublereal) lwkopt;
return 0;
/* End of DSYEV */
} /* dsyev_ */
#undef a_ref