LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ zlatrs()

subroutine zlatrs ( character uplo,
character trans,
character diag,
character normin,
integer n,
complex*16, dimension( lda, * ) a,
integer lda,
complex*16, dimension( * ) x,
double precision scale,
double precision, dimension( * ) cnorm,
integer info )

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

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

Purpose:
!>
!> ZLATRS 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.  Here A is an upper or lower
!> triangular matrix, 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
!> ZTRSV 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]A
!>          A is COMPLEX*16 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 COMPLEX*16 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 DOUBLE PRECISION
!>          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 DOUBLE PRECISION 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.
Further Details:
!>
!>  A rough bound on x is computed; if that is less than overflow, ZTRSV
!>  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 ZTRSV 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 ZTRSV if 1/M(n) and 1/G(n) are both greater
!>  than max(underflow, 1/overflow).
!> 

Definition at line 235 of file zlatrs.f.

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