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

◆ clarhs()

subroutine clarhs ( character*3  path,
character  xtype,
character  uplo,
character  trans,
integer  m,
integer  n,
integer  kl,
integer  ku,
integer  nrhs,
complex, dimension( lda, * )  a,
integer  lda,
complex, dimension( ldx, * )  x,
integer  ldx,
complex, dimension( ldb, * )  b,
integer  ldb,
integer, dimension( 4 )  iseed,
integer  info 
)

CLARHS

Purpose:
 CLARHS 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) = A, A**T or A**H, depending on TRANS.
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  (No transpose)
          = 'T':  B := A**T * X  (Transpose)
          = 'C':  B := A**H * X  (Conjugate transpose)
[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 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 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 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
          CLATMS).  Modified on exit.
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -i, the i-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 206 of file clarhs.f.

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 A( LDA, * ), B( LDB, * ), X( LDX, * )
221* ..
222*
223* =====================================================================
224*
225* .. Parameters ..
226 COMPLEX ONE, ZERO
227 parameter( one = ( 1.0e+0, 0.0e+0 ),
228 $ zero = ( 0.0e+0, 0.0e+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 cgbmv, cgemm, chbmv, chemm, chpmv, clacpy,
243 $ ctrmm, xerbla
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, 'Complex 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.
272 $ .NOT.( 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( 'CLARHS', -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 clarnv( 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 cgemm( 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 chemm( '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 csymm( '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 cgbmv( 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 chbmv( 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 csbmv( 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 chpmv( 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 cspmv( 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 clacpy( '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 ctrmm( '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 clacpy( '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 ctpmv( 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 clacpy( '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 ctbmv( 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( 'CLARHS', -info )
435 END IF
436*
437 RETURN
438*
439* End of CLARHS
440*
subroutine csbmv(uplo, n, k, alpha, a, lda, x, incx, beta, y, incy)
CSBMV
Definition csbmv.f:152
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine cgbmv(trans, m, n, kl, ku, alpha, a, lda, x, incx, beta, y, incy)
CGBMV
Definition cgbmv.f:190
subroutine cgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
CGEMM
Definition cgemm.f:188
subroutine chbmv(uplo, n, k, alpha, a, lda, x, incx, beta, y, incy)
CHBMV
Definition chbmv.f:187
subroutine chemm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
CHEMM
Definition chemm.f:191
subroutine csymm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
CSYMM
Definition csymm.f:189
subroutine cspmv(uplo, n, alpha, ap, x, incx, beta, y, incy)
CSPMV computes a matrix-vector product for complex vectors using a complex symmetric packed matrix
Definition cspmv.f:151
subroutine chpmv(uplo, n, alpha, ap, x, incx, beta, y, incy)
CHPMV
Definition chpmv.f:149
subroutine clacpy(uplo, m, n, a, lda, b, ldb)
CLACPY copies all or part of one two-dimensional array to another.
Definition clacpy.f:103
subroutine clarnv(idist, iseed, n, x)
CLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition clarnv.f:99
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
logical function lsamen(n, ca, cb)
LSAMEN
Definition lsamen.f:74
subroutine ctbmv(uplo, trans, diag, n, k, a, lda, x, incx)
CTBMV
Definition ctbmv.f:186
subroutine ctpmv(uplo, trans, diag, n, ap, x, incx)
CTPMV
Definition ctpmv.f:142
subroutine ctrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
CTRMM
Definition ctrmm.f:177
Here is the call graph for this function:
Here is the caller graph for this function: