LAPACK 3.12.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.
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
                    principal minor of order i of B is not positive.
                    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 202 of file sspgvd.f.

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