Error Bounds for the Generalized Symmetric Definite Eigenproblem



next up previous contents index
Next: Further Details: Error Up: Accuracy and Stability Previous: Further Details: Error

Error Bounds for the Generalized Symmetric Definite Eigenproblem

 

There are three types of problems to consider.   In all cases A and B are real symmetric (or complex Hermitian) and B is positive definite. These decompositions are computed for real symmetric matrices by the driver routines xSYGV, xSPGV and (for type 1 only) xSBGV (see subsection 2.2.5.1). These decompositions are computed for complex Hermitian matrices by the driver routines xHEGV, xHPGV and (for type 1 only) xHBGV (see subsection 2.2.5.1). In each of the following three decompositions, is real and diagonal with diagonal entries , and the columns of Z are linearly independent vectors. The are called eigenvalues and the are eigenvectors.               

  1. . The eigendecomposition may be written and (or and if A and B are complex). This may also be written .

  2. . The eigendecomposition may be written and ( and if A and B are complex). This may also be written .

  3. . The eigendecomposition may be written and ( and if A and B are complex). This may also be written .

The approximate error bounds for the computed eigenvalues are

The approximate error bounds   for the computed eigenvectors , which bound the acute angles between the computed eigenvectors and true eigenvectors , are    

These bounds are computed differently, depending on which of the above three problems are to be solved. The following code fragments show how.

  1. First we consider error bounds for problem 1.

       EPSMCH = SLAMCH( 'E' )
    *  Solve the eigenproblem A - lambda B (ITYPE = 1)
       ITYPE = 1
    *  Compute the norms of A and B
       ANORM = SLANSY( '1', UPLO, N, A, LDA, WORK )
       BNORM = SLANSY( '1', UPLO, N, B, LDB, WORK )
    *  The eigenvalues are returned in W
    *  The eigenvectors are returned in A
       CALL SSYGV( ITYPE, 'V', UPLO, N, A, LDA, B, LDB, W,
      $            WORK, LWORK, INFO )
       IF( INFO.GT.0 .AND. INFO.LE.N ) THEN
             PRINT *,'SSYGV did not converge'
       ELSE IF( INFO.GT.N ) THEN
          PRINT *,'B not positive definite'
       ELSE IF ( N.GT.0 ) THEN
    *     Get reciprocal condition number RCONDB of Cholesky
    *      factor of B
          CALL STRCON( '1', UPLO, 'N', N, B, LDB, RCONDB,
         $             WORK, IWORK, INFO )
          RCONDB = MAX( RCONDB, EPSMCH )
          CALL SDISNA( 'Eigenvectors', N, N, W, RCONDZ, INFO )
          DO 10 I = 1, N
          EERRBD( I ) = ( EPSMCH / RCONDB**2 ) *
         $              ( ANORM / BNORM + ABS( W(I) ) )
             ZERRBD( I ) = ( EPSMCH / RCONDB**3 ) *
    $                   ( ( ANORM / BNORM ) / RCONDZ(I) +
    $                   ( ABS( W(I) ) / RCONDZ(I) ) * RCONDB )
    10    CONTINUE
       END IF
    
     

    For example, if ,

    then ANORM = 120231, BNORM = 120, and RCOND = .8326, and the approximate eigenvalues, approximate error bounds, and true errors are

    This code fragment cannot be adapted to use xSBGV (or xHBGV), because xSBGV does not return a conventional Cholesky factor in B, but rather a ``split'' Choleksy factorization (performed by xPBSTF). A future LAPACK release will include error bounds for xSBGV.

  2. Problem types 2 and 3 have the same error bounds. We illustrate only type 2.    

       EPSMCH = SLAMCH( 'E' )
    *  Solve the eigenproblem A*B - lambda I (ITYPE = 2)
       ITYPE = 2
    *  Compute the norms of A and B
       ANORM = SLANSY( '1', UPLO, N, A, LDA, WORK )
       BNORM = SLANSY( '1', UPLO, N, B, LDB, WORK )
    *  The eigenvalues are returned in W
    *  The eigenvectors are returned in A
       CALL SSYGV( ITYPE, 'V', UPLO, N, A, LDA, B, LDB, W,
    $              WORK, LWORK, INFO )
       IF( INFO.GT.0 .AND. INFO.LE.N ) THEN
          PRINT *,'SSYGV did not converge'
       ELSE IF( INFO.GT.N ) THEN
          PRINT *,'B not positive definite'
       ELSE IF ( N.GT.0 ) THEN
    *     Get reciprocal condition number RCONDB of Cholesky
    *      factor of B
          CALL STRCON( '1', UPLO, 'N', N, B, LDB, RCONDB,
    $                  WORK, IWORK, INFO )
          RCONDB = MAX( RCONDB, EPSMCH )
          CALL SDISNA( 'Eigenvectors', N, N, W, RCONDZ, INFO )
          DO 10 I = 1, N
             EERRBD(I) = ( ANORM * BNORM ) * EPSMCH + 
      $                  ( EPSMCH / RCONDB**2 ) * ABS( W(I) )
             ZERRBD(I) = ( EPSMCH / RCONDB ) * ( ( ANORM *
      $                    BNORM ) / RCONDZ(I) + 1.0 / RCONDB )
    10       CONTINUE
          END IF
    
     

    For the same A and B as above, the approximate eigenvalues, approximate error bounds, and true errors are





next up previous contents index
Next: Further Details: Error Up: Accuracy and Stability Previous: Further Details: Error




Tue Nov 29 14:03:33 EST 1994