LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ssyt01_aa.f
Go to the documentation of this file.
1*> \brief \b SSYT01_AA
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 SSYT01_AA( UPLO, N, A, LDA, AFAC, LDAFAC, IPIV,
12* C, LDC, RWORK, RESID )
13*
14* .. Scalar Arguments ..
15* CHARACTER UPLO
16* INTEGER LDA, LDAFAC, LDC, N
17* REAL RESID
18* ..
19* .. Array Arguments ..
20* INTEGER IPIV( * )
21* REAL A( LDA, * ), AFAC( LDAFAC, * ), C( LDC, * ),
22* $ RWORK( * )
23* ..
24*
25*
26*> \par Purpose:
27* =============
28*>
29*> \verbatim
30*>
31*> SSYT01_AA reconstructs a symmetric indefinite matrix A from its
32*> block L*D*L' or U*D*U' factorization and computes the residual
33*> norm( C - A ) / ( N * norm(A) * EPS ),
34*> where C is the reconstructed matrix and EPS is the machine epsilon.
35*> \endverbatim
36*
37* Arguments:
38* ==========
39*
40*> \param[in] UPLO
41*> \verbatim
42*> UPLO is CHARACTER*1
43*> Specifies whether the upper or lower triangular part of the
44*> symmetric matrix A is stored:
45*> = 'U': Upper triangular
46*> = 'L': Lower triangular
47*> \endverbatim
48*>
49*> \param[in] N
50*> \verbatim
51*> N is INTEGER
52*> The number of rows and columns of the matrix A. N >= 0.
53*> \endverbatim
54*>
55*> \param[in] A
56*> \verbatim
57*> A is REAL array, dimension (LDA,N)
58*> The original symmetric matrix A.
59*> \endverbatim
60*>
61*> \param[in] LDA
62*> \verbatim
63*> LDA is INTEGER
64*> The leading dimension of the array A. LDA >= max(1,N)
65*> \endverbatim
66*>
67*> \param[in] AFAC
68*> \verbatim
69*> AFAC is REAL array, dimension (LDAFAC,N)
70*> The factored form of the matrix A. AFAC contains the block
71*> diagonal matrix D and the multipliers used to obtain the
72*> factor L or U from the block L*D*L' or U*D*U' factorization
73*> as computed by SSYTRF.
74*> \endverbatim
75*>
76*> \param[in] LDAFAC
77*> \verbatim
78*> LDAFAC is INTEGER
79*> The leading dimension of the array AFAC. LDAFAC >= max(1,N).
80*> \endverbatim
81*>
82*> \param[in] IPIV
83*> \verbatim
84*> IPIV is INTEGER array, dimension (N)
85*> The pivot indices from SSYTRF.
86*> \endverbatim
87*>
88*> \param[out] C
89*> \verbatim
90*> C is REAL array, dimension (LDC,N)
91*> \endverbatim
92*>
93*> \param[in] LDC
94*> \verbatim
95*> LDC is INTEGER
96*> The leading dimension of the array C. LDC >= max(1,N).
97*> \endverbatim
98*>
99*> \param[out] RWORK
100*> \verbatim
101*> RWORK is REAL array, dimension (N)
102*> \endverbatim
103*>
104*> \param[out] RESID
105*> \verbatim
106*> RESID is REAL
107*> If UPLO = 'L', norm(L*D*L' - A) / ( N * norm(A) * EPS )
108*> If UPLO = 'U', norm(U*D*U' - A) / ( N * norm(A) * EPS )
109*> \endverbatim
110*
111* Authors:
112* ========
113*
114*> \author Univ. of Tennessee
115*> \author Univ. of California Berkeley
116*> \author Univ. of Colorado Denver
117*> \author NAG Ltd.
118*
119*> \ingroup real_lin
120*
121* =====================================================================
122 SUBROUTINE ssyt01_aa( UPLO, N, A, LDA, AFAC, LDAFAC, IPIV, C,
123 $ LDC, RWORK, RESID )
124*
125* -- LAPACK test routine --
126* -- LAPACK is a software package provided by Univ. of Tennessee, --
127* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
128*
129* .. Scalar Arguments ..
130 CHARACTER UPLO
131 INTEGER LDA, LDAFAC, LDC, N
132 REAL RESID
133* ..
134* .. Array Arguments ..
135 INTEGER IPIV( * )
136 REAL A( LDA, * ), AFAC( LDAFAC, * ), C( LDC, * ),
137 $ rwork( * )
138* ..
139*
140* =====================================================================
141*
142* .. Parameters ..
143 REAL ZERO, ONE
144 parameter( zero = 0.0e+0, one = 1.0e+0 )
145* ..
146* .. Local Scalars ..
147 INTEGER I, J
148 REAL ANORM, EPS
149* ..
150* .. External Functions ..
151 LOGICAL LSAME
152 REAL SLAMCH, SLANSY
153 EXTERNAL lsame, slamch, slansy
154* ..
155* .. External Subroutines ..
156 EXTERNAL slaset, slavsy, sswap, strmm, slacpy
157* ..
158* .. Intrinsic Functions ..
159 INTRINSIC real
160* ..
161* .. Executable Statements ..
162*
163* Quick exit if N = 0.
164*
165 IF( n.LE.0 ) THEN
166 resid = zero
167 RETURN
168 END IF
169*
170* Determine EPS and the norm of A.
171*
172 eps = slamch( 'Epsilon' )
173 anorm = slansy( '1', uplo, n, a, lda, rwork )
174*
175* Initialize C to the tridiagonal matrix T.
176*
177 CALL slaset( 'Full', n, n, zero, zero, c, ldc )
178 CALL slacpy( 'F', 1, n, afac( 1, 1 ), ldafac+1, c( 1, 1 ), ldc+1 )
179 IF( n.GT.1 ) THEN
180 IF( lsame( uplo, 'U' ) ) THEN
181 CALL slacpy( 'F', 1, n-1, afac( 1, 2 ), ldafac+1, c( 1, 2 ),
182 $ ldc+1 )
183 CALL slacpy( 'F', 1, n-1, afac( 1, 2 ), ldafac+1, c( 2, 1 ),
184 $ ldc+1 )
185 ELSE
186 CALL slacpy( 'F', 1, n-1, afac( 2, 1 ), ldafac+1, c( 1, 2 ),
187 $ ldc+1 )
188 CALL slacpy( 'F', 1, n-1, afac( 2, 1 ), ldafac+1, c( 2, 1 ),
189 $ ldc+1 )
190 ENDIF
191*
192* Call STRMM to form the product U' * D (or L * D ).
193*
194 IF( lsame( uplo, 'U' ) ) THEN
195 CALL strmm( 'Left', uplo, 'Transpose', 'Unit', n-1, n,
196 $ one, afac( 1, 2 ), ldafac, c( 2, 1 ), ldc )
197 ELSE
198 CALL strmm( 'Left', uplo, 'No transpose', 'Unit', n-1, n,
199 $ one, afac( 2, 1 ), ldafac, c( 2, 1 ), ldc )
200 END IF
201*
202* Call STRMM again to multiply by U (or L ).
203*
204 IF( lsame( uplo, 'U' ) ) THEN
205 CALL strmm( 'Right', uplo, 'No transpose', 'Unit', n, n-1,
206 $ one, afac( 1, 2 ), ldafac, c( 1, 2 ), ldc )
207 ELSE
208 CALL strmm( 'Right', uplo, 'Transpose', 'Unit', n, n-1,
209 $ one, afac( 2, 1 ), ldafac, c( 1, 2 ), ldc )
210 END IF
211 ENDIF
212*
213* Apply symmetric pivots
214*
215 DO j = n, 1, -1
216 i = ipiv( j )
217 IF( i.NE.j )
218 $ CALL sswap( n, c( j, 1 ), ldc, c( i, 1 ), ldc )
219 END DO
220 DO j = n, 1, -1
221 i = ipiv( j )
222 IF( i.NE.j )
223 $ CALL sswap( n, c( 1, j ), 1, c( 1, i ), 1 )
224 END DO
225*
226*
227* Compute the difference C - A .
228*
229 IF( lsame( uplo, 'U' ) ) THEN
230 DO j = 1, n
231 DO i = 1, j
232 c( i, j ) = c( i, j ) - a( i, j )
233 END DO
234 END DO
235 ELSE
236 DO j = 1, n
237 DO i = j, n
238 c( i, j ) = c( i, j ) - a( i, j )
239 END DO
240 END DO
241 END IF
242*
243* Compute norm( C - A ) / ( N * norm(A) * EPS )
244*
245 resid = slansy( '1', uplo, n, c, ldc, rwork )
246*
247 IF( anorm.LE.zero ) THEN
248 IF( resid.NE.zero )
249 $ resid = one / eps
250 ELSE
251 resid = ( ( resid / real( n ) ) / anorm ) / eps
252 END IF
253*
254 RETURN
255*
256* End of SSYT01_AA
257*
258 END
subroutine slacpy(uplo, m, n, a, lda, b, ldb)
SLACPY copies all or part of one two-dimensional array to another.
Definition slacpy.f:103
subroutine slaset(uplo, m, n, alpha, beta, a, lda)
SLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition slaset.f:110
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82
subroutine strmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
STRMM
Definition strmm.f:177
subroutine slavsy(uplo, trans, diag, n, nrhs, a, lda, ipiv, b, ldb, info)
SLAVSY
Definition slavsy.f:155
subroutine ssyt01_aa(uplo, n, a, lda, afac, ldafac, ipiv, c, ldc, rwork, resid)
SSYT01_AA
Definition ssyt01_aa.f:124