LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
dlasyf.f
Go to the documentation of this file.
1 *> \brief \b DLASYF computes a partial factorization of a real symmetric matrix, using the 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 + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/dlasyf.f">
11 *> [TGZ]</a>
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/dlasyf.f">
13 *> [ZIP]</a>
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/dlasyf.f">
15 *> [TXT]</a>
16 *> \endhtmlonly
17 *
18 * Definition:
19 * ===========
20 *
21 * SUBROUTINE DLASYF( UPLO, N, NB, KB, A, LDA, IPIV, W, LDW, INFO )
22 *
23 * .. Scalar Arguments ..
24 * CHARACTER UPLO
25 * INTEGER INFO, KB, LDA, LDW, N, NB
26 * ..
27 * .. Array Arguments ..
28 * INTEGER IPIV( * )
29 * DOUBLE PRECISION A( LDA, * ), W( LDW, * )
30 * ..
31 *
32 *
33 *> \par Purpose:
34 * =============
35 *>
36 *> \verbatim
37 *>
38 *> DLASYF computes a partial factorization of a real symmetric matrix A
39 *> using the Bunch-Kaufman diagonal pivoting method. The partial
40 *> 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 is an auxiliary routine called by DSYTRF. It uses blocked code
52 *> (calling Level 3 BLAS) to update the submatrix A11 (if UPLO = 'U') or
53 *> 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. If UPLO = 'U', the leading
93 *> n-by-n upper triangular part of A contains the upper
94 *> triangular part of the matrix A, and the strictly lower
95 *> triangular part of A is not referenced. If UPLO = 'L', the
96 *> leading n-by-n lower triangular part of A contains the lower
97 *> triangular part of the matrix A, and the strictly upper
98 *> triangular part of A is not referenced.
99 *> On exit, A contains details of the partial factorization.
100 *> \endverbatim
101 *>
102 *> \param[in] LDA
103 *> \verbatim
104 *> LDA is INTEGER
105 *> The leading dimension of the array A. LDA >= max(1,N).
106 *> \endverbatim
107 *>
108 *> \param[out] IPIV
109 *> \verbatim
110 *> IPIV is INTEGER array, dimension (N)
111 *> Details of the interchanges and the block structure of D.
112 *> If UPLO = 'U', only the last KB elements of IPIV are set;
113 *> if UPLO = 'L', only the first KB elements are set.
114 *>
115 *> If IPIV(k) > 0, then rows and columns k and IPIV(k) were
116 *> interchanged and D(k,k) is a 1-by-1 diagonal block.
117 *> If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0, then rows and
118 *> columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k)
119 *> is a 2-by-2 diagonal block. If UPLO = 'L' and IPIV(k) =
120 *> IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were
121 *> interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
122 *> \endverbatim
123 *>
124 *> \param[out] W
125 *> \verbatim
126 *> W is DOUBLE PRECISION array, dimension (LDW,NB)
127 *> \endverbatim
128 *>
129 *> \param[in] LDW
130 *> \verbatim
131 *> LDW is INTEGER
132 *> The leading dimension of the array W. LDW >= max(1,N).
133 *> \endverbatim
134 *>
135 *> \param[out] INFO
136 *> \verbatim
137 *> INFO is INTEGER
138 *> = 0: successful exit
139 *> > 0: if INFO = k, D(k,k) is exactly zero. The factorization
140 *> has been completed, but the block diagonal matrix D is
141 *> exactly singular.
142 *> \endverbatim
143 *
144 * Authors:
145 * ========
146 *
147 *> \author Univ. of Tennessee
148 *> \author Univ. of California Berkeley
149 *> \author Univ. of Colorado Denver
150 *> \author NAG Ltd.
151 *
152 *> \date September 2012
153 *
154 *> \ingroup doubleSYcomputational
155 *
156 * =====================================================================
157  SUBROUTINE dlasyf( UPLO, N, NB, KB, A, LDA, IPIV, W, LDW, INFO )
158 *
159 * -- LAPACK computational routine (version 3.4.2) --
160 * -- LAPACK is a software package provided by Univ. of Tennessee, --
161 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
162 * September 2012
163 *
164 * .. Scalar Arguments ..
165  CHARACTER uplo
166  INTEGER info, kb, lda, ldw, n, nb
167 * ..
168 * .. Array Arguments ..
169  INTEGER ipiv( * )
170  DOUBLE PRECISION a( lda, * ), w( ldw, * )
171 * ..
172 *
173 * =====================================================================
174 *
175 * .. Parameters ..
176  DOUBLE PRECISION zero, one
177  parameter( zero = 0.0d+0, one = 1.0d+0 )
178  DOUBLE PRECISION eight, sevten
179  parameter( eight = 8.0d+0, sevten = 17.0d+0 )
180 * ..
181 * .. Local Scalars ..
182  INTEGER imax, j, jb, jj, jmax, jp, k, kk, kkw, kp,
183  $ kstep, kw
184  DOUBLE PRECISION absakk, alpha, colmax, d11, d21, d22, r1,
185  $ rowmax, t
186 * ..
187 * .. External Functions ..
188  LOGICAL lsame
189  INTEGER idamax
190  EXTERNAL lsame, idamax
191 * ..
192 * .. External Subroutines ..
193  EXTERNAL dcopy, dgemm, dgemv, dscal, dswap
194 * ..
195 * .. Intrinsic Functions ..
196  INTRINSIC abs, max, min, sqrt
197 * ..
198 * .. Executable Statements ..
199 *
200  info = 0
201 *
202 * Initialize ALPHA for use in choosing pivot block size.
203 *
204  alpha = ( one+sqrt( sevten ) ) / eight
205 *
206  IF( lsame( uplo, 'U' ) ) THEN
207 *
208 * Factorize the trailing columns of A using the upper triangle
209 * of A and working backwards, and compute the matrix W = U12*D
210 * for use in updating A11
211 *
212 * K is the main loop index, decreasing from N in steps of 1 or 2
213 *
214 * KW is the column of W which corresponds to column K of A
215 *
216  k = n
217  10 continue
218  kw = nb + k - n
219 *
220 * Exit from loop
221 *
222  IF( ( k.LE.n-nb+1 .AND. nb.LT.n ) .OR. k.LT.1 )
223  $ go to 30
224 *
225 * Copy column K of A to column KW of W and update it
226 *
227  CALL dcopy( k, a( 1, k ), 1, w( 1, kw ), 1 )
228  IF( k.LT.n )
229  $ CALL dgemv( 'No transpose', k, n-k, -one, a( 1, k+1 ), lda,
230  $ w( k, kw+1 ), ldw, one, w( 1, kw ), 1 )
231 *
232  kstep = 1
233 *
234 * Determine rows and columns to be interchanged and whether
235 * a 1-by-1 or 2-by-2 pivot block will be used
236 *
237  absakk = abs( w( k, kw ) )
238 *
239 * IMAX is the row-index of the largest off-diagonal element in
240 * column K, and COLMAX is its absolute value
241 *
242  IF( k.GT.1 ) THEN
243  imax = idamax( k-1, w( 1, kw ), 1 )
244  colmax = abs( w( imax, kw ) )
245  ELSE
246  colmax = zero
247  END IF
248 *
249  IF( max( absakk, colmax ).EQ.zero ) THEN
250 *
251 * Column K is zero: set INFO and continue
252 *
253  IF( info.EQ.0 )
254  $ info = k
255  kp = k
256  ELSE
257  IF( absakk.GE.alpha*colmax ) THEN
258 *
259 * no interchange, use 1-by-1 pivot block
260 *
261  kp = k
262  ELSE
263 *
264 * Copy column IMAX to column KW-1 of W and update it
265 *
266  CALL dcopy( imax, a( 1, imax ), 1, w( 1, kw-1 ), 1 )
267  CALL dcopy( k-imax, a( imax, imax+1 ), lda,
268  $ w( imax+1, kw-1 ), 1 )
269  IF( k.LT.n )
270  $ CALL dgemv( 'No transpose', k, n-k, -one, a( 1, k+1 ),
271  $ lda, w( imax, kw+1 ), ldw, one,
272  $ w( 1, kw-1 ), 1 )
273 *
274 * JMAX is the column-index of the largest off-diagonal
275 * element in row IMAX, and ROWMAX is its absolute value
276 *
277  jmax = imax + idamax( k-imax, w( imax+1, kw-1 ), 1 )
278  rowmax = abs( w( jmax, kw-1 ) )
279  IF( imax.GT.1 ) THEN
280  jmax = idamax( imax-1, w( 1, kw-1 ), 1 )
281  rowmax = max( rowmax, abs( w( jmax, kw-1 ) ) )
282  END IF
283 *
284  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
285 *
286 * no interchange, use 1-by-1 pivot block
287 *
288  kp = k
289  ELSE IF( abs( w( imax, kw-1 ) ).GE.alpha*rowmax ) THEN
290 *
291 * interchange rows and columns K and IMAX, use 1-by-1
292 * pivot block
293 *
294  kp = imax
295 *
296 * copy column KW-1 of W to column KW
297 *
298  CALL dcopy( k, w( 1, kw-1 ), 1, w( 1, kw ), 1 )
299  ELSE
300 *
301 * interchange rows and columns K-1 and IMAX, use 2-by-2
302 * pivot block
303 *
304  kp = imax
305  kstep = 2
306  END IF
307  END IF
308 *
309  kk = k - kstep + 1
310  kkw = nb + kk - n
311 *
312 * Updated column KP is already stored in column KKW of W
313 *
314  IF( kp.NE.kk ) THEN
315 *
316 * Copy non-updated column KK to column KP
317 *
318  a( kp, k ) = a( kk, k )
319  CALL dcopy( k-1-kp, a( kp+1, kk ), 1, a( kp, kp+1 ),
320  $ lda )
321  CALL dcopy( kp, a( 1, kk ), 1, a( 1, kp ), 1 )
322 *
323 * Interchange rows KK and KP in last KK columns of A and W
324 *
325  CALL dswap( n-kk+1, a( kk, kk ), lda, a( kp, kk ), lda )
326  CALL dswap( n-kk+1, w( kk, kkw ), ldw, w( kp, kkw ),
327  $ ldw )
328  END IF
329 *
330  IF( kstep.EQ.1 ) THEN
331 *
332 * 1-by-1 pivot block D(k): column KW of W now holds
333 *
334 * W(k) = U(k)*D(k)
335 *
336 * where U(k) is the k-th column of U
337 *
338 * Store U(k) in column k of A
339 *
340  CALL dcopy( k, w( 1, kw ), 1, a( 1, k ), 1 )
341  r1 = one / a( k, k )
342  CALL dscal( k-1, r1, a( 1, k ), 1 )
343  ELSE
344 *
345 * 2-by-2 pivot block D(k): columns KW and KW-1 of W now
346 * hold
347 *
348 * ( W(k-1) W(k) ) = ( U(k-1) U(k) )*D(k)
349 *
350 * where U(k) and U(k-1) are the k-th and (k-1)-th columns
351 * of U
352 *
353  IF( k.GT.2 ) THEN
354 *
355 * Store U(k) and U(k-1) in columns k and k-1 of A
356 *
357  d21 = w( k-1, kw )
358  d11 = w( k, kw ) / d21
359  d22 = w( k-1, kw-1 ) / d21
360  t = one / ( d11*d22-one )
361  d21 = t / d21
362  DO 20 j = 1, k - 2
363  a( j, k-1 ) = d21*( d11*w( j, kw-1 )-w( j, kw ) )
364  a( j, k ) = d21*( d22*w( j, kw )-w( j, kw-1 ) )
365  20 continue
366  END IF
367 *
368 * Copy D(k) to A
369 *
370  a( k-1, k-1 ) = w( k-1, kw-1 )
371  a( k-1, k ) = w( k-1, kw )
372  a( k, k ) = w( k, kw )
373  END IF
374  END IF
375 *
376 * Store details of the interchanges in IPIV
377 *
378  IF( kstep.EQ.1 ) THEN
379  ipiv( k ) = kp
380  ELSE
381  ipiv( k ) = -kp
382  ipiv( k-1 ) = -kp
383  END IF
384 *
385 * Decrease K and return to the start of the main loop
386 *
387  k = k - kstep
388  go to 10
389 *
390  30 continue
391 *
392 * Update the upper triangle of A11 (= A(1:k,1:k)) as
393 *
394 * A11 := A11 - U12*D*U12**T = A11 - U12*W**T
395 *
396 * computing blocks of NB columns at a time
397 *
398  DO 50 j = ( ( k-1 ) / nb )*nb + 1, 1, -nb
399  jb = min( nb, k-j+1 )
400 *
401 * Update the upper triangle of the diagonal block
402 *
403  DO 40 jj = j, j + jb - 1
404  CALL dgemv( 'No transpose', jj-j+1, n-k, -one,
405  $ a( j, k+1 ), lda, w( jj, kw+1 ), ldw, one,
406  $ a( j, jj ), 1 )
407  40 continue
408 *
409 * Update the rectangular superdiagonal block
410 *
411  CALL dgemm( 'No transpose', 'Transpose', j-1, jb, n-k, -one,
412  $ a( 1, k+1 ), lda, w( j, kw+1 ), ldw, one,
413  $ a( 1, j ), lda )
414  50 continue
415 *
416 * Put U12 in standard form by partially undoing the interchanges
417 * in columns k+1:n
418 *
419  j = k + 1
420  60 continue
421  jj = j
422  jp = ipiv( j )
423  IF( jp.LT.0 ) THEN
424  jp = -jp
425  j = j + 1
426  END IF
427  j = j + 1
428  IF( jp.NE.jj .AND. j.LE.n )
429  $ CALL dswap( n-j+1, a( jp, j ), lda, a( jj, j ), lda )
430  IF( j.LE.n )
431  $ go to 60
432 *
433 * Set KB to the number of columns factorized
434 *
435  kb = n - k
436 *
437  ELSE
438 *
439 * Factorize the leading columns of A using the lower triangle
440 * of A and working forwards, and compute the matrix W = L21*D
441 * for use in updating A22
442 *
443 * K is the main loop index, increasing from 1 in steps of 1 or 2
444 *
445  k = 1
446  70 continue
447 *
448 * Exit from loop
449 *
450  IF( ( k.GE.nb .AND. nb.LT.n ) .OR. k.GT.n )
451  $ go to 90
452 *
453 * Copy column K of A to column K of W and update it
454 *
455  CALL dcopy( n-k+1, a( k, k ), 1, w( k, k ), 1 )
456  CALL dgemv( 'No transpose', n-k+1, k-1, -one, a( k, 1 ), lda,
457  $ w( k, 1 ), ldw, one, w( k, k ), 1 )
458 *
459  kstep = 1
460 *
461 * Determine rows and columns to be interchanged and whether
462 * a 1-by-1 or 2-by-2 pivot block will be used
463 *
464  absakk = abs( w( k, k ) )
465 *
466 * IMAX is the row-index of the largest off-diagonal element in
467 * column K, and COLMAX is its absolute value
468 *
469  IF( k.LT.n ) THEN
470  imax = k + idamax( n-k, w( k+1, k ), 1 )
471  colmax = abs( w( imax, k ) )
472  ELSE
473  colmax = zero
474  END IF
475 *
476  IF( max( absakk, colmax ).EQ.zero ) THEN
477 *
478 * Column K is zero: set INFO and continue
479 *
480  IF( info.EQ.0 )
481  $ info = k
482  kp = k
483  ELSE
484  IF( absakk.GE.alpha*colmax ) THEN
485 *
486 * no interchange, use 1-by-1 pivot block
487 *
488  kp = k
489  ELSE
490 *
491 * Copy column IMAX to column K+1 of W and update it
492 *
493  CALL dcopy( imax-k, a( imax, k ), lda, w( k, k+1 ), 1 )
494  CALL dcopy( n-imax+1, a( imax, imax ), 1, w( imax, k+1 ),
495  $ 1 )
496  CALL dgemv( 'No transpose', n-k+1, k-1, -one, a( k, 1 ),
497  $ lda, w( imax, 1 ), ldw, one, w( k, k+1 ), 1 )
498 *
499 * JMAX is the column-index of the largest off-diagonal
500 * element in row IMAX, and ROWMAX is its absolute value
501 *
502  jmax = k - 1 + idamax( imax-k, w( k, k+1 ), 1 )
503  rowmax = abs( w( jmax, k+1 ) )
504  IF( imax.LT.n ) THEN
505  jmax = imax + idamax( n-imax, w( imax+1, k+1 ), 1 )
506  rowmax = max( rowmax, abs( w( jmax, k+1 ) ) )
507  END IF
508 *
509  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
510 *
511 * no interchange, use 1-by-1 pivot block
512 *
513  kp = k
514  ELSE IF( abs( w( imax, k+1 ) ).GE.alpha*rowmax ) THEN
515 *
516 * interchange rows and columns K and IMAX, use 1-by-1
517 * pivot block
518 *
519  kp = imax
520 *
521 * copy column K+1 of W to column K
522 *
523  CALL dcopy( n-k+1, w( k, k+1 ), 1, w( k, k ), 1 )
524  ELSE
525 *
526 * interchange rows and columns K+1 and IMAX, use 2-by-2
527 * pivot block
528 *
529  kp = imax
530  kstep = 2
531  END IF
532  END IF
533 *
534  kk = k + kstep - 1
535 *
536 * Updated column KP is already stored in column KK of W
537 *
538  IF( kp.NE.kk ) THEN
539 *
540 * Copy non-updated column KK to column KP
541 *
542  a( kp, k ) = a( kk, k )
543  CALL dcopy( kp-k-1, a( k+1, kk ), 1, a( kp, k+1 ), lda )
544  CALL dcopy( n-kp+1, a( kp, kk ), 1, a( kp, kp ), 1 )
545 *
546 * Interchange rows KK and KP in first KK columns of A and W
547 *
548  CALL dswap( kk, a( kk, 1 ), lda, a( kp, 1 ), lda )
549  CALL dswap( kk, w( kk, 1 ), ldw, w( kp, 1 ), ldw )
550  END IF
551 *
552  IF( kstep.EQ.1 ) THEN
553 *
554 * 1-by-1 pivot block D(k): column k of W now holds
555 *
556 * W(k) = L(k)*D(k)
557 *
558 * where L(k) is the k-th column of L
559 *
560 * Store L(k) in column k of A
561 *
562  CALL dcopy( n-k+1, w( k, k ), 1, a( k, k ), 1 )
563  IF( k.LT.n ) THEN
564  r1 = one / a( k, k )
565  CALL dscal( n-k, r1, a( k+1, k ), 1 )
566  END IF
567  ELSE
568 *
569 * 2-by-2 pivot block D(k): columns k and k+1 of W now hold
570 *
571 * ( W(k) W(k+1) ) = ( L(k) L(k+1) )*D(k)
572 *
573 * where L(k) and L(k+1) are the k-th and (k+1)-th columns
574 * of L
575 *
576  IF( k.LT.n-1 ) THEN
577 *
578 * Store L(k) and L(k+1) in columns k and k+1 of A
579 *
580  d21 = w( k+1, k )
581  d11 = w( k+1, k+1 ) / d21
582  d22 = w( k, k ) / d21
583  t = one / ( d11*d22-one )
584  d21 = t / d21
585  DO 80 j = k + 2, n
586  a( j, k ) = d21*( d11*w( j, k )-w( j, k+1 ) )
587  a( j, k+1 ) = d21*( d22*w( j, k+1 )-w( j, k ) )
588  80 continue
589  END IF
590 *
591 * Copy D(k) to A
592 *
593  a( k, k ) = w( k, k )
594  a( k+1, k ) = w( k+1, k )
595  a( k+1, k+1 ) = w( k+1, k+1 )
596  END IF
597  END IF
598 *
599 * Store details of the interchanges in IPIV
600 *
601  IF( kstep.EQ.1 ) THEN
602  ipiv( k ) = kp
603  ELSE
604  ipiv( k ) = -kp
605  ipiv( k+1 ) = -kp
606  END IF
607 *
608 * Increase K and return to the start of the main loop
609 *
610  k = k + kstep
611  go to 70
612 *
613  90 continue
614 *
615 * Update the lower triangle of A22 (= A(k:n,k:n)) as
616 *
617 * A22 := A22 - L21*D*L21**T = A22 - L21*W**T
618 *
619 * computing blocks of NB columns at a time
620 *
621  DO 110 j = k, n, nb
622  jb = min( nb, n-j+1 )
623 *
624 * Update the lower triangle of the diagonal block
625 *
626  DO 100 jj = j, j + jb - 1
627  CALL dgemv( 'No transpose', j+jb-jj, k-1, -one,
628  $ a( jj, 1 ), lda, w( jj, 1 ), ldw, one,
629  $ a( jj, jj ), 1 )
630  100 continue
631 *
632 * Update the rectangular subdiagonal block
633 *
634  IF( j+jb.LE.n )
635  $ CALL dgemm( 'No transpose', 'Transpose', n-j-jb+1, jb,
636  $ k-1, -one, a( j+jb, 1 ), lda, w( j, 1 ), ldw,
637  $ one, a( j+jb, j ), lda )
638  110 continue
639 *
640 * Put L21 in standard form by partially undoing the interchanges
641 * in columns 1:k-1
642 *
643  j = k - 1
644  120 continue
645  jj = j
646  jp = ipiv( j )
647  IF( jp.LT.0 ) THEN
648  jp = -jp
649  j = j - 1
650  END IF
651  j = j - 1
652  IF( jp.NE.jj .AND. j.GE.1 )
653  $ CALL dswap( j, a( jp, 1 ), lda, a( jj, 1 ), lda )
654  IF( j.GE.1 )
655  $ go to 120
656 *
657 * Set KB to the number of columns factorized
658 *
659  kb = k - 1
660 *
661  END IF
662  return
663 *
664 * End of DLASYF
665 *
666  END