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