LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zdrvls ( logical, dimension( * )  DOTYPE,
integer  NM,
integer, dimension( * )  MVAL,
integer  NN,
integer, dimension( * )  NVAL,
integer  NNS,
integer, dimension( * )  NSVAL,
integer  NNB,
integer, dimension( * )  NBVAL,
integer, dimension( * )  NXVAL,
double precision  THRESH,
logical  TSTERR,
complex*16, dimension( * )  A,
complex*16, dimension( * )  COPYA,
complex*16, dimension( * )  B,
complex*16, dimension( * )  COPYB,
complex*16, dimension( * )  C,
double precision, dimension( * )  S,
double precision, dimension( * )  COPYS,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

ZDRVLS

Purpose:
 ZDRVLS tests the least squares driver routines ZGELS, CGELSS, ZGELSY
 and CGELSD.
Parameters
[in]DOTYPE
          DOTYPE is LOGICAL array, dimension (NTYPES)
          The matrix types to be used for testing.  Matrices of type j
          (for 1 <= j <= NTYPES) are used for testing if DOTYPE(j) =
          .TRUE.; if DOTYPE(j) = .FALSE., then type j is not used.
          The matrix of type j is generated as follows:
          j=1: A = U*D*V where U and V are random unitary matrices
               and D has random entries (> 0.1) taken from a uniform
               distribution (0,1). A is full rank.
          j=2: The same of 1, but A is scaled up.
          j=3: The same of 1, but A is scaled down.
          j=4: A = U*D*V where U and V are random unitary matrices
               and D has 3*min(M,N)/4 random entries (> 0.1) taken
               from a uniform distribution (0,1) and the remaining
               entries set to 0. A is rank-deficient.
          j=5: The same of 4, but A is scaled up.
          j=6: The same of 5, but A is scaled down.
[in]NM
          NM is INTEGER
          The number of values of M contained in the vector MVAL.
[in]MVAL
          MVAL is INTEGER array, dimension (NM)
          The values of the matrix row dimension M.
[in]NN
          NN is INTEGER
          The number of values of N contained in the vector NVAL.
[in]NVAL
          NVAL is INTEGER array, dimension (NN)
          The values of the matrix column dimension N.
[in]NNB
          NNB is INTEGER
          The number of values of NB and NX contained in the
          vectors NBVAL and NXVAL.  The blocking parameters are used
          in pairs (NB,NX).
[in]NBVAL
          NBVAL is INTEGER array, dimension (NNB)
          The values of the blocksize NB.
[in]NXVAL
          NXVAL is INTEGER array, dimension (NNB)
          The values of the crossover point NX.
[in]NNS
          NNS is INTEGER
          The number of values of NRHS contained in the vector NSVAL.
[in]NSVAL
          NSVAL is INTEGER array, dimension (NNS)
          The values of the number of right hand sides NRHS.
[in]THRESH
          THRESH is DOUBLE PRECISION
          The threshold value for the test ratios.  A result is
          included in the output file if RESULT >= THRESH.  To have
          every test ratio printed, use THRESH = 0.
[in]TSTERR
          TSTERR is LOGICAL
          Flag that indicates whether error exits are to be tested.
[out]A
          A is COMPLEX*16 array, dimension (MMAX*NMAX)
          where MMAX is the maximum value of M in MVAL and NMAX is the
          maximum value of N in NVAL.
[out]COPYA
          COPYA is COMPLEX*16 array, dimension (MMAX*NMAX)
[out]B
          B is COMPLEX*16 array, dimension (MMAX*NSMAX)
          where MMAX is the maximum value of M in MVAL and NSMAX is the
          maximum value of NRHS in NSVAL.
[out]COPYB
          COPYB is COMPLEX*16 array, dimension (MMAX*NSMAX)
[out]C
          C is COMPLEX*16 array, dimension (MMAX*NSMAX)
[out]S
          S is DOUBLE PRECISION array, dimension
                      (min(MMAX,NMAX))
[out]COPYS
          COPYS is DOUBLE PRECISION array, dimension
                      (min(MMAX,NMAX))
[out]WORK
          WORK is COMPLEX*16 array, dimension
                      (MMAX*NMAX + 4*NMAX + MMAX).
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (5*NMAX-1)
[out]IWORK
          IWORK is INTEGER array, dimension (15*NMAX)
[in]NOUT
          NOUT is INTEGER
          The unit number for output.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2015

Definition at line 211 of file zdrvls.f.

211 *
212 * -- LAPACK test routine (version 3.6.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 2015
216 *
217 * .. Scalar Arguments ..
218  LOGICAL tsterr
219  INTEGER nm, nn, nnb, nns, nout
220  DOUBLE PRECISION thresh
221 * ..
222 * .. Array Arguments ..
223  LOGICAL dotype( * )
224  INTEGER iwork( * ), mval( * ), nbval( * ), nsval( * ),
225  $ nval( * ), nxval( * )
226  DOUBLE PRECISION copys( * ), rwork( * ), s( * )
227  COMPLEX*16 a( * ), b( * ), c( * ), copya( * ), copyb( * ),
228  $ work( * )
229 * ..
230 *
231 * =====================================================================
232 *
233 * .. Parameters ..
234  INTEGER ntests
235  parameter ( ntests = 14 )
236  INTEGER smlsiz
237  parameter ( smlsiz = 25 )
238  DOUBLE PRECISION one, zero
239  parameter ( one = 1.0d+0, zero = 0.0d+0 )
240  COMPLEX*16 cone, czero
241  parameter ( cone = ( 1.0d+0, 0.0d+0 ),
242  $ czero = ( 0.0d+0, 0.0d+0 ) )
243 * ..
244 * .. Local Scalars ..
245  CHARACTER trans
246  CHARACTER*3 path
247  INTEGER crank, i, im, in, inb, info, ins, irank,
248  $ iscale, itran, itype, j, k, lda, ldb, ldwork,
249  $ lwlsy, lwork, m, mnmin, n, nb, ncols, nerrs,
250  $ nfail, nrhs, nrows, nrun, rank
251  DOUBLE PRECISION eps, norma, normb, rcond
252 * ..
253 * .. Local Arrays ..
254  INTEGER iseed( 4 ), iseedy( 4 )
255  DOUBLE PRECISION result( ntests )
256 * ..
257 * .. External Functions ..
258  DOUBLE PRECISION dasum, dlamch, zqrt12, zqrt14, zqrt17
259  EXTERNAL dasum, dlamch, zqrt12, zqrt14, zqrt17
260 * ..
261 * .. External Subroutines ..
262  EXTERNAL alaerh, alahd, alasvm, daxpy, dlasrt, xlaenv,
265  $ zqrt16
266 * ..
267 * .. Intrinsic Functions ..
268  INTRINSIC dble, max, min, sqrt
269 * ..
270 * .. Scalars in Common ..
271  LOGICAL lerr, ok
272  CHARACTER*32 srnamt
273  INTEGER infot, iounit
274 * ..
275 * .. Common blocks ..
276  COMMON / infoc / infot, iounit, ok, lerr
277  COMMON / srnamc / srnamt
278 * ..
279 * .. Data statements ..
280  DATA iseedy / 1988, 1989, 1990, 1991 /
281 * ..
282 * .. Executable Statements ..
283 *
284 * Initialize constants and the random number seed.
285 *
286  path( 1: 1 ) = 'Zomplex precision'
287  path( 2: 3 ) = 'LS'
288  nrun = 0
289  nfail = 0
290  nerrs = 0
291  DO 10 i = 1, 4
292  iseed( i ) = iseedy( i )
293  10 CONTINUE
294  eps = dlamch( 'Epsilon' )
295 *
296 * Threshold for rank estimation
297 *
298  rcond = sqrt( eps ) - ( sqrt( eps )-eps ) / 2
299 *
300 * Test the error exits
301 *
302  CALL xlaenv( 9, smlsiz )
303  IF( tsterr )
304  $ CALL zerrls( path, nout )
305 *
306 * Print the header if NM = 0 or NN = 0 and THRESH = 0.
307 *
308  IF( ( nm.EQ.0 .OR. nn.EQ.0 ) .AND. thresh.EQ.zero )
309  $ CALL alahd( nout, path )
310  infot = 0
311 *
312  DO 140 im = 1, nm
313  m = mval( im )
314  lda = max( 1, m )
315 *
316  DO 130 in = 1, nn
317  n = nval( in )
318  mnmin = min( m, n )
319  ldb = max( 1, m, n )
320 *
321  DO 120 ins = 1, nns
322  nrhs = nsval( ins )
323  lwork = max( 1, ( m+nrhs )*( n+2 ), ( n+nrhs )*( m+2 ),
324  $ m*n+4*mnmin+max( m, n ), 2*n+m )
325 *
326  DO 110 irank = 1, 2
327  DO 100 iscale = 1, 3
328  itype = ( irank-1 )*3 + iscale
329  IF( .NOT.dotype( itype ) )
330  $ GO TO 100
331 *
332  IF( irank.EQ.1 ) THEN
333 *
334 * Test ZGELS
335 *
336 * Generate a matrix of scaling type ISCALE
337 *
338  CALL zqrt13( iscale, m, n, copya, lda, norma,
339  $ iseed )
340  DO 40 inb = 1, nnb
341  nb = nbval( inb )
342  CALL xlaenv( 1, nb )
343  CALL xlaenv( 3, nxval( inb ) )
344 *
345  DO 30 itran = 1, 2
346  IF( itran.EQ.1 ) THEN
347  trans = 'N'
348  nrows = m
349  ncols = n
350  ELSE
351  trans = 'C'
352  nrows = n
353  ncols = m
354  END IF
355  ldwork = max( 1, ncols )
356 *
357 * Set up a consistent rhs
358 *
359  IF( ncols.GT.0 ) THEN
360  CALL zlarnv( 2, iseed, ncols*nrhs,
361  $ work )
362  CALL zdscal( ncols*nrhs,
363  $ one / dble( ncols ), work,
364  $ 1 )
365  END IF
366  CALL zgemm( trans, 'No transpose', nrows,
367  $ nrhs, ncols, cone, copya, lda,
368  $ work, ldwork, czero, b, ldb )
369  CALL zlacpy( 'Full', nrows, nrhs, b, ldb,
370  $ copyb, ldb )
371 *
372 * Solve LS or overdetermined system
373 *
374  IF( m.GT.0 .AND. n.GT.0 ) THEN
375  CALL zlacpy( 'Full', m, n, copya, lda,
376  $ a, lda )
377  CALL zlacpy( 'Full', nrows, nrhs,
378  $ copyb, ldb, b, ldb )
379  END IF
380  srnamt = 'ZGELS '
381  CALL zgels( trans, m, n, nrhs, a, lda, b,
382  $ ldb, work, lwork, info )
383 *
384  IF( info.NE.0 )
385  $ CALL alaerh( path, 'ZGELS ', info, 0,
386  $ trans, m, n, nrhs, -1, nb,
387  $ itype, nfail, nerrs,
388  $ nout )
389 *
390 * Check correctness of results
391 *
392  ldwork = max( 1, nrows )
393  IF( nrows.GT.0 .AND. nrhs.GT.0 )
394  $ CALL zlacpy( 'Full', nrows, nrhs,
395  $ copyb, ldb, c, ldb )
396  CALL zqrt16( trans, m, n, nrhs, copya,
397  $ lda, b, ldb, c, ldb, rwork,
398  $ result( 1 ) )
399 *
400  IF( ( itran.EQ.1 .AND. m.GE.n ) .OR.
401  $ ( itran.EQ.2 .AND. m.LT.n ) ) THEN
402 *
403 * Solving LS system
404 *
405  result( 2 ) = zqrt17( trans, 1, m, n,
406  $ nrhs, copya, lda, b, ldb,
407  $ copyb, ldb, c, work,
408  $ lwork )
409  ELSE
410 *
411 * Solving overdetermined system
412 *
413  result( 2 ) = zqrt14( trans, m, n,
414  $ nrhs, copya, lda, b, ldb,
415  $ work, lwork )
416  END IF
417 *
418 * Print information about the tests that
419 * did not pass the threshold.
420 *
421  DO 20 k = 1, 2
422  IF( result( k ).GE.thresh ) THEN
423  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
424  $ CALL alahd( nout, path )
425  WRITE( nout, fmt = 9999 )trans, m,
426  $ n, nrhs, nb, itype, k,
427  $ result( k )
428  nfail = nfail + 1
429  END IF
430  20 CONTINUE
431  nrun = nrun + 2
432  30 CONTINUE
433  40 CONTINUE
434  END IF
435 *
436 * Generate a matrix of scaling type ISCALE and rank
437 * type IRANK.
438 *
439  CALL zqrt15( iscale, irank, m, n, nrhs, copya, lda,
440  $ copyb, ldb, copys, rank, norma, normb,
441  $ iseed, work, lwork )
442 *
443 * workspace used: MAX(M+MIN(M,N),NRHS*MIN(M,N),2*N+M)
444 *
445  ldwork = max( 1, m )
446 *
447 * Loop for testing different block sizes.
448 *
449  DO 90 inb = 1, nnb
450  nb = nbval( inb )
451  CALL xlaenv( 1, nb )
452  CALL xlaenv( 3, nxval( inb ) )
453 *
454 * Test ZGELSY
455 *
456 * ZGELSY: Compute the minimum-norm solution
457 * X to min( norm( A * X - B ) )
458 * using the rank-revealing orthogonal
459 * factorization.
460 *
461  CALL zlacpy( 'Full', m, n, copya, lda, a, lda )
462  CALL zlacpy( 'Full', m, nrhs, copyb, ldb, b,
463  $ ldb )
464 *
465 * Initialize vector IWORK.
466 *
467  DO 70 j = 1, n
468  iwork( j ) = 0
469  70 CONTINUE
470 *
471 * Set LWLSY to the adequate value.
472 *
473  lwlsy = mnmin + max( 2*mnmin, nb*( n+1 ),
474  $ mnmin+nb*nrhs )
475  lwlsy = max( 1, lwlsy )
476 *
477  srnamt = 'ZGELSY'
478  CALL zgelsy( m, n, nrhs, a, lda, b, ldb, iwork,
479  $ rcond, crank, work, lwlsy, rwork,
480  $ info )
481  IF( info.NE.0 )
482  $ CALL alaerh( path, 'ZGELSY', info, 0, ' ', m,
483  $ n, nrhs, -1, nb, itype, nfail,
484  $ nerrs, nout )
485 *
486 * workspace used: 2*MNMIN+NB*NB+NB*MAX(N,NRHS)
487 *
488 * Test 3: Compute relative error in svd
489 * workspace: M*N + 4*MIN(M,N) + MAX(M,N)
490 *
491  result( 3 ) = zqrt12( crank, crank, a, lda,
492  $ copys, work, lwork, rwork )
493 *
494 * Test 4: Compute error in solution
495 * workspace: M*NRHS + M
496 *
497  CALL zlacpy( 'Full', m, nrhs, copyb, ldb, work,
498  $ ldwork )
499  CALL zqrt16( 'No transpose', m, n, nrhs, copya,
500  $ lda, b, ldb, work, ldwork, rwork,
501  $ result( 4 ) )
502 *
503 * Test 5: Check norm of r'*A
504 * workspace: NRHS*(M+N)
505 *
506  result( 5 ) = zero
507  IF( m.GT.crank )
508  $ result( 5 ) = zqrt17( 'No transpose', 1, m,
509  $ n, nrhs, copya, lda, b, ldb,
510  $ copyb, ldb, c, work, lwork )
511 *
512 * Test 6: Check if x is in the rowspace of A
513 * workspace: (M+NRHS)*(N+2)
514 *
515  result( 6 ) = zero
516 *
517  IF( n.GT.crank )
518  $ result( 6 ) = zqrt14( 'No transpose', m, n,
519  $ nrhs, copya, lda, b, ldb,
520  $ work, lwork )
521 *
522 * Test ZGELSS
523 *
524 * ZGELSS: Compute the minimum-norm solution
525 * X to min( norm( A * X - B ) )
526 * using the SVD.
527 *
528  CALL zlacpy( 'Full', m, n, copya, lda, a, lda )
529  CALL zlacpy( 'Full', m, nrhs, copyb, ldb, b,
530  $ ldb )
531  srnamt = 'ZGELSS'
532  CALL zgelss( m, n, nrhs, a, lda, b, ldb, s,
533  $ rcond, crank, work, lwork, rwork,
534  $ info )
535 *
536  IF( info.NE.0 )
537  $ CALL alaerh( path, 'ZGELSS', info, 0, ' ', m,
538  $ n, nrhs, -1, nb, itype, nfail,
539  $ nerrs, nout )
540 *
541 * workspace used: 3*min(m,n) +
542 * max(2*min(m,n),nrhs,max(m,n))
543 *
544 * Test 7: Compute relative error in svd
545 *
546  IF( rank.GT.0 ) THEN
547  CALL daxpy( mnmin, -one, copys, 1, s, 1 )
548  result( 7 ) = dasum( mnmin, s, 1 ) /
549  $ dasum( mnmin, copys, 1 ) /
550  $ ( eps*dble( mnmin ) )
551  ELSE
552  result( 7 ) = zero
553  END IF
554 *
555 * Test 8: Compute error in solution
556 *
557  CALL zlacpy( 'Full', m, nrhs, copyb, ldb, work,
558  $ ldwork )
559  CALL zqrt16( 'No transpose', m, n, nrhs, copya,
560  $ lda, b, ldb, work, ldwork, rwork,
561  $ result( 8 ) )
562 *
563 * Test 9: Check norm of r'*A
564 *
565  result( 9 ) = zero
566  IF( m.GT.crank )
567  $ result( 9 ) = zqrt17( 'No transpose', 1, m,
568  $ n, nrhs, copya, lda, b, ldb,
569  $ copyb, ldb, c, work, lwork )
570 *
571 * Test 10: Check if x is in the rowspace of A
572 *
573  result( 10 ) = zero
574  IF( n.GT.crank )
575  $ result( 10 ) = zqrt14( 'No transpose', m, n,
576  $ nrhs, copya, lda, b, ldb,
577  $ work, lwork )
578 *
579 * Test ZGELSD
580 *
581 * ZGELSD: Compute the minimum-norm solution X
582 * to min( norm( A * X - B ) ) using a
583 * divide and conquer SVD.
584 *
585  CALL xlaenv( 9, 25 )
586 *
587  CALL zlacpy( 'Full', m, n, copya, lda, a, lda )
588  CALL zlacpy( 'Full', m, nrhs, copyb, ldb, b,
589  $ ldb )
590 *
591  srnamt = 'ZGELSD'
592  CALL zgelsd( m, n, nrhs, a, lda, b, ldb, s,
593  $ rcond, crank, work, lwork, rwork,
594  $ iwork, info )
595  IF( info.NE.0 )
596  $ CALL alaerh( path, 'ZGELSD', info, 0, ' ', m,
597  $ n, nrhs, -1, nb, itype, nfail,
598  $ nerrs, nout )
599 *
600 * Test 11: Compute relative error in svd
601 *
602  IF( rank.GT.0 ) THEN
603  CALL daxpy( mnmin, -one, copys, 1, s, 1 )
604  result( 11 ) = dasum( mnmin, s, 1 ) /
605  $ dasum( mnmin, copys, 1 ) /
606  $ ( eps*dble( mnmin ) )
607  ELSE
608  result( 11 ) = zero
609  END IF
610 *
611 * Test 12: Compute error in solution
612 *
613  CALL zlacpy( 'Full', m, nrhs, copyb, ldb, work,
614  $ ldwork )
615  CALL zqrt16( 'No transpose', m, n, nrhs, copya,
616  $ lda, b, ldb, work, ldwork, rwork,
617  $ result( 12 ) )
618 *
619 * Test 13: Check norm of r'*A
620 *
621  result( 13 ) = zero
622  IF( m.GT.crank )
623  $ result( 13 ) = zqrt17( 'No transpose', 1, m,
624  $ n, nrhs, copya, lda, b, ldb,
625  $ copyb, ldb, c, work, lwork )
626 *
627 * Test 14: Check if x is in the rowspace of A
628 *
629  result( 14 ) = zero
630  IF( n.GT.crank )
631  $ result( 14 ) = zqrt14( 'No transpose', m, n,
632  $ nrhs, copya, lda, b, ldb,
633  $ work, lwork )
634 *
635 * Print information about the tests that did not
636 * pass the threshold.
637 *
638  DO 80 k = 3, ntests
639  IF( result( k ).GE.thresh ) THEN
640  IF( nfail.EQ.0 .AND. nerrs.EQ.0 )
641  $ CALL alahd( nout, path )
642  WRITE( nout, fmt = 9998 )m, n, nrhs, nb,
643  $ itype, k, result( k )
644  nfail = nfail + 1
645  END IF
646  80 CONTINUE
647  nrun = nrun + 12
648 *
649  90 CONTINUE
650  100 CONTINUE
651  110 CONTINUE
652  120 CONTINUE
653  130 CONTINUE
654  140 CONTINUE
655 *
656 * Print a summary of the results.
657 *
658  CALL alasvm( path, nout, nfail, nrun, nerrs )
659 *
660  9999 FORMAT( ' TRANS=''', a1, ''', M=', i5, ', N=', i5, ', NRHS=', i4,
661  $ ', NB=', i4, ', type', i2, ', test(', i2, ')=', g12.5 )
662  9998 FORMAT( ' M=', i5, ', N=', i5, ', NRHS=', i4, ', NB=', i4,
663  $ ', type', i2, ', test(', i2, ')=', g12.5 )
664  RETURN
665 *
666 * End of ZDRVLS
667 *
subroutine alasvm(TYPE, NOUT, NFAIL, NRUN, NERRS)
ALASVM
Definition: alasvm.f:75
subroutine dlasrt(ID, N, D, INFO)
DLASRT sorts numbers in increasing or decreasing order.
Definition: dlasrt.f:90
subroutine alahd(IOUNIT, PATH)
ALAHD
Definition: alahd.f:95
subroutine zqrt16(TRANS, M, N, NRHS, A, LDA, X, LDX, B, LDB, RWORK, RESID)
ZQRT16
Definition: zqrt16.f:135
subroutine zgels(TRANS, M, N, NRHS, A, LDA, B, LDB, WORK, LWORK, INFO)
ZGELS solves overdetermined or underdetermined systems for GE matrices
Definition: zgels.f:184
subroutine alaerh(PATH, SUBNAM, INFO, INFOE, OPTS, M, N, KL, KU, N5, IMAT, NFAIL, NERRS, NOUT)
ALAERH
Definition: alaerh.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:105
subroutine zlarnv(IDIST, ISEED, N, X)
ZLARNV returns a vector of random numbers from a uniform or normal distribution.
Definition: zlarnv.f:101
subroutine zgelsy(M, N, NRHS, A, LDA, B, LDB, JPVT, RCOND, RANK, WORK, LWORK, RWORK, INFO)
ZGELSY solves overdetermined or underdetermined systems for GE matrices
Definition: zgelsy.f:212
double precision function zqrt12(M, N, A, LDA, S, WORK, LWORK, RWORK)
ZQRT12
Definition: zqrt12.f:99
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine zgelss(M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, INFO)
ZGELSS solves overdetermined or underdetermined systems for GE matrices
Definition: zgelss.f:180
double precision function zqrt14(TRANS, M, N, NRHS, A, LDA, X, LDX, WORK, LWORK)
ZQRT14
Definition: zqrt14.f:118
subroutine zgelsd(M, N, NRHS, A, LDA, B, LDB, S, RCOND, RANK, WORK, LWORK, RWORK, IWORK, INFO)
ZGELSD computes the minimum-norm solution to a linear least squares problem for GE matrices ...
Definition: zgelsd.f:227
subroutine zgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
ZGEMM
Definition: zgemm.f:189
subroutine daxpy(N, DA, DX, INCX, DY, INCY)
DAXPY
Definition: daxpy.f:54
subroutine xlaenv(ISPEC, NVALUE)
XLAENV
Definition: xlaenv.f:83
subroutine zqrt13(SCALE, M, N, A, LDA, NORMA, ISEED)
ZQRT13
Definition: zqrt13.f:93
double precision function dasum(N, DX, INCX)
DASUM
Definition: dasum.f:53
subroutine zerrls(PATH, NUNIT)
ZERRLS
Definition: zerrls.f:57
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:54
double precision function zqrt17(TRANS, IRESID, M, N, NRHS, A, LDA, X, LDX, B, LDB, C, WORK, LWORK)
ZQRT17
Definition: zqrt17.f:152
subroutine zqrt15(SCALE, RKSEL, M, N, NRHS, A, LDA, B, LDB, S, RANK, NORMA, NORMB, ISEED, WORK, LWORK)
ZQRT15
Definition: zqrt15.f:151

Here is the call graph for this function:

Here is the caller graph for this function: