LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zlarhs.f
Go to the documentation of this file.
1*> \brief \b ZLARHS
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 ZLARHS( PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS,
12* A, LDA, X, LDX, B, LDB, ISEED, INFO )
13*
14* .. Scalar Arguments ..
15* CHARACTER TRANS, UPLO, XTYPE
16* CHARACTER*3 PATH
17* INTEGER INFO, KL, KU, LDA, LDB, LDX, M, N, NRHS
18* ..
19* .. Array Arguments ..
20* INTEGER ISEED( 4 )
21* COMPLEX*16 A( LDA, * ), B( LDB, * ), X( LDX, * )
22* ..
23*
24*
25*> \par Purpose:
26* =============
27*>
28*> \verbatim
29*>
30*> ZLARHS chooses a set of NRHS random solution vectors and sets
31*> up the right hand sides for the linear system
32*> op(A) * X = B,
33*> where op(A) = A, A**T, or A**H, depending on TRANS.
34*> \endverbatim
35*
36* Arguments:
37* ==========
38*
39*> \param[in] PATH
40*> \verbatim
41*> PATH is CHARACTER*3
42*> The type of the complex matrix A. PATH may be given in any
43*> combination of upper and lower case. Valid paths include
44*> xGE: General m x n matrix
45*> xGB: General banded matrix
46*> xPO: Hermitian positive definite, 2-D storage
47*> xPP: Hermitian positive definite packed
48*> xPB: Hermitian positive definite banded
49*> xHE: Hermitian indefinite, 2-D storage
50*> xHP: Hermitian indefinite packed
51*> xHB: Hermitian indefinite banded
52*> xSY: Symmetric indefinite, 2-D storage
53*> xSP: Symmetric indefinite packed
54*> xSB: Symmetric indefinite banded
55*> xTR: Triangular
56*> xTP: Triangular packed
57*> xTB: Triangular banded
58*> xQR: General m x n matrix
59*> xLQ: General m x n matrix
60*> xQL: General m x n matrix
61*> xRQ: General m x n matrix
62*> where the leading character indicates the precision.
63*> \endverbatim
64*>
65*> \param[in] XTYPE
66*> \verbatim
67*> XTYPE is CHARACTER*1
68*> Specifies how the exact solution X will be determined:
69*> = 'N': New solution; generate a random X.
70*> = 'C': Computed; use value of X on entry.
71*> \endverbatim
72*>
73*> \param[in] UPLO
74*> \verbatim
75*> UPLO is CHARACTER*1
76*> Used only if A is symmetric or triangular; specifies whether
77*> the upper or lower triangular part of the matrix A is stored.
78*> = 'U': Upper triangular
79*> = 'L': Lower triangular
80*> \endverbatim
81*>
82*> \param[in] TRANS
83*> \verbatim
84*> TRANS is CHARACTER*1
85*> Used only if A is nonsymmetric; specifies the operation
86*> applied to the matrix A.
87*> = 'N': B := A * X (No transpose)
88*> = 'T': B := A**T * X (Transpose)
89*> = 'C': B := A**H * X (Conjugate transpose)
90*> \endverbatim
91*>
92*> \param[in] M
93*> \verbatim
94*> M is INTEGER
95*> The number of rows of the matrix A. M >= 0.
96*> \endverbatim
97*>
98*> \param[in] N
99*> \verbatim
100*> N is INTEGER
101*> The number of columns of the matrix A. N >= 0.
102*> \endverbatim
103*>
104*> \param[in] KL
105*> \verbatim
106*> KL is INTEGER
107*> Used only if A is a band matrix; specifies the number of
108*> subdiagonals of A if A is a general band matrix or if A is
109*> symmetric or triangular and UPLO = 'L'; specifies the number
110*> of superdiagonals of A if A is symmetric or triangular and
111*> UPLO = 'U'. 0 <= KL <= M-1.
112*> \endverbatim
113*>
114*> \param[in] KU
115*> \verbatim
116*> KU is INTEGER
117*> Used only if A is a general band matrix or if A is
118*> triangular.
119*>
120*> If PATH = xGB, specifies the number of superdiagonals of A,
121*> and 0 <= KU <= N-1.
122*>
123*> If PATH = xTR, xTP, or xTB, specifies whether or not the
124*> matrix has unit diagonal:
125*> = 1: matrix has non-unit diagonal (default)
126*> = 2: matrix has unit diagonal
127*> \endverbatim
128*>
129*> \param[in] NRHS
130*> \verbatim
131*> NRHS is INTEGER
132*> The number of right hand side vectors in the system A*X = B.
133*> \endverbatim
134*>
135*> \param[in] A
136*> \verbatim
137*> A is COMPLEX*16 array, dimension (LDA,N)
138*> The test matrix whose type is given by PATH.
139*> \endverbatim
140*>
141*> \param[in] LDA
142*> \verbatim
143*> LDA is INTEGER
144*> The leading dimension of the array A.
145*> If PATH = xGB, LDA >= KL+KU+1.
146*> If PATH = xPB, xSB, xHB, or xTB, LDA >= KL+1.
147*> Otherwise, LDA >= max(1,M).
148*> \endverbatim
149*>
150*> \param[in,out] X
151*> \verbatim
152*> X is or output) COMPLEX*16 array, dimension (LDX,NRHS)
153*> On entry, if XTYPE = 'C' (for 'Computed'), then X contains
154*> the exact solution to the system of linear equations.
155*> On exit, if XTYPE = 'N' (for 'New'), then X is initialized
156*> with random values.
157*> \endverbatim
158*>
159*> \param[in] LDX
160*> \verbatim
161*> LDX is INTEGER
162*> The leading dimension of the array X. If TRANS = 'N',
163*> LDX >= max(1,N); if TRANS = 'T', LDX >= max(1,M).
164*> \endverbatim
165*>
166*> \param[out] B
167*> \verbatim
168*> B is COMPLEX*16 array, dimension (LDB,NRHS)
169*> The right hand side vector(s) for the system of equations,
170*> computed from B = op(A) * X, where op(A) is determined by
171*> TRANS.
172*> \endverbatim
173*>
174*> \param[in] LDB
175*> \verbatim
176*> LDB is INTEGER
177*> The leading dimension of the array B. If TRANS = 'N',
178*> LDB >= max(1,M); if TRANS = 'T', LDB >= max(1,N).
179*> \endverbatim
180*>
181*> \param[in,out] ISEED
182*> \verbatim
183*> ISEED is INTEGER array, dimension (4)
184*> The seed vector for the random number generator (used in
185*> ZLATMS). Modified on exit.
186*> \endverbatim
187*>
188*> \param[out] INFO
189*> \verbatim
190*> INFO is INTEGER
191*> = 0: successful exit
192*> < 0: if INFO = -i, the i-th argument had an illegal value
193*> \endverbatim
194*
195* Authors:
196* ========
197*
198*> \author Univ. of Tennessee
199*> \author Univ. of California Berkeley
200*> \author Univ. of Colorado Denver
201*> \author NAG Ltd.
202*
203*> \ingroup complex16_eig
204*
205* =====================================================================
206 SUBROUTINE zlarhs( PATH, XTYPE, UPLO, TRANS, M, N, KL, KU, NRHS,
207 $ A, LDA, X, LDX, B, LDB, ISEED, INFO )
208*
209* -- LAPACK test routine --
210* -- LAPACK is a software package provided by Univ. of Tennessee, --
211* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
212*
213* .. Scalar Arguments ..
214 CHARACTER TRANS, UPLO, XTYPE
215 CHARACTER*3 PATH
216 INTEGER INFO, KL, KU, LDA, LDB, LDX, M, N, NRHS
217* ..
218* .. Array Arguments ..
219 INTEGER ISEED( 4 )
220 COMPLEX*16 A( LDA, * ), B( LDB, * ), X( LDX, * )
221* ..
222*
223* =====================================================================
224*
225* .. Parameters ..
226 COMPLEX*16 ONE, ZERO
227 parameter( one = ( 1.0d+0, 0.0d+0 ),
228 $ zero = ( 0.0d+0, 0.0d+0 ) )
229* ..
230* .. Local Scalars ..
231 LOGICAL BAND, GEN, NOTRAN, QRS, SYM, TRAN, TRI
232 CHARACTER C1, DIAG
233 CHARACTER*2 C2
234 INTEGER J, MB, NX
235* ..
236* .. External Functions ..
237 LOGICAL LSAME, LSAMEN
238 EXTERNAL lsame, lsamen
239* ..
240* .. External Subroutines ..
241 EXTERNAL xerbla, zgbmv, zgemm, zhbmv, zhemm, zhpmv,
243 $ ztpmv, ztrmm
244* ..
245* .. Intrinsic Functions ..
246 INTRINSIC max
247* ..
248* .. Executable Statements ..
249*
250* Test the input parameters.
251*
252 info = 0
253 c1 = path( 1: 1 )
254 c2 = path( 2: 3 )
255 tran = lsame( trans, 'T' ) .OR. lsame( trans, 'C' )
256 notran = .NOT.tran
257 gen = lsame( path( 2: 2 ), 'G' )
258 qrs = lsame( path( 2: 2 ), 'Q' ) .OR. lsame( path( 3: 3 ), 'Q' )
259 sym = lsame( path( 2: 2 ), 'P' ) .OR.
260 $ lsame( path( 2: 2 ), 'S' ) .OR. lsame( path( 2: 2 ), 'H' )
261 tri = lsame( path( 2: 2 ), 'T' )
262 band = lsame( path( 3: 3 ), 'B' )
263 IF( .NOT.lsame( c1, 'Zomplex precision' ) ) THEN
264 info = -1
265 ELSE IF( .NOT.( lsame( xtype, 'N' ) .OR. lsame( xtype, 'C' ) ) )
266 $ THEN
267 info = -2
268 ELSE IF( ( sym .OR. tri ) .AND. .NOT.
269 $ ( lsame( uplo, 'U' ) .OR. lsame( uplo, 'L' ) ) ) THEN
270 info = -3
271 ELSE IF( ( gen .OR. qrs ) .AND. .NOT.
272 $ ( tran .OR. lsame( trans, 'N' ) ) ) THEN
273 info = -4
274 ELSE IF( m.LT.0 ) THEN
275 info = -5
276 ELSE IF( n.LT.0 ) THEN
277 info = -6
278 ELSE IF( band .AND. kl.LT.0 ) THEN
279 info = -7
280 ELSE IF( band .AND. ku.LT.0 ) THEN
281 info = -8
282 ELSE IF( nrhs.LT.0 ) THEN
283 info = -9
284 ELSE IF( ( .NOT.band .AND. lda.LT.max( 1, m ) ) .OR.
285 $ ( band .AND. ( sym .OR. tri ) .AND. lda.LT.kl+1 ) .OR.
286 $ ( band .AND. gen .AND. lda.LT.kl+ku+1 ) ) THEN
287 info = -11
288 ELSE IF( ( notran .AND. ldx.LT.max( 1, n ) ) .OR.
289 $ ( tran .AND. ldx.LT.max( 1, m ) ) ) THEN
290 info = -13
291 ELSE IF( ( notran .AND. ldb.LT.max( 1, m ) ) .OR.
292 $ ( tran .AND. ldb.LT.max( 1, n ) ) ) THEN
293 info = -15
294 END IF
295 IF( info.NE.0 ) THEN
296 CALL xerbla( 'ZLARHS', -info )
297 RETURN
298 END IF
299*
300* Initialize X to NRHS random vectors unless XTYPE = 'C'.
301*
302 IF( tran ) THEN
303 nx = m
304 mb = n
305 ELSE
306 nx = n
307 mb = m
308 END IF
309 IF( .NOT.lsame( xtype, 'C' ) ) THEN
310 DO 10 j = 1, nrhs
311 CALL zlarnv( 2, iseed, n, x( 1, j ) )
312 10 CONTINUE
313 END IF
314*
315* Multiply X by op(A) using an appropriate
316* matrix multiply routine.
317*
318 IF( lsamen( 2, c2, 'GE' ) .OR. lsamen( 2, c2, 'QR' ) .OR.
319 $ lsamen( 2, c2, 'LQ' ) .OR. lsamen( 2, c2, 'QL' ) .OR.
320 $ lsamen( 2, c2, 'RQ' ) ) THEN
321*
322* General matrix
323*
324 CALL zgemm( trans, 'N', mb, nrhs, nx, one, a, lda, x, ldx,
325 $ zero, b, ldb )
326*
327 ELSE IF( lsamen( 2, c2, 'PO' ) .OR. lsamen( 2, c2, 'HE' ) ) THEN
328*
329* Hermitian matrix, 2-D storage
330*
331 CALL zhemm( 'Left', uplo, n, nrhs, one, a, lda, x, ldx, zero,
332 $ b, ldb )
333*
334 ELSE IF( lsamen( 2, c2, 'SY' ) ) THEN
335*
336* Symmetric matrix, 2-D storage
337*
338 CALL zsymm( 'Left', uplo, n, nrhs, one, a, lda, x, ldx, zero,
339 $ b, ldb )
340*
341 ELSE IF( lsamen( 2, c2, 'GB' ) ) THEN
342*
343* General matrix, band storage
344*
345 DO 20 j = 1, nrhs
346 CALL zgbmv( trans, m, n, kl, ku, one, a, lda, x( 1, j ), 1,
347 $ zero, b( 1, j ), 1 )
348 20 CONTINUE
349*
350 ELSE IF( lsamen( 2, c2, 'PB' ) .OR. lsamen( 2, c2, 'HB' ) ) THEN
351*
352* Hermitian matrix, band storage
353*
354 DO 30 j = 1, nrhs
355 CALL zhbmv( uplo, n, kl, one, a, lda, x( 1, j ), 1, zero,
356 $ b( 1, j ), 1 )
357 30 CONTINUE
358*
359 ELSE IF( lsamen( 2, c2, 'SB' ) ) THEN
360*
361* Symmetric matrix, band storage
362*
363 DO 40 j = 1, nrhs
364 CALL zsbmv( uplo, n, kl, one, a, lda, x( 1, j ), 1, zero,
365 $ b( 1, j ), 1 )
366 40 CONTINUE
367*
368 ELSE IF( lsamen( 2, c2, 'PP' ) .OR. lsamen( 2, c2, 'HP' ) ) THEN
369*
370* Hermitian matrix, packed storage
371*
372 DO 50 j = 1, nrhs
373 CALL zhpmv( uplo, n, one, a, x( 1, j ), 1, zero, b( 1, j ),
374 $ 1 )
375 50 CONTINUE
376*
377 ELSE IF( lsamen( 2, c2, 'SP' ) ) THEN
378*
379* Symmetric matrix, packed storage
380*
381 DO 60 j = 1, nrhs
382 CALL zspmv( uplo, n, one, a, x( 1, j ), 1, zero, b( 1, j ),
383 $ 1 )
384 60 CONTINUE
385*
386 ELSE IF( lsamen( 2, c2, 'TR' ) ) THEN
387*
388* Triangular matrix. Note that for triangular matrices,
389* KU = 1 => non-unit triangular
390* KU = 2 => unit triangular
391*
392 CALL zlacpy( 'Full', n, nrhs, x, ldx, b, ldb )
393 IF( ku.EQ.2 ) THEN
394 diag = 'U'
395 ELSE
396 diag = 'N'
397 END IF
398 CALL ztrmm( 'Left', uplo, trans, diag, n, nrhs, one, a, lda, b,
399 $ ldb )
400*
401 ELSE IF( lsamen( 2, c2, 'TP' ) ) THEN
402*
403* Triangular matrix, packed storage
404*
405 CALL zlacpy( 'Full', n, nrhs, x, ldx, b, ldb )
406 IF( ku.EQ.2 ) THEN
407 diag = 'U'
408 ELSE
409 diag = 'N'
410 END IF
411 DO 70 j = 1, nrhs
412 CALL ztpmv( uplo, trans, diag, n, a, b( 1, j ), 1 )
413 70 CONTINUE
414*
415 ELSE IF( lsamen( 2, c2, 'TB' ) ) THEN
416*
417* Triangular matrix, banded storage
418*
419 CALL zlacpy( 'Full', n, nrhs, x, ldx, b, ldb )
420 IF( ku.EQ.2 ) THEN
421 diag = 'U'
422 ELSE
423 diag = 'N'
424 END IF
425 DO 80 j = 1, nrhs
426 CALL ztbmv( uplo, trans, diag, n, kl, a, lda, b( 1, j ), 1 )
427 80 CONTINUE
428*
429 ELSE
430*
431* If none of the above, set INFO = -1 and return
432*
433 info = -1
434 CALL xerbla( 'ZLARHS', -info )
435 END IF
436*
437 RETURN
438*
439* End of ZLARHS
440*
441 END
subroutine zlarhs(path, xtype, uplo, trans, m, n, kl, ku, nrhs, a, lda, x, ldx, b, ldb, iseed, info)
ZLARHS
Definition zlarhs.f:208
subroutine zsbmv(uplo, n, k, alpha, a, lda, x, incx, beta, y, incy)
ZSBMV
Definition zsbmv.f:152
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zgbmv(trans, m, n, kl, ku, alpha, a, lda, x, incx, beta, y, incy)
ZGBMV
Definition zgbmv.f:190
subroutine zgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
ZGEMM
Definition zgemm.f:188
subroutine zhbmv(uplo, n, k, alpha, a, lda, x, incx, beta, y, incy)
ZHBMV
Definition zhbmv.f:187
subroutine zhemm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
ZHEMM
Definition zhemm.f:191
subroutine zsymm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
ZSYMM
Definition zsymm.f:189
subroutine zspmv(uplo, n, alpha, ap, x, incx, beta, y, incy)
ZSPMV computes a matrix-vector product for complex vectors using a complex symmetric packed matrix
Definition zspmv.f:151
subroutine zhpmv(uplo, n, alpha, ap, x, incx, beta, y, incy)
ZHPMV
Definition zhpmv.f:149
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
subroutine zlarnv(idist, iseed, n, x)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition zlarnv.f:99
subroutine ztbmv(uplo, trans, diag, n, k, a, lda, x, incx)
ZTBMV
Definition ztbmv.f:186
subroutine ztpmv(uplo, trans, diag, n, ap, x, incx)
ZTPMV
Definition ztpmv.f:142
subroutine ztrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
ZTRMM
Definition ztrmm.f:177