LAPACK  3.7.0
LAPACK: Linear Algebra PACKage
ssytf2_rk.f
Go to the documentation of this file.
1 *> \brief \b SSYTF2_RK computes the factorization of a real symmetric indefinite matrix using the bounded Bunch-Kaufman (rook) diagonal pivoting method (BLAS2 unblocked algorithm).
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 *> \htmlonly
9 *> Download SSYTF2_RK + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/ssytf2_rk.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/ssytf2_rk.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/ssytf2_rk.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE SSYTF2_RK( UPLO, N, A, LDA, E, IPIV, INFO )
22 *
23 * .. Scalar Arguments ..
24 * CHARACTER UPLO
25 * INTEGER INFO, LDA, N
26 * ..
27 * .. Array Arguments ..
28 * INTEGER IPIV( * )
29 * REAL A( LDA, * ), E ( * )
30 * ..
31 *
32 *
33 *> \par Purpose:
34 * =============
35 *>
36 *> \verbatim
37 *> SSYTF2_RK computes the factorization of a real symmetric matrix A
38 *> using the bounded Bunch-Kaufman (rook) diagonal pivoting method:
39 *>
40 *> A = P*U*D*(U**T)*(P**T) or A = P*L*D*(L**T)*(P**T),
41 *>
42 *> where U (or L) is unit upper (or lower) triangular matrix,
43 *> U**T (or L**T) is the transpose of U (or L), P is a permutation
44 *> matrix, P**T is the transpose of P, and D is symmetric and block
45 *> diagonal with 1-by-1 and 2-by-2 diagonal blocks.
46 *>
47 *> This is the unblocked version of the algorithm, calling Level 2 BLAS.
48 *> For more information see Further Details section.
49 *> \endverbatim
50 *
51 * Arguments:
52 * ==========
53 *
54 *> \param[in] UPLO
55 *> \verbatim
56 *> UPLO is CHARACTER*1
57 *> Specifies whether the upper or lower triangular part of the
58 *> symmetric matrix A is stored:
59 *> = 'U': Upper triangular
60 *> = 'L': Lower triangular
61 *> \endverbatim
62 *>
63 *> \param[in] N
64 *> \verbatim
65 *> N is INTEGER
66 *> The order of the matrix A. N >= 0.
67 *> \endverbatim
68 *>
69 *> \param[in,out] A
70 *> \verbatim
71 *> A is REAL array, dimension (LDA,N)
72 *> On entry, the symmetric matrix A.
73 *> If UPLO = 'U': the leading N-by-N upper triangular part
74 *> of A contains the upper triangular part of the matrix A,
75 *> and the strictly lower triangular part of A is not
76 *> referenced.
77 *>
78 *> If UPLO = 'L': the leading N-by-N lower triangular part
79 *> of A contains the lower triangular part of the matrix A,
80 *> and the strictly upper triangular part of A is not
81 *> referenced.
82 *>
83 *> On exit, contains:
84 *> a) ONLY diagonal elements of the symmetric block diagonal
85 *> matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
86 *> (superdiagonal (or subdiagonal) elements of D
87 *> are stored on exit in array E), and
88 *> b) If UPLO = 'U': factor U in the superdiagonal part of A.
89 *> If UPLO = 'L': factor L in the subdiagonal part of A.
90 *> \endverbatim
91 *>
92 *> \param[in] LDA
93 *> \verbatim
94 *> LDA is INTEGER
95 *> The leading dimension of the array A. LDA >= max(1,N).
96 *> \endverbatim
97 *>
98 *> \param[out] E
99 *> \verbatim
100 *> E is REAL array, dimension (N)
101 *> On exit, contains the superdiagonal (or subdiagonal)
102 *> elements of the symmetric block diagonal matrix D
103 *> with 1-by-1 or 2-by-2 diagonal blocks, where
104 *> If UPLO = 'U': E(i) = D(i-1,i), i=2:N, E(1) is set to 0;
105 *> If UPLO = 'L': E(i) = D(i+1,i), i=1:N-1, E(N) is set to 0.
106 *>
107 *> NOTE: For 1-by-1 diagonal block D(k), where
108 *> 1 <= k <= N, the element E(k) is set to 0 in both
109 *> UPLO = 'U' or UPLO = 'L' cases.
110 *> \endverbatim
111 *>
112 *> \param[out] IPIV
113 *> \verbatim
114 *> IPIV is INTEGER array, dimension (N)
115 *> IPIV describes the permutation matrix P in the factorization
116 *> of matrix A as follows. The absolute value of IPIV(k)
117 *> represents the index of row and column that were
118 *> interchanged with the k-th row and column. The value of UPLO
119 *> describes the order in which the interchanges were applied.
120 *> Also, the sign of IPIV represents the block structure of
121 *> the symmetric block diagonal matrix D with 1-by-1 or 2-by-2
122 *> diagonal blocks which correspond to 1 or 2 interchanges
123 *> at each factorization step. For more info see Further
124 *> Details section.
125 *>
126 *> If UPLO = 'U',
127 *> ( in factorization order, k decreases from N to 1 ):
128 *> a) A single positive entry IPIV(k) > 0 means:
129 *> D(k,k) is a 1-by-1 diagonal block.
130 *> If IPIV(k) != k, rows and columns k and IPIV(k) were
131 *> interchanged in the matrix A(1:N,1:N);
132 *> If IPIV(k) = k, no interchange occurred.
133 *>
134 *> b) A pair of consecutive negative entries
135 *> IPIV(k) < 0 and IPIV(k-1) < 0 means:
136 *> D(k-1:k,k-1:k) is a 2-by-2 diagonal block.
137 *> (NOTE: negative entries in IPIV appear ONLY in pairs).
138 *> 1) If -IPIV(k) != k, rows and columns
139 *> k and -IPIV(k) were interchanged
140 *> in the matrix A(1:N,1:N).
141 *> If -IPIV(k) = k, no interchange occurred.
142 *> 2) If -IPIV(k-1) != k-1, rows and columns
143 *> k-1 and -IPIV(k-1) were interchanged
144 *> in the matrix A(1:N,1:N).
145 *> If -IPIV(k-1) = k-1, no interchange occurred.
146 *>
147 *> c) In both cases a) and b), always ABS( IPIV(k) ) <= k.
148 *>
149 *> d) NOTE: Any entry IPIV(k) is always NONZERO on output.
150 *>
151 *> If UPLO = 'L',
152 *> ( in factorization order, k increases from 1 to N ):
153 *> a) A single positive entry IPIV(k) > 0 means:
154 *> D(k,k) is a 1-by-1 diagonal block.
155 *> If IPIV(k) != k, rows and columns k and IPIV(k) were
156 *> interchanged in the matrix A(1:N,1:N).
157 *> If IPIV(k) = k, no interchange occurred.
158 *>
159 *> b) A pair of consecutive negative entries
160 *> IPIV(k) < 0 and IPIV(k+1) < 0 means:
161 *> D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
162 *> (NOTE: negative entries in IPIV appear ONLY in pairs).
163 *> 1) If -IPIV(k) != k, rows and columns
164 *> k and -IPIV(k) were interchanged
165 *> in the matrix A(1:N,1:N).
166 *> If -IPIV(k) = k, no interchange occurred.
167 *> 2) If -IPIV(k+1) != k+1, rows and columns
168 *> k-1 and -IPIV(k-1) were interchanged
169 *> in the matrix A(1:N,1:N).
170 *> If -IPIV(k+1) = k+1, no interchange occurred.
171 *>
172 *> c) In both cases a) and b), always ABS( IPIV(k) ) >= k.
173 *>
174 *> d) NOTE: Any entry IPIV(k) is always NONZERO on output.
175 *> \endverbatim
176 *>
177 *> \param[out] INFO
178 *> \verbatim
179 *> INFO is INTEGER
180 *> = 0: successful exit
181 *>
182 *> < 0: If INFO = -k, the k-th argument had an illegal value
183 *>
184 *> > 0: If INFO = k, the matrix A is singular, because:
185 *> If UPLO = 'U': column k in the upper
186 *> triangular part of A contains all zeros.
187 *> If UPLO = 'L': column k in the lower
188 *> triangular part of A contains all zeros.
189 *>
190 *> Therefore D(k,k) is exactly zero, and superdiagonal
191 *> elements of column k of U (or subdiagonal elements of
192 *> column k of L ) are all zeros. The factorization has
193 *> been completed, but the block diagonal matrix D is
194 *> exactly singular, and division by zero will occur if
195 *> it is used to solve a system of equations.
196 *>
197 *> NOTE: INFO only stores the first occurrence of
198 *> a singularity, any subsequent occurrence of singularity
199 *> is not stored in INFO even though the factorization
200 *> always completes.
201 *> \endverbatim
202 *
203 * Authors:
204 * ========
205 *
206 *> \author Univ. of Tennessee
207 *> \author Univ. of California Berkeley
208 *> \author Univ. of Colorado Denver
209 *> \author NAG Ltd.
210 *
211 *> \date December 2016
212 *
213 *> \ingroup singleSYcomputational
214 *
215 *> \par Further Details:
216 * =====================
217 *>
218 *> \verbatim
219 *> TODO: put further details
220 *> \endverbatim
221 *
222 *> \par Contributors:
223 * ==================
224 *>
225 *> \verbatim
226 *>
227 *> December 2016, Igor Kozachenko,
228 *> Computer Science Division,
229 *> University of California, Berkeley
230 *>
231 *> September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
232 *> School of Mathematics,
233 *> University of Manchester
234 *>
235 *> 01-01-96 - Based on modifications by
236 *> J. Lewis, Boeing Computer Services Company
237 *> A. Petitet, Computer Science Dept.,
238 *> Univ. of Tenn., Knoxville abd , USA
239 *> \endverbatim
240 *
241 * =====================================================================
242  SUBROUTINE ssytf2_rk( UPLO, N, A, LDA, E, IPIV, INFO )
243 *
244 * -- LAPACK computational routine (version 3.7.0) --
245 * -- LAPACK is a software package provided by Univ. of Tennessee, --
246 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
247 * December 2016
248 *
249 * .. Scalar Arguments ..
250  CHARACTER UPLO
251  INTEGER INFO, LDA, N
252 * ..
253 * .. Array Arguments ..
254  INTEGER IPIV( * )
255  REAL A( lda, * ), E( * )
256 * ..
257 *
258 * =====================================================================
259 *
260 * .. Parameters ..
261  REAL ZERO, ONE
262  parameter ( zero = 0.0e+0, one = 1.0e+0 )
263  REAL EIGHT, SEVTEN
264  parameter ( eight = 8.0e+0, sevten = 17.0e+0 )
265 * ..
266 * .. Local Scalars ..
267  LOGICAL UPPER, DONE
268  INTEGER I, IMAX, J, JMAX, ITEMP, K, KK, KP, KSTEP,
269  $ p, ii
270  REAL ABSAKK, ALPHA, COLMAX, D11, D12, D21, D22,
271  $ rowmax, stemp, t, wk, wkm1, wkp1, sfmin
272 * ..
273 * .. External Functions ..
274  LOGICAL LSAME
275  INTEGER ISAMAX
276  REAL SLAMCH
277  EXTERNAL lsame, isamax, slamch
278 * ..
279 * .. External Subroutines ..
280  EXTERNAL sscal, sswap, ssyr, xerbla
281 * ..
282 * .. Intrinsic Functions ..
283  INTRINSIC abs, max, sqrt
284 * ..
285 * .. Executable Statements ..
286 *
287 * Test the input parameters.
288 *
289  info = 0
290  upper = lsame( uplo, 'U' )
291  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
292  info = -1
293  ELSE IF( n.LT.0 ) THEN
294  info = -2
295  ELSE IF( lda.LT.max( 1, n ) ) THEN
296  info = -4
297  END IF
298  IF( info.NE.0 ) THEN
299  CALL xerbla( 'SSYTF2_RK', -info )
300  RETURN
301  END IF
302 *
303 * Initialize ALPHA for use in choosing pivot block size.
304 *
305  alpha = ( one+sqrt( sevten ) ) / eight
306 *
307 * Compute machine safe minimum
308 *
309  sfmin = slamch( 'S' )
310 *
311  IF( upper ) THEN
312 *
313 * Factorize A as U*D*U**T using the upper triangle of A
314 *
315 * Initilize the first entry of array E, where superdiagonal
316 * elements of D are stored
317 *
318  e( 1 ) = zero
319 *
320 * K is the main loop index, decreasing from N to 1 in steps of
321 * 1 or 2
322 *
323  k = n
324  10 CONTINUE
325 *
326 * If K < 1, exit from loop
327 *
328  IF( k.LT.1 )
329  $ GO TO 34
330  kstep = 1
331  p = k
332 *
333 * Determine rows and columns to be interchanged and whether
334 * a 1-by-1 or 2-by-2 pivot block will be used
335 *
336  absakk = abs( a( k, k ) )
337 *
338 * IMAX is the row-index of the largest off-diagonal element in
339 * column K, and COLMAX is its absolute value.
340 * Determine both COLMAX and IMAX.
341 *
342  IF( k.GT.1 ) THEN
343  imax = isamax( k-1, a( 1, k ), 1 )
344  colmax = abs( a( imax, k ) )
345  ELSE
346  colmax = zero
347  END IF
348 *
349  IF( (max( absakk, colmax ).EQ.zero) ) THEN
350 *
351 * Column K is zero or underflow: set INFO and continue
352 *
353  IF( info.EQ.0 )
354  $ info = k
355  kp = k
356 *
357 * Set E( K ) to zero
358 *
359  IF( k.GT.1 )
360  $ e( k ) = zero
361 *
362  ELSE
363 *
364 * Test for interchange
365 *
366 * Equivalent to testing for (used to handle NaN and Inf)
367 * ABSAKK.GE.ALPHA*COLMAX
368 *
369  IF( .NOT.( absakk.LT.alpha*colmax ) ) THEN
370 *
371 * no interchange,
372 * use 1-by-1 pivot block
373 *
374  kp = k
375  ELSE
376 *
377  done = .false.
378 *
379 * Loop until pivot found
380 *
381  12 CONTINUE
382 *
383 * Begin pivot search loop body
384 *
385 * JMAX is the column-index of the largest off-diagonal
386 * element in row IMAX, and ROWMAX is its absolute value.
387 * Determine both ROWMAX and JMAX.
388 *
389  IF( imax.NE.k ) THEN
390  jmax = imax + isamax( k-imax, a( imax, imax+1 ),
391  $ lda )
392  rowmax = abs( a( imax, jmax ) )
393  ELSE
394  rowmax = zero
395  END IF
396 *
397  IF( imax.GT.1 ) THEN
398  itemp = isamax( imax-1, a( 1, imax ), 1 )
399  stemp = abs( a( itemp, imax ) )
400  IF( stemp.GT.rowmax ) THEN
401  rowmax = stemp
402  jmax = itemp
403  END IF
404  END IF
405 *
406 * Equivalent to testing for (used to handle NaN and Inf)
407 * ABS( A( IMAX, IMAX ) ).GE.ALPHA*ROWMAX
408 *
409  IF( .NOT.( abs( a( imax, imax ) ).LT.alpha*rowmax ) )
410  $ THEN
411 *
412 * interchange rows and columns K and IMAX,
413 * use 1-by-1 pivot block
414 *
415  kp = imax
416  done = .true.
417 *
418 * Equivalent to testing for ROWMAX .EQ. COLMAX,
419 * used to handle NaN and Inf
420 *
421  ELSE IF( ( p.EQ.jmax ).OR.( rowmax.LE.colmax ) ) THEN
422 *
423 * interchange rows and columns K+1 and IMAX,
424 * use 2-by-2 pivot block
425 *
426  kp = imax
427  kstep = 2
428  done = .true.
429  ELSE
430 *
431 * Pivot NOT found, set variables and repeat
432 *
433  p = imax
434  colmax = rowmax
435  imax = jmax
436  END IF
437 *
438 * End pivot search loop body
439 *
440  IF( .NOT. done ) GOTO 12
441 *
442  END IF
443 *
444 * Swap TWO rows and TWO columns
445 *
446 * First swap
447 *
448  IF( ( kstep.EQ.2 ) .AND. ( p.NE.k ) ) THEN
449 *
450 * Interchange rows and column K and P in the leading
451 * submatrix A(1:k,1:k) if we have a 2-by-2 pivot
452 *
453  IF( p.GT.1 )
454  $ CALL sswap( p-1, a( 1, k ), 1, a( 1, p ), 1 )
455  IF( p.LT.(k-1) )
456  $ CALL sswap( k-p-1, a( p+1, k ), 1, a( p, p+1 ),
457  $ lda )
458  t = a( k, k )
459  a( k, k ) = a( p, p )
460  a( p, p ) = t
461 *
462 * Convert upper triangle of A into U form by applying
463 * the interchanges in columns k+1:N.
464 *
465  IF( k.LT.n )
466  $ CALL sswap( n-k, a( k, k+1 ), lda, a( p, k+1 ), lda )
467 *
468  END IF
469 *
470 * Second swap
471 *
472  kk = k - kstep + 1
473  IF( kp.NE.kk ) THEN
474 *
475 * Interchange rows and columns KK and KP in the leading
476 * submatrix A(1:k,1:k)
477 *
478  IF( kp.GT.1 )
479  $ CALL sswap( kp-1, a( 1, kk ), 1, a( 1, kp ), 1 )
480  IF( ( kk.GT.1 ) .AND. ( kp.LT.(kk-1) ) )
481  $ CALL sswap( kk-kp-1, a( kp+1, kk ), 1, a( kp, kp+1 ),
482  $ lda )
483  t = a( kk, kk )
484  a( kk, kk ) = a( kp, kp )
485  a( kp, kp ) = t
486  IF( kstep.EQ.2 ) THEN
487  t = a( k-1, k )
488  a( k-1, k ) = a( kp, k )
489  a( kp, k ) = t
490  END IF
491 *
492 * Convert upper triangle of A into U form by applying
493 * the interchanges in columns k+1:N.
494 *
495  IF( k.LT.n )
496  $ CALL sswap( n-k, a( kk, k+1 ), lda, a( kp, k+1 ),
497  $ lda )
498 *
499  END IF
500 *
501 * Update the leading submatrix
502 *
503  IF( kstep.EQ.1 ) THEN
504 *
505 * 1-by-1 pivot block D(k): column k now holds
506 *
507 * W(k) = U(k)*D(k)
508 *
509 * where U(k) is the k-th column of U
510 *
511  IF( k.GT.1 ) THEN
512 *
513 * Perform a rank-1 update of A(1:k-1,1:k-1) and
514 * store U(k) in column k
515 *
516  IF( abs( a( k, k ) ).GE.sfmin ) THEN
517 *
518 * Perform a rank-1 update of A(1:k-1,1:k-1) as
519 * A := A - U(k)*D(k)*U(k)**T
520 * = A - W(k)*1/D(k)*W(k)**T
521 *
522  d11 = one / a( k, k )
523  CALL ssyr( uplo, k-1, -d11, a( 1, k ), 1, a, lda )
524 *
525 * Store U(k) in column k
526 *
527  CALL sscal( k-1, d11, a( 1, k ), 1 )
528  ELSE
529 *
530 * Store L(k) in column K
531 *
532  d11 = a( k, k )
533  DO 16 ii = 1, k - 1
534  a( ii, k ) = a( ii, k ) / d11
535  16 CONTINUE
536 *
537 * Perform a rank-1 update of A(k+1:n,k+1:n) as
538 * A := A - U(k)*D(k)*U(k)**T
539 * = A - W(k)*(1/D(k))*W(k)**T
540 * = A - (W(k)/D(k))*(D(k))*(W(k)/D(K))**T
541 *
542  CALL ssyr( uplo, k-1, -d11, a( 1, k ), 1, a, lda )
543  END IF
544 *
545 * Store the superdiagonal element of D in array E
546 *
547  e( k ) = zero
548 *
549  END IF
550 *
551  ELSE
552 *
553 * 2-by-2 pivot block D(k): columns k and k-1 now hold
554 *
555 * ( W(k-1) W(k) ) = ( U(k-1) U(k) )*D(k)
556 *
557 * where U(k) and U(k-1) are the k-th and (k-1)-th columns
558 * of U
559 *
560 * Perform a rank-2 update of A(1:k-2,1:k-2) as
561 *
562 * A := A - ( U(k-1) U(k) )*D(k)*( U(k-1) U(k) )**T
563 * = A - ( ( A(k-1)A(k) )*inv(D(k)) ) * ( A(k-1)A(k) )**T
564 *
565 * and store L(k) and L(k+1) in columns k and k+1
566 *
567  IF( k.GT.2 ) THEN
568 *
569  d12 = a( k-1, k )
570  d22 = a( k-1, k-1 ) / d12
571  d11 = a( k, k ) / d12
572  t = one / ( d11*d22-one )
573 *
574  DO 30 j = k - 2, 1, -1
575 *
576  wkm1 = t*( d11*a( j, k-1 )-a( j, k ) )
577  wk = t*( d22*a( j, k )-a( j, k-1 ) )
578 *
579  DO 20 i = j, 1, -1
580  a( i, j ) = a( i, j ) - (a( i, k ) / d12 )*wk -
581  $ ( a( i, k-1 ) / d12 )*wkm1
582  20 CONTINUE
583 *
584 * Store U(k) and U(k-1) in cols k and k-1 for row J
585 *
586  a( j, k ) = wk / d12
587  a( j, k-1 ) = wkm1 / d12
588 *
589  30 CONTINUE
590 *
591  END IF
592 *
593 * Copy superdiagonal elements of D(K) to E(K) and
594 * ZERO out superdiagonal entry of A
595 *
596  e( k ) = a( k-1, k )
597  e( k-1 ) = zero
598  a( k-1, k ) = zero
599 *
600  END IF
601 *
602 * End column K is nonsingular
603 *
604  END IF
605 *
606 * Store details of the interchanges in IPIV
607 *
608  IF( kstep.EQ.1 ) THEN
609  ipiv( k ) = kp
610  ELSE
611  ipiv( k ) = -p
612  ipiv( k-1 ) = -kp
613  END IF
614 *
615 * Decrease K and return to the start of the main loop
616 *
617  k = k - kstep
618  GO TO 10
619 *
620  34 CONTINUE
621 *
622  ELSE
623 *
624 * Factorize A as L*D*L**T using the lower triangle of A
625 *
626 * Initilize the unused last entry of the subdiagonal array E.
627 *
628  e( n ) = zero
629 *
630 * K is the main loop index, increasing from 1 to N in steps of
631 * 1 or 2
632 *
633  k = 1
634  40 CONTINUE
635 *
636 * If K > N, exit from loop
637 *
638  IF( k.GT.n )
639  $ GO TO 64
640  kstep = 1
641  p = k
642 *
643 * Determine rows and columns to be interchanged and whether
644 * a 1-by-1 or 2-by-2 pivot block will be used
645 *
646  absakk = abs( a( k, k ) )
647 *
648 * IMAX is the row-index of the largest off-diagonal element in
649 * column K, and COLMAX is its absolute value.
650 * Determine both COLMAX and IMAX.
651 *
652  IF( k.LT.n ) THEN
653  imax = k + isamax( n-k, a( k+1, k ), 1 )
654  colmax = abs( a( imax, k ) )
655  ELSE
656  colmax = zero
657  END IF
658 *
659  IF( ( max( absakk, colmax ).EQ.zero ) ) THEN
660 *
661 * Column K is zero or underflow: set INFO and continue
662 *
663  IF( info.EQ.0 )
664  $ info = k
665  kp = k
666 *
667 * Set E( K ) to zero
668 *
669  IF( k.LT.n )
670  $ e( k ) = zero
671 *
672  ELSE
673 *
674 * Test for interchange
675 *
676 * Equivalent to testing for (used to handle NaN and Inf)
677 * ABSAKK.GE.ALPHA*COLMAX
678 *
679  IF( .NOT.( absakk.LT.alpha*colmax ) ) THEN
680 *
681 * no interchange, use 1-by-1 pivot block
682 *
683  kp = k
684 *
685  ELSE
686 *
687  done = .false.
688 *
689 * Loop until pivot found
690 *
691  42 CONTINUE
692 *
693 * Begin pivot search loop body
694 *
695 * JMAX is the column-index of the largest off-diagonal
696 * element in row IMAX, and ROWMAX is its absolute value.
697 * Determine both ROWMAX and JMAX.
698 *
699  IF( imax.NE.k ) THEN
700  jmax = k - 1 + isamax( imax-k, a( imax, k ), lda )
701  rowmax = abs( a( imax, jmax ) )
702  ELSE
703  rowmax = zero
704  END IF
705 *
706  IF( imax.LT.n ) THEN
707  itemp = imax + isamax( n-imax, a( imax+1, imax ),
708  $ 1 )
709  stemp = abs( a( itemp, imax ) )
710  IF( stemp.GT.rowmax ) THEN
711  rowmax = stemp
712  jmax = itemp
713  END IF
714  END IF
715 *
716 * Equivalent to testing for (used to handle NaN and Inf)
717 * ABS( A( IMAX, IMAX ) ).GE.ALPHA*ROWMAX
718 *
719  IF( .NOT.( abs( a( imax, imax ) ).LT.alpha*rowmax ) )
720  $ THEN
721 *
722 * interchange rows and columns K and IMAX,
723 * use 1-by-1 pivot block
724 *
725  kp = imax
726  done = .true.
727 *
728 * Equivalent to testing for ROWMAX .EQ. COLMAX,
729 * used to handle NaN and Inf
730 *
731  ELSE IF( ( p.EQ.jmax ).OR.( rowmax.LE.colmax ) ) THEN
732 *
733 * interchange rows and columns K+1 and IMAX,
734 * use 2-by-2 pivot block
735 *
736  kp = imax
737  kstep = 2
738  done = .true.
739  ELSE
740 *
741 * Pivot NOT found, set variables and repeat
742 *
743  p = imax
744  colmax = rowmax
745  imax = jmax
746  END IF
747 *
748 * End pivot search loop body
749 *
750  IF( .NOT. done ) GOTO 42
751 *
752  END IF
753 *
754 * Swap TWO rows and TWO columns
755 *
756 * First swap
757 *
758  IF( ( kstep.EQ.2 ) .AND. ( p.NE.k ) ) THEN
759 *
760 * Interchange rows and column K and P in the trailing
761 * submatrix A(k:n,k:n) if we have a 2-by-2 pivot
762 *
763  IF( p.LT.n )
764  $ CALL sswap( n-p, a( p+1, k ), 1, a( p+1, p ), 1 )
765  IF( p.GT.(k+1) )
766  $ CALL sswap( p-k-1, a( k+1, k ), 1, a( p, k+1 ), lda )
767  t = a( k, k )
768  a( k, k ) = a( p, p )
769  a( p, p ) = t
770 *
771 * Convert lower triangle of A into L form by applying
772 * the interchanges in columns 1:k-1.
773 *
774  IF ( k.GT.1 )
775  $ CALL sswap( k-1, a( k, 1 ), lda, a( p, 1 ), lda )
776 *
777  END IF
778 *
779 * Second swap
780 *
781  kk = k + kstep - 1
782  IF( kp.NE.kk ) THEN
783 *
784 * Interchange rows and columns KK and KP in the trailing
785 * submatrix A(k:n,k:n)
786 *
787  IF( kp.LT.n )
788  $ CALL sswap( n-kp, a( kp+1, kk ), 1, a( kp+1, kp ), 1 )
789  IF( ( kk.LT.n ) .AND. ( kp.GT.(kk+1) ) )
790  $ CALL sswap( kp-kk-1, a( kk+1, kk ), 1, a( kp, kk+1 ),
791  $ lda )
792  t = a( kk, kk )
793  a( kk, kk ) = a( kp, kp )
794  a( kp, kp ) = t
795  IF( kstep.EQ.2 ) THEN
796  t = a( k+1, k )
797  a( k+1, k ) = a( kp, k )
798  a( kp, k ) = t
799  END IF
800 *
801 * Convert lower triangle of A into L form by applying
802 * the interchanges in columns 1:k-1.
803 *
804  IF ( k.GT.1 )
805  $ CALL sswap( k-1, a( kk, 1 ), lda, a( kp, 1 ), lda )
806 *
807  END IF
808 *
809 * Update the trailing submatrix
810 *
811  IF( kstep.EQ.1 ) THEN
812 *
813 * 1-by-1 pivot block D(k): column k now holds
814 *
815 * W(k) = L(k)*D(k)
816 *
817 * where L(k) is the k-th column of L
818 *
819  IF( k.LT.n ) THEN
820 *
821 * Perform a rank-1 update of A(k+1:n,k+1:n) and
822 * store L(k) in column k
823 *
824  IF( abs( a( k, k ) ).GE.sfmin ) THEN
825 *
826 * Perform a rank-1 update of A(k+1:n,k+1:n) as
827 * A := A - L(k)*D(k)*L(k)**T
828 * = A - W(k)*(1/D(k))*W(k)**T
829 *
830  d11 = one / a( k, k )
831  CALL ssyr( uplo, n-k, -d11, a( k+1, k ), 1,
832  $ a( k+1, k+1 ), lda )
833 *
834 * Store L(k) in column k
835 *
836  CALL sscal( n-k, d11, a( k+1, k ), 1 )
837  ELSE
838 *
839 * Store L(k) in column k
840 *
841  d11 = a( k, k )
842  DO 46 ii = k + 1, n
843  a( ii, k ) = a( ii, k ) / d11
844  46 CONTINUE
845 *
846 * Perform a rank-1 update of A(k+1:n,k+1:n) as
847 * A := A - L(k)*D(k)*L(k)**T
848 * = A - W(k)*(1/D(k))*W(k)**T
849 * = A - (W(k)/D(k))*(D(k))*(W(k)/D(K))**T
850 *
851  CALL ssyr( uplo, n-k, -d11, a( k+1, k ), 1,
852  $ a( k+1, k+1 ), lda )
853  END IF
854 *
855 * Store the subdiagonal element of D in array E
856 *
857  e( k ) = zero
858 *
859  END IF
860 *
861  ELSE
862 *
863 * 2-by-2 pivot block D(k): columns k and k+1 now hold
864 *
865 * ( W(k) W(k+1) ) = ( L(k) L(k+1) )*D(k)
866 *
867 * where L(k) and L(k+1) are the k-th and (k+1)-th columns
868 * of L
869 *
870 *
871 * Perform a rank-2 update of A(k+2:n,k+2:n) as
872 *
873 * A := A - ( L(k) L(k+1) ) * D(k) * ( L(k) L(k+1) )**T
874 * = A - ( ( A(k)A(k+1) )*inv(D(k) ) * ( A(k)A(k+1) )**T
875 *
876 * and store L(k) and L(k+1) in columns k and k+1
877 *
878  IF( k.LT.n-1 ) THEN
879 *
880  d21 = a( k+1, k )
881  d11 = a( k+1, k+1 ) / d21
882  d22 = a( k, k ) / d21
883  t = one / ( d11*d22-one )
884 *
885  DO 60 j = k + 2, n
886 *
887 * Compute D21 * ( W(k)W(k+1) ) * inv(D(k)) for row J
888 *
889  wk = t*( d11*a( j, k )-a( j, k+1 ) )
890  wkp1 = t*( d22*a( j, k+1 )-a( j, k ) )
891 *
892 * Perform a rank-2 update of A(k+2:n,k+2:n)
893 *
894  DO 50 i = j, n
895  a( i, j ) = a( i, j ) - ( a( i, k ) / d21 )*wk -
896  $ ( a( i, k+1 ) / d21 )*wkp1
897  50 CONTINUE
898 *
899 * Store L(k) and L(k+1) in cols k and k+1 for row J
900 *
901  a( j, k ) = wk / d21
902  a( j, k+1 ) = wkp1 / d21
903 *
904  60 CONTINUE
905 *
906  END IF
907 *
908 * Copy subdiagonal elements of D(K) to E(K) and
909 * ZERO out subdiagonal entry of A
910 *
911  e( k ) = a( k+1, k )
912  e( k+1 ) = zero
913  a( k+1, k ) = zero
914 *
915  END IF
916 *
917 * End column K is nonsingular
918 *
919  END IF
920 *
921 * Store details of the interchanges in IPIV
922 *
923  IF( kstep.EQ.1 ) THEN
924  ipiv( k ) = kp
925  ELSE
926  ipiv( k ) = -p
927  ipiv( k+1 ) = -kp
928  END IF
929 *
930 * Increase K and return to the start of the main loop
931 *
932  k = k + kstep
933  GO TO 40
934 *
935  64 CONTINUE
936 *
937  END IF
938 *
939  RETURN
940 *
941 * End of SSYTF2_RK
942 *
943  END
subroutine ssytf2_rk(UPLO, N, A, LDA, E, IPIV, INFO)
SSYTF2_RK computes the factorization of a real symmetric indefinite matrix using the bounded Bunch-Ka...
Definition: ssytf2_rk.f:243
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine sscal(N, SA, SX, INCX)
SSCAL
Definition: sscal.f:55
subroutine sswap(N, SX, INCX, SY, INCY)
SSWAP
Definition: sswap.f:53
subroutine ssyr(UPLO, N, ALPHA, X, INCX, A, LDA)
SSYR
Definition: ssyr.f:134