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

◆ dlatps()

subroutine dlatps ( character uplo,
character trans,
character diag,
character normin,
integer n,
double precision, dimension( * ) ap,
double precision, dimension( * ) x,
double precision scale,
double precision, dimension( * ) cnorm,
integer info )

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

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

Purpose:
!>
!> DLATPS 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
!> DTPSV 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 DOUBLE PRECISION 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 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, DTPSV
!>  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 DTPSV 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 DTPSV if 1/M(n) and 1/G(n) are both greater
!>  than max(underflow, 1/overflow).
!> 

Definition at line 225 of file dlatps.f.

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