LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zqrt12.f
Go to the documentation of this file.
1*> \brief \b ZQRT12
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8* Definition:
9* ===========
10*
11* DOUBLE PRECISION FUNCTION ZQRT12( M, N, A, LDA, S, WORK, LWORK,
12* RWORK )
13*
14* .. Scalar Arguments ..
15* INTEGER LDA, LWORK, M, N
16* ..
17* .. Array Arguments ..
18* DOUBLE PRECISION RWORK( * ), S( * )
19* COMPLEX*16 A( LDA, * ), WORK( LWORK )
20* ..
21*
22*
23*> \par Purpose:
24* =============
25*>
26*> \verbatim
27*>
28*> ZQRT12 computes the singular values `svlues' of the upper trapezoid
29*> of A(1:M,1:N) and returns the ratio
30*>
31*> || svlues - s||/(||s||*eps*max(M,N))
32*> \endverbatim
33*
34* Arguments:
35* ==========
36*
37*> \param[in] M
38*> \verbatim
39*> M is INTEGER
40*> The number of rows of the matrix A.
41*> \endverbatim
42*>
43*> \param[in] N
44*> \verbatim
45*> N is INTEGER
46*> The number of columns of the matrix A.
47*> \endverbatim
48*>
49*> \param[in] A
50*> \verbatim
51*> A is COMPLEX*16 array, dimension (LDA,N)
52*> The M-by-N matrix A. Only the upper trapezoid is referenced.
53*> \endverbatim
54*>
55*> \param[in] LDA
56*> \verbatim
57*> LDA is INTEGER
58*> The leading dimension of the array A.
59*> \endverbatim
60*>
61*> \param[in] S
62*> \verbatim
63*> S is DOUBLE PRECISION array, dimension (min(M,N))
64*> The singular values of the matrix A.
65*> \endverbatim
66*>
67*> \param[out] WORK
68*> \verbatim
69*> WORK is COMPLEX*16 array, dimension (LWORK)
70*> \endverbatim
71*>
72*> \param[in] LWORK
73*> \verbatim
74*> LWORK is INTEGER
75*> The length of the array WORK. LWORK >= M*N + 2*min(M,N) +
76*> max(M,N).
77*> \endverbatim
78*>
79*> \param[out] RWORK
80*> \verbatim
81*> RWORK is DOUBLE PRECISION array, dimension (2*min(M,N))
82*> \endverbatim
83*
84* Authors:
85* ========
86*
87*> \author Univ. of Tennessee
88*> \author Univ. of California Berkeley
89*> \author Univ. of Colorado Denver
90*> \author NAG Ltd.
91*
92*> \ingroup complex16_lin
93*
94* =====================================================================
95 DOUBLE PRECISION FUNCTION zqrt12( M, N, A, LDA, S, WORK, LWORK,
96 $ RWORK )
97*
98* -- LAPACK test routine --
99* -- LAPACK is a software package provided by Univ. of Tennessee, --
100* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
101*
102* .. Scalar Arguments ..
103 INTEGER lda, lwork, m, n
104* ..
105* .. Array Arguments ..
106 DOUBLE PRECISION rwork( * ), s( * )
107 COMPLEX*16 a( lda, * ), work( lwork )
108* ..
109*
110* =====================================================================
111*
112* .. Parameters ..
113 DOUBLE PRECISION zero, one
114 parameter( zero = 0.0d0, one = 1.0d0 )
115* ..
116* .. Local Scalars ..
117 INTEGER i, info, iscl, j, mn
118 DOUBLE PRECISION anrm, bignum, nrmsvl, smlnum
119* ..
120* .. Local Arrays ..
121 DOUBLE PRECISION dummy( 1 )
122* ..
123* .. External Functions ..
124 DOUBLE PRECISION dasum, dlamch, dnrm2, zlange
125 EXTERNAL dasum, dlamch, dnrm2, zlange
126* ..
127* .. External Subroutines ..
128 EXTERNAL daxpy, dbdsqr, dlascl, xerbla, zgebd2, zlascl,
129 $ zlaset
130* ..
131* .. Intrinsic Functions ..
132 INTRINSIC dble, dcmplx, max, min
133* ..
134* .. Executable Statements ..
135*
136 zqrt12 = zero
137*
138* Test that enough workspace is supplied
139*
140 IF( lwork.LT.m*n+2*min( m, n )+max( m, n ) ) THEN
141 CALL xerbla( 'ZQRT12', 7 )
142 RETURN
143 END IF
144*
145* Quick return if possible
146*
147 mn = min( m, n )
148 IF( mn.LE.zero )
149 $ RETURN
150*
151 nrmsvl = dnrm2( mn, s, 1 )
152*
153* Copy upper triangle of A into work
154*
155 CALL zlaset( 'Full', m, n, dcmplx( zero ), dcmplx( zero ), work,
156 $ m )
157 DO j = 1, n
158 DO i = 1, min( j, m )
159 work( ( j-1 )*m+i ) = a( i, j )
160 END DO
161 END DO
162*
163* Get machine parameters
164*
165 smlnum = dlamch( 'S' ) / dlamch( 'P' )
166 bignum = one / smlnum
167*
168* Scale work if max entry outside range [SMLNUM,BIGNUM]
169*
170 anrm = zlange( 'M', m, n, work, m, dummy )
171 iscl = 0
172 IF( anrm.GT.zero .AND. anrm.LT.smlnum ) THEN
173*
174* Scale matrix norm up to SMLNUM
175*
176 CALL zlascl( 'G', 0, 0, anrm, smlnum, m, n, work, m, info )
177 iscl = 1
178 ELSE IF( anrm.GT.bignum ) THEN
179*
180* Scale matrix norm down to BIGNUM
181*
182 CALL zlascl( 'G', 0, 0, anrm, bignum, m, n, work, m, info )
183 iscl = 1
184 END IF
185*
186 IF( anrm.NE.zero ) THEN
187*
188* Compute SVD of work
189*
190 CALL zgebd2( m, n, work, m, rwork( 1 ), rwork( mn+1 ),
191 $ work( m*n+1 ), work( m*n+mn+1 ),
192 $ work( m*n+2*mn+1 ), info )
193 CALL dbdsqr( 'Upper', mn, 0, 0, 0, rwork( 1 ), rwork( mn+1 ),
194 $ dummy, mn, dummy, 1, dummy, mn, rwork( 2*mn+1 ),
195 $ info )
196*
197 IF( iscl.EQ.1 ) THEN
198 IF( anrm.GT.bignum ) THEN
199 CALL dlascl( 'G', 0, 0, bignum, anrm, mn, 1, rwork( 1 ),
200 $ mn, info )
201 END IF
202 IF( anrm.LT.smlnum ) THEN
203 CALL dlascl( 'G', 0, 0, smlnum, anrm, mn, 1, rwork( 1 ),
204 $ mn, info )
205 END IF
206 END IF
207*
208 ELSE
209*
210 DO i = 1, mn
211 rwork( i ) = zero
212 END DO
213 END IF
214*
215* Compare s and singular values of work
216*
217 CALL daxpy( mn, -one, s, 1, rwork( 1 ), 1 )
218 zqrt12 = dasum( mn, rwork( 1 ), 1 ) /
219 $ ( dlamch( 'Epsilon' )*dble( max( m, n ) ) )
220*
221 IF( nrmsvl.NE.zero )
222 $ zqrt12 = zqrt12 / nrmsvl
223*
224 RETURN
225*
226* End of ZQRT12
227*
228 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
double precision function dasum(n, dx, incx)
DASUM
Definition dasum.f:71
subroutine daxpy(n, da, dx, incx, dy, incy)
DAXPY
Definition daxpy.f:89
subroutine dbdsqr(uplo, n, ncvt, nru, ncc, d, e, vt, ldvt, u, ldu, c, ldc, work, info)
DBDSQR
Definition dbdsqr.f:241
subroutine zgebd2(m, n, a, lda, d, e, tauq, taup, work, info)
ZGEBD2 reduces a general matrix to bidiagonal form using an unblocked algorithm.
Definition zgebd2.f:189
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
subroutine zlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
ZLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition zlascl.f:143
subroutine dlascl(type, kl, ku, cfrom, cto, m, n, a, lda, info)
DLASCL multiplies a general rectangular matrix by a real scalar defined as cto/cfrom.
Definition dlascl.f:143
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
real(wp) function dnrm2(n, x, incx)
DNRM2
Definition dnrm2.f90:89
double precision function zqrt12(m, n, a, lda, s, work, lwork, rwork)
ZQRT12
Definition zqrt12.f:97