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

DLARHS

Purpose:
 DLARHS 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 or A' (transpose of A).
Parameters
[in]PATH
          PATH is CHARACTER*3
          The type of the real matrix A.  PATH may be given in any
          combination of upper and lower case.  Valid types include
             xGE:  General m x n matrix
             xGB:  General banded matrix
             xPO:  Symmetric positive definite, 2-D storage
             xPP:  Symmetric positive definite packed
             xPB:  Symmetric positive definite 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
          Specifies whether the upper or lower triangular part of the
          matrix A is stored, if A is symmetric.
          = 'U':  Upper triangular
          = 'L':  Lower triangular
[in]TRANS
          TRANS is CHARACTER*1
          Specifies the operation applied to the matrix A.
          = 'N':  System is  A * x = b
          = 'T':  System is  A'* x = b
          = 'C':  System is  A'* x = b
[in]M
          M is INTEGER
          The number or 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 DOUBLE PRECISION 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) DOUBLE PRECISION 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 DOUBLE PRECISION 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
          DLATMS).  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.
Date
November 2011

Definition at line 206 of file dlarhs.f.

206 *
207 * -- LAPACK test routine (version 3.4.0) --
208 * -- LAPACK is a software package provided by Univ. of Tennessee, --
209 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
210 * November 2011
211 *
212 * .. Scalar Arguments ..
213  CHARACTER trans, uplo, xtype
214  CHARACTER*3 path
215  INTEGER info, kl, ku, lda, ldb, ldx, m, n, nrhs
216 * ..
217 * .. Array Arguments ..
218  INTEGER iseed( 4 )
219  DOUBLE PRECISION a( lda, * ), b( ldb, * ), x( ldx, * )
220 * ..
221 *
222 * =====================================================================
223 *
224 * .. Parameters ..
225  DOUBLE PRECISION one, zero
226  parameter ( one = 1.0d+0, zero = 0.0d+0 )
227 * ..
228 * .. Local Scalars ..
229  LOGICAL band, gen, notran, qrs, sym, tran, tri
230  CHARACTER c1, diag
231  CHARACTER*2 c2
232  INTEGER j, mb, nx
233 * ..
234 * .. External Functions ..
235  LOGICAL lsame, lsamen
236  EXTERNAL lsame, lsamen
237 * ..
238 * .. External Subroutines ..
239  EXTERNAL dgbmv, dgemm, dlacpy, dlarnv, dsbmv, dspmv,
241 * ..
242 * .. Intrinsic Functions ..
243  INTRINSIC max
244 * ..
245 * .. Executable Statements ..
246 *
247 * Test the input parameters.
248 *
249  info = 0
250  c1 = path( 1: 1 )
251  c2 = path( 2: 3 )
252  tran = lsame( trans, 'T' ) .OR. lsame( trans, 'C' )
253  notran = .NOT.tran
254  gen = lsame( path( 2: 2 ), 'G' )
255  qrs = lsame( path( 2: 2 ), 'Q' ) .OR. lsame( path( 3: 3 ), 'Q' )
256  sym = lsame( path( 2: 2 ), 'P' ) .OR. lsame( path( 2: 2 ), 'S' )
257  tri = lsame( path( 2: 2 ), 'T' )
258  band = lsame( path( 3: 3 ), 'B' )
259  IF( .NOT.lsame( c1, 'Double precision' ) ) THEN
260  info = -1
261  ELSE IF( .NOT.( lsame( xtype, 'N' ) .OR. lsame( xtype, 'C' ) ) )
262  $ THEN
263  info = -2
264  ELSE IF( ( sym .OR. tri ) .AND. .NOT.
265  $ ( lsame( uplo, 'U' ) .OR. lsame( uplo, 'L' ) ) ) THEN
266  info = -3
267  ELSE IF( ( gen .OR. qrs ) .AND. .NOT.
268  $ ( tran .OR. lsame( trans, 'N' ) ) ) THEN
269  info = -4
270  ELSE IF( m.LT.0 ) THEN
271  info = -5
272  ELSE IF( n.LT.0 ) THEN
273  info = -6
274  ELSE IF( band .AND. kl.LT.0 ) THEN
275  info = -7
276  ELSE IF( band .AND. ku.LT.0 ) THEN
277  info = -8
278  ELSE IF( nrhs.LT.0 ) THEN
279  info = -9
280  ELSE IF( ( .NOT.band .AND. lda.LT.max( 1, m ) ) .OR.
281  $ ( band .AND. ( sym .OR. tri ) .AND. lda.LT.kl+1 ) .OR.
282  $ ( band .AND. gen .AND. lda.LT.kl+ku+1 ) ) THEN
283  info = -11
284  ELSE IF( ( notran .AND. ldx.LT.max( 1, n ) ) .OR.
285  $ ( tran .AND. ldx.LT.max( 1, m ) ) ) THEN
286  info = -13
287  ELSE IF( ( notran .AND. ldb.LT.max( 1, m ) ) .OR.
288  $ ( tran .AND. ldb.LT.max( 1, n ) ) ) THEN
289  info = -15
290  END IF
291  IF( info.NE.0 ) THEN
292  CALL xerbla( 'DLARHS', -info )
293  RETURN
294  END IF
295 *
296 * Initialize X to NRHS random vectors unless XTYPE = 'C'.
297 *
298  IF( tran ) THEN
299  nx = m
300  mb = n
301  ELSE
302  nx = n
303  mb = m
304  END IF
305  IF( .NOT.lsame( xtype, 'C' ) ) THEN
306  DO 10 j = 1, nrhs
307  CALL dlarnv( 2, iseed, n, x( 1, j ) )
308  10 CONTINUE
309  END IF
310 *
311 * Multiply X by op( A ) using an appropriate
312 * matrix multiply routine.
313 *
314  IF( lsamen( 2, c2, 'GE' ) .OR. lsamen( 2, c2, 'QR' ) .OR.
315  $ lsamen( 2, c2, 'LQ' ) .OR. lsamen( 2, c2, 'QL' ) .OR.
316  $ lsamen( 2, c2, 'RQ' ) ) THEN
317 *
318 * General matrix
319 *
320  CALL dgemm( trans, 'N', mb, nrhs, nx, one, a, lda, x, ldx,
321  $ zero, b, ldb )
322 *
323  ELSE IF( lsamen( 2, c2, 'PO' ) .OR. lsamen( 2, c2, 'SY' ) ) THEN
324 *
325 * Symmetric matrix, 2-D storage
326 *
327  CALL dsymm( 'Left', uplo, n, nrhs, one, a, lda, x, ldx, zero,
328  $ b, ldb )
329 *
330  ELSE IF( lsamen( 2, c2, 'GB' ) ) THEN
331 *
332 * General matrix, band storage
333 *
334  DO 20 j = 1, nrhs
335  CALL dgbmv( trans, mb, nx, kl, ku, one, a, lda, x( 1, j ),
336  $ 1, zero, b( 1, j ), 1 )
337  20 CONTINUE
338 *
339  ELSE IF( lsamen( 2, c2, 'PB' ) ) THEN
340 *
341 * Symmetric matrix, band storage
342 *
343  DO 30 j = 1, nrhs
344  CALL dsbmv( uplo, n, kl, one, a, lda, x( 1, j ), 1, zero,
345  $ b( 1, j ), 1 )
346  30 CONTINUE
347 *
348  ELSE IF( lsamen( 2, c2, 'PP' ) .OR. lsamen( 2, c2, 'SP' ) ) THEN
349 *
350 * Symmetric matrix, packed storage
351 *
352  DO 40 j = 1, nrhs
353  CALL dspmv( uplo, n, one, a, x( 1, j ), 1, zero, b( 1, j ),
354  $ 1 )
355  40 CONTINUE
356 *
357  ELSE IF( lsamen( 2, c2, 'TR' ) ) THEN
358 *
359 * Triangular matrix. Note that for triangular matrices,
360 * KU = 1 => non-unit triangular
361 * KU = 2 => unit triangular
362 *
363  CALL dlacpy( 'Full', n, nrhs, x, ldx, b, ldb )
364  IF( ku.EQ.2 ) THEN
365  diag = 'U'
366  ELSE
367  diag = 'N'
368  END IF
369  CALL dtrmm( 'Left', uplo, trans, diag, n, nrhs, one, a, lda, b,
370  $ ldb )
371 *
372  ELSE IF( lsamen( 2, c2, 'TP' ) ) THEN
373 *
374 * Triangular matrix, packed storage
375 *
376  CALL dlacpy( 'Full', n, nrhs, x, ldx, b, ldb )
377  IF( ku.EQ.2 ) THEN
378  diag = 'U'
379  ELSE
380  diag = 'N'
381  END IF
382  DO 50 j = 1, nrhs
383  CALL dtpmv( uplo, trans, diag, n, a, b( 1, j ), 1 )
384  50 CONTINUE
385 *
386  ELSE IF( lsamen( 2, c2, 'TB' ) ) THEN
387 *
388 * Triangular matrix, banded storage
389 *
390  CALL dlacpy( 'Full', n, nrhs, x, ldx, b, ldb )
391  IF( ku.EQ.2 ) THEN
392  diag = 'U'
393  ELSE
394  diag = 'N'
395  END IF
396  DO 60 j = 1, nrhs
397  CALL dtbmv( uplo, trans, diag, n, kl, a, lda, b( 1, j ), 1 )
398  60 CONTINUE
399 *
400  ELSE
401 *
402 * If PATH is none of the above, return with an error code.
403 *
404  info = -1
405  CALL xerbla( 'DLARHS', -info )
406  END IF
407 *
408  RETURN
409 *
410 * End of DLARHS
411 *
logical function lsamen(N, CA, CB)
LSAMEN
Definition: lsamen.f:76
subroutine dsymm(SIDE, UPLO, M, N, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DSYMM
Definition: dsymm.f:191
subroutine dtrmm(SIDE, UPLO, TRANSA, DIAG, M, N, ALPHA, A, LDA, B, LDB)
DTRMM
Definition: dtrmm.f:179
subroutine dlacpy(UPLO, M, N, A, LDA, B, LDB)
DLACPY copies all or part of one two-dimensional array to another.
Definition: dlacpy.f:105
subroutine dspmv(UPLO, N, ALPHA, AP, X, INCX, BETA, Y, INCY)
DSPMV
Definition: dspmv.f:149
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dgbmv(TRANS, M, N, KL, KU, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGBMV
Definition: dgbmv.f:187
subroutine dsbmv(UPLO, N, K, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DSBMV
Definition: dsbmv.f:186
subroutine dtbmv(UPLO, TRANS, DIAG, N, K, A, LDA, X, INCX)
DTBMV
Definition: dtbmv.f:188
subroutine dtpmv(UPLO, TRANS, DIAG, N, AP, X, INCX)
DTPMV
Definition: dtpmv.f:144
subroutine dlarnv(IDIST, ISEED, N, X)
DLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: dlarnv.f:99
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: