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

◆ dlarhs()

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) = A or A**T, depending on TRANS.
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
          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 = Transpose)
[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.

Definition at line 203 of file dlarhs.f.

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