LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine slatrs ( character  UPLO,
character  TRANS,
character  DIAG,
character  NORMIN,
integer  N,
real, dimension( lda, * )  A,
integer  LDA,
real, dimension( * )  X,
real  SCALE,
real, dimension( * )  CNORM,
integer  INFO 
)

SLATRS solves a triangular system of equations with the scale factor set to prevent overflow.

Download SLATRS + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 SLATRS solves one of the triangular systems

    A *x = s*b  or  A**T*x = s*b

 with scaling to prevent overflow.  Here A is an upper or lower
 triangular matrix, A**T denotes the transpose of A, x and b are
 n-element vectors, and s is a scaling factor, usually less than
 or equal to 1, chosen so that the components of x will be less than
 the overflow threshold.  If the unscaled problem will not cause
 overflow, the Level 2 BLAS routine STRSV is called.  If the matrix A
 is singular (A(j,j) = 0 for some j), then s is set to 0 and a
 non-trivial solution to A*x = 0 is returned.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          Specifies whether the matrix A is upper or lower triangular.
          = 'U':  Upper triangular
          = 'L':  Lower triangular
[in]TRANS
          TRANS is CHARACTER*1
          Specifies the operation applied to A.
          = 'N':  Solve A * x = s*b  (No transpose)
          = 'T':  Solve A**T* x = s*b  (Transpose)
          = 'C':  Solve A**T* x = s*b  (Conjugate transpose = Transpose)
[in]DIAG
          DIAG is CHARACTER*1
          Specifies whether or not the matrix A is unit triangular.
          = 'N':  Non-unit triangular
          = 'U':  Unit triangular
[in]NORMIN
          NORMIN is CHARACTER*1
          Specifies whether CNORM has been set or not.
          = 'Y':  CNORM contains the column norms on entry
          = 'N':  CNORM is not set on entry.  On exit, the norms will
                  be computed and stored in CNORM.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]A
          A is REAL array, dimension (LDA,N)
          The triangular matrix A.  If UPLO = 'U', the leading n by n
          upper triangular part of the array A contains the upper
          triangular matrix, and the strictly lower triangular part of
          A is not referenced.  If UPLO = 'L', the leading n by n lower
          triangular part of the array A contains the lower triangular
          matrix, and the strictly upper triangular part of A is not
          referenced.  If DIAG = 'U', the diagonal elements of A are
          also not referenced and are assumed to be 1.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max (1,N).
[in,out]X
          X is REAL array, dimension (N)
          On entry, the right hand side b of the triangular system.
          On exit, X is overwritten by the solution vector x.
[out]SCALE
          SCALE is REAL
          The scaling factor s for the triangular system
             A * x = s*b  or  A**T* x = s*b.
          If SCALE = 0, the matrix A is singular or badly scaled, and
          the vector x is an exact or approximate solution to A*x = 0.
[in,out]CNORM
          CNORM is REAL array, dimension (N)

          If NORMIN = 'Y', CNORM is an input argument and CNORM(j)
          contains the norm of the off-diagonal part of the j-th column
          of A.  If TRANS = 'N', CNORM(j) must be greater than or equal
          to the infinity-norm, and if TRANS = 'T' or 'C', CNORM(j)
          must be greater than or equal to the 1-norm.

          If NORMIN = 'N', CNORM is an output argument and CNORM(j)
          returns the 1-norm of the offdiagonal part of the j-th column
          of A.
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -k, the k-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012
Further Details:
  A rough bound on x is computed; if that is less than overflow, STRSV
  is called, otherwise, specific code is used which checks for possible
  overflow or divide-by-zero at every operation.

  A columnwise scheme is used for solving A*x = b.  The basic algorithm
  if A is lower triangular is

       x[1:n] := b[1:n]
       for j = 1, ..., n
            x(j) := x(j) / A(j,j)
            x[j+1:n] := x[j+1:n] - x(j) * A[j+1:n,j]
       end

  Define bounds on the components of x after j iterations of the loop:
     M(j) = bound on x[1:j]
     G(j) = bound on x[j+1:n]
  Initially, let M(0) = 0 and G(0) = max{x(i), i=1,...,n}.

  Then for iteration j+1 we have
     M(j+1) <= G(j) / | A(j+1,j+1) |
     G(j+1) <= G(j) + M(j+1) * | A[j+2:n,j+1] |
            <= G(j) ( 1 + CNORM(j+1) / | A(j+1,j+1) | )

  where CNORM(j+1) is greater than or equal to the infinity-norm of
  column j+1 of A, not counting the diagonal.  Hence

     G(j) <= G(0) product ( 1 + CNORM(i) / | A(i,i) | )
                  1<=i<=j
  and

     |x(j)| <= ( G(0) / |A(j,j)| ) product ( 1 + CNORM(i) / |A(i,i)| )
                                   1<=i< j

  Since |x(j)| <= M(j), we use the Level 2 BLAS routine STRSV if the
  reciprocal of the largest M(j), j=1,..,n, is larger than
  max(underflow, 1/overflow).

  The bound on x(j) is also used to determine when a step in the
  columnwise method can be performed without fear of overflow.  If
  the computed bound is greater than a large constant, x is scaled to
  prevent overflow, but if the bound overflows, x is set to 0, x(j) to
  1, and scale to 0, and a non-trivial solution to A*x = 0 is found.

  Similarly, a row-wise scheme is used to solve A**T*x = b.  The basic
  algorithm for A upper triangular is

       for j = 1, ..., n
            x(j) := ( b(j) - A[1:j-1,j]**T * x[1:j-1] ) / A(j,j)
       end

  We simultaneously compute two bounds
       G(j) = bound on ( b(i) - A[1:i-1,i]**T * x[1:i-1] ), 1<=i<=j
       M(j) = bound on x(i), 1<=i<=j

  The initial values are G(0) = 0, M(0) = max{b(i), i=1,..,n}, and we
  add the constraint G(j) >= G(j-1) and M(j) >= M(j-1) for j >= 1.
  Then the bound on x(j) is

       M(j) <= M(j-1) * ( 1 + CNORM(j) ) / | A(j,j) |

            <= M(0) * product ( ( 1 + CNORM(i) ) / |A(i,i)| )
                      1<=i<=j

  and we can safely call STRSV if 1/M(n) and 1/G(n) are both greater
  than max(underflow, 1/overflow).

Definition at line 240 of file slatrs.f.

240 *
241 * -- LAPACK auxiliary routine (version 3.4.2) --
242 * -- LAPACK is a software package provided by Univ. of Tennessee, --
243 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
244 * September 2012
245 *
246 * .. Scalar Arguments ..
247  CHARACTER diag, normin, trans, uplo
248  INTEGER info, lda, n
249  REAL scale
250 * ..
251 * .. Array Arguments ..
252  REAL a( lda, * ), cnorm( * ), x( * )
253 * ..
254 *
255 * =====================================================================
256 *
257 * .. Parameters ..
258  REAL zero, half, one
259  parameter ( zero = 0.0e+0, half = 0.5e+0, one = 1.0e+0 )
260 * ..
261 * .. Local Scalars ..
262  LOGICAL notran, nounit, upper
263  INTEGER i, imax, j, jfirst, jinc, jlast
264  REAL bignum, grow, rec, smlnum, sumj, tjj, tjjs,
265  $ tmax, tscal, uscal, xbnd, xj, xmax
266 * ..
267 * .. External Functions ..
268  LOGICAL lsame
269  INTEGER isamax
270  REAL sasum, sdot, slamch
271  EXTERNAL lsame, isamax, sasum, sdot, slamch
272 * ..
273 * .. External Subroutines ..
274  EXTERNAL saxpy, sscal, strsv, xerbla
275 * ..
276 * .. Intrinsic Functions ..
277  INTRINSIC abs, max, min
278 * ..
279 * .. Executable Statements ..
280 *
281  info = 0
282  upper = lsame( uplo, 'U' )
283  notran = lsame( trans, 'N' )
284  nounit = lsame( diag, 'N' )
285 *
286 * Test the input parameters.
287 *
288  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
289  info = -1
290  ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
291  $ lsame( trans, 'C' ) ) THEN
292  info = -2
293  ELSE IF( .NOT.nounit .AND. .NOT.lsame( diag, 'U' ) ) THEN
294  info = -3
295  ELSE IF( .NOT.lsame( normin, 'Y' ) .AND. .NOT.
296  $ lsame( normin, 'N' ) ) THEN
297  info = -4
298  ELSE IF( n.LT.0 ) THEN
299  info = -5
300  ELSE IF( lda.LT.max( 1, n ) ) THEN
301  info = -7
302  END IF
303  IF( info.NE.0 ) THEN
304  CALL xerbla( 'SLATRS', -info )
305  RETURN
306  END IF
307 *
308 * Quick return if possible
309 *
310  IF( n.EQ.0 )
311  $ RETURN
312 *
313 * Determine machine dependent parameters to control overflow.
314 *
315  smlnum = slamch( 'Safe minimum' ) / slamch( 'Precision' )
316  bignum = one / smlnum
317  scale = one
318 *
319  IF( lsame( normin, 'N' ) ) THEN
320 *
321 * Compute the 1-norm of each column, not including the diagonal.
322 *
323  IF( upper ) THEN
324 *
325 * A is upper triangular.
326 *
327  DO 10 j = 1, n
328  cnorm( j ) = sasum( j-1, a( 1, j ), 1 )
329  10 CONTINUE
330  ELSE
331 *
332 * A is lower triangular.
333 *
334  DO 20 j = 1, n - 1
335  cnorm( j ) = sasum( n-j, a( j+1, j ), 1 )
336  20 CONTINUE
337  cnorm( n ) = zero
338  END IF
339  END IF
340 *
341 * Scale the column norms by TSCAL if the maximum element in CNORM is
342 * greater than BIGNUM.
343 *
344  imax = isamax( n, cnorm, 1 )
345  tmax = cnorm( imax )
346  IF( tmax.LE.bignum ) THEN
347  tscal = one
348  ELSE
349  tscal = one / ( smlnum*tmax )
350  CALL sscal( n, tscal, cnorm, 1 )
351  END IF
352 *
353 * Compute a bound on the computed solution vector to see if the
354 * Level 2 BLAS routine STRSV can be used.
355 *
356  j = isamax( n, x, 1 )
357  xmax = abs( x( j ) )
358  xbnd = xmax
359  IF( notran ) THEN
360 *
361 * Compute the growth in A * x = b.
362 *
363  IF( upper ) THEN
364  jfirst = n
365  jlast = 1
366  jinc = -1
367  ELSE
368  jfirst = 1
369  jlast = n
370  jinc = 1
371  END IF
372 *
373  IF( tscal.NE.one ) THEN
374  grow = zero
375  GO TO 50
376  END IF
377 *
378  IF( nounit ) THEN
379 *
380 * A is non-unit triangular.
381 *
382 * Compute GROW = 1/G(j) and XBND = 1/M(j).
383 * Initially, G(0) = max{x(i), i=1,...,n}.
384 *
385  grow = one / max( xbnd, smlnum )
386  xbnd = grow
387  DO 30 j = jfirst, jlast, jinc
388 *
389 * Exit the loop if the growth factor is too small.
390 *
391  IF( grow.LE.smlnum )
392  $ GO TO 50
393 *
394 * M(j) = G(j-1) / abs(A(j,j))
395 *
396  tjj = abs( a( j, j ) )
397  xbnd = min( xbnd, min( one, tjj )*grow )
398  IF( tjj+cnorm( j ).GE.smlnum ) THEN
399 *
400 * G(j) = G(j-1)*( 1 + CNORM(j) / abs(A(j,j)) )
401 *
402  grow = grow*( tjj / ( tjj+cnorm( j ) ) )
403  ELSE
404 *
405 * G(j) could overflow, set GROW to 0.
406 *
407  grow = zero
408  END IF
409  30 CONTINUE
410  grow = xbnd
411  ELSE
412 *
413 * A is unit triangular.
414 *
415 * Compute GROW = 1/G(j), where G(0) = max{x(i), i=1,...,n}.
416 *
417  grow = min( one, one / max( xbnd, smlnum ) )
418  DO 40 j = jfirst, jlast, jinc
419 *
420 * Exit the loop if the growth factor is too small.
421 *
422  IF( grow.LE.smlnum )
423  $ GO TO 50
424 *
425 * G(j) = G(j-1)*( 1 + CNORM(j) )
426 *
427  grow = grow*( one / ( one+cnorm( j ) ) )
428  40 CONTINUE
429  END IF
430  50 CONTINUE
431 *
432  ELSE
433 *
434 * Compute the growth in A**T * x = b.
435 *
436  IF( upper ) THEN
437  jfirst = 1
438  jlast = n
439  jinc = 1
440  ELSE
441  jfirst = n
442  jlast = 1
443  jinc = -1
444  END IF
445 *
446  IF( tscal.NE.one ) THEN
447  grow = zero
448  GO TO 80
449  END IF
450 *
451  IF( nounit ) THEN
452 *
453 * A is non-unit triangular.
454 *
455 * Compute GROW = 1/G(j) and XBND = 1/M(j).
456 * Initially, M(0) = max{x(i), i=1,...,n}.
457 *
458  grow = one / max( xbnd, smlnum )
459  xbnd = grow
460  DO 60 j = jfirst, jlast, jinc
461 *
462 * Exit the loop if the growth factor is too small.
463 *
464  IF( grow.LE.smlnum )
465  $ GO TO 80
466 *
467 * G(j) = max( G(j-1), M(j-1)*( 1 + CNORM(j) ) )
468 *
469  xj = one + cnorm( j )
470  grow = min( grow, xbnd / xj )
471 *
472 * M(j) = M(j-1)*( 1 + CNORM(j) ) / abs(A(j,j))
473 *
474  tjj = abs( a( j, j ) )
475  IF( xj.GT.tjj )
476  $ xbnd = xbnd*( tjj / xj )
477  60 CONTINUE
478  grow = min( grow, xbnd )
479  ELSE
480 *
481 * A is unit triangular.
482 *
483 * Compute GROW = 1/G(j), where G(0) = max{x(i), i=1,...,n}.
484 *
485  grow = min( one, one / max( xbnd, smlnum ) )
486  DO 70 j = jfirst, jlast, jinc
487 *
488 * Exit the loop if the growth factor is too small.
489 *
490  IF( grow.LE.smlnum )
491  $ GO TO 80
492 *
493 * G(j) = ( 1 + CNORM(j) )*G(j-1)
494 *
495  xj = one + cnorm( j )
496  grow = grow / xj
497  70 CONTINUE
498  END IF
499  80 CONTINUE
500  END IF
501 *
502  IF( ( grow*tscal ).GT.smlnum ) THEN
503 *
504 * Use the Level 2 BLAS solve if the reciprocal of the bound on
505 * elements of X is not too small.
506 *
507  CALL strsv( uplo, trans, diag, n, a, lda, x, 1 )
508  ELSE
509 *
510 * Use a Level 1 BLAS solve, scaling intermediate results.
511 *
512  IF( xmax.GT.bignum ) THEN
513 *
514 * Scale X so that its components are less than or equal to
515 * BIGNUM in absolute value.
516 *
517  scale = bignum / xmax
518  CALL sscal( n, scale, x, 1 )
519  xmax = bignum
520  END IF
521 *
522  IF( notran ) THEN
523 *
524 * Solve A * x = b
525 *
526  DO 100 j = jfirst, jlast, jinc
527 *
528 * Compute x(j) = b(j) / A(j,j), scaling x if necessary.
529 *
530  xj = abs( x( j ) )
531  IF( nounit ) THEN
532  tjjs = a( j, j )*tscal
533  ELSE
534  tjjs = tscal
535  IF( tscal.EQ.one )
536  $ GO TO 95
537  END IF
538  tjj = abs( tjjs )
539  IF( tjj.GT.smlnum ) THEN
540 *
541 * abs(A(j,j)) > SMLNUM:
542 *
543  IF( tjj.LT.one ) THEN
544  IF( xj.GT.tjj*bignum ) THEN
545 *
546 * Scale x by 1/b(j).
547 *
548  rec = one / xj
549  CALL sscal( n, rec, x, 1 )
550  scale = scale*rec
551  xmax = xmax*rec
552  END IF
553  END IF
554  x( j ) = x( j ) / tjjs
555  xj = abs( x( j ) )
556  ELSE IF( tjj.GT.zero ) THEN
557 *
558 * 0 < abs(A(j,j)) <= SMLNUM:
559 *
560  IF( xj.GT.tjj*bignum ) THEN
561 *
562 * Scale x by (1/abs(x(j)))*abs(A(j,j))*BIGNUM
563 * to avoid overflow when dividing by A(j,j).
564 *
565  rec = ( tjj*bignum ) / xj
566  IF( cnorm( j ).GT.one ) THEN
567 *
568 * Scale by 1/CNORM(j) to avoid overflow when
569 * multiplying x(j) times column j.
570 *
571  rec = rec / cnorm( j )
572  END IF
573  CALL sscal( n, rec, x, 1 )
574  scale = scale*rec
575  xmax = xmax*rec
576  END IF
577  x( j ) = x( j ) / tjjs
578  xj = abs( x( j ) )
579  ELSE
580 *
581 * A(j,j) = 0: Set x(1:n) = 0, x(j) = 1, and
582 * scale = 0, and compute a solution to A*x = 0.
583 *
584  DO 90 i = 1, n
585  x( i ) = zero
586  90 CONTINUE
587  x( j ) = one
588  xj = one
589  scale = zero
590  xmax = zero
591  END IF
592  95 CONTINUE
593 *
594 * Scale x if necessary to avoid overflow when adding a
595 * multiple of column j of A.
596 *
597  IF( xj.GT.one ) THEN
598  rec = one / xj
599  IF( cnorm( j ).GT.( bignum-xmax )*rec ) THEN
600 *
601 * Scale x by 1/(2*abs(x(j))).
602 *
603  rec = rec*half
604  CALL sscal( n, rec, x, 1 )
605  scale = scale*rec
606  END IF
607  ELSE IF( xj*cnorm( j ).GT.( bignum-xmax ) ) THEN
608 *
609 * Scale x by 1/2.
610 *
611  CALL sscal( n, half, x, 1 )
612  scale = scale*half
613  END IF
614 *
615  IF( upper ) THEN
616  IF( j.GT.1 ) THEN
617 *
618 * Compute the update
619 * x(1:j-1) := x(1:j-1) - x(j) * A(1:j-1,j)
620 *
621  CALL saxpy( j-1, -x( j )*tscal, a( 1, j ), 1, x,
622  $ 1 )
623  i = isamax( j-1, x, 1 )
624  xmax = abs( x( i ) )
625  END IF
626  ELSE
627  IF( j.LT.n ) THEN
628 *
629 * Compute the update
630 * x(j+1:n) := x(j+1:n) - x(j) * A(j+1:n,j)
631 *
632  CALL saxpy( n-j, -x( j )*tscal, a( j+1, j ), 1,
633  $ x( j+1 ), 1 )
634  i = j + isamax( n-j, x( j+1 ), 1 )
635  xmax = abs( x( i ) )
636  END IF
637  END IF
638  100 CONTINUE
639 *
640  ELSE
641 *
642 * Solve A**T * x = b
643 *
644  DO 140 j = jfirst, jlast, jinc
645 *
646 * Compute x(j) = b(j) - sum A(k,j)*x(k).
647 * k<>j
648 *
649  xj = abs( x( j ) )
650  uscal = tscal
651  rec = one / max( xmax, one )
652  IF( cnorm( j ).GT.( bignum-xj )*rec ) THEN
653 *
654 * If x(j) could overflow, scale x by 1/(2*XMAX).
655 *
656  rec = rec*half
657  IF( nounit ) THEN
658  tjjs = a( j, j )*tscal
659  ELSE
660  tjjs = tscal
661  END IF
662  tjj = abs( tjjs )
663  IF( tjj.GT.one ) THEN
664 *
665 * Divide by A(j,j) when scaling x if A(j,j) > 1.
666 *
667  rec = min( one, rec*tjj )
668  uscal = uscal / tjjs
669  END IF
670  IF( rec.LT.one ) THEN
671  CALL sscal( n, rec, x, 1 )
672  scale = scale*rec
673  xmax = xmax*rec
674  END IF
675  END IF
676 *
677  sumj = zero
678  IF( uscal.EQ.one ) THEN
679 *
680 * If the scaling needed for A in the dot product is 1,
681 * call SDOT to perform the dot product.
682 *
683  IF( upper ) THEN
684  sumj = sdot( j-1, a( 1, j ), 1, x, 1 )
685  ELSE IF( j.LT.n ) THEN
686  sumj = sdot( n-j, a( j+1, j ), 1, x( j+1 ), 1 )
687  END IF
688  ELSE
689 *
690 * Otherwise, use in-line code for the dot product.
691 *
692  IF( upper ) THEN
693  DO 110 i = 1, j - 1
694  sumj = sumj + ( a( i, j )*uscal )*x( i )
695  110 CONTINUE
696  ELSE IF( j.LT.n ) THEN
697  DO 120 i = j + 1, n
698  sumj = sumj + ( a( i, j )*uscal )*x( i )
699  120 CONTINUE
700  END IF
701  END IF
702 *
703  IF( uscal.EQ.tscal ) THEN
704 *
705 * Compute x(j) := ( x(j) - sumj ) / A(j,j) if 1/A(j,j)
706 * was not used to scale the dotproduct.
707 *
708  x( j ) = x( j ) - sumj
709  xj = abs( x( j ) )
710  IF( nounit ) THEN
711  tjjs = a( j, j )*tscal
712  ELSE
713  tjjs = tscal
714  IF( tscal.EQ.one )
715  $ GO TO 135
716  END IF
717 *
718 * Compute x(j) = x(j) / A(j,j), scaling if necessary.
719 *
720  tjj = abs( tjjs )
721  IF( tjj.GT.smlnum ) THEN
722 *
723 * abs(A(j,j)) > SMLNUM:
724 *
725  IF( tjj.LT.one ) THEN
726  IF( xj.GT.tjj*bignum ) THEN
727 *
728 * Scale X by 1/abs(x(j)).
729 *
730  rec = one / xj
731  CALL sscal( n, rec, x, 1 )
732  scale = scale*rec
733  xmax = xmax*rec
734  END IF
735  END IF
736  x( j ) = x( j ) / tjjs
737  ELSE IF( tjj.GT.zero ) THEN
738 *
739 * 0 < abs(A(j,j)) <= SMLNUM:
740 *
741  IF( xj.GT.tjj*bignum ) THEN
742 *
743 * Scale x by (1/abs(x(j)))*abs(A(j,j))*BIGNUM.
744 *
745  rec = ( tjj*bignum ) / xj
746  CALL sscal( n, rec, x, 1 )
747  scale = scale*rec
748  xmax = xmax*rec
749  END IF
750  x( j ) = x( j ) / tjjs
751  ELSE
752 *
753 * A(j,j) = 0: Set x(1:n) = 0, x(j) = 1, and
754 * scale = 0, and compute a solution to A**T*x = 0.
755 *
756  DO 130 i = 1, n
757  x( i ) = zero
758  130 CONTINUE
759  x( j ) = one
760  scale = zero
761  xmax = zero
762  END IF
763  135 CONTINUE
764  ELSE
765 *
766 * Compute x(j) := x(j) / A(j,j) - sumj if the dot
767 * product has already been divided by 1/A(j,j).
768 *
769  x( j ) = x( j ) / tjjs - sumj
770  END IF
771  xmax = max( xmax, abs( x( j ) ) )
772  140 CONTINUE
773  END IF
774  scale = scale / tscal
775  END IF
776 *
777 * Scale the column norms by 1/TSCAL for return.
778 *
779  IF( tscal.NE.one ) THEN
780  CALL sscal( n, one / tscal, cnorm, 1 )
781  END IF
782 *
783  RETURN
784 *
785 * End of SLATRS
786 *
real function sdot(N, SX, INCX, SY, INCY)
SDOT
Definition: sdot.f:53
integer function isamax(N, SX, INCX)
ISAMAX
Definition: isamax.f:53
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
real function sasum(N, SX, INCX)
SASUM
Definition: sasum.f:54
subroutine strsv(UPLO, TRANS, DIAG, N, A, LDA, X, INCX)
STRSV
Definition: strsv.f:151
subroutine saxpy(N, SA, SX, INCX, SY, INCY)
SAXPY
Definition: saxpy.f:54
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:55
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
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: