LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zlarhs ( character*3  PATH,
character  XTYPE,
character  UPLO,
character  TRANS,
integer  M,
integer  N,
integer  KL,
integer  KU,
integer  NRHS,
complex*16, dimension( lda, * )  A,
integer  LDA,
complex*16, dimension( ldx, * )  X,
integer  LDX,
complex*16, dimension( ldb, * )  B,
integer  LDB,
integer, dimension( 4 )  ISEED,
integer  INFO 
)

ZLARHS

Purpose:
 ZLARHS chooses a set of NRHS random solution vectors and sets
 up the right hand sides for the linear system
    op( A ) * X = B,
 where op( A ) may be A, A**T (transpose of A), or A**H (conjugate
 transpose of A).
Parameters
[in]PATH
          PATH is CHARACTER*3
          The type of the complex matrix A.  PATH may be given in any
          combination of upper and lower case.  Valid paths include
             xGE:  General m x n matrix
             xGB:  General banded matrix
             xPO:  Hermitian positive definite, 2-D storage
             xPP:  Hermitian positive definite packed
             xPB:  Hermitian positive definite banded
             xHE:  Hermitian indefinite, 2-D storage
             xHP:  Hermitian indefinite packed
             xHB:  Hermitian indefinite banded
             xSY:  Symmetric indefinite, 2-D storage
             xSP:  Symmetric indefinite packed
             xSB:  Symmetric indefinite banded
             xTR:  Triangular
             xTP:  Triangular packed
             xTB:  Triangular banded
             xQR:  General m x n matrix
             xLQ:  General m x n matrix
             xQL:  General m x n matrix
             xRQ:  General m x n matrix
          where the leading character indicates the precision.
[in]XTYPE
          XTYPE is CHARACTER*1
          Specifies how the exact solution X will be determined:
          = 'N':  New solution; generate a random X.
          = 'C':  Computed; use value of X on entry.
[in]UPLO
          UPLO is CHARACTER*1
          Used only if A is symmetric or triangular; specifies whether
          the upper or lower triangular part of the matrix A is stored.
          = 'U':  Upper triangular
          = 'L':  Lower triangular
[in]TRANS
          TRANS is CHARACTER*1
          Used only if A is nonsymmetric; specifies the operation
          applied to the matrix A.
          = 'N':  B := A    * X
          = 'T':  B := A**T * X
          = 'C':  B := A**H * X
[in]M
          M is INTEGER
          The number of rows of the matrix A.  M >= 0.
[in]N
          N is INTEGER
          The number of columns of the matrix A.  N >= 0.
[in]KL
          KL is INTEGER
          Used only if A is a band matrix; specifies the number of
          subdiagonals of A if A is a general band matrix or if A is
          symmetric or triangular and UPLO = 'L'; specifies the number
          of superdiagonals of A if A is symmetric or triangular and
          UPLO = 'U'.  0 <= KL <= M-1.
[in]KU
          KU is INTEGER
          Used only if A is a general band matrix or if A is
          triangular.

          If PATH = xGB, specifies the number of superdiagonals of A,
          and 0 <= KU <= N-1.

          If PATH = xTR, xTP, or xTB, specifies whether or not the
          matrix has unit diagonal:
          = 1:  matrix has non-unit diagonal (default)
          = 2:  matrix has unit diagonal
[in]NRHS
          NRHS is INTEGER
          The number of right hand side vectors in the system A*X = B.
[in]A
          A is COMPLEX*16 array, dimension (LDA,N)
          The test matrix whose type is given by PATH.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.
          If PATH = xGB, LDA >= KL+KU+1.
          If PATH = xPB, xSB, xHB, or xTB, LDA >= KL+1.
          Otherwise, LDA >= max(1,M).
[in,out]X
          X is or output) COMPLEX*16 array, dimension (LDX,NRHS)
          On entry, if XTYPE = 'C' (for 'Computed'), then X contains
          the exact solution to the system of linear equations.
          On exit, if XTYPE = 'N' (for 'New'), then X is initialized
          with random values.
[in]LDX
          LDX is INTEGER
          The leading dimension of the array X.  If TRANS = 'N',
          LDX >= max(1,N); if TRANS = 'T', LDX >= max(1,M).
[out]B
          B is COMPLEX*16 array, dimension (LDB,NRHS)
          The right hand side vector(s) for the system of equations,
          computed from B = op(A) * X, where op(A) is determined by
          TRANS.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  If TRANS = 'N',
          LDB >= max(1,M); if TRANS = 'T', LDB >= max(1,N).
[in,out]ISEED
          ISEED is INTEGER array, dimension (4)
          The seed vector for the random number generator (used in
          ZLATMS).  Modified on exit.
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -k, the k-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011

Definition at line 211 of file zlarhs.f.

211 *
212 * -- LAPACK test routine (version 3.4.0) --
213 * -- LAPACK is a software package provided by Univ. of Tennessee, --
214 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
215 * November 2011
216 *
217 * .. Scalar Arguments ..
218  CHARACTER trans, uplo, xtype
219  CHARACTER*3 path
220  INTEGER info, kl, ku, lda, ldb, ldx, m, n, nrhs
221 * ..
222 * .. Array Arguments ..
223  INTEGER iseed( 4 )
224  COMPLEX*16 a( lda, * ), b( ldb, * ), x( ldx, * )
225 * ..
226 *
227 * =====================================================================
228 *
229 * .. Parameters ..
230  COMPLEX*16 one, zero
231  parameter ( one = ( 1.0d+0, 0.0d+0 ),
232  $ zero = ( 0.0d+0, 0.0d+0 ) )
233 * ..
234 * .. Local Scalars ..
235  LOGICAL band, gen, notran, qrs, sym, tran, tri
236  CHARACTER c1, diag
237  CHARACTER*2 c2
238  INTEGER j, mb, nx
239 * ..
240 * .. External Functions ..
241  LOGICAL lsame, lsamen
242  EXTERNAL lsame, lsamen
243 * ..
244 * .. External Subroutines ..
245  EXTERNAL xerbla, zgbmv, zgemm, zhbmv, zhemm, zhpmv,
247  $ ztpmv, ztrmm
248 * ..
249 * .. Intrinsic Functions ..
250  INTRINSIC max
251 * ..
252 * .. Executable Statements ..
253 *
254 * Test the input parameters.
255 *
256  info = 0
257  c1 = path( 1: 1 )
258  c2 = path( 2: 3 )
259  tran = lsame( trans, 'T' ) .OR. lsame( trans, 'C' )
260  notran = .NOT.tran
261  gen = lsame( path( 2: 2 ), 'G' )
262  qrs = lsame( path( 2: 2 ), 'Q' ) .OR. lsame( path( 3: 3 ), 'Q' )
263  sym = lsame( path( 2: 2 ), 'P' ) .OR.
264  $ lsame( path( 2: 2 ), 'S' ) .OR. lsame( path( 2: 2 ), 'H' )
265  tri = lsame( path( 2: 2 ), 'T' )
266  band = lsame( path( 3: 3 ), 'B' )
267  IF( .NOT.lsame( c1, 'Zomplex precision' ) ) THEN
268  info = -1
269  ELSE IF( .NOT.( lsame( xtype, 'N' ) .OR. lsame( xtype, 'C' ) ) )
270  $ THEN
271  info = -2
272  ELSE IF( ( sym .OR. tri ) .AND. .NOT.
273  $ ( lsame( uplo, 'U' ) .OR. lsame( uplo, 'L' ) ) ) THEN
274  info = -3
275  ELSE IF( ( gen .OR. qrs ) .AND. .NOT.
276  $ ( tran .OR. lsame( trans, 'N' ) ) ) THEN
277  info = -4
278  ELSE IF( m.LT.0 ) THEN
279  info = -5
280  ELSE IF( n.LT.0 ) THEN
281  info = -6
282  ELSE IF( band .AND. kl.LT.0 ) THEN
283  info = -7
284  ELSE IF( band .AND. ku.LT.0 ) THEN
285  info = -8
286  ELSE IF( nrhs.LT.0 ) THEN
287  info = -9
288  ELSE IF( ( .NOT.band .AND. lda.LT.max( 1, m ) ) .OR.
289  $ ( band .AND. ( sym .OR. tri ) .AND. lda.LT.kl+1 ) .OR.
290  $ ( band .AND. gen .AND. lda.LT.kl+ku+1 ) ) THEN
291  info = -11
292  ELSE IF( ( notran .AND. ldx.LT.max( 1, n ) ) .OR.
293  $ ( tran .AND. ldx.LT.max( 1, m ) ) ) THEN
294  info = -13
295  ELSE IF( ( notran .AND. ldb.LT.max( 1, m ) ) .OR.
296  $ ( tran .AND. ldb.LT.max( 1, n ) ) ) THEN
297  info = -15
298  END IF
299  IF( info.NE.0 ) THEN
300  CALL xerbla( 'ZLARHS', -info )
301  RETURN
302  END IF
303 *
304 * Initialize X to NRHS random vectors unless XTYPE = 'C'.
305 *
306  IF( tran ) THEN
307  nx = m
308  mb = n
309  ELSE
310  nx = n
311  mb = m
312  END IF
313  IF( .NOT.lsame( xtype, 'C' ) ) THEN
314  DO 10 j = 1, nrhs
315  CALL zlarnv( 2, iseed, n, x( 1, j ) )
316  10 CONTINUE
317  END IF
318 *
319 * Multiply X by op( A ) using an appropriate
320 * matrix multiply routine.
321 *
322  IF( lsamen( 2, c2, 'GE' ) .OR. lsamen( 2, c2, 'QR' ) .OR.
323  $ lsamen( 2, c2, 'LQ' ) .OR. lsamen( 2, c2, 'QL' ) .OR.
324  $ lsamen( 2, c2, 'RQ' ) ) THEN
325 *
326 * General matrix
327 *
328  CALL zgemm( trans, 'N', mb, nrhs, nx, one, a, lda, x, ldx,
329  $ zero, b, ldb )
330 *
331  ELSE IF( lsamen( 2, c2, 'PO' ) .OR. lsamen( 2, c2, 'HE' ) ) THEN
332 *
333 * Hermitian matrix, 2-D storage
334 *
335  CALL zhemm( 'Left', uplo, n, nrhs, one, a, lda, x, ldx, zero,
336  $ b, ldb )
337 *
338  ELSE IF( lsamen( 2, c2, 'SY' ) ) THEN
339 *
340 * Symmetric matrix, 2-D storage
341 *
342  CALL zsymm( 'Left', uplo, n, nrhs, one, a, lda, x, ldx, zero,
343  $ b, ldb )
344 *
345  ELSE IF( lsamen( 2, c2, 'GB' ) ) THEN
346 *
347 * General matrix, band storage
348 *
349  DO 20 j = 1, nrhs
350  CALL zgbmv( trans, m, n, kl, ku, one, a, lda, x( 1, j ), 1,
351  $ zero, b( 1, j ), 1 )
352  20 CONTINUE
353 *
354  ELSE IF( lsamen( 2, c2, 'PB' ) .OR. lsamen( 2, c2, 'HB' ) ) THEN
355 *
356 * Hermitian matrix, band storage
357 *
358  DO 30 j = 1, nrhs
359  CALL zhbmv( uplo, n, kl, one, a, lda, x( 1, j ), 1, zero,
360  $ b( 1, j ), 1 )
361  30 CONTINUE
362 *
363  ELSE IF( lsamen( 2, c2, 'SB' ) ) THEN
364 *
365 * Symmetric matrix, band storage
366 *
367  DO 40 j = 1, nrhs
368  CALL zsbmv( uplo, n, kl, one, a, lda, x( 1, j ), 1, zero,
369  $ b( 1, j ), 1 )
370  40 CONTINUE
371 *
372  ELSE IF( lsamen( 2, c2, 'PP' ) .OR. lsamen( 2, c2, 'HP' ) ) THEN
373 *
374 * Hermitian matrix, packed storage
375 *
376  DO 50 j = 1, nrhs
377  CALL zhpmv( uplo, n, one, a, x( 1, j ), 1, zero, b( 1, j ),
378  $ 1 )
379  50 CONTINUE
380 *
381  ELSE IF( lsamen( 2, c2, 'SP' ) ) THEN
382 *
383 * Symmetric matrix, packed storage
384 *
385  DO 60 j = 1, nrhs
386  CALL zspmv( uplo, n, one, a, x( 1, j ), 1, zero, b( 1, j ),
387  $ 1 )
388  60 CONTINUE
389 *
390  ELSE IF( lsamen( 2, c2, 'TR' ) ) THEN
391 *
392 * Triangular matrix. Note that for triangular matrices,
393 * KU = 1 => non-unit triangular
394 * KU = 2 => unit triangular
395 *
396  CALL zlacpy( 'Full', n, nrhs, x, ldx, b, ldb )
397  IF( ku.EQ.2 ) THEN
398  diag = 'U'
399  ELSE
400  diag = 'N'
401  END IF
402  CALL ztrmm( 'Left', uplo, trans, diag, n, nrhs, one, a, lda, b,
403  $ ldb )
404 *
405  ELSE IF( lsamen( 2, c2, 'TP' ) ) THEN
406 *
407 * Triangular matrix, packed storage
408 *
409  CALL zlacpy( 'Full', n, nrhs, x, ldx, b, ldb )
410  IF( ku.EQ.2 ) THEN
411  diag = 'U'
412  ELSE
413  diag = 'N'
414  END IF
415  DO 70 j = 1, nrhs
416  CALL ztpmv( uplo, trans, diag, n, a, b( 1, j ), 1 )
417  70 CONTINUE
418 *
419  ELSE IF( lsamen( 2, c2, 'TB' ) ) THEN
420 *
421 * Triangular matrix, banded storage
422 *
423  CALL zlacpy( 'Full', n, nrhs, x, ldx, b, ldb )
424  IF( ku.EQ.2 ) THEN
425  diag = 'U'
426  ELSE
427  diag = 'N'
428  END IF
429  DO 80 j = 1, nrhs
430  CALL ztbmv( uplo, trans, diag, n, kl, a, lda, b( 1, j ), 1 )
431  80 CONTINUE
432 *
433  ELSE
434 *
435 * If none of the above, set INFO = -1 and return
436 *
437  info = -1
438  CALL xerbla( 'ZLARHS', -info )
439  END IF
440 *
441  RETURN
442 *
443 * End of ZLARHS
444 *
subroutine zlacpy(UPLO, M, N, A, LDA, B, LDB)
ZLACPY copies all or part of one two-dimensional array to another.
Definition: zlacpy.f:105
subroutine zlarnv(IDIST, ISEED, N, X)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: zlarnv.f:101
logical function lsamen(N, CA, CB)
LSAMEN
Definition: lsamen.f:76
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:153
subroutine zgbmv(TRANS, M, N, KL, KU, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGBMV
Definition: zgbmv.f:189
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:189
subroutine zsymm(SIDE, UPLO, M, N, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZSYMM
Definition: zsymm.f:191
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ztbmv(UPLO, TRANS, DIAG, N, K, A, LDA, X, INCX)
ZTBMV
Definition: ztbmv.f:188
subroutine ztrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
ZTRMM
Definition: ztrmm.f:179
subroutine zhemm(SIDE, UPLO, M, N, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZHEMM
Definition: zhemm.f:193
subroutine ztpmv(UPLO, TRANS, DIAG, N, AP, X, INCX)
ZTPMV
Definition: ztpmv.f:144
subroutine zhbmv(UPLO, N, K, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZHBMV
Definition: zhbmv.f:189
subroutine zsbmv(UPLO, N, K, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZSBMV
Definition: zsbmv.f:154
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zhpmv(UPLO, N, ALPHA, AP, X, INCX, BETA, Y, INCY)
ZHPMV
Definition: zhpmv.f:151

Here is the call graph for this function: