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

◆ dlatbs()

subroutine dlatbs ( character uplo,
character trans,
character diag,
character normin,
integer n,
integer kd,
double precision, dimension( ldab, * ) ab,
integer ldab,
double precision, dimension( * ) x,
double precision scale,
double precision, dimension( * ) cnorm,
integer info )

DLATBS solves a triangular banded system of equations.

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

Purpose:
!>
!> DLATBS 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 band matrix.  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 DTBSV 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]KD
!>          KD is INTEGER
!>          The number of subdiagonals or superdiagonals in the
!>          triangular matrix A.  KD >= 0.
!> 
[in]AB
!>          AB is DOUBLE PRECISION array, dimension (LDAB,N)
!>          The upper or lower triangular band matrix A, stored in the
!>          first KD+1 rows of the array. The j-th column of A is stored
!>          in the j-th column of the array AB as follows:
!>          if UPLO = 'U', AB(kd+1+i-j,j) = A(i,j) for max(1,j-kd)<=i<=j;
!>          if UPLO = 'L', AB(1+i-j,j)    = A(i,j) for j<=i<=min(n,j+kd).
!> 
[in]LDAB
!>          LDAB is INTEGER
!>          The leading dimension of the array AB.  LDAB >= KD+1.
!> 
[in,out]X
!>          X is DOUBLE PRECISION 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  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 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, DTBSV
!>  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 DTBSV 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 DTBSV if 1/M(n) and 1/G(n) are both greater
!>  than max(underflow, 1/overflow).
!> 

Definition at line 238 of file dlatbs.f.

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