LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dspgv ( 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  INFO 
)

DSPGV

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

Purpose:
 DSPGV 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.
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 (3*N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
          > 0:  DPPTRF or DSPEV returned an error code:
             <= N:  if INFO = i, DSPEV 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.
Date
November 2015

Definition at line 163 of file dspgv.f.

163 *
164 * -- LAPACK driver routine (version 3.6.0) --
165 * -- LAPACK is a software package provided by Univ. of Tennessee, --
166 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
167 * November 2015
168 *
169 * .. Scalar Arguments ..
170  CHARACTER jobz, uplo
171  INTEGER info, itype, ldz, n
172 * ..
173 * .. Array Arguments ..
174  DOUBLE PRECISION ap( * ), bp( * ), w( * ), work( * ),
175  $ z( ldz, * )
176 * ..
177 *
178 * =====================================================================
179 *
180 * .. Local Scalars ..
181  LOGICAL upper, wantz
182  CHARACTER trans
183  INTEGER j, neig
184 * ..
185 * .. External Functions ..
186  LOGICAL lsame
187  EXTERNAL lsame
188 * ..
189 * .. External Subroutines ..
190  EXTERNAL dpptrf, dspev, dspgst, dtpmv, dtpsv, xerbla
191 * ..
192 * .. Executable Statements ..
193 *
194 * Test the input parameters.
195 *
196  wantz = lsame( jobz, 'V' )
197  upper = lsame( uplo, 'U' )
198 *
199  info = 0
200  IF( itype.LT.1 .OR. itype.GT.3 ) THEN
201  info = -1
202  ELSE IF( .NOT.( wantz .OR. lsame( jobz, 'N' ) ) ) THEN
203  info = -2
204  ELSE IF( .NOT.( upper .OR. lsame( uplo, 'L' ) ) ) THEN
205  info = -3
206  ELSE IF( n.LT.0 ) THEN
207  info = -4
208  ELSE IF( ldz.LT.1 .OR. ( wantz .AND. ldz.LT.n ) ) THEN
209  info = -9
210  END IF
211  IF( info.NE.0 ) THEN
212  CALL xerbla( 'DSPGV ', -info )
213  RETURN
214  END IF
215 *
216 * Quick return if possible
217 *
218  IF( n.EQ.0 )
219  $ RETURN
220 *
221 * Form a Cholesky factorization of B.
222 *
223  CALL dpptrf( uplo, n, bp, info )
224  IF( info.NE.0 ) THEN
225  info = n + info
226  RETURN
227  END IF
228 *
229 * Transform problem to standard eigenvalue problem and solve.
230 *
231  CALL dspgst( itype, uplo, n, ap, bp, info )
232  CALL dspev( jobz, uplo, n, ap, w, z, ldz, work, info )
233 *
234  IF( wantz ) THEN
235 *
236 * Backtransform eigenvectors to the original problem.
237 *
238  neig = n
239  IF( info.GT.0 )
240  $ neig = info - 1
241  IF( itype.EQ.1 .OR. itype.EQ.2 ) THEN
242 *
243 * For A*x=(lambda)*B*x and A*B*x=(lambda)*x;
244 * backtransform eigenvectors: x = inv(L)**T*y or inv(U)*y
245 *
246  IF( upper ) THEN
247  trans = 'N'
248  ELSE
249  trans = 'T'
250  END IF
251 *
252  DO 10 j = 1, neig
253  CALL dtpsv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
254  $ 1 )
255  10 CONTINUE
256 *
257  ELSE IF( itype.EQ.3 ) THEN
258 *
259 * For B*A*x=(lambda)*x;
260 * backtransform eigenvectors: x = L*y or U**T*y
261 *
262  IF( upper ) THEN
263  trans = 'T'
264  ELSE
265  trans = 'N'
266  END IF
267 *
268  DO 20 j = 1, neig
269  CALL dtpmv( uplo, trans, 'Non-unit', n, bp, z( 1, j ),
270  $ 1 )
271  20 CONTINUE
272  END IF
273  END IF
274  RETURN
275 *
276 * End of DSPGV
277 *
subroutine dspgst(ITYPE, UPLO, N, AP, BP, INFO)
DSPGST
Definition: dspgst.f:115
subroutine dspev(JOBZ, UPLO, N, AP, W, Z, LDZ, WORK, INFO)
DSPEV computes the eigenvalues and, optionally, the left and/or right eigenvectors for OTHER matrice...
Definition: dspev.f:132
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dpptrf(UPLO, N, AP, INFO)
DPPTRF
Definition: dpptrf.f:121
subroutine dtpsv(UPLO, TRANS, DIAG, N, AP, X, INCX)
DTPSV
Definition: dtpsv.f:146
subroutine dtpmv(UPLO, TRANS, DIAG, N, AP, X, INCX)
DTPMV
Definition: dtpmv.f:144
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: