File: node144.html

package info (click to toggle)
scalapack-doc 1.5-11
  • links: PTS
  • area: main
  • in suites: bullseye, buster, stretch
  • size: 10,336 kB
  • ctags: 4,931
  • sloc: makefile: 47; sh: 18
file content (269 lines) | stat: -rw-r--r-- 18,790 bytes parent folder | download | duplicates (4)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
<!DOCTYPE HTML PUBLIC "-//IETF//DTD HTML 2.0//EN">
<!--Converted with LaTeX2HTML 96.1-h (September 30, 1996) by Nikos Drakos (nikos@cbl.leeds.ac.uk), CBLU, University of Leeds -->
<HTML>
<HEAD>
<TITLE>Error Bounds for the Generalized Symmetric Definite
Eigenproblem</TITLE>
<META NAME="description" CONTENT="Error Bounds for the Generalized Symmetric Definite
Eigenproblem">
<META NAME="keywords" CONTENT="slug">
<META NAME="resource-type" CONTENT="document">
<META NAME="distribution" CONTENT="global">
<LINK REL=STYLESHEET HREF="slug.css">
</HEAD>
<BODY LANG="EN" >
 <A NAME="tex2html4007" HREF="node145.html"><IMG WIDTH=37 HEIGHT=24 ALIGN=BOTTOM ALT="next" SRC="http://www.netlib.org/utk/icons/next_motif.gif"></A> <A NAME="tex2html4005" HREF="node132.html"><IMG WIDTH=26 HEIGHT=24 ALIGN=BOTTOM ALT="up" SRC="http://www.netlib.org/utk/icons/up_motif.gif"></A> <A NAME="tex2html4001" HREF="node143.html"><IMG WIDTH=63 HEIGHT=24 ALIGN=BOTTOM ALT="previous" SRC="http://www.netlib.org/utk/icons/previous_motif.gif"></A> <A NAME="tex2html4009" HREF="node1.html"><IMG WIDTH=65 HEIGHT=24 ALIGN=BOTTOM ALT="contents" SRC="http://www.netlib.org/utk/icons/contents_motif.gif"></A> <A NAME="tex2html4010" HREF="node190.html"><IMG WIDTH=43 HEIGHT=24 ALIGN=BOTTOM ALT="index" SRC="http://www.netlib.org/utk/icons/index_motif.gif"></A> <BR>
<B> Next:</B> <A NAME="tex2html4008" HREF="node145.html">Troubleshooting</A>
<B>Up:</B> <A NAME="tex2html4006" HREF="node132.html">Accuracy and Stability</A>
<B> Previous:</B> <A NAME="tex2html4002" HREF="node143.html">Further Details:  Error Bounds </A>
<BR> <P>
<H1><A NAME="SECTION04690000000000000000">Error Bounds for the Generalized Symmetric Definite
Eigenproblem</A></H1>
<A NAME="secgendef">&#160;</A>
<P>
Three types of problems must be considered. 
<A NAME="6003">&#160;</A> 
In all cases <I>A</I> and <I>B</I>
are real symmetric (or complex Hermitian) and <I>B</I> is positive definite.
These decompositions are computed for real symmetric matrices
by the driver routines
PxSYGVX (see section <A HREF="node49.html#secGSEP">3.2.4</A>).
These decompositions are computed for complex Hermitian matrices
by the driver routines
PxHEGVX (see subsection <A HREF="node49.html#secGSEP">3.2.4</A>).
In each of the following three decompositions, 
<IMG WIDTH=12 HEIGHT=13 ALIGN=BOTTOM ALT="tex2html_wrap_inline12784" SRC="img74.gif"> is real and diagonal with diagonal entries
<IMG WIDTH=101 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline19079" SRC="img701.gif">, and
the columns <IMG WIDTH=11 HEIGHT=17 ALIGN=MIDDLE ALT="tex2html_wrap_inline18630" SRC="img638.gif"> of <I>Z</I> are linearly independent vectors.
The <IMG WIDTH=13 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline18626" SRC="img637.gif"> are called
<B>eigenvalues</B> and the <IMG WIDTH=11 HEIGHT=17 ALIGN=MIDDLE ALT="tex2html_wrap_inline18630" SRC="img638.gif"> are 
<B>eigenvectors</B>.<A NAME="6008">&#160;</A><A NAME="6009">&#160;</A>
<A NAME="6010">&#160;</A><A NAME="6011">&#160;</A>
<P>
<OL>
<LI> <IMG WIDTH=57 HEIGHT=22 ALIGN=MIDDLE ALT="tex2html_wrap_inline18087" SRC="img550.gif">.
The eigendecomposition may be written <IMG WIDTH=84 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline12895" SRC="img86.gif"> and
<IMG WIDTH=82 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline12903" SRC="img88.gif"> (or <IMG WIDTH=86 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline12913" SRC="img90.gif"> and <IMG WIDTH=84 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline12921" SRC="img92.gif">
if <I>A</I> and <I>B</I> are complex).
This may also be written <IMG WIDTH=91 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline19103" SRC="img702.gif">.
<LI> <IMG WIDTH=67 HEIGHT=22 ALIGN=MIDDLE ALT="tex2html_wrap_inline18089" SRC="img551.gif">.
The eigendecomposition may be written <IMG WIDTH=113 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline19107" SRC="img703.gif"> and
<IMG WIDTH=82 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline12903" SRC="img88.gif"> (<IMG WIDTH=115 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline19111" SRC="img704.gif"> and <IMG WIDTH=84 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline12921" SRC="img92.gif"> if <I>A</I>
and <I>B</I> are complex).
This may also be written <IMG WIDTH=91 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline19119" SRC="img705.gif">.
<LI> <IMG WIDTH=67 HEIGHT=22 ALIGN=MIDDLE ALT="tex2html_wrap_inline18091" SRC="img552.gif">.
The eigendecomposition may be written <IMG WIDTH=84 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline12895" SRC="img86.gif">
and <IMG WIDTH=100 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline12905" SRC="img89.gif"> (<IMG WIDTH=86 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline12913" SRC="img90.gif"> and <IMG WIDTH=102 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline12923" SRC="img93.gif"> if <I>A</I>
and <I>B</I> are complex).
This may also be written <IMG WIDTH=91 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline19135" SRC="img706.gif">.
<P>
</OL>
<P>
The approximate error bounds
for the computed eigenvalues <IMG WIDTH=101 HEIGHT=33 ALIGN=MIDDLE ALT="tex2html_wrap_inline18636" SRC="img640.gif">
are
<BR><IMG WIDTH=333 HEIGHT=22 ALIGN=BOTTOM ALT="displaymath19053" SRC="img707.gif"><BR>
The approximate error 
bounds<A NAME="6025">&#160;</A><A NAME="6026">&#160;</A>  
for the computed eigenvectors <IMG WIDTH=11 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline18638" SRC="img642.gif">,
which bound the acute angles between the computed eigenvectors and true
eigenvectors <IMG WIDTH=11 HEIGHT=17 ALIGN=MIDDLE ALT="tex2html_wrap_inline18630" SRC="img638.gif">, are
<A NAME="6028">&#160;</A>
<A NAME="6029">&#160;</A>
<BR><IMG WIDTH=331 HEIGHT=17 ALIGN=BOTTOM ALT="displaymath18597" SRC="img643.gif"><BR>
These bounds are computed differently, depending on which of the above three
problems are to be solved. The following code fragments show how.
<P>
<OL>
<LI>
First we consider error bounds for problem 1.
<P>
<PRE>      EPSMCH = PSLAMCH( ICTXT, 'E' )
      UNFL = PSLAMCH( ICTXT, 'U' )
*     Solve the eigenproblem A - lambda B (ITYPE = 1)
      ITYPE = 1
*     Compute the norms of A and B
      ANORM = PSLANSY( '1', UPLO, N, A, IA, JA, DESCA, WORK )
      BNORM = PSLANSY( '1', UPLO, N, B, IB, JB, DESCB, WORK )
*     The eigenvalues are returned in W
*     The eigenvectors are returned in A
      SUBROUTINE PSSYGVX( ITYPE, 'V', 'A', UPLO, N, A, IA, JA,
     $                    DESCA, B, IB, JB, DESCB, VL, VU, IL, IU,
     $                    UNFL, M, NZ, W, -1.0, Z, IZ, JZ, DESCZ,
     $                    WORK, LWORK, IWORK, LIWORK, IFAIL, ICLUSTR,
     $                    GAP, INFO )
      IF( INFO.GT.0 ) THEN
         PRINT *,'PSSYGVX did not converge, or B not positive definite'
      ELSE IF( N.GT.0 ) THEN
*        Get reciprocal condition number RCONDB of Cholesky factor of B
         CALL PSTRCON( '1', UPLO, 'N', N, B, IB, JB, DESCB, RCONDB,
     $                 WORK, LWORK, IWORK, LIWORK, 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</PRE>
<A NAME="6033">&#160;</A>
<P>
For example, if <IMG WIDTH=315 HEIGHT=30 ALIGN=MIDDLE ALT="tex2html_wrap_inline18242" SRC="img571.gif">,
<BR><IMG WIDTH=502 HEIGHT=67 ALIGN=BOTTOM ALT="displaymath19055" SRC="img708.gif"><BR>
then <IMG WIDTH=119 HEIGHT=13 ALIGN=BOTTOM ALT="tex2html_wrap_inline19145" SRC="img709.gif">,
<IMG WIDTH=94 HEIGHT=13 ALIGN=BOTTOM ALT="tex2html_wrap_inline19147" SRC="img710.gif">, and
<IMG WIDTH=117 HEIGHT=13 ALIGN=BOTTOM ALT="tex2html_wrap_inline19149" SRC="img711.gif">, and
the approximate eigenvalues, approximate error bounds,
and true errors are
shown below.
<P>
<BR><IMG WIDTH=533 HEIGHT=90 ALIGN=BOTTOM ALT="tabular6044" SRC="img712.gif"><BR>

<LI>
Problem types 2 and 3 have the same error bounds. We illustrate
only type 2. 
<A NAME="6064">&#160;</A>  
<A NAME="6065">&#160;</A>
<P>
<PRE>      EPSMCH = PSLAMCH( ICTXT, 'E' )
*     Solve the eigenproblem A*B - lambda I (ITYPE = 2)
      ITYPE = 2
*     Compute the norms of A and B
      ANORM = PSLANSY( '1', UPLO, N, A, IA, JA, DESCA, WORK )
      BNORM = PSLANSY( '1', UPLO, N, B, IB, JB, DESCB, WORK )
*     The eigenvalues are returned in W
*     The eigenvectors are returned in A
      SUBROUTINE PSSYGVX( ITYPE, 'V', 'A', UPLO, N, A, IA, JA,
     $                    DESCA, B, IB, JB, DESCB, VL, VU, IL, IU,
     $                    UNFL, M, NZ, W, -1.0, Z, IZ, JZ, DESCZ,
     $                    WORK, LWORK, IWORK, LIWORK, IFAIL, ICLUSTR,
     $                    GAP, INFO )
      IF( INFO.GT.0 .AND. INFO.LE.N ) THEN
         PRINT *,'PSSYGVX 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 PSTRCON( '1', UPLO, 'N', N, B, IB, JB, DESCB, RCONDB,
     $                 WORK, LWORK, IWORK, LIWORK, 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</PRE>
<A NAME="6066">&#160;</A>
<P>
For the same <I>A</I> and <I>B</I> as above, the approximate eigenvalues, 
approximate error bounds, and true errors are
shown below.
<P>
<BR><IMG WIDTH=571 HEIGHT=90 ALIGN=BOTTOM ALT="tabular6068" SRC="img713.gif"><BR>
</OL>
<P>
<B>Further Details:  Error Bounds for the Generalized Symmetric Definite
Eigenproblem</B><A NAME="secGSEPFurtherDetails">&#160;</A>
<P>
The error analysis of the driver routine PxSYGVX, or PxHEGVX in the complex case
<A NAME="6086">&#160;</A><A NAME="6087">&#160;</A><A NAME="6088">&#160;</A><A NAME="6089">&#160;</A>
(see section <A HREF="node49.html#secGSEP">3.2.4</A>),
goes as follows. 
In all cases 
<IMG WIDTH=171 HEIGHT=28 ALIGN=MIDDLE ALT="tex2html_wrap_inline18712" SRC="img656.gif"> is 
the <B>absolute gap</B><A NAME="6093">&#160;</A><A NAME="6094">&#160;</A>
between <IMG WIDTH=13 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline18626" SRC="img637.gif"> and the nearest other eigenvalue.
<A NAME="6095">&#160;</A>
<P>
<OL>
<LI> <IMG WIDTH=57 HEIGHT=22 ALIGN=MIDDLE ALT="tex2html_wrap_inline18087" SRC="img550.gif">. 
The computed eigenvalues <IMG WIDTH=13 HEIGHT=33 ALIGN=MIDDLE ALT="tex2html_wrap_inline19245" SRC="img714.gif"> can differ
from true eigenvalues <IMG WIDTH=13 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline18626" SRC="img637.gif"> by at most about
<BR><IMG WIDTH=478 HEIGHT=22 ALIGN=BOTTOM ALT="displaymath19056" SRC="img715.gif"><BR>
<P>
<A NAME="6102">&#160;</A>  
<A NAME="6103">&#160;</A>  
The angular difference between the computed eigenvector 
<IMG WIDTH=11 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline18638" SRC="img642.gif"> and a true eigenvector <IMG WIDTH=11 HEIGHT=17 ALIGN=MIDDLE ALT="tex2html_wrap_inline18630" SRC="img638.gif"> is
<BR><IMG WIDTH=498 HEIGHT=44 ALIGN=BOTTOM ALT="displaymath19057" SRC="img716.gif"><BR>
<LI> <IMG WIDTH=67 HEIGHT=22 ALIGN=MIDDLE ALT="tex2html_wrap_inline18089" SRC="img551.gif"> or <IMG WIDTH=67 HEIGHT=22 ALIGN=MIDDLE ALT="tex2html_wrap_inline18091" SRC="img552.gif">.
The computed eigenvalues <IMG WIDTH=13 HEIGHT=33 ALIGN=MIDDLE ALT="tex2html_wrap_inline19245" SRC="img714.gif"> can differ
from true eigenvalues <IMG WIDTH=13 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline18626" SRC="img637.gif"> by at most about
<BR><IMG WIDTH=469 HEIGHT=22 ALIGN=BOTTOM ALT="displaymath19058" SRC="img717.gif"><BR>
<P>
<A NAME="6115">&#160;</A>  
<A NAME="6116">&#160;</A>
<P>
The angular difference between the computed eigenvector 
<IMG WIDTH=11 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline18638" SRC="img642.gif"> and a true eigenvector <IMG WIDTH=11 HEIGHT=17 ALIGN=MIDDLE ALT="tex2html_wrap_inline18630" SRC="img638.gif"> is
<BR><IMG WIDTH=490 HEIGHT=48 ALIGN=BOTTOM ALT="displaymath19059" SRC="img718.gif"><BR></OL>
<P>
The code fragments above replace <I>p</I>(<I>n</I>) by 1 and make sure
neither <TT>RCONDB</TT> nor <TT>RCONDZ</TT> is so small as to cause
overflow when used as divisors in the expressions for error bounds.
<P>
These error bounds are large when <I>B</I> is ill-conditioned with respect to
inversion (<IMG WIDTH=42 HEIGHT=26 ALIGN=MIDDLE ALT="tex2html_wrap_inline19277" SRC="img719.gif"> is large). Often, the eigenvalues
and eigenvectors are much better conditioned than indicated here.
We mention two ways to get tighter bounds.
The first way is effective when the diagonal entries of <I>B</I> differ
widely in magnitude:<A NAME="tex2html1464" HREF="footnode.html#6125"><IMG  ALIGN=BOTTOM ALT="gif" SRC="http://www.netlib.org/utk/icons/foot_motif.gif"></A>
<OL>
<LI> <IMG WIDTH=57 HEIGHT=22 ALIGN=MIDDLE ALT="tex2html_wrap_inline18087" SRC="img550.gif">. Let <IMG WIDTH=199 HEIGHT=39 ALIGN=MIDDLE ALT="tex2html_wrap_inline19283" SRC="img720.gif"> 
be a diagonal matrix.
Then replace <I>B</I> by <I>DBD</I> and <I>A</I> by <I>DAD</I> in the above bounds.
<LI> <IMG WIDTH=67 HEIGHT=22 ALIGN=MIDDLE ALT="tex2html_wrap_inline18089" SRC="img551.gif"> or <IMG WIDTH=67 HEIGHT=22 ALIGN=MIDDLE ALT="tex2html_wrap_inline18091" SRC="img552.gif">. 
Let <IMG WIDTH=199 HEIGHT=39 ALIGN=MIDDLE ALT="tex2html_wrap_inline19283" SRC="img720.gif"> 
be a diagonal matrix.
Then replace <I>B</I> by <I>DBD</I> and <I>A</I> by <IMG WIDTH=77 HEIGHT=15 ALIGN=BOTTOM ALT="tex2html_wrap_inline19305" SRC="img721.gif"> in the above bounds. 
</OL>
<P>
The second way to get tighter bounds does not actually supply guaranteed 
bounds, but its estimates are often better in practice. 
It is not guaranteed because it assumes the algorithm is backward stable,
which is not necessarily true when <I>B</I> is ill-conditioned.
<A NAME="6138">&#160;</A>
<A NAME="6139">&#160;</A>
It estimates the <B>chordal distance</B> between a 
true eigenvalue <IMG WIDTH=13 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline18626" SRC="img637.gif"> and a computed eigenvalue <IMG WIDTH=13 HEIGHT=33 ALIGN=MIDDLE ALT="tex2html_wrap_inline19245" SRC="img714.gif">:
<A NAME="6142">&#160;</A> 
<BR><IMG WIDTH=367 HEIGHT=56 ALIGN=BOTTOM ALT="displaymath19060" SRC="img722.gif"><BR>
To interpret this measure, we write <IMG WIDTH=74 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline19313" SRC="img723.gif">
and <IMG WIDTH=74 HEIGHT=33 ALIGN=MIDDLE ALT="tex2html_wrap_inline19315" SRC="img724.gif">. Then
<IMG WIDTH=171 HEIGHT=34 ALIGN=MIDDLE ALT="tex2html_wrap_inline19317" SRC="img725.gif">.
In other words, if <IMG WIDTH=13 HEIGHT=33 ALIGN=MIDDLE ALT="tex2html_wrap_inline19245" SRC="img714.gif"> represents the one-dimensional subspace
<IMG WIDTH=11 HEIGHT=16 ALIGN=BOTTOM ALT="tex2html_wrap_inline17614" SRC="img477.gif"> consisting of the line through the origin with slope <IMG WIDTH=13 HEIGHT=33 ALIGN=MIDDLE ALT="tex2html_wrap_inline19245" SRC="img714.gif">,
and <IMG WIDTH=13 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline18626" SRC="img637.gif"> represents the analogous subspace <IMG WIDTH=11 HEIGHT=12 ALIGN=BOTTOM ALT="tex2html_wrap_inline17602" SRC="img475.gif">, then
<IMG WIDTH=61 HEIGHT=32 ALIGN=MIDDLE ALT="tex2html_wrap_inline19329" SRC="img726.gif">
is the sine of the acute angle <IMG WIDTH=51 HEIGHT=32 ALIGN=MIDDLE ALT="tex2html_wrap_inline17870" SRC="img519.gif"> between these 
subspaces. 
<A NAME="6157">&#160;</A>
<A NAME="6158">&#160;</A>
Thus, <IMG WIDTH=10 HEIGHT=17 ALIGN=MIDDLE ALT="tex2html_wrap_inline19333" SRC="img727.gif"> is bounded by one and is small when both arguments are 
large.<A NAME="tex2html1470" HREF="footnode.html#8108"><IMG  ALIGN=BOTTOM ALT="gif" SRC="http://www.netlib.org/utk/icons/foot_motif.gif"></A>
It applies only to the first problem, <IMG WIDTH=57 HEIGHT=22 ALIGN=MIDDLE ALT="tex2html_wrap_inline18087" SRC="img550.gif">:
<P>
<BLOCKQUOTE> Suppose a computed eigenvalue <IMG WIDTH=13 HEIGHT=33 ALIGN=MIDDLE ALT="tex2html_wrap_inline19245" SRC="img714.gif"> of <IMG WIDTH=57 HEIGHT=22 ALIGN=MIDDLE ALT="tex2html_wrap_inline18087" SRC="img550.gif"> is
the exact eigenvalue of a perturbed problem <IMG WIDTH=153 HEIGHT=26 ALIGN=MIDDLE ALT="tex2html_wrap_inline19345" SRC="img728.gif">.
Let <IMG WIDTH=14 HEIGHT=17 ALIGN=MIDDLE ALT="tex2html_wrap_inline17852" SRC="img518.gif"> be the unit eigenvector (<IMG WIDTH=68 HEIGHT=26 ALIGN=MIDDLE ALT="tex2html_wrap_inline19349" SRC="img729.gif">) for the exact 
eigenvalue <IMG WIDTH=13 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline18626" SRC="img637.gif">. 
Then if <IMG WIDTH=27 HEIGHT=26 ALIGN=MIDDLE ALT="tex2html_wrap_inline19353" SRC="img730.gif"> is small compared with 
<IMG WIDTH=26 HEIGHT=26 ALIGN=MIDDLE ALT="tex2html_wrap_inline19355" SRC="img731.gif">, and if <IMG WIDTH=27 HEIGHT=26 ALIGN=MIDDLE ALT="tex2html_wrap_inline19357" SRC="img732.gif"> is small compared with <IMG WIDTH=28 HEIGHT=26 ALIGN=MIDDLE ALT="tex2html_wrap_inline19359" SRC="img733.gif">, we have
<BR><IMG WIDTH=393 HEIGHT=52 ALIGN=BOTTOM ALT="displaymath19061" SRC="img734.gif"><BR>
Thus <IMG WIDTH=199 HEIGHT=41 ALIGN=MIDDLE ALT="tex2html_wrap_inline19363" SRC="img735.gif"> is a condition number for
eigenvalue <IMG WIDTH=13 HEIGHT=25 ALIGN=MIDDLE ALT="tex2html_wrap_inline18626" SRC="img637.gif">.
<A NAME="6166">&#160;</A>
</BLOCKQUOTE><HR><A NAME="tex2html4007" HREF="node145.html"><IMG WIDTH=37 HEIGHT=24 ALIGN=BOTTOM ALT="next" SRC="http://www.netlib.org/utk/icons/next_motif.gif"></A> <A NAME="tex2html4005" HREF="node132.html"><IMG WIDTH=26 HEIGHT=24 ALIGN=BOTTOM ALT="up" SRC="http://www.netlib.org/utk/icons/up_motif.gif"></A> <A NAME="tex2html4001" HREF="node143.html"><IMG WIDTH=63 HEIGHT=24 ALIGN=BOTTOM ALT="previous" SRC="http://www.netlib.org/utk/icons/previous_motif.gif"></A> <A NAME="tex2html4009" HREF="node1.html"><IMG WIDTH=65 HEIGHT=24 ALIGN=BOTTOM ALT="contents" SRC="http://www.netlib.org/utk/icons/contents_motif.gif"></A> <A NAME="tex2html4010" HREF="node190.html"><IMG WIDTH=43 HEIGHT=24 ALIGN=BOTTOM ALT="index" SRC="http://www.netlib.org/utk/icons/index_motif.gif"></A> <BR>
<B> Next:</B> <A NAME="tex2html4008" HREF="node145.html">Troubleshooting</A>
<B>Up:</B> <A NAME="tex2html4006" HREF="node132.html">Accuracy and Stability</A>
<B> Previous:</B> <A NAME="tex2html4002" HREF="node143.html">Further Details:  Error Bounds </A>
<P><ADDRESS>
<I>Susan Blackford <BR>
Tue May 13 09:21:01 EDT 1997</I>
</ADDRESS>
</BODY>
</HTML>