LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
dgtt01.f
Go to the documentation of this file.
1 *> \brief \b DGTT01
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 * Definition:
9 * ===========
10 *
11 * SUBROUTINE DGTT01( N, DL, D, DU, DLF, DF, DUF, DU2, IPIV, WORK,
12 * LDWORK, RWORK, RESID )
13 *
14 * .. Scalar Arguments ..
15 * INTEGER LDWORK, N
16 * DOUBLE PRECISION RESID
17 * ..
18 * .. Array Arguments ..
19 * INTEGER IPIV( * )
20 * DOUBLE PRECISION D( * ), DF( * ), DL( * ), DLF( * ), DU( * ),
21 * $ DU2( * ), DUF( * ), RWORK( * ),
22 * $ WORK( LDWORK, * )
23 * ..
24 *
25 *
26 *> \par Purpose:
27 * =============
28 *>
29 *> \verbatim
30 *>
31 *> DGTT01 reconstructs a tridiagonal matrix A from its LU factorization
32 *> and computes the residual
33 *> norm(L*U - A) / ( norm(A) * EPS ),
34 *> where EPS is the machine epsilon.
35 *> \endverbatim
36 *
37 * Arguments:
38 * ==========
39 *
40 *> \param[in] N
41 *> \verbatim
42 *> N is INTEGTER
43 *> The order of the matrix A. N >= 0.
44 *> \endverbatim
45 *>
46 *> \param[in] DL
47 *> \verbatim
48 *> DL is DOUBLE PRECISION array, dimension (N-1)
49 *> The (n-1) sub-diagonal elements of A.
50 *> \endverbatim
51 *>
52 *> \param[in] D
53 *> \verbatim
54 *> D is DOUBLE PRECISION array, dimension (N)
55 *> The diagonal elements of A.
56 *> \endverbatim
57 *>
58 *> \param[in] DU
59 *> \verbatim
60 *> DU is DOUBLE PRECISION array, dimension (N-1)
61 *> The (n-1) super-diagonal elements of A.
62 *> \endverbatim
63 *>
64 *> \param[in] DLF
65 *> \verbatim
66 *> DLF is DOUBLE PRECISION array, dimension (N-1)
67 *> The (n-1) multipliers that define the matrix L from the
68 *> LU factorization of A.
69 *> \endverbatim
70 *>
71 *> \param[in] DF
72 *> \verbatim
73 *> DF is DOUBLE PRECISION array, dimension (N)
74 *> The n diagonal elements of the upper triangular matrix U from
75 *> the LU factorization of A.
76 *> \endverbatim
77 *>
78 *> \param[in] DUF
79 *> \verbatim
80 *> DUF is DOUBLE PRECISION array, dimension (N-1)
81 *> The (n-1) elements of the first super-diagonal of U.
82 *> \endverbatim
83 *>
84 *> \param[in] DU2
85 *> \verbatim
86 *> DU2 is DOUBLE PRECISION array, dimension (N-2)
87 *> The (n-2) elements of the second super-diagonal of U.
88 *> \endverbatim
89 *>
90 *> \param[in] IPIV
91 *> \verbatim
92 *> IPIV is INTEGER array, dimension (N)
93 *> The pivot indices; for 1 <= i <= n, row i of the matrix was
94 *> interchanged with row IPIV(i). IPIV(i) will always be either
95 *> i or i+1; IPIV(i) = i indicates a row interchange was not
96 *> required.
97 *> \endverbatim
98 *>
99 *> \param[out] WORK
100 *> \verbatim
101 *> WORK is DOUBLE PRECISION array, dimension (LDWORK,N)
102 *> \endverbatim
103 *>
104 *> \param[in] LDWORK
105 *> \verbatim
106 *> LDWORK is INTEGER
107 *> The leading dimension of the array WORK. LDWORK >= max(1,N).
108 *> \endverbatim
109 *>
110 *> \param[out] RWORK
111 *> \verbatim
112 *> RWORK is DOUBLE PRECISION array, dimension (N)
113 *> \endverbatim
114 *>
115 *> \param[out] RESID
116 *> \verbatim
117 *> RESID is DOUBLE PRECISION
118 *> The scaled residual: norm(L*U - A) / (norm(A) * EPS)
119 *> \endverbatim
120 *
121 * Authors:
122 * ========
123 *
124 *> \author Univ. of Tennessee
125 *> \author Univ. of California Berkeley
126 *> \author Univ. of Colorado Denver
127 *> \author NAG Ltd.
128 *
129 *> \date November 2011
130 *
131 *> \ingroup double_lin
132 *
133 * =====================================================================
134  SUBROUTINE dgtt01( N, DL, D, DU, DLF, DF, DUF, DU2, IPIV, WORK,
135  $ ldwork, rwork, resid )
136 *
137 * -- LAPACK test routine (version 3.4.0) --
138 * -- LAPACK is a software package provided by Univ. of Tennessee, --
139 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
140 * November 2011
141 *
142 * .. Scalar Arguments ..
143  INTEGER LDWORK, N
144  DOUBLE PRECISION RESID
145 * ..
146 * .. Array Arguments ..
147  INTEGER IPIV( * )
148  DOUBLE PRECISION D( * ), DF( * ), DL( * ), DLF( * ), DU( * ),
149  $ du2( * ), duf( * ), rwork( * ),
150  $ work( ldwork, * )
151 * ..
152 *
153 * =====================================================================
154 *
155 * .. Parameters ..
156  DOUBLE PRECISION ONE, ZERO
157  parameter ( one = 1.0d+0, zero = 0.0d+0 )
158 * ..
159 * .. Local Scalars ..
160  INTEGER I, IP, J, LASTJ
161  DOUBLE PRECISION ANORM, EPS, LI
162 * ..
163 * .. External Functions ..
164  DOUBLE PRECISION DLAMCH, DLANGT, DLANHS
165  EXTERNAL dlamch, dlangt, dlanhs
166 * ..
167 * .. Intrinsic Functions ..
168  INTRINSIC min
169 * ..
170 * .. External Subroutines ..
171  EXTERNAL daxpy, dswap
172 * ..
173 * .. Executable Statements ..
174 *
175 * Quick return if possible
176 *
177  IF( n.LE.0 ) THEN
178  resid = zero
179  RETURN
180  END IF
181 *
182  eps = dlamch( 'Epsilon' )
183 *
184 * Copy the matrix U to WORK.
185 *
186  DO 20 j = 1, n
187  DO 10 i = 1, n
188  work( i, j ) = zero
189  10 CONTINUE
190  20 CONTINUE
191  DO 30 i = 1, n
192  IF( i.EQ.1 ) THEN
193  work( i, i ) = df( i )
194  IF( n.GE.2 )
195  $ work( i, i+1 ) = duf( i )
196  IF( n.GE.3 )
197  $ work( i, i+2 ) = du2( i )
198  ELSE IF( i.EQ.n ) THEN
199  work( i, i ) = df( i )
200  ELSE
201  work( i, i ) = df( i )
202  work( i, i+1 ) = duf( i )
203  IF( i.LT.n-1 )
204  $ work( i, i+2 ) = du2( i )
205  END IF
206  30 CONTINUE
207 *
208 * Multiply on the left by L.
209 *
210  lastj = n
211  DO 40 i = n - 1, 1, -1
212  li = dlf( i )
213  CALL daxpy( lastj-i+1, li, work( i, i ), ldwork,
214  $ work( i+1, i ), ldwork )
215  ip = ipiv( i )
216  IF( ip.EQ.i ) THEN
217  lastj = min( i+2, n )
218  ELSE
219  CALL dswap( lastj-i+1, work( i, i ), ldwork, work( i+1, i ),
220  $ ldwork )
221  END IF
222  40 CONTINUE
223 *
224 * Subtract the matrix A.
225 *
226  work( 1, 1 ) = work( 1, 1 ) - d( 1 )
227  IF( n.GT.1 ) THEN
228  work( 1, 2 ) = work( 1, 2 ) - du( 1 )
229  work( n, n-1 ) = work( n, n-1 ) - dl( n-1 )
230  work( n, n ) = work( n, n ) - d( n )
231  DO 50 i = 2, n - 1
232  work( i, i-1 ) = work( i, i-1 ) - dl( i-1 )
233  work( i, i ) = work( i, i ) - d( i )
234  work( i, i+1 ) = work( i, i+1 ) - du( i )
235  50 CONTINUE
236  END IF
237 *
238 * Compute the 1-norm of the tridiagonal matrix A.
239 *
240  anorm = dlangt( '1', n, dl, d, du )
241 *
242 * Compute the 1-norm of WORK, which is only guaranteed to be
243 * upper Hessenberg.
244 *
245  resid = dlanhs( '1', n, work, ldwork, rwork )
246 *
247 * Compute norm(L*U - A) / (norm(A) * EPS)
248 *
249  IF( anorm.LE.zero ) THEN
250  IF( resid.NE.zero )
251  $ resid = one / eps
252  ELSE
253  resid = ( resid / anorm ) / eps
254  END IF
255 *
256  RETURN
257 *
258 * End of DGTT01
259 *
260  END
subroutine daxpy(N, DA, DX, INCX, DY, INCY)
DAXPY
Definition: daxpy.f:54
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
subroutine dgtt01(N, DL, D, DU, DLF, DF, DUF, DU2, IPIV, WORK, LDWORK, RWORK, RESID)
DGTT01
Definition: dgtt01.f:136