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

◆ zbdt05()

subroutine zbdt05 ( integer  m,
integer  n,
complex*16, dimension( lda, * )  a,
integer  lda,
double precision, dimension( * )  s,
integer  ns,
complex*16, dimension( * )  u,
integer  ldu,
complex*16, dimension( ldvt, * )  vt,
integer  ldvt,
complex*16, dimension( * )  work,
double precision  resid 
)

ZBDT05

Purpose:
 ZBDT05 reconstructs a bidiagonal matrix B from its (partial) SVD:
    S = U' * B * V
 where U and V are orthogonal matrices and S is diagonal.

 The test ratio to test the singular value decomposition is
    RESID = norm( S - U' * B * V ) / ( n * norm(B) * EPS )
 where VT = V' and EPS is the machine precision.
Parameters
[in]M
          M is INTEGER
          The number of rows of the matrices A and U.
[in]N
          N is INTEGER
          The number of columns of the matrices A and VT.
[in]A
          A is COMPLEX*16 array, dimension (LDA,N)
          The m by n matrix A.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,M).
[in]S
          S is DOUBLE PRECISION array, dimension (NS)
          The singular values from the (partial) SVD of B, sorted in
          decreasing order.
[in]NS
          NS is INTEGER
          The number of singular values/vectors from the (partial)
          SVD of B.
[in]U
          U is COMPLEX*16 array, dimension (LDU,NS)
          The n by ns orthogonal matrix U in S = U' * B * V.
[in]LDU
          LDU is INTEGER
          The leading dimension of the array U.  LDU >= max(1,N)
[in]VT
          VT is COMPLEX*16 array, dimension (LDVT,N)
          The n by ns orthogonal matrix V in S = U' * B * V.
[in]LDVT
          LDVT is INTEGER
          The leading dimension of the array VT.
[out]WORK
          WORK is COMPLEX*16 array, dimension (M,N)
[out]RESID
          RESID is DOUBLE PRECISION
          The test ratio:  norm(S - U' * A * V) / ( n * norm(A) * EPS )
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 123 of file zbdt05.f.

125*
126* -- LAPACK test routine --
127* -- LAPACK is a software package provided by Univ. of Tennessee, --
128* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
129*
130* .. Scalar Arguments ..
131 INTEGER LDA, LDU, LDVT, M, N, NS
132 DOUBLE PRECISION RESID
133* ..
134* .. Array Arguments ..
135 DOUBLE PRECISION S( * )
136 COMPLEX*16 A( LDA, * ), U( * ), VT( LDVT, * ), WORK( * )
137* ..
138*
139* ======================================================================
140*
141* .. Parameters ..
142 DOUBLE PRECISION ZERO, ONE
143 parameter( zero = 0.0d+0, one = 1.0d+0 )
144 COMPLEX*16 CZERO, CONE
145 parameter( czero = ( 0.0d+0, 0.0d+0 ),
146 $ cone = ( 1.0d+0, 0.0d+0 ) )
147* ..
148* .. Local Scalars ..
149 INTEGER I, J
150 DOUBLE PRECISION ANORM, EPS
151* ..
152* .. Local Arrays ..
153 DOUBLE PRECISION DUM( 1 )
154* ..
155* .. External Functions ..
156 LOGICAL LSAME
157 INTEGER IDAMAX
158 DOUBLE PRECISION DASUM, DZASUM, DLAMCH, ZLANGE
159 EXTERNAL lsame, idamax, dasum, dzasum, dlamch, zlange
160* ..
161* .. External Subroutines ..
162 EXTERNAL zgemm
163* ..
164* .. Intrinsic Functions ..
165 INTRINSIC abs, dble, max, min
166* ..
167* .. Executable Statements ..
168*
169* Quick return if possible.
170*
171 resid = zero
172 IF( min( m, n ).LE.0 .OR. ns.LE.0 )
173 $ RETURN
174*
175 eps = dlamch( 'Precision' )
176 anorm = zlange( 'M', m, n, a, lda, dum )
177*
178* Compute U' * A * V.
179*
180 CALL zgemm( 'N', 'C', m, ns, n, cone, a, lda, vt,
181 $ ldvt, czero, work( 1+ns*ns ), m )
182 CALL zgemm( 'C', 'N', ns, ns, m, -cone, u, ldu, work( 1+ns*ns ),
183 $ m, czero, work, ns )
184*
185* norm(S - U' * B * V)
186*
187 j = 0
188 DO 10 i = 1, ns
189 work( j+i ) = work( j+i ) + dcmplx( s( i ), zero )
190 resid = max( resid, dzasum( ns, work( j+1 ), 1 ) )
191 j = j + ns
192 10 CONTINUE
193*
194 IF( anorm.LE.zero ) THEN
195 IF( resid.NE.zero )
196 $ resid = one / eps
197 ELSE
198 IF( anorm.GE.resid ) THEN
199 resid = ( resid / anorm ) / ( dble( n )*eps )
200 ELSE
201 IF( anorm.LT.one ) THEN
202 resid = ( min( resid, dble( n )*anorm ) / anorm ) /
203 $ ( dble( n )*eps )
204 ELSE
205 resid = min( resid / anorm, dble( n ) ) /
206 $ ( dble( n )*eps )
207 END IF
208 END IF
209 END IF
210*
211 RETURN
212*
213* End of ZBDT05
214*
double precision function dasum(n, dx, incx)
DASUM
Definition dasum.f:71
double precision function dzasum(n, zx, incx)
DZASUM
Definition dzasum.f:72
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:188
integer function idamax(n, dx, incx)
IDAMAX
Definition idamax.f:71
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
double precision function zlange(norm, m, n, a, lda, work)
ZLANGE returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlange.f:115
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: