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

◆ dspgvd()

subroutine dspgvd ( integer  itype,
character  jobz,
character  uplo,
integer  n,
double precision, dimension( * )  ap,
double precision, dimension( * )  bp,
double precision, dimension( * )  w,
double precision, dimension( ldz, * )  z,
integer  ldz,
double precision, dimension( * )  work,
integer  lwork,
integer, dimension( * )  iwork,
integer  liwork,
integer  info 
)

DSPGVD

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

Purpose:
 DSPGVD 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 DOUBLE PRECISION 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 DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (N)
          If INFO = 0, the eigenvalues in ascending order.
[out]Z
          Z is DOUBLE PRECISION 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 DOUBLE PRECISION 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:  DPPTRF or DSPEVD returned an error code:
             <= N:  if INFO = i, DSPEVD 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 dspgvd.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 DOUBLE PRECISION 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 EXTERNAL lsame
229* ..
230* .. External Subroutines ..
231 EXTERNAL dpptrf, dspevd, dspgst, dtpmv, dtpsv, xerbla
232* ..
233* .. Intrinsic Functions ..
234 INTRINSIC dble, max
235* ..
236* .. Executable Statements ..
237*
238* Test the input parameters.
239*
240 wantz = lsame( jobz, 'V' )
241 upper = lsame( uplo, 'U' )
242 lquery = ( lwork.EQ.-1 .OR. liwork.EQ.-1 )
243*
244 info = 0
245 IF( itype.LT.1 .OR. itype.GT.3 ) THEN
246 info = -1
247 ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
248 info = -2
249 ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
250 info = -3
251 ELSE IF( n.LT.0 ) THEN
252 info = -4
253 ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
254 info = -9
255 END IF
256*
257 IF( info.EQ.0 ) THEN
258 IF( n.LE.1 ) THEN
259 liwmin = 1
260 lwmin = 1
261 ELSE
262 IF( wantz ) THEN
263 liwmin = 3 + 5*n
264 lwmin = 1 + 6*n + 2*n**2
265 ELSE
266 liwmin = 1
267 lwmin = 2*n
268 END IF
269 END IF
270 work( 1 ) = lwmin
271 iwork( 1 ) = liwmin
272 IF( lwork.LT.lwmin .AND. .NOT.lquery ) THEN
273 info = -11
274 ELSE IF( liwork.LT.liwmin .AND. .NOT.lquery ) THEN
275 info = -13
276 END IF
277 END IF
278*
279 IF( info.NE.0 ) THEN
280 CALL xerbla( 'DSPGVD', -info )
281 RETURN
282 ELSE IF( lquery ) THEN
283 RETURN
284 END IF
285*
286* Quick return if possible
287*
288 IF( n.EQ.0 )
289 $ RETURN
290*
291* Form a Cholesky factorization of BP.
292*
293 CALL dpptrf( uplo, n, bp, info )
294 IF( info.NE.0 ) THEN
295 info = n + info
296 RETURN
297 END IF
298*
299* Transform problem to standard eigenvalue problem and solve.
300*
301 CALL dspgst( itype, uplo, n, ap, bp, info )
302 CALL dspevd( jobz, uplo, n, ap, w, z, ldz, work, lwork, iwork,
303 $ liwork, info )
304 lwmin = int( max( dble( lwmin ), dble( work( 1 ) ) ) )
305 liwmin = int( max( dble( liwmin ), dble( iwork( 1 ) ) ) )
306*
307 IF( wantz ) THEN
308*
309* Backtransform eigenvectors to the original problem.
310*
311 neig = n
312 IF( info.GT.0 )
313 $ neig = info - 1
314 IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
315*
316* For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
317* backtransform eigenvectors: x = inv(L)**T *y or inv(U)*y
318*
319 IF( upper ) THEN
320 trans = 'N'
321 ELSE
322 trans = 'T'
323 END IF
324*
325 DO 10 j = 1, neig
326 CALL dtpsv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
327 $ 1 )
328 10 CONTINUE
329*
330 ELSE IF( itype.EQ.3 ) THEN
331*
332* For B*A*x=(lambda)*x;
333* backtransform eigenvectors: x = L*y or U**T *y
334*
335 IF( upper ) THEN
336 trans = 'T'
337 ELSE
338 trans = 'N'
339 END IF
340*
341 DO 20 j = 1, neig
342 CALL dtpmv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
343 $ 1 )
344 20 CONTINUE
345 END IF
346 END IF
347*
348 work( 1 ) = lwmin
349 iwork( 1 ) = liwmin
350*
351 RETURN
352*
353* End of DSPGVD
354*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dspevd(jobz, uplo, n, ap, w, z, ldz, work, lwork, iwork, liwork, info)
DSPEVD computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition dspevd.f:172
subroutine dspgst(itype, uplo, n, ap, bp, info)
DSPGST
Definition dspgst.f:113
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dpptrf(uplo, n, ap, info)
DPPTRF
Definition dpptrf.f:119
subroutine dtpmv(uplo, trans, diag, n, ap, x, incx)
DTPMV
Definition dtpmv.f:142
subroutine dtpsv(uplo, trans, diag, n, ap, x, incx)
DTPSV
Definition dtpsv.f:144
Here is the call graph for this function:
Here is the caller graph for this function: