LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ sspgvd()

subroutine sspgvd ( integer  ITYPE,
character  JOBZ,
character  UPLO,
integer  N,
real, dimension( * )  AP,
real, dimension( * )  BP,
real, dimension( * )  W,
real, dimension( ldz, * )  Z,
integer  LDZ,
real, dimension( * )  WORK,
integer  LWORK,
integer, dimension( * )  IWORK,
integer  LIWORK,
integer  INFO 
)

SSPGVD

Download SSPGVD + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 SSPGVD computes all the eigenvalues, and optionally, the eigenvectors
 of a real generalized symmetric-definite eigenproblem, of the form
 A*x=(lambda)*B*x,  A*Bx=(lambda)*x,  or B*A*x=(lambda)*x.  Here A and
 B are assumed to be symmetric, stored in packed format, and B is also
 positive definite.
 If eigenvectors are desired, it uses a divide and conquer algorithm.

 The divide and conquer algorithm makes very mild assumptions about
 floating point arithmetic. It will work on machines with a guard
 digit in add/subtract, or on those binary machines without guard
 digits which subtract like the Cray X-MP, Cray Y-MP, Cray C-90, or
 Cray-2. It could conceivably fail on hexadecimal or decimal machines
 without guard digits, but we know of none.
Parameters
[in]ITYPE
          ITYPE is INTEGER
          Specifies the problem type to be solved:
          = 1:  A*x = (lambda)*B*x
          = 2:  A*B*x = (lambda)*x
          = 3:  B*A*x = (lambda)*x
[in]JOBZ
          JOBZ is CHARACTER*1
          = 'N':  Compute eigenvalues only;
          = 'V':  Compute eigenvalues and eigenvectors.
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangles of A and B are stored;
          = 'L':  Lower triangles of A and B are stored.
[in]N
          N is INTEGER
          The order of the matrices A and B.  N >= 0.
[in,out]AP
          AP is REAL array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the symmetric matrix
          A, packed columnwise in a linear array.  The j-th column of A
          is stored in the array AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2*n-j)/2) = A(i,j) for j<=i<=n.

          On exit, the contents of AP are destroyed.
[in,out]BP
          BP is REAL array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the symmetric matrix
          B, packed columnwise in a linear array.  The j-th column of B
          is stored in the array BP as follows:
          if UPLO = 'U', BP(i + (j-1)*j/2) = B(i,j) for 1<=i<=j;
          if UPLO = 'L', BP(i + (j-1)*(2*n-j)/2) = B(i,j) for j<=i<=n.

          On exit, the triangular factor U or L from the Cholesky
          factorization B = U**T*U or B = L*L**T, in the same storage
          format as B.
[out]W
          W is REAL array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is REAL array, dimension (LDZ, N)
          If JOBZ = 'V', then if INFO = 0, Z contains the matrix Z of
          eigenvectors.  The eigenvectors are normalized as follows:
          if ITYPE = 1 or 2, Z**T*B*Z = I;
          if ITYPE = 3, Z**T*inv(B)*Z = I.
          If JOBZ = 'N', then Z is not referenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= 1, and if
          JOBZ = 'V', LDZ >= max(1,N).
[out]WORK
          WORK is REAL array, dimension (MAX(1,LWORK))
          On exit, if INFO = 0, WORK(1) returns the required LWORK.
[in]LWORK
          LWORK is INTEGER
          The dimension of the array WORK.
          If N <= 1,               LWORK >= 1.
          If JOBZ = 'N' and N > 1, LWORK >= 2*N.
          If JOBZ = 'V' and N > 1, LWORK >= 1 + 6*N + 2*N**2.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the required sizes of the WORK and IWORK
          arrays, returns these values as the first entries of the WORK
          and IWORK arrays, and no error message related to LWORK or
          LIWORK is issued by XERBLA.
[out]IWORK
          IWORK is INTEGER array, dimension (MAX(1,LIWORK))
          On exit, if INFO = 0, IWORK(1) returns the required LIWORK.
[in]LIWORK
          LIWORK is INTEGER
          The dimension of the array IWORK.
          If JOBZ  = 'N' or N <= 1, LIWORK >= 1.
          If JOBZ  = 'V' and N > 1, LIWORK >= 3 + 5*N.

          If LIWORK = -1, then a workspace query is assumed; the
          routine only calculates the required sizes of the WORK and
          IWORK arrays, returns these values as the first entries of
          the WORK and IWORK arrays, and no error message related to
          LWORK or LIWORK is issued by XERBLA.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  SPPTRF or SSPEVD returned an error code:
             <= N:  if INFO = i, SSPEVD failed to converge;
                    i off-diagonal elements of an intermediate
                    tridiagonal form did not converge to zero;
             > N:   if INFO = N + i, for 1 <= i <= N, then the leading
                    minor of order i of B is not positive definite.
                    The factorization of B could not be completed and
                    no eigenvalues or eigenvectors were computed.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Contributors:
Mark Fahey, Department of Mathematics, Univ. of Kentucky, USA

Definition at line 208 of file sspgvd.f.

210*
211* -- LAPACK driver routine --
212* -- LAPACK is a software package provided by Univ. of Tennessee, --
213* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
214*
215* .. Scalar Arguments ..
216 CHARACTER JOBZ, UPLO
217 INTEGER INFO, ITYPE, LDZ, LIWORK, LWORK, N
218* ..
219* .. Array Arguments ..
220 INTEGER IWORK( * )
221 REAL AP( * ), BP( * ), W( * ), WORK( * ),
222 $ Z( LDZ, * )
223* ..
224*
225* =====================================================================
226*
227* .. Local Scalars ..
228 LOGICAL LQUERY, UPPER, WANTZ
229 CHARACTER TRANS
230 INTEGER J, LIWMIN, LWMIN, NEIG
231* ..
232* .. External Functions ..
233 LOGICAL LSAME
234 EXTERNAL lsame
235* ..
236* .. External Subroutines ..
237 EXTERNAL spptrf, sspevd, sspgst, stpmv, stpsv, xerbla
238* ..
239* .. Intrinsic Functions ..
240 INTRINSIC max, real
241* ..
242* .. Executable Statements ..
243*
244* Test the input parameters.
245*
246 wantz = lsame( jobz, 'V' )
247 upper = lsame( uplo, 'U' )
248 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
249*
250 info = 0
251 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
252 info = -1
253 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
254 info = -2
255 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
256 info = -3
257 ELSE IF( n.LT.0 ) THEN
258 info = -4
259 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
260 info = -9
261 END IF
262*
263 IF( info.EQ.0 ) THEN
264 IF( n.LE.1 ) THEN
265 liwmin = 1
266 lwmin = 1
267 ELSE
268 IF( wantz ) THEN
269 liwmin = 3 + 5*n
270 lwmin = 1 + 6*n + 2*n**2
271 ELSE
272 liwmin = 1
273 lwmin = 2*n
274 END IF
275 END IF
276 work( 1 ) = lwmin
277 iwork( 1 ) = liwmin
278 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
279 info = -11
280 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
281 info = -13
282 END IF
283 END IF
284*
285 IF( info.NE.0 ) THEN
286 CALL xerbla( 'SSPGVD', -info )
287 RETURN
288 ELSE IF( lquery ) THEN
289 RETURN
290 END IF
291*
292* Quick return if possible
293*
294 IF( n.EQ.0 )
295 $ RETURN
296*
297* Form a Cholesky factorization of BP.
298*
299 CALL spptrf( uplo, n, bp, info )
300 IF( info.NE.0 ) THEN
301 info = n + info
302 RETURN
303 END IF
304*
305* Transform problem to standard eigenvalue problem and solve.
306*
307 CALL sspgst( itype, uplo, n, ap, bp, info )
308 CALL sspevd( jobz, uplo, n, ap, w, z, ldz, work, lwork, iwork,
309 $ liwork, info )
310 lwmin = int( max( real( lwmin ), real( work( 1 ) ) ) )
311 liwmin = int( max( real( liwmin ), real( iwork( 1 ) ) ) )
312*
313 IF( wantz ) THEN
314*
315* Backtransform eigenvectors to the original problem.
316*
317 neig = n
318 IF( info.GT.0 )
319 $ neig = info - 1
320 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
321*
322* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
323* backtransform eigenvectors: x = inv(L)**T *y or inv(U)*y
324*
325 IF( upper ) THEN
326 trans = 'N'
327 ELSE
328 trans = 'T'
329 END IF
330*
331 DO 10 j = 1, neig
332 CALL stpsv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
333 $ 1 )
334 10 CONTINUE
335*
336 ELSE IF( itype.EQ.3 ) THEN
337*
338* For B*A*x=(lambda)*x;
339* backtransform eigenvectors: x = L*y or U**T *y
340*
341 IF( upper ) THEN
342 trans = 'T'
343 ELSE
344 trans = 'N'
345 END IF
346*
347 DO 20 j = 1, neig
348 CALL stpmv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
349 $ 1 )
350 20 CONTINUE
351 END IF
352 END IF
353*
354 work( 1 ) = lwmin
355 iwork( 1 ) = liwmin
356*
357 RETURN
358*
359* End of SSPGVD
360*
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:53
subroutine sspgst(ITYPE, UPLO, N, AP, BP, INFO)
SSPGST
Definition: sspgst.f:113
subroutine spptrf(UPLO, N, AP, INFO)
SPPTRF
Definition: spptrf.f:119
subroutine sspevd(JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, LWORK, IWORK, LIWORK, INFO)
SSPEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: sspevd.f:178
subroutine stpmv(UPLO, TRANS, DIAG, N, AP, X, INCX)
STPMV
Definition: stpmv.f:142
subroutine stpsv(UPLO, TRANS, DIAG, N, AP, X, INCX)
STPSV
Definition: stpsv.f:144
Here is the call graph for this function:
Here is the caller graph for this function: