LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zhpt01()

subroutine zhpt01 ( character  uplo,
integer  n,
complex*16, dimension( * )  a,
complex*16, dimension( * )  afac,
integer, dimension( * )  ipiv,
complex*16, dimension( ldc, * )  c,
integer  ldc,
double precision, dimension( * )  rwork,
double precision  resid 
)

ZHPT01

Purpose:
 ZHPT01 reconstructs a Hermitian indefinite packed matrix A from its
 block L*D*L' or U*D*U' factorization and computes the residual
    norm( C - A ) / ( N * norm(A) * EPS ),
 where C is the reconstructed matrix, EPS is the machine epsilon,
 L' is the conjugate transpose of L, and U' is the conjugate transpose
 of U.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          Specifies whether the upper or lower triangular part of the
          Hermitian matrix A is stored:
          = 'U':  Upper triangular
          = 'L':  Lower triangular
[in]N
          N is INTEGER
          The number of rows and columns of the matrix A.  N >= 0.
[in]A
          A is COMPLEX*16 array, dimension (N*(N+1)/2)
          The original Hermitian matrix A, stored as a packed
          triangular matrix.
[in]AFAC
          AFAC is COMPLEX*16 array, dimension (N*(N+1)/2)
          The factored form of the matrix A, stored as a packed
          triangular matrix.  AFAC contains the block diagonal matrix D
          and the multipliers used to obtain the factor L or U from the
          block L*D*L' or U*D*U' factorization as computed by ZHPTRF.
[in]IPIV
          IPIV is INTEGER array, dimension (N)
          The pivot indices from ZHPTRF.
[out]C
          C is COMPLEX*16 array, dimension (LDC,N)
[in]LDC
          LDC is INTEGER
          The leading dimension of the array C.  LDC >= max(1,N).
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (N)
[out]RESID
          RESID is DOUBLE PRECISION
          If UPLO = 'L', norm(L*D*L' - A) / ( N * norm(A) * EPS )
          If UPLO = 'U', norm(U*D*U' - A) / ( N * norm(A) * EPS )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 112 of file zhpt01.f.

113*
114* -- LAPACK test routine --
115* -- LAPACK is a software package provided by Univ. of Tennessee, --
116* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
117*
118* .. Scalar Arguments ..
119 CHARACTER UPLO
120 INTEGER LDC, N
121 DOUBLE PRECISION RESID
122* ..
123* .. Array Arguments ..
124 INTEGER IPIV( * )
125 DOUBLE PRECISION RWORK( * )
126 COMPLEX*16 A( * ), AFAC( * ), C( LDC, * )
127* ..
128*
129* =====================================================================
130*
131* .. Parameters ..
132 DOUBLE PRECISION ZERO, ONE
133 parameter( zero = 0.0d+0, one = 1.0d+0 )
134 COMPLEX*16 CZERO, CONE
135 parameter( czero = ( 0.0d+0, 0.0d+0 ),
136 $ cone = ( 1.0d+0, 0.0d+0 ) )
137* ..
138* .. Local Scalars ..
139 INTEGER I, INFO, J, JC
140 DOUBLE PRECISION ANORM, EPS
141* ..
142* .. External Functions ..
143 LOGICAL LSAME
144 DOUBLE PRECISION DLAMCH, ZLANHE, ZLANHP
145 EXTERNAL lsame, dlamch, zlanhe, zlanhp
146* ..
147* .. External Subroutines ..
148 EXTERNAL zlaset, zlavhp
149* ..
150* .. Intrinsic Functions ..
151 INTRINSIC dble, dimag
152* ..
153* .. Executable Statements ..
154*
155* Quick exit if N = 0.
156*
157 IF( n.LE.0 ) THEN
158 resid = zero
159 RETURN
160 END IF
161*
162* Determine EPS and the norm of A.
163*
164 eps = dlamch( 'Epsilon' )
165 anorm = zlanhp( '1', uplo, n, a, rwork )
166*
167* Check the imaginary parts of the diagonal elements and return with
168* an error code if any are nonzero.
169*
170 jc = 1
171 IF( lsame( uplo, 'U' ) ) THEN
172 DO 10 j = 1, n
173 IF( dimag( afac( jc ) ).NE.zero ) THEN
174 resid = one / eps
175 RETURN
176 END IF
177 jc = jc + j + 1
178 10 CONTINUE
179 ELSE
180 DO 20 j = 1, n
181 IF( dimag( afac( jc ) ).NE.zero ) THEN
182 resid = one / eps
183 RETURN
184 END IF
185 jc = jc + n - j + 1
186 20 CONTINUE
187 END IF
188*
189* Initialize C to the identity matrix.
190*
191 CALL zlaset( 'Full', n, n, czero, cone, c, ldc )
192*
193* Call ZLAVHP to form the product D * U' (or D * L' ).
194*
195 CALL zlavhp( uplo, 'Conjugate', 'Non-unit', n, n, afac, ipiv, c,
196 $ ldc, info )
197*
198* Call ZLAVHP again to multiply by U ( or L ).
199*
200 CALL zlavhp( uplo, 'No transpose', 'Unit', n, n, afac, ipiv, c,
201 $ ldc, info )
202*
203* Compute the difference C - A .
204*
205 IF( lsame( uplo, 'U' ) ) THEN
206 jc = 0
207 DO 40 j = 1, n
208 DO 30 i = 1, j - 1
209 c( i, j ) = c( i, j ) - a( jc+i )
210 30 CONTINUE
211 c( j, j ) = c( j, j ) - dble( a( jc+j ) )
212 jc = jc + j
213 40 CONTINUE
214 ELSE
215 jc = 1
216 DO 60 j = 1, n
217 c( j, j ) = c( j, j ) - dble( a( jc ) )
218 DO 50 i = j + 1, n
219 c( i, j ) = c( i, j ) - a( jc+i-j )
220 50 CONTINUE
221 jc = jc + n - j + 1
222 60 CONTINUE
223 END IF
224*
225* Compute norm( C - A ) / ( N * norm(A) * EPS )
226*
227 resid = zlanhe( '1', uplo, n, c, ldc, rwork )
228*
229 IF( anorm.LE.zero ) THEN
230 IF( resid.NE.zero )
231 $ resid = one / eps
232 ELSE
233 resid = ( ( resid / dble( n ) ) / anorm ) / eps
234 END IF
235*
236 RETURN
237*
238* End of ZHPT01
239*
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlanhe(norm, uplo, n, a, lda, work)
ZLANHE returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanhe.f:124
double precision function zlanhp(norm, uplo, n, ap, work)
ZLANHP returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanhp.f:117
subroutine zlaset(uplo, m, n, alpha, beta, a, lda)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition zlaset.f:106
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine zlavhp(uplo, trans, diag, n, nrhs, a, ipiv, b, ldb, info)
ZLAVHP
Definition zlavhp.f:131
Here is the call graph for this function:
Here is the caller graph for this function: