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

SLATPS solves a triangular system of equations with the matrix held in packed storage.

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

Purpose:
 SLATPS solves one of the triangular systems

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

 with scaling to prevent overflow, where A is an upper or lower
 triangular matrix stored in packed form.  Here 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
 STPSV 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]AP
          AP is REAL array, dimension (N*(N+1)/2)
          The upper or lower triangular matrix A, packed columnwise in
          a linear array.  The j-th column of A is stored in the array
          AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=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, STPSV
  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 STPSV 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 STPSV if 1/M(n) and 1/G(n) are both greater
  than max(underflow, 1/overflow).

Definition at line 231 of file slatps.f.

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