LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dsytf2_rook ( character  UPLO,
integer  N,
double precision, dimension( lda, * )  A,
integer  LDA,
integer, dimension( * )  IPIV,
integer  INFO 
)

DSYTF2_ROOK computes the factorization of a real symmetric indefinite matrix using the bounded Bunch-Kaufman ("rook") diagonal pivoting method (unblocked algorithm).

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

Purpose:
 DSYTF2_ROOK computes the factorization of a real symmetric matrix A
 using the bounded Bunch-Kaufman ("rook") diagonal pivoting method:

    A = U*D*U**T  or  A = L*D*L**T

 where U (or L) is a product of permutation and unit upper (lower)
 triangular matrices, U**T is the transpose of U, and D is symmetric and
 block diagonal with 1-by-1 and 2-by-2 diagonal blocks.

 This is the unblocked version of the algorithm, calling Level 2 BLAS.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          Specifies whether the upper or lower triangular part of the
          symmetric matrix A is stored:
          = 'U':  Upper triangular
          = 'L':  Lower triangular
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]A
          A is DOUBLE PRECISION array, dimension (LDA,N)
          On entry, the symmetric matrix A.  If UPLO = 'U', the leading
          n-by-n upper triangular part of A contains the upper
          triangular part of the matrix A, and the strictly lower
          triangular part of A is not referenced.  If UPLO = 'L', the
          leading n-by-n lower triangular part of A contains the lower
          triangular part of the matrix A, and the strictly upper
          triangular part of A is not referenced.

          On exit, the block diagonal matrix D and the multipliers used
          to obtain the factor U or L (see below for further details).
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]IPIV
          IPIV is INTEGER array, dimension (N)
          Details of the interchanges and the block structure of D.

          If UPLO = 'U':
             If IPIV(k) > 0, then rows and columns k and IPIV(k)
             were interchanged and D(k,k) is a 1-by-1 diagonal block.

             If IPIV(k) < 0 and IPIV(k-1) < 0, then rows and
             columns k and -IPIV(k) were interchanged and rows and
             columns k-1 and -IPIV(k-1) were inerchaged,
             D(k-1:k,k-1:k) is a 2-by-2 diagonal block.

          If UPLO = 'L':
             If IPIV(k) > 0, then rows and columns k and IPIV(k)
             were interchanged and D(k,k) is a 1-by-1 diagonal block.

             If IPIV(k) < 0 and IPIV(k+1) < 0, then rows and
             columns k and -IPIV(k) were interchanged and rows and
             columns k+1 and -IPIV(k+1) were inerchaged,
             D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -k, the k-th argument had an illegal value
          > 0: if INFO = k, D(k,k) is exactly zero.  The factorization
               has been completed, but the block diagonal matrix D is
               exactly singular, and division by zero will occur if it
               is used to solve a system of equations.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2013
Further Details:
  If UPLO = 'U', then A = U*D*U**T, where
     U = P(n)*U(n)* ... *P(k)U(k)* ...,
  i.e., U is a product of terms P(k)*U(k), where k decreases from n to
  1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1
  and 2-by-2 diagonal blocks D(k).  P(k) is a permutation matrix as
  defined by IPIV(k), and U(k) is a unit upper triangular matrix, such
  that if the diagonal block D(k) is of order s (s = 1 or 2), then

             (   I    v    0   )   k-s
     U(k) =  (   0    I    0   )   s
             (   0    0    I   )   n-k
                k-s   s   n-k

  If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k).
  If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k),
  and A(k,k), and v overwrites A(1:k-2,k-1:k).

  If UPLO = 'L', then A = L*D*L**T, where
     L = P(1)*L(1)* ... *P(k)*L(k)* ...,
  i.e., L is a product of terms P(k)*L(k), where k increases from 1 to
  n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1
  and 2-by-2 diagonal blocks D(k).  P(k) is a permutation matrix as
  defined by IPIV(k), and L(k) is a unit lower triangular matrix, such
  that if the diagonal block D(k) is of order s (s = 1 or 2), then

             (   I    0     0   )  k-1
     L(k) =  (   0    I     0   )  s
             (   0    v     I   )  n-k-s+1
                k-1   s  n-k-s+1

  If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k).
  If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k),
  and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
Contributors:
  November 2013,     Igor Kozachenko,
                  Computer Science Division,
                  University of California, Berkeley

  September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
                  School of Mathematics,
                  University of Manchester

  01-01-96 - Based on modifications by
    J. Lewis, Boeing Computer Services Company
    A. Petitet, Computer Science Dept., Univ. of Tenn., Knoxville abd , USA

Definition at line 196 of file dsytf2_rook.f.

196 *
197 * -- LAPACK computational routine (version 3.5.0) --
198 * -- LAPACK is a software package provided by Univ. of Tennessee, --
199 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
200 * November 2013
201 *
202 * .. Scalar Arguments ..
203  CHARACTER uplo
204  INTEGER info, lda, n
205 * ..
206 * .. Array Arguments ..
207  INTEGER ipiv( * )
208  DOUBLE PRECISION a( lda, * )
209 * ..
210 *
211 * =====================================================================
212 *
213 * .. Parameters ..
214  DOUBLE PRECISION zero, one
215  parameter ( zero = 0.0d+0, one = 1.0d+0 )
216  DOUBLE PRECISION eight, sevten
217  parameter ( eight = 8.0d+0, sevten = 17.0d+0 )
218 * ..
219 * .. Local Scalars ..
220  LOGICAL upper, done
221  INTEGER i, imax, j, jmax, itemp, k, kk, kp, kstep,
222  $ p, ii
223  DOUBLE PRECISION absakk, alpha, colmax, d11, d12, d21, d22,
224  $ rowmax, dtemp, t, wk, wkm1, wkp1, sfmin
225 * ..
226 * .. External Functions ..
227  LOGICAL lsame
228  INTEGER idamax
229  DOUBLE PRECISION dlamch
230  EXTERNAL lsame, idamax, dlamch
231 * ..
232 * .. External Subroutines ..
233  EXTERNAL dscal, dswap, dsyr, xerbla
234 * ..
235 * .. Intrinsic Functions ..
236  INTRINSIC abs, max, sqrt
237 * ..
238 * .. Executable Statements ..
239 *
240 * Test the input parameters.
241 *
242  info = 0
243  upper = lsame( uplo, 'U' )
244  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
245  info = -1
246  ELSE IF( n.LT.0 ) THEN
247  info = -2
248  ELSE IF( lda.LT.max( 1, n ) ) THEN
249  info = -4
250  END IF
251  IF( info.NE.0 ) THEN
252  CALL xerbla( 'DSYTF2_ROOK', -info )
253  RETURN
254  END IF
255 *
256 * Initialize ALPHA for use in choosing pivot block size.
257 *
258  alpha = ( one+sqrt( sevten ) ) / eight
259 *
260 * Compute machine safe minimum
261 *
262  sfmin = dlamch( 'S' )
263 *
264  IF( upper ) THEN
265 *
266 * Factorize A as U*D*U**T using the upper triangle of A
267 *
268 * K is the main loop index, decreasing from N to 1 in steps of
269 * 1 or 2
270 *
271  k = n
272  10 CONTINUE
273 *
274 * If K < 1, exit from loop
275 *
276  IF( k.LT.1 )
277  $ GO TO 70
278  kstep = 1
279  p = k
280 *
281 * Determine rows and columns to be interchanged and whether
282 * a 1-by-1 or 2-by-2 pivot block will be used
283 *
284  absakk = abs( a( k, k ) )
285 *
286 * IMAX is the row-index of the largest off-diagonal element in
287 * column K, and COLMAX is its absolute value.
288 * Determine both COLMAX and IMAX.
289 *
290  IF( k.GT.1 ) THEN
291  imax = idamax( k-1, a( 1, k ), 1 )
292  colmax = abs( a( imax, k ) )
293  ELSE
294  colmax = zero
295  END IF
296 *
297  IF( (max( absakk, colmax ).EQ.zero) ) THEN
298 *
299 * Column K is zero or underflow: set INFO and continue
300 *
301  IF( info.EQ.0 )
302  $ info = k
303  kp = k
304  ELSE
305 *
306 * Test for interchange
307 *
308 * Equivalent to testing for (used to handle NaN and Inf)
309 * ABSAKK.GE.ALPHA*COLMAX
310 *
311  IF( .NOT.( absakk.LT.alpha*colmax ) ) THEN
312 *
313 * no interchange,
314 * use 1-by-1 pivot block
315 *
316  kp = k
317  ELSE
318 *
319  done = .false.
320 *
321 * Loop until pivot found
322 *
323  12 CONTINUE
324 *
325 * Begin pivot search loop body
326 *
327 * JMAX is the column-index of the largest off-diagonal
328 * element in row IMAX, and ROWMAX is its absolute value.
329 * Determine both ROWMAX and JMAX.
330 *
331  IF( imax.NE.k ) THEN
332  jmax = imax + idamax( k-imax, a( imax, imax+1 ),
333  $ lda )
334  rowmax = abs( a( imax, jmax ) )
335  ELSE
336  rowmax = zero
337  END IF
338 *
339  IF( imax.GT.1 ) THEN
340  itemp = idamax( imax-1, a( 1, imax ), 1 )
341  dtemp = abs( a( itemp, imax ) )
342  IF( dtemp.GT.rowmax ) THEN
343  rowmax = dtemp
344  jmax = itemp
345  END IF
346  END IF
347 *
348 * Equivalent to testing for (used to handle NaN and Inf)
349 * ABS( A( IMAX, IMAX ) ).GE.ALPHA*ROWMAX
350 *
351  IF( .NOT.( abs( a( imax, imax ) ).LT.alpha*rowmax ) )
352  $ THEN
353 *
354 * interchange rows and columns K and IMAX,
355 * use 1-by-1 pivot block
356 *
357  kp = imax
358  done = .true.
359 *
360 * Equivalent to testing for ROWMAX .EQ. COLMAX,
361 * used to handle NaN and Inf
362 *
363  ELSE IF( ( p.EQ.jmax ).OR.( rowmax.LE.colmax ) ) THEN
364 *
365 * interchange rows and columns K+1 and IMAX,
366 * use 2-by-2 pivot block
367 *
368  kp = imax
369  kstep = 2
370  done = .true.
371  ELSE
372 *
373 * Pivot NOT found, set variables and repeat
374 *
375  p = imax
376  colmax = rowmax
377  imax = jmax
378  END IF
379 *
380 * End pivot search loop body
381 *
382  IF( .NOT. done ) GOTO 12
383 *
384  END IF
385 *
386 * Swap TWO rows and TWO columns
387 *
388 * First swap
389 *
390  IF( ( kstep.EQ.2 ) .AND. ( p.NE.k ) ) THEN
391 *
392 * Interchange rows and column K and P in the leading
393 * submatrix A(1:k,1:k) if we have a 2-by-2 pivot
394 *
395  IF( p.GT.1 )
396  $ CALL dswap( p-1, a( 1, k ), 1, a( 1, p ), 1 )
397  IF( p.LT.(k-1) )
398  $ CALL dswap( k-p-1, a( p+1, k ), 1, a( p, p+1 ),
399  $ lda )
400  t = a( k, k )
401  a( k, k ) = a( p, p )
402  a( p, p ) = t
403  END IF
404 *
405 * Second swap
406 *
407  kk = k - kstep + 1
408  IF( kp.NE.kk ) THEN
409 *
410 * Interchange rows and columns KK and KP in the leading
411 * submatrix A(1:k,1:k)
412 *
413  IF( kp.GT.1 )
414  $ CALL dswap( kp-1, a( 1, kk ), 1, a( 1, kp ), 1 )
415  IF( ( kk.GT.1 ) .AND. ( kp.LT.(kk-1) ) )
416  $ CALL dswap( kk-kp-1, a( kp+1, kk ), 1, a( kp, kp+1 ),
417  $ lda )
418  t = a( kk, kk )
419  a( kk, kk ) = a( kp, kp )
420  a( kp, kp ) = t
421  IF( kstep.EQ.2 ) THEN
422  t = a( k-1, k )
423  a( k-1, k ) = a( kp, k )
424  a( kp, k ) = t
425  END IF
426  END IF
427 *
428 * Update the leading submatrix
429 *
430  IF( kstep.EQ.1 ) THEN
431 *
432 * 1-by-1 pivot block D(k): column k now holds
433 *
434 * W(k) = U(k)*D(k)
435 *
436 * where U(k) is the k-th column of U
437 *
438  IF( k.GT.1 ) THEN
439 *
440 * Perform a rank-1 update of A(1:k-1,1:k-1) and
441 * store U(k) in column k
442 *
443  IF( abs( a( k, k ) ).GE.sfmin ) THEN
444 *
445 * Perform a rank-1 update of A(1:k-1,1:k-1) as
446 * A := A - U(k)*D(k)*U(k)**T
447 * = A - W(k)*1/D(k)*W(k)**T
448 *
449  d11 = one / a( k, k )
450  CALL dsyr( uplo, k-1, -d11, a( 1, k ), 1, a, lda )
451 *
452 * Store U(k) in column k
453 *
454  CALL dscal( k-1, d11, a( 1, k ), 1 )
455  ELSE
456 *
457 * Store L(k) in column K
458 *
459  d11 = a( k, k )
460  DO 16 ii = 1, k - 1
461  a( ii, k ) = a( ii, k ) / d11
462  16 CONTINUE
463 *
464 * Perform a rank-1 update of A(k+1:n,k+1:n) as
465 * A := A - U(k)*D(k)*U(k)**T
466 * = A - W(k)*(1/D(k))*W(k)**T
467 * = A - (W(k)/D(k))*(D(k))*(W(k)/D(K))**T
468 *
469  CALL dsyr( uplo, k-1, -d11, a( 1, k ), 1, a, lda )
470  END IF
471  END IF
472 *
473  ELSE
474 *
475 * 2-by-2 pivot block D(k): columns k and k-1 now hold
476 *
477 * ( W(k-1) W(k) ) = ( U(k-1) U(k) )*D(k)
478 *
479 * where U(k) and U(k-1) are the k-th and (k-1)-th columns
480 * of U
481 *
482 * Perform a rank-2 update of A(1:k-2,1:k-2) as
483 *
484 * A := A - ( U(k-1) U(k) )*D(k)*( U(k-1) U(k) )**T
485 * = A - ( ( A(k-1)A(k) )*inv(D(k)) ) * ( A(k-1)A(k) )**T
486 *
487 * and store L(k) and L(k+1) in columns k and k+1
488 *
489  IF( k.GT.2 ) THEN
490 *
491  d12 = a( k-1, k )
492  d22 = a( k-1, k-1 ) / d12
493  d11 = a( k, k ) / d12
494  t = one / ( d11*d22-one )
495 *
496  DO 30 j = k - 2, 1, -1
497 *
498  wkm1 = t*( d11*a( j, k-1 )-a( j, k ) )
499  wk = t*( d22*a( j, k )-a( j, k-1 ) )
500 *
501  DO 20 i = j, 1, -1
502  a( i, j ) = a( i, j ) - (a( i, k ) / d12 )*wk -
503  $ ( a( i, k-1 ) / d12 )*wkm1
504  20 CONTINUE
505 *
506 * Store U(k) and U(k-1) in cols k and k-1 for row J
507 *
508  a( j, k ) = wk / d12
509  a( j, k-1 ) = wkm1 / d12
510 *
511  30 CONTINUE
512 *
513  END IF
514 *
515  END IF
516  END IF
517 *
518 * Store details of the interchanges in IPIV
519 *
520  IF( kstep.EQ.1 ) THEN
521  ipiv( k ) = kp
522  ELSE
523  ipiv( k ) = -p
524  ipiv( k-1 ) = -kp
525  END IF
526 *
527 * Decrease K and return to the start of the main loop
528 *
529  k = k - kstep
530  GO TO 10
531 *
532  ELSE
533 *
534 * Factorize A as L*D*L**T using the lower triangle of A
535 *
536 * K is the main loop index, increasing from 1 to N in steps of
537 * 1 or 2
538 *
539  k = 1
540  40 CONTINUE
541 *
542 * If K > N, exit from loop
543 *
544  IF( k.GT.n )
545  $ GO TO 70
546  kstep = 1
547  p = k
548 *
549 * Determine rows and columns to be interchanged and whether
550 * a 1-by-1 or 2-by-2 pivot block will be used
551 *
552  absakk = abs( a( k, k ) )
553 *
554 * IMAX is the row-index of the largest off-diagonal element in
555 * column K, and COLMAX is its absolute value.
556 * Determine both COLMAX and IMAX.
557 *
558  IF( k.LT.n ) THEN
559  imax = k + idamax( n-k, a( k+1, k ), 1 )
560  colmax = abs( a( imax, k ) )
561  ELSE
562  colmax = zero
563  END IF
564 *
565  IF( ( max( absakk, colmax ).EQ.zero ) ) THEN
566 *
567 * Column K is zero or underflow: set INFO and continue
568 *
569  IF( info.EQ.0 )
570  $ info = k
571  kp = k
572  ELSE
573 *
574 * Test for interchange
575 *
576 * Equivalent to testing for (used to handle NaN and Inf)
577 * ABSAKK.GE.ALPHA*COLMAX
578 *
579  IF( .NOT.( absakk.LT.alpha*colmax ) ) THEN
580 *
581 * no interchange, use 1-by-1 pivot block
582 *
583  kp = k
584  ELSE
585 *
586  done = .false.
587 *
588 * Loop until pivot found
589 *
590  42 CONTINUE
591 *
592 * Begin pivot search loop body
593 *
594 * JMAX is the column-index of the largest off-diagonal
595 * element in row IMAX, and ROWMAX is its absolute value.
596 * Determine both ROWMAX and JMAX.
597 *
598  IF( imax.NE.k ) THEN
599  jmax = k - 1 + idamax( imax-k, a( imax, k ), lda )
600  rowmax = abs( a( imax, jmax ) )
601  ELSE
602  rowmax = zero
603  END IF
604 *
605  IF( imax.LT.n ) THEN
606  itemp = imax + idamax( n-imax, a( imax+1, imax ),
607  $ 1 )
608  dtemp = abs( a( itemp, imax ) )
609  IF( dtemp.GT.rowmax ) THEN
610  rowmax = dtemp
611  jmax = itemp
612  END IF
613  END IF
614 *
615 * Equivalent to testing for (used to handle NaN and Inf)
616 * ABS( A( IMAX, IMAX ) ).GE.ALPHA*ROWMAX
617 *
618  IF( .NOT.( abs( a( imax, imax ) ).LT.alpha*rowmax ) )
619  $ THEN
620 *
621 * interchange rows and columns K and IMAX,
622 * use 1-by-1 pivot block
623 *
624  kp = imax
625  done = .true.
626 *
627 * Equivalent to testing for ROWMAX .EQ. COLMAX,
628 * used to handle NaN and Inf
629 *
630  ELSE IF( ( p.EQ.jmax ).OR.( rowmax.LE.colmax ) ) THEN
631 *
632 * interchange rows and columns K+1 and IMAX,
633 * use 2-by-2 pivot block
634 *
635  kp = imax
636  kstep = 2
637  done = .true.
638  ELSE
639 *
640 * Pivot NOT found, set variables and repeat
641 *
642  p = imax
643  colmax = rowmax
644  imax = jmax
645  END IF
646 *
647 * End pivot search loop body
648 *
649  IF( .NOT. done ) GOTO 42
650 *
651  END IF
652 *
653 * Swap TWO rows and TWO columns
654 *
655 * First swap
656 *
657  IF( ( kstep.EQ.2 ) .AND. ( p.NE.k ) ) THEN
658 *
659 * Interchange rows and column K and P in the trailing
660 * submatrix A(k:n,k:n) if we have a 2-by-2 pivot
661 *
662  IF( p.LT.n )
663  $ CALL dswap( n-p, a( p+1, k ), 1, a( p+1, p ), 1 )
664  IF( p.GT.(k+1) )
665  $ CALL dswap( p-k-1, a( k+1, k ), 1, a( p, k+1 ), lda )
666  t = a( k, k )
667  a( k, k ) = a( p, p )
668  a( p, p ) = t
669  END IF
670 *
671 * Second swap
672 *
673  kk = k + kstep - 1
674  IF( kp.NE.kk ) THEN
675 *
676 * Interchange rows and columns KK and KP in the trailing
677 * submatrix A(k:n,k:n)
678 *
679  IF( kp.LT.n )
680  $ CALL dswap( n-kp, a( kp+1, kk ), 1, a( kp+1, kp ), 1 )
681  IF( ( kk.LT.n ) .AND. ( kp.GT.(kk+1) ) )
682  $ CALL dswap( kp-kk-1, a( kk+1, kk ), 1, a( kp, kk+1 ),
683  $ lda )
684  t = a( kk, kk )
685  a( kk, kk ) = a( kp, kp )
686  a( kp, kp ) = t
687  IF( kstep.EQ.2 ) THEN
688  t = a( k+1, k )
689  a( k+1, k ) = a( kp, k )
690  a( kp, k ) = t
691  END IF
692  END IF
693 *
694 * Update the trailing submatrix
695 *
696  IF( kstep.EQ.1 ) THEN
697 *
698 * 1-by-1 pivot block D(k): column k now holds
699 *
700 * W(k) = L(k)*D(k)
701 *
702 * where L(k) is the k-th column of L
703 *
704  IF( k.LT.n ) THEN
705 *
706 * Perform a rank-1 update of A(k+1:n,k+1:n) and
707 * store L(k) in column k
708 *
709  IF( abs( a( k, k ) ).GE.sfmin ) THEN
710 *
711 * Perform a rank-1 update of A(k+1:n,k+1:n) as
712 * A := A - L(k)*D(k)*L(k)**T
713 * = A - W(k)*(1/D(k))*W(k)**T
714 *
715  d11 = one / a( k, k )
716  CALL dsyr( uplo, n-k, -d11, a( k+1, k ), 1,
717  $ a( k+1, k+1 ), lda )
718 *
719 * Store L(k) in column k
720 *
721  CALL dscal( n-k, d11, a( k+1, k ), 1 )
722  ELSE
723 *
724 * Store L(k) in column k
725 *
726  d11 = a( k, k )
727  DO 46 ii = k + 1, n
728  a( ii, k ) = a( ii, k ) / d11
729  46 CONTINUE
730 *
731 * Perform a rank-1 update of A(k+1:n,k+1:n) as
732 * A := A - L(k)*D(k)*L(k)**T
733 * = A - W(k)*(1/D(k))*W(k)**T
734 * = A - (W(k)/D(k))*(D(k))*(W(k)/D(K))**T
735 *
736  CALL dsyr( uplo, n-k, -d11, a( k+1, k ), 1,
737  $ a( k+1, k+1 ), lda )
738  END IF
739  END IF
740 *
741  ELSE
742 *
743 * 2-by-2 pivot block D(k): columns k and k+1 now hold
744 *
745 * ( W(k) W(k+1) ) = ( L(k) L(k+1) )*D(k)
746 *
747 * where L(k) and L(k+1) are the k-th and (k+1)-th columns
748 * of L
749 *
750 *
751 * Perform a rank-2 update of A(k+2:n,k+2:n) as
752 *
753 * A := A - ( L(k) L(k+1) ) * D(k) * ( L(k) L(k+1) )**T
754 * = A - ( ( A(k)A(k+1) )*inv(D(k) ) * ( A(k)A(k+1) )**T
755 *
756 * and store L(k) and L(k+1) in columns k and k+1
757 *
758  IF( k.LT.n-1 ) THEN
759 *
760  d21 = a( k+1, k )
761  d11 = a( k+1, k+1 ) / d21
762  d22 = a( k, k ) / d21
763  t = one / ( d11*d22-one )
764 *
765  DO 60 j = k + 2, n
766 *
767 * Compute D21 * ( W(k)W(k+1) ) * inv(D(k)) for row J
768 *
769  wk = t*( d11*a( j, k )-a( j, k+1 ) )
770  wkp1 = t*( d22*a( j, k+1 )-a( j, k ) )
771 *
772 * Perform a rank-2 update of A(k+2:n,k+2:n)
773 *
774  DO 50 i = j, n
775  a( i, j ) = a( i, j ) - ( a( i, k ) / d21 )*wk -
776  $ ( a( i, k+1 ) / d21 )*wkp1
777  50 CONTINUE
778 *
779 * Store L(k) and L(k+1) in cols k and k+1 for row J
780 *
781  a( j, k ) = wk / d21
782  a( j, k+1 ) = wkp1 / d21
783 *
784  60 CONTINUE
785 *
786  END IF
787 *
788  END IF
789  END IF
790 *
791 * Store details of the interchanges in IPIV
792 *
793  IF( kstep.EQ.1 ) THEN
794  ipiv( k ) = kp
795  ELSE
796  ipiv( k ) = -p
797  ipiv( k+1 ) = -kp
798  END IF
799 *
800 * Increase K and return to the start of the main loop
801 *
802  k = k + kstep
803  GO TO 40
804 *
805  END IF
806 *
807  70 CONTINUE
808 *
809  RETURN
810 *
811 * End of DSYTF2_ROOK
812 *
integer function idamax(N, DX, INCX)
IDAMAX
Definition: idamax.f:53
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dsyr(UPLO, N, ALPHA, X, INCX, A, LDA)
DSYR
Definition: dsyr.f:134
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
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: