LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
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 *> \date December 2016
120 *
121 *
122 *> \ingroup real_lin
123 *
124 * =====================================================================
125  SUBROUTINE ssyt01_aa( UPLO, N, A, LDA, AFAC, LDAFAC, IPIV, C,
126  $ ldc, rwork, resid )
127 *
128 * -- LAPACK test routine (version 3.7.0) --
129 * -- LAPACK is a software package provided by Univ. of Tennessee, --
130 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
131 * December 2016
132 *
133 * .. Scalar Arguments ..
134  CHARACTER UPLO
135  INTEGER LDA, LDAFAC, LDC, N
136  REAL RESID
137 * ..
138 * .. Array Arguments ..
139  INTEGER IPIV( * )
140  REAL A( lda, * ), AFAC( ldafac, * ), C( ldc, * ),
141  $ rwork( * )
142 * ..
143 *
144 * =====================================================================
145 *
146 * .. Parameters ..
147  REAL ZERO, ONE
148  parameter ( zero = 0.0e+0, one = 1.0e+0 )
149 * ..
150 * .. Local Scalars ..
151  INTEGER I, J
152  REAL ANORM, EPS
153 * ..
154 * .. External Functions ..
155  LOGICAL LSAME
156  REAL SLAMCH, SLANSY
157  EXTERNAL lsame, slamch, slansy
158 * ..
159 * .. External Subroutines ..
160  EXTERNAL slaset, slavsy
161 * ..
162 * .. Intrinsic Functions ..
163  INTRINSIC dble
164 * ..
165 * .. Executable Statements ..
166 *
167 * Quick exit if N = 0.
168 *
169  IF( n.LE.0 ) THEN
170  resid = zero
171  RETURN
172  END IF
173 *
174 * Determine EPS and the norm of A.
175 *
176  eps = slamch( 'Epsilon' )
177  anorm = slansy( '1', uplo, n, a, lda, rwork )
178 *
179 * Initialize C to the tridiagonal matrix T.
180 *
181  CALL slaset( 'Full', n, n, zero, zero, c, ldc )
182  CALL slacpy( 'F', 1, n, afac( 1, 1 ), ldafac+1, c( 1, 1 ), ldc+1 )
183  IF( n.GT.1 ) THEN
184  IF( lsame( uplo, 'U' ) ) THEN
185  CALL slacpy( 'F', 1, n-1, afac( 1, 2 ), ldafac+1, c( 1, 2 ),
186  $ ldc+1 )
187  CALL slacpy( 'F', 1, n-1, afac( 1, 2 ), ldafac+1, c( 2, 1 ),
188  $ ldc+1 )
189  ELSE
190  CALL slacpy( 'F', 1, n-1, afac( 2, 1 ), ldafac+1, c( 1, 2 ),
191  $ ldc+1 )
192  CALL slacpy( 'F', 1, n-1, afac( 2, 1 ), ldafac+1, c( 2, 1 ),
193  $ ldc+1 )
194  ENDIF
195 *
196 * Call STRMM to form the product U' * D (or L * D ).
197 *
198  IF( lsame( uplo, 'U' ) ) THEN
199  CALL strmm( 'Left', uplo, 'Transpose', 'Unit', n-1, n,
200  $ one, afac( 1, 2 ), ldafac, c( 2, 1 ), ldc )
201  ELSE
202  CALL strmm( 'Left', uplo, 'No transpose', 'Unit', n-1, n,
203  $ one, afac( 2, 1 ), ldafac, c( 2, 1 ), ldc )
204  END IF
205 *
206 * Call STRMM again to multiply by U (or L ).
207 *
208  IF( lsame( uplo, 'U' ) ) THEN
209  CALL strmm( 'Right', uplo, 'No transpose', 'Unit', n, n-1,
210  $ one, afac( 1, 2 ), ldafac, c( 1, 2 ), ldc )
211  ELSE
212  CALL strmm( 'Right', uplo, 'Transpose', 'Unit', n, n-1,
213  $ one, afac( 2, 1 ), ldafac, c( 1, 2 ), ldc )
214  END IF
215  ENDIF
216 *
217 * Apply symmetric pivots
218 *
219  DO j = n, 1, -1
220  i = ipiv( j )
221  IF( i.NE.j )
222  $ CALL sswap( n, c( j, 1 ), ldc, c( i, 1 ), ldc )
223  END DO
224  DO j = n, 1, -1
225  i = ipiv( j )
226  IF( i.NE.j )
227  $ CALL sswap( n, c( 1, j ), 1, c( 1, i ), 1 )
228  END DO
229 *
230 *
231 * Compute the difference C - A .
232 *
233  IF( lsame( uplo, 'U' ) ) THEN
234  DO j = 1, n
235  DO i = 1, j
236  c( i, j ) = c( i, j ) - a( i, j )
237  END DO
238  END DO
239  ELSE
240  DO j = 1, n
241  DO i = j, n
242  c( i, j ) = c( i, j ) - a( i, j )
243  END DO
244  END DO
245  END IF
246 *
247 * Compute norm( C - A ) / ( N * norm(A) * EPS )
248 *
249  resid = slansy( '1', uplo, n, c, ldc, rwork )
250 *
251  IF( anorm.LE.zero ) THEN
252  IF( resid.NE.zero )
253  $ resid = one / eps
254  ELSE
255  resid = ( ( resid / dble( n ) ) / anorm ) / eps
256  END IF
257 *
258  RETURN
259 *
260 * End of SSYT01_AA
261 *
262  END
subroutine ssyt01_aa(UPLO, N, A, LDA, AFAC, LDAFAC, IPIV, C, LDC, RWORK, RESID)
SSYT01_AA
Definition: ssyt01_aa.f:127
subroutine slavsy(UPLO, TRANS, DIAG, N, NRHS, A, LDA, IPIV, B, LDB, INFO)
SLAVSY
Definition: slavsy.f:157
subroutine strmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
STRMM
Definition: strmm.f:179
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:112
subroutine sswap(N, SX, INCX, SY, INCY)
SSWAP
Definition: sswap.f:53
subroutine slacpy(UPLO, M, N, A, LDA, B, LDB)
SLACPY copies all or part of one two-dimensional array to another.
Definition: slacpy.f:105