LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine sdrvls ( 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,
real  THRESH,
logical  TSTERR,
real, dimension( * )  A,
real, dimension( * )  COPYA,
real, dimension( * )  B,
real, dimension( * )  COPYB,
real, dimension( * )  C,
real, dimension( * )  S,
real, dimension( * )  COPYS,
real, dimension( * )  WORK,
integer, dimension( * )  IWORK,
integer  NOUT 
)

SDRVLS

Purpose:
 SDRVLS tests the least squares driver routines SGELS, SGELSS, SGELSY
 and SGELSD.
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 orthogonal 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 orthogonal 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]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]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]THRESH
          THRESH is REAL
          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 REAL 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 REAL array, dimension (MMAX*NMAX)
[out]B
          B is REAL 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 REAL array, dimension (MMAX*NSMAX)
[out]C
          C is REAL array, dimension (MMAX*NSMAX)
[out]S
          S is REAL array, dimension
                      (min(MMAX,NMAX))
[out]COPYS
          COPYS is REAL array, dimension
                      (min(MMAX,NMAX))
[out]WORK
          WORK is REAL array,
                      dimension (MMAX*NMAX + 4*NMAX + MMAX).
[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 205 of file sdrvls.f.

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

Here is the call graph for this function:

Here is the caller graph for this function: