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

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

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

Purpose:
 CLATPS solves one of the triangular systems

    A * x = s*b,  A**T * x = s*b,  or  A**H * 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, A**H denotes the conjugate 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 CTPSV 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**H * x = s*b  (Conjugate 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 COMPLEX 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 COMPLEX 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,  A**T * x = s*b,  or  A**H * 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, CTPSV
  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 CTPSV 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  or
  A**H *x = b.  The basic algorithm for A upper triangular is

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

  We simultaneously compute two bounds
       G(j) = bound on ( b(i) - A[1:i-1,i]' * 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 CTPSV if 1/M(n) and 1/G(n) are both greater
  than max(underflow, 1/overflow).

Definition at line 233 of file clatps.f.

233 *
234 * -- LAPACK auxiliary routine (version 3.4.2) --
235 * -- LAPACK is a software package provided by Univ. of Tennessee, --
236 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
237 * September 2012
238 *
239 * .. Scalar Arguments ..
240  CHARACTER diag, normin, trans, uplo
241  INTEGER info, n
242  REAL scale
243 * ..
244 * .. Array Arguments ..
245  REAL cnorm( * )
246  COMPLEX ap( * ), x( * )
247 * ..
248 *
249 * =====================================================================
250 *
251 * .. Parameters ..
252  REAL zero, half, one, two
253  parameter ( zero = 0.0e+0, half = 0.5e+0, one = 1.0e+0,
254  $ two = 2.0e+0 )
255 * ..
256 * .. Local Scalars ..
257  LOGICAL notran, nounit, upper
258  INTEGER i, imax, ip, j, jfirst, jinc, jlast, jlen
259  REAL bignum, grow, rec, smlnum, tjj, tmax, tscal,
260  $ xbnd, xj, xmax
261  COMPLEX csumj, tjjs, uscal, zdum
262 * ..
263 * .. External Functions ..
264  LOGICAL lsame
265  INTEGER icamax, isamax
266  REAL scasum, slamch
267  COMPLEX cdotc, cdotu, cladiv
268  EXTERNAL lsame, icamax, isamax, scasum, slamch, cdotc,
269  $ cdotu, cladiv
270 * ..
271 * .. External Subroutines ..
272  EXTERNAL caxpy, csscal, ctpsv, slabad, sscal, xerbla
273 * ..
274 * .. Intrinsic Functions ..
275  INTRINSIC abs, aimag, cmplx, conjg, max, min, real
276 * ..
277 * .. Statement Functions ..
278  REAL cabs1, cabs2
279 * ..
280 * .. Statement Function definitions ..
281  cabs1( zdum ) = abs( REAL( ZDUM ) ) + abs( aimag( zdum ) )
282  cabs2( zdum ) = abs( REAL( ZDUM ) / 2. ) +
283  $ abs( aimag( zdum ) / 2. )
284 * ..
285 * .. Executable Statements ..
286 *
287  info = 0
288  upper = lsame( uplo, 'U' )
289  notran = lsame( trans, 'N' )
290  nounit = lsame( diag, 'N' )
291 *
292 * Test the input parameters.
293 *
294  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
295  info = -1
296  ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
297  $ lsame( trans, 'C' ) ) THEN
298  info = -2
299  ELSE IF( .NOT.nounit .AND. .NOT.lsame( diag, 'U' ) ) THEN
300  info = -3
301  ELSE IF( .NOT.lsame( normin, 'Y' ) .AND. .NOT.
302  $ lsame( normin, 'N' ) ) THEN
303  info = -4
304  ELSE IF( n.LT.0 ) THEN
305  info = -5
306  END IF
307  IF( info.NE.0 ) THEN
308  CALL xerbla( 'CLATPS', -info )
309  RETURN
310  END IF
311 *
312 * Quick return if possible
313 *
314  IF( n.EQ.0 )
315  $ RETURN
316 *
317 * Determine machine dependent parameters to control overflow.
318 *
319  smlnum = slamch( 'Safe minimum' )
320  bignum = one / smlnum
321  CALL slabad( smlnum, bignum )
322  smlnum = smlnum / slamch( 'Precision' )
323  bignum = one / smlnum
324  scale = one
325 *
326  IF( lsame( normin, 'N' ) ) THEN
327 *
328 * Compute the 1-norm of each column, not including the diagonal.
329 *
330  IF( upper ) THEN
331 *
332 * A is upper triangular.
333 *
334  ip = 1
335  DO 10 j = 1, n
336  cnorm( j ) = scasum( j-1, ap( ip ), 1 )
337  ip = ip + j
338  10 CONTINUE
339  ELSE
340 *
341 * A is lower triangular.
342 *
343  ip = 1
344  DO 20 j = 1, n - 1
345  cnorm( j ) = scasum( n-j, ap( ip+1 ), 1 )
346  ip = ip + n - j + 1
347  20 CONTINUE
348  cnorm( n ) = zero
349  END IF
350  END IF
351 *
352 * Scale the column norms by TSCAL if the maximum element in CNORM is
353 * greater than BIGNUM/2.
354 *
355  imax = isamax( n, cnorm, 1 )
356  tmax = cnorm( imax )
357  IF( tmax.LE.bignum*half ) THEN
358  tscal = one
359  ELSE
360  tscal = half / ( smlnum*tmax )
361  CALL sscal( n, tscal, cnorm, 1 )
362  END IF
363 *
364 * Compute a bound on the computed solution vector to see if the
365 * Level 2 BLAS routine CTPSV can be used.
366 *
367  xmax = zero
368  DO 30 j = 1, n
369  xmax = max( xmax, cabs2( x( j ) ) )
370  30 CONTINUE
371  xbnd = xmax
372  IF( notran ) THEN
373 *
374 * Compute the growth in A * x = b.
375 *
376  IF( upper ) THEN
377  jfirst = n
378  jlast = 1
379  jinc = -1
380  ELSE
381  jfirst = 1
382  jlast = n
383  jinc = 1
384  END IF
385 *
386  IF( tscal.NE.one ) THEN
387  grow = zero
388  GO TO 60
389  END IF
390 *
391  IF( nounit ) THEN
392 *
393 * A is non-unit triangular.
394 *
395 * Compute GROW = 1/G(j) and XBND = 1/M(j).
396 * Initially, G(0) = max{x(i), i=1,...,n}.
397 *
398  grow = half / max( xbnd, smlnum )
399  xbnd = grow
400  ip = jfirst*( jfirst+1 ) / 2
401  jlen = n
402  DO 40 j = jfirst, jlast, jinc
403 *
404 * Exit the loop if the growth factor is too small.
405 *
406  IF( grow.LE.smlnum )
407  $ GO TO 60
408 *
409  tjjs = ap( ip )
410  tjj = cabs1( tjjs )
411 *
412  IF( tjj.GE.smlnum ) THEN
413 *
414 * M(j) = G(j-1) / abs(A(j,j))
415 *
416  xbnd = min( xbnd, min( one, tjj )*grow )
417  ELSE
418 *
419 * M(j) could overflow, set XBND to 0.
420 *
421  xbnd = zero
422  END IF
423 *
424  IF( tjj+cnorm( j ).GE.smlnum ) THEN
425 *
426 * G(j) = G(j-1)*( 1 + CNORM(j) / abs(A(j,j)) )
427 *
428  grow = grow*( tjj / ( tjj+cnorm( j ) ) )
429  ELSE
430 *
431 * G(j) could overflow, set GROW to 0.
432 *
433  grow = zero
434  END IF
435  ip = ip + jinc*jlen
436  jlen = jlen - 1
437  40 CONTINUE
438  grow = xbnd
439  ELSE
440 *
441 * A is unit triangular.
442 *
443 * Compute GROW = 1/G(j), where G(0) = max{x(i), i=1,...,n}.
444 *
445  grow = min( one, half / max( xbnd, smlnum ) )
446  DO 50 j = jfirst, jlast, jinc
447 *
448 * Exit the loop if the growth factor is too small.
449 *
450  IF( grow.LE.smlnum )
451  $ GO TO 60
452 *
453 * G(j) = G(j-1)*( 1 + CNORM(j) )
454 *
455  grow = grow*( one / ( one+cnorm( j ) ) )
456  50 CONTINUE
457  END IF
458  60 CONTINUE
459 *
460  ELSE
461 *
462 * Compute the growth in A**T * x = b or A**H * x = b.
463 *
464  IF( upper ) THEN
465  jfirst = 1
466  jlast = n
467  jinc = 1
468  ELSE
469  jfirst = n
470  jlast = 1
471  jinc = -1
472  END IF
473 *
474  IF( tscal.NE.one ) THEN
475  grow = zero
476  GO TO 90
477  END IF
478 *
479  IF( nounit ) THEN
480 *
481 * A is non-unit triangular.
482 *
483 * Compute GROW = 1/G(j) and XBND = 1/M(j).
484 * Initially, M(0) = max{x(i), i=1,...,n}.
485 *
486  grow = half / max( xbnd, smlnum )
487  xbnd = grow
488  ip = jfirst*( jfirst+1 ) / 2
489  jlen = 1
490  DO 70 j = jfirst, jlast, jinc
491 *
492 * Exit the loop if the growth factor is too small.
493 *
494  IF( grow.LE.smlnum )
495  $ GO TO 90
496 *
497 * G(j) = max( G(j-1), M(j-1)*( 1 + CNORM(j) ) )
498 *
499  xj = one + cnorm( j )
500  grow = min( grow, xbnd / xj )
501 *
502  tjjs = ap( ip )
503  tjj = cabs1( tjjs )
504 *
505  IF( tjj.GE.smlnum ) THEN
506 *
507 * M(j) = M(j-1)*( 1 + CNORM(j) ) / abs(A(j,j))
508 *
509  IF( xj.GT.tjj )
510  $ xbnd = xbnd*( tjj / xj )
511  ELSE
512 *
513 * M(j) could overflow, set XBND to 0.
514 *
515  xbnd = zero
516  END IF
517  jlen = jlen + 1
518  ip = ip + jinc*jlen
519  70 CONTINUE
520  grow = min( grow, xbnd )
521  ELSE
522 *
523 * A is unit triangular.
524 *
525 * Compute GROW = 1/G(j), where G(0) = max{x(i), i=1,...,n}.
526 *
527  grow = min( one, half / max( xbnd, smlnum ) )
528  DO 80 j = jfirst, jlast, jinc
529 *
530 * Exit the loop if the growth factor is too small.
531 *
532  IF( grow.LE.smlnum )
533  $ GO TO 90
534 *
535 * G(j) = ( 1 + CNORM(j) )*G(j-1)
536 *
537  xj = one + cnorm( j )
538  grow = grow / xj
539  80 CONTINUE
540  END IF
541  90 CONTINUE
542  END IF
543 *
544  IF( ( grow*tscal ).GT.smlnum ) THEN
545 *
546 * Use the Level 2 BLAS solve if the reciprocal of the bound on
547 * elements of X is not too small.
548 *
549  CALL ctpsv( uplo, trans, diag, n, ap, x, 1 )
550  ELSE
551 *
552 * Use a Level 1 BLAS solve, scaling intermediate results.
553 *
554  IF( xmax.GT.bignum*half ) THEN
555 *
556 * Scale X so that its components are less than or equal to
557 * BIGNUM in absolute value.
558 *
559  scale = ( bignum*half ) / xmax
560  CALL csscal( n, scale, x, 1 )
561  xmax = bignum
562  ELSE
563  xmax = xmax*two
564  END IF
565 *
566  IF( notran ) THEN
567 *
568 * Solve A * x = b
569 *
570  ip = jfirst*( jfirst+1 ) / 2
571  DO 110 j = jfirst, jlast, jinc
572 *
573 * Compute x(j) = b(j) / A(j,j), scaling x if necessary.
574 *
575  xj = cabs1( x( j ) )
576  IF( nounit ) THEN
577  tjjs = ap( ip )*tscal
578  ELSE
579  tjjs = tscal
580  IF( tscal.EQ.one )
581  $ GO TO 105
582  END IF
583  tjj = cabs1( tjjs )
584  IF( tjj.GT.smlnum ) THEN
585 *
586 * abs(A(j,j)) > SMLNUM:
587 *
588  IF( tjj.LT.one ) THEN
589  IF( xj.GT.tjj*bignum ) THEN
590 *
591 * Scale x by 1/b(j).
592 *
593  rec = one / xj
594  CALL csscal( n, rec, x, 1 )
595  scale = scale*rec
596  xmax = xmax*rec
597  END IF
598  END IF
599  x( j ) = cladiv( x( j ), tjjs )
600  xj = cabs1( x( j ) )
601  ELSE IF( tjj.GT.zero ) THEN
602 *
603 * 0 < abs(A(j,j)) <= SMLNUM:
604 *
605  IF( xj.GT.tjj*bignum ) THEN
606 *
607 * Scale x by (1/abs(x(j)))*abs(A(j,j))*BIGNUM
608 * to avoid overflow when dividing by A(j,j).
609 *
610  rec = ( tjj*bignum ) / xj
611  IF( cnorm( j ).GT.one ) THEN
612 *
613 * Scale by 1/CNORM(j) to avoid overflow when
614 * multiplying x(j) times column j.
615 *
616  rec = rec / cnorm( j )
617  END IF
618  CALL csscal( n, rec, x, 1 )
619  scale = scale*rec
620  xmax = xmax*rec
621  END IF
622  x( j ) = cladiv( x( j ), tjjs )
623  xj = cabs1( x( j ) )
624  ELSE
625 *
626 * A(j,j) = 0: Set x(1:n) = 0, x(j) = 1, and
627 * scale = 0, and compute a solution to A*x = 0.
628 *
629  DO 100 i = 1, n
630  x( i ) = zero
631  100 CONTINUE
632  x( j ) = one
633  xj = one
634  scale = zero
635  xmax = zero
636  END IF
637  105 CONTINUE
638 *
639 * Scale x if necessary to avoid overflow when adding a
640 * multiple of column j of A.
641 *
642  IF( xj.GT.one ) THEN
643  rec = one / xj
644  IF( cnorm( j ).GT.( bignum-xmax )*rec ) THEN
645 *
646 * Scale x by 1/(2*abs(x(j))).
647 *
648  rec = rec*half
649  CALL csscal( n, rec, x, 1 )
650  scale = scale*rec
651  END IF
652  ELSE IF( xj*cnorm( j ).GT.( bignum-xmax ) ) THEN
653 *
654 * Scale x by 1/2.
655 *
656  CALL csscal( n, half, x, 1 )
657  scale = scale*half
658  END IF
659 *
660  IF( upper ) THEN
661  IF( j.GT.1 ) THEN
662 *
663 * Compute the update
664 * x(1:j-1) := x(1:j-1) - x(j) * A(1:j-1,j)
665 *
666  CALL caxpy( j-1, -x( j )*tscal, ap( ip-j+1 ), 1, x,
667  $ 1 )
668  i = icamax( j-1, x, 1 )
669  xmax = cabs1( x( i ) )
670  END IF
671  ip = ip - j
672  ELSE
673  IF( j.LT.n ) THEN
674 *
675 * Compute the update
676 * x(j+1:n) := x(j+1:n) - x(j) * A(j+1:n,j)
677 *
678  CALL caxpy( n-j, -x( j )*tscal, ap( ip+1 ), 1,
679  $ x( j+1 ), 1 )
680  i = j + icamax( n-j, x( j+1 ), 1 )
681  xmax = cabs1( x( i ) )
682  END IF
683  ip = ip + n - j + 1
684  END IF
685  110 CONTINUE
686 *
687  ELSE IF( lsame( trans, 'T' ) ) THEN
688 *
689 * Solve A**T * x = b
690 *
691  ip = jfirst*( jfirst+1 ) / 2
692  jlen = 1
693  DO 150 j = jfirst, jlast, jinc
694 *
695 * Compute x(j) = b(j) - sum A(k,j)*x(k).
696 * k<>j
697 *
698  xj = cabs1( x( j ) )
699  uscal = tscal
700  rec = one / max( xmax, one )
701  IF( cnorm( j ).GT.( bignum-xj )*rec ) THEN
702 *
703 * If x(j) could overflow, scale x by 1/(2*XMAX).
704 *
705  rec = rec*half
706  IF( nounit ) THEN
707  tjjs = ap( ip )*tscal
708  ELSE
709  tjjs = tscal
710  END IF
711  tjj = cabs1( tjjs )
712  IF( tjj.GT.one ) THEN
713 *
714 * Divide by A(j,j) when scaling x if A(j,j) > 1.
715 *
716  rec = min( one, rec*tjj )
717  uscal = cladiv( uscal, tjjs )
718  END IF
719  IF( rec.LT.one ) THEN
720  CALL csscal( n, rec, x, 1 )
721  scale = scale*rec
722  xmax = xmax*rec
723  END IF
724  END IF
725 *
726  csumj = zero
727  IF( uscal.EQ.cmplx( one ) ) THEN
728 *
729 * If the scaling needed for A in the dot product is 1,
730 * call CDOTU to perform the dot product.
731 *
732  IF( upper ) THEN
733  csumj = cdotu( j-1, ap( ip-j+1 ), 1, x, 1 )
734  ELSE IF( j.LT.n ) THEN
735  csumj = cdotu( n-j, ap( ip+1 ), 1, x( j+1 ), 1 )
736  END IF
737  ELSE
738 *
739 * Otherwise, use in-line code for the dot product.
740 *
741  IF( upper ) THEN
742  DO 120 i = 1, j - 1
743  csumj = csumj + ( ap( ip-j+i )*uscal )*x( i )
744  120 CONTINUE
745  ELSE IF( j.LT.n ) THEN
746  DO 130 i = 1, n - j
747  csumj = csumj + ( ap( ip+i )*uscal )*x( j+i )
748  130 CONTINUE
749  END IF
750  END IF
751 *
752  IF( uscal.EQ.cmplx( tscal ) ) THEN
753 *
754 * Compute x(j) := ( x(j) - CSUMJ ) / A(j,j) if 1/A(j,j)
755 * was not used to scale the dotproduct.
756 *
757  x( j ) = x( j ) - csumj
758  xj = cabs1( x( j ) )
759  IF( nounit ) THEN
760 *
761 * Compute x(j) = x(j) / A(j,j), scaling if necessary.
762 *
763  tjjs = ap( ip )*tscal
764  ELSE
765  tjjs = tscal
766  IF( tscal.EQ.one )
767  $ GO TO 145
768  END IF
769  tjj = cabs1( tjjs )
770  IF( tjj.GT.smlnum ) THEN
771 *
772 * abs(A(j,j)) > SMLNUM:
773 *
774  IF( tjj.LT.one ) THEN
775  IF( xj.GT.tjj*bignum ) THEN
776 *
777 * Scale X by 1/abs(x(j)).
778 *
779  rec = one / xj
780  CALL csscal( n, rec, x, 1 )
781  scale = scale*rec
782  xmax = xmax*rec
783  END IF
784  END IF
785  x( j ) = cladiv( x( j ), tjjs )
786  ELSE IF( tjj.GT.zero ) THEN
787 *
788 * 0 < abs(A(j,j)) <= SMLNUM:
789 *
790  IF( xj.GT.tjj*bignum ) THEN
791 *
792 * Scale x by (1/abs(x(j)))*abs(A(j,j))*BIGNUM.
793 *
794  rec = ( tjj*bignum ) / xj
795  CALL csscal( n, rec, x, 1 )
796  scale = scale*rec
797  xmax = xmax*rec
798  END IF
799  x( j ) = cladiv( x( j ), tjjs )
800  ELSE
801 *
802 * A(j,j) = 0: Set x(1:n) = 0, x(j) = 1, and
803 * scale = 0 and compute a solution to A**T *x = 0.
804 *
805  DO 140 i = 1, n
806  x( i ) = zero
807  140 CONTINUE
808  x( j ) = one
809  scale = zero
810  xmax = zero
811  END IF
812  145 CONTINUE
813  ELSE
814 *
815 * Compute x(j) := x(j) / A(j,j) - CSUMJ if the dot
816 * product has already been divided by 1/A(j,j).
817 *
818  x( j ) = cladiv( x( j ), tjjs ) - csumj
819  END IF
820  xmax = max( xmax, cabs1( x( j ) ) )
821  jlen = jlen + 1
822  ip = ip + jinc*jlen
823  150 CONTINUE
824 *
825  ELSE
826 *
827 * Solve A**H * x = b
828 *
829  ip = jfirst*( jfirst+1 ) / 2
830  jlen = 1
831  DO 190 j = jfirst, jlast, jinc
832 *
833 * Compute x(j) = b(j) - sum A(k,j)*x(k).
834 * k<>j
835 *
836  xj = cabs1( x( j ) )
837  uscal = tscal
838  rec = one / max( xmax, one )
839  IF( cnorm( j ).GT.( bignum-xj )*rec ) THEN
840 *
841 * If x(j) could overflow, scale x by 1/(2*XMAX).
842 *
843  rec = rec*half
844  IF( nounit ) THEN
845  tjjs = conjg( ap( ip ) )*tscal
846  ELSE
847  tjjs = tscal
848  END IF
849  tjj = cabs1( tjjs )
850  IF( tjj.GT.one ) THEN
851 *
852 * Divide by A(j,j) when scaling x if A(j,j) > 1.
853 *
854  rec = min( one, rec*tjj )
855  uscal = cladiv( uscal, tjjs )
856  END IF
857  IF( rec.LT.one ) THEN
858  CALL csscal( n, rec, x, 1 )
859  scale = scale*rec
860  xmax = xmax*rec
861  END IF
862  END IF
863 *
864  csumj = zero
865  IF( uscal.EQ.cmplx( one ) ) THEN
866 *
867 * If the scaling needed for A in the dot product is 1,
868 * call CDOTC to perform the dot product.
869 *
870  IF( upper ) THEN
871  csumj = cdotc( j-1, ap( ip-j+1 ), 1, x, 1 )
872  ELSE IF( j.LT.n ) THEN
873  csumj = cdotc( n-j, ap( ip+1 ), 1, x( j+1 ), 1 )
874  END IF
875  ELSE
876 *
877 * Otherwise, use in-line code for the dot product.
878 *
879  IF( upper ) THEN
880  DO 160 i = 1, j - 1
881  csumj = csumj + ( conjg( ap( ip-j+i ) )*uscal )*
882  $ x( i )
883  160 CONTINUE
884  ELSE IF( j.LT.n ) THEN
885  DO 170 i = 1, n - j
886  csumj = csumj + ( conjg( ap( ip+i ) )*uscal )*
887  $ x( j+i )
888  170 CONTINUE
889  END IF
890  END IF
891 *
892  IF( uscal.EQ.cmplx( tscal ) ) THEN
893 *
894 * Compute x(j) := ( x(j) - CSUMJ ) / A(j,j) if 1/A(j,j)
895 * was not used to scale the dotproduct.
896 *
897  x( j ) = x( j ) - csumj
898  xj = cabs1( x( j ) )
899  IF( nounit ) THEN
900 *
901 * Compute x(j) = x(j) / A(j,j), scaling if necessary.
902 *
903  tjjs = conjg( ap( ip ) )*tscal
904  ELSE
905  tjjs = tscal
906  IF( tscal.EQ.one )
907  $ GO TO 185
908  END IF
909  tjj = cabs1( tjjs )
910  IF( tjj.GT.smlnum ) THEN
911 *
912 * abs(A(j,j)) > SMLNUM:
913 *
914  IF( tjj.LT.one ) THEN
915  IF( xj.GT.tjj*bignum ) THEN
916 *
917 * Scale X by 1/abs(x(j)).
918 *
919  rec = one / xj
920  CALL csscal( n, rec, x, 1 )
921  scale = scale*rec
922  xmax = xmax*rec
923  END IF
924  END IF
925  x( j ) = cladiv( x( j ), tjjs )
926  ELSE IF( tjj.GT.zero ) THEN
927 *
928 * 0 < abs(A(j,j)) <= SMLNUM:
929 *
930  IF( xj.GT.tjj*bignum ) THEN
931 *
932 * Scale x by (1/abs(x(j)))*abs(A(j,j))*BIGNUM.
933 *
934  rec = ( tjj*bignum ) / xj
935  CALL csscal( n, rec, x, 1 )
936  scale = scale*rec
937  xmax = xmax*rec
938  END IF
939  x( j ) = cladiv( x( j ), tjjs )
940  ELSE
941 *
942 * A(j,j) = 0: Set x(1:n) = 0, x(j) = 1, and
943 * scale = 0 and compute a solution to A**H *x = 0.
944 *
945  DO 180 i = 1, n
946  x( i ) = zero
947  180 CONTINUE
948  x( j ) = one
949  scale = zero
950  xmax = zero
951  END IF
952  185 CONTINUE
953  ELSE
954 *
955 * Compute x(j) := x(j) / A(j,j) - CSUMJ if the dot
956 * product has already been divided by 1/A(j,j).
957 *
958  x( j ) = cladiv( x( j ), tjjs ) - csumj
959  END IF
960  xmax = max( xmax, cabs1( x( j ) ) )
961  jlen = jlen + 1
962  ip = ip + jinc*jlen
963  190 CONTINUE
964  END IF
965  scale = scale / tscal
966  END IF
967 *
968 * Scale the column norms by 1/TSCAL for return.
969 *
970  IF( tscal.NE.one ) THEN
971  CALL sscal( n, one / tscal, cnorm, 1 )
972  END IF
973 *
974  RETURN
975 *
976 * End of CLATPS
977 *
real function scasum(N, CX, INCX)
SCASUM
Definition: scasum.f:54
integer function isamax(N, SX, INCX)
ISAMAX
Definition: isamax.f:53
subroutine slabad(SMALL, LARGE)
SLABAD
Definition: slabad.f:76
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
complex function cdotu(N, CX, INCX, CY, INCY)
CDOTU
Definition: cdotu.f:54
complex function cdotc(N, CX, INCX, CY, INCY)
CDOTC
Definition: cdotc.f:54
integer function icamax(N, CX, INCX)
ICAMAX
Definition: icamax.f:53
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:55
subroutine ctpsv(UPLO, TRANS, DIAG, N, AP, X, INCX)
CTPSV
Definition: ctpsv.f:146
complex function cladiv(X, Y)
CLADIV performs complex division in real arithmetic, avoiding unnecessary overflow.
Definition: cladiv.f:66
real function slamch(CMACH)
SLAMCH
Definition: slamch.f:69
subroutine caxpy(N, CA, CX, INCX, CY, INCY)
CAXPY
Definition: caxpy.f:53
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine csscal(N, SA, CX, INCX)
CSSCAL
Definition: csscal.f:54

Here is the call graph for this function:

Here is the caller graph for this function: