LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zqrt15.f
Go to the documentation of this file.
1*> \brief \b ZQRT15
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 ZQRT15( SCALE, RKSEL, M, N, NRHS, A, LDA, B, LDB, S,
12* RANK, NORMA, NORMB, ISEED, WORK, LWORK )
13*
14* .. Scalar Arguments ..
15* INTEGER LDA, LDB, LWORK, M, N, NRHS, RANK, RKSEL, SCALE
16* DOUBLE PRECISION NORMA, NORMB
17* ..
18* .. Array Arguments ..
19* INTEGER ISEED( 4 )
20* DOUBLE PRECISION S( * )
21* COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( LWORK )
22* ..
23*
24*
25*> \par Purpose:
26* =============
27*>
28*> \verbatim
29*>
30*> ZQRT15 generates a matrix with full or deficient rank and of various
31*> norms.
32*> \endverbatim
33*
34* Arguments:
35* ==========
36*
37*> \param[in] SCALE
38*> \verbatim
39*> SCALE is INTEGER
40*> SCALE = 1: normally scaled matrix
41*> SCALE = 2: matrix scaled up
42*> SCALE = 3: matrix scaled down
43*> \endverbatim
44*>
45*> \param[in] RKSEL
46*> \verbatim
47*> RKSEL is INTEGER
48*> RKSEL = 1: full rank matrix
49*> RKSEL = 2: rank-deficient matrix
50*> \endverbatim
51*>
52*> \param[in] M
53*> \verbatim
54*> M is INTEGER
55*> The number of rows of the matrix A.
56*> \endverbatim
57*>
58*> \param[in] N
59*> \verbatim
60*> N is INTEGER
61*> The number of columns of A.
62*> \endverbatim
63*>
64*> \param[in] NRHS
65*> \verbatim
66*> NRHS is INTEGER
67*> The number of columns of B.
68*> \endverbatim
69*>
70*> \param[out] A
71*> \verbatim
72*> A is COMPLEX*16 array, dimension (LDA,N)
73*> The M-by-N matrix A.
74*> \endverbatim
75*>
76*> \param[in] LDA
77*> \verbatim
78*> LDA is INTEGER
79*> The leading dimension of the array A.
80*> \endverbatim
81*>
82*> \param[out] B
83*> \verbatim
84*> B is COMPLEX*16 array, dimension (LDB, NRHS)
85*> A matrix that is in the range space of matrix A.
86*> \endverbatim
87*>
88*> \param[in] LDB
89*> \verbatim
90*> LDB is INTEGER
91*> The leading dimension of the array B.
92*> \endverbatim
93*>
94*> \param[out] S
95*> \verbatim
96*> S is DOUBLE PRECISION array, dimension MIN(M,N)
97*> Singular values of A.
98*> \endverbatim
99*>
100*> \param[out] RANK
101*> \verbatim
102*> RANK is INTEGER
103*> number of nonzero singular values of A.
104*> \endverbatim
105*>
106*> \param[out] NORMA
107*> \verbatim
108*> NORMA is DOUBLE PRECISION
109*> one-norm norm of A.
110*> \endverbatim
111*>
112*> \param[out] NORMB
113*> \verbatim
114*> NORMB is DOUBLE PRECISION
115*> one-norm norm of B.
116*> \endverbatim
117*>
118*> \param[in,out] ISEED
119*> \verbatim
120*> ISEED is integer array, dimension (4)
121*> seed for random number generator.
122*> \endverbatim
123*>
124*> \param[out] WORK
125*> \verbatim
126*> WORK is COMPLEX*16 array, dimension (LWORK)
127*> \endverbatim
128*>
129*> \param[in] LWORK
130*> \verbatim
131*> LWORK is INTEGER
132*> length of work space required.
133*> LWORK >= MAX(M+MIN(M,N),NRHS*MIN(M,N),2*N+M)
134*> \endverbatim
135*
136* Authors:
137* ========
138*
139*> \author Univ. of Tennessee
140*> \author Univ. of California Berkeley
141*> \author Univ. of Colorado Denver
142*> \author NAG Ltd.
143*
144*> \ingroup complex16_lin
145*
146* =====================================================================
147 SUBROUTINE zqrt15( SCALE, RKSEL, M, N, NRHS, A, LDA, B, LDB, S,
148 $ RANK, NORMA, NORMB, ISEED, WORK, LWORK )
149*
150* -- LAPACK test routine --
151* -- LAPACK is a software package provided by Univ. of Tennessee, --
152* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
153*
154* .. Scalar Arguments ..
155 INTEGER LDA, LDB, LWORK, M, N, NRHS, RANK, RKSEL, SCALE
156 DOUBLE PRECISION NORMA, NORMB
157* ..
158* .. Array Arguments ..
159 INTEGER ISEED( 4 )
160 DOUBLE PRECISION S( * )
161 COMPLEX*16 A( LDA, * ), B( LDB, * ), WORK( LWORK )
162* ..
163*
164* =====================================================================
165*
166* .. Parameters ..
167 DOUBLE PRECISION ZERO, ONE, TWO, SVMIN
168 parameter( zero = 0.0d+0, one = 1.0d+0, two = 2.0d+0,
169 $ svmin = 0.1d+0 )
170 COMPLEX*16 CZERO, CONE
171 parameter( czero = ( 0.0d+0, 0.0d+0 ),
172 $ cone = ( 1.0d+0, 0.0d+0 ) )
173* ..
174* .. Local Scalars ..
175 INTEGER INFO, J, MN
176 DOUBLE PRECISION BIGNUM, EPS, SMLNUM, TEMP
177* ..
178* .. Local Arrays ..
179 DOUBLE PRECISION DUMMY( 1 )
180* ..
181* .. External Functions ..
182 DOUBLE PRECISION DASUM, DLAMCH, DLARND, DZNRM2, ZLANGE
183 EXTERNAL dasum, dlamch, dlarnd, dznrm2, zlange
184* ..
185* .. External Subroutines ..
186 EXTERNAL dlabad, dlaord, dlascl, xerbla, zdscal, zgemm,
188* ..
189* .. Intrinsic Functions ..
190 INTRINSIC abs, dcmplx, max, min
191* ..
192* .. Executable Statements ..
193*
194 mn = min( m, n )
195 IF( lwork.LT.max( m+mn, mn*nrhs, 2*n+m ) ) THEN
196 CALL xerbla( 'ZQRT15', 16 )
197 RETURN
198 END IF
199*
200 smlnum = dlamch( 'Safe minimum' )
201 bignum = one / smlnum
202 CALL dlabad( smlnum, bignum )
203 eps = dlamch( 'Epsilon' )
204 smlnum = ( smlnum / eps ) / eps
205 bignum = one / smlnum
206*
207* Determine rank and (unscaled) singular values
208*
209 IF( rksel.EQ.1 ) THEN
210 rank = mn
211 ELSE IF( rksel.EQ.2 ) THEN
212 rank = ( 3*mn ) / 4
213 DO 10 j = rank + 1, mn
214 s( j ) = zero
215 10 CONTINUE
216 ELSE
217 CALL xerbla( 'ZQRT15', 2 )
218 END IF
219*
220 IF( rank.GT.0 ) THEN
221*
222* Nontrivial case
223*
224 s( 1 ) = one
225 DO 30 j = 2, rank
226 20 CONTINUE
227 temp = dlarnd( 1, iseed )
228 IF( temp.GT.svmin ) THEN
229 s( j ) = abs( temp )
230 ELSE
231 GO TO 20
232 END IF
233 30 CONTINUE
234 CALL dlaord( 'Decreasing', rank, s, 1 )
235*
236* Generate 'rank' columns of a random orthogonal matrix in A
237*
238 CALL zlarnv( 2, iseed, m, work )
239 CALL zdscal( m, one / dznrm2( m, work, 1 ), work, 1 )
240 CALL zlaset( 'Full', m, rank, czero, cone, a, lda )
241 CALL zlarf( 'Left', m, rank, work, 1, dcmplx( two ), a, lda,
242 $ work( m+1 ) )
243*
244* workspace used: m+mn
245*
246* Generate consistent rhs in the range space of A
247*
248 CALL zlarnv( 2, iseed, rank*nrhs, work )
249 CALL zgemm( 'No transpose', 'No transpose', m, nrhs, rank,
250 $ cone, a, lda, work, rank, czero, b, ldb )
251*
252* work space used: <= mn *nrhs
253*
254* generate (unscaled) matrix A
255*
256 DO 40 j = 1, rank
257 CALL zdscal( m, s( j ), a( 1, j ), 1 )
258 40 CONTINUE
259 IF( rank.LT.n )
260 $ CALL zlaset( 'Full', m, n-rank, czero, czero,
261 $ a( 1, rank+1 ), lda )
262 CALL zlaror( 'Right', 'No initialization', m, n, a, lda, iseed,
263 $ work, info )
264*
265 ELSE
266*
267* work space used 2*n+m
268*
269* Generate null matrix and rhs
270*
271 DO 50 j = 1, mn
272 s( j ) = zero
273 50 CONTINUE
274 CALL zlaset( 'Full', m, n, czero, czero, a, lda )
275 CALL zlaset( 'Full', m, nrhs, czero, czero, b, ldb )
276*
277 END IF
278*
279* Scale the matrix
280*
281 IF( scale.NE.1 ) THEN
282 norma = zlange( 'Max', m, n, a, lda, dummy )
283 IF( norma.NE.zero ) THEN
284 IF( scale.EQ.2 ) THEN
285*
286* matrix scaled up
287*
288 CALL zlascl( 'General', 0, 0, norma, bignum, m, n, a,
289 $ lda, info )
290 CALL dlascl( 'General', 0, 0, norma, bignum, mn, 1, s,
291 $ mn, info )
292 CALL zlascl( 'General', 0, 0, norma, bignum, m, nrhs, b,
293 $ ldb, info )
294 ELSE IF( scale.EQ.3 ) THEN
295*
296* matrix scaled down
297*
298 CALL zlascl( 'General', 0, 0, norma, smlnum, m, n, a,
299 $ lda, info )
300 CALL dlascl( 'General', 0, 0, norma, smlnum, mn, 1, s,
301 $ mn, info )
302 CALL zlascl( 'General', 0, 0, norma, smlnum, m, nrhs, b,
303 $ ldb, info )
304 ELSE
305 CALL xerbla( 'ZQRT15', 1 )
306 RETURN
307 END IF
308 END IF
309 END IF
310*
311 norma = dasum( mn, s, 1 )
312 normb = zlange( 'One-norm', m, nrhs, b, ldb, dummy )
313*
314 RETURN
315*
316* End of ZQRT15
317*
318 END
subroutine dlabad(SMALL, LARGE)
DLABAD
Definition: dlabad.f:74
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 xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:78
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:187
subroutine zqrt15(SCALE, RKSEL, M, N, NRHS, A, LDA, B, LDB, S, RANK, NORMA, NORMB, ISEED, WORK, LWORK)
ZQRT15
Definition: zqrt15.f:149
subroutine zlaror(SIDE, INIT, M, N, A, LDA, ISEED, X, INFO)
ZLAROR
Definition: zlaror.f:158
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 zlarnv(IDIST, ISEED, N, X)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: zlarnv.f:99
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
subroutine zlarf(SIDE, M, N, V, INCV, TAU, C, LDC, WORK)
ZLARF applies an elementary reflector to a general rectangular matrix.
Definition: zlarf.f:128
subroutine dlaord(JOB, N, X, INCX)
DLAORD
Definition: dlaord.f:73