LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
dlasyf.f
Go to the documentation of this file.
1 *> \brief \b DLASYF computes a partial factorization of a real symmetric matrix using the Bunch-Kaufman 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 *>
113 *> If UPLO = 'U':
114 *> Only the last KB elements of IPIV are set.
115 *>
116 *> If IPIV(k) > 0, then rows and columns k and IPIV(k) were
117 *> interchanged and D(k,k) is a 1-by-1 diagonal block.
118 *>
119 *> If IPIV(k) = IPIV(k-1) < 0, then rows and columns
120 *> k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k)
121 *> is a 2-by-2 diagonal block.
122 *>
123 *> If UPLO = 'L':
124 *> Only the first KB elements of IPIV are set.
125 *>
126 *> If IPIV(k) > 0, then rows and columns k and IPIV(k) were
127 *> interchanged and D(k,k) is a 1-by-1 diagonal block.
128 *>
129 *> If IPIV(k) = IPIV(k+1) < 0, then rows and columns
130 *> k+1 and -IPIV(k) were interchanged and D(k:k+1,k:k+1)
131 *> is a 2-by-2 diagonal block.
132 *> \endverbatim
133 *>
134 *> \param[out] W
135 *> \verbatim
136 *> W is DOUBLE PRECISION array, dimension (LDW,NB)
137 *> \endverbatim
138 *>
139 *> \param[in] LDW
140 *> \verbatim
141 *> LDW is INTEGER
142 *> The leading dimension of the array W. LDW >= max(1,N).
143 *> \endverbatim
144 *>
145 *> \param[out] INFO
146 *> \verbatim
147 *> INFO is INTEGER
148 *> = 0: successful exit
149 *> > 0: if INFO = k, D(k,k) is exactly zero. The factorization
150 *> has been completed, but the block diagonal matrix D is
151 *> exactly singular.
152 *> \endverbatim
153 *
154 * Authors:
155 * ========
156 *
157 *> \author Univ. of Tennessee
158 *> \author Univ. of California Berkeley
159 *> \author Univ. of Colorado Denver
160 *> \author NAG Ltd.
161 *
162 *> \date November 2013
163 *
164 *> \ingroup doubleSYcomputational
165 *
166 *> \par Contributors:
167 * ==================
168 *>
169 *> \verbatim
170 *>
171 *> November 2013, Igor Kozachenko,
172 *> Computer Science Division,
173 *> University of California, Berkeley
174 *> \endverbatim
175 *
176 * =====================================================================
177  SUBROUTINE dlasyf( UPLO, N, NB, KB, A, LDA, IPIV, W, LDW, INFO )
178 *
179 * -- LAPACK computational routine (version 3.5.0) --
180 * -- LAPACK is a software package provided by Univ. of Tennessee, --
181 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
182 * November 2013
183 *
184 * .. Scalar Arguments ..
185  CHARACTER UPLO
186  INTEGER INFO, KB, LDA, LDW, N, NB
187 * ..
188 * .. Array Arguments ..
189  INTEGER IPIV( * )
190  DOUBLE PRECISION A( lda, * ), W( ldw, * )
191 * ..
192 *
193 * =====================================================================
194 *
195 * .. Parameters ..
196  DOUBLE PRECISION ZERO, ONE
197  parameter ( zero = 0.0d+0, one = 1.0d+0 )
198  DOUBLE PRECISION EIGHT, SEVTEN
199  parameter ( eight = 8.0d+0, sevten = 17.0d+0 )
200 * ..
201 * .. Local Scalars ..
202  INTEGER IMAX, J, JB, JJ, JMAX, JP, K, KK, KKW, KP,
203  $ kstep, kw
204  DOUBLE PRECISION ABSAKK, ALPHA, COLMAX, D11, D21, D22, R1,
205  $ rowmax, t
206 * ..
207 * .. External Functions ..
208  LOGICAL LSAME
209  INTEGER IDAMAX
210  EXTERNAL lsame, idamax
211 * ..
212 * .. External Subroutines ..
213  EXTERNAL dcopy, dgemm, dgemv, dscal, dswap
214 * ..
215 * .. Intrinsic Functions ..
216  INTRINSIC abs, max, min, sqrt
217 * ..
218 * .. Executable Statements ..
219 *
220  info = 0
221 *
222 * Initialize ALPHA for use in choosing pivot block size.
223 *
224  alpha = ( one+sqrt( sevten ) ) / eight
225 *
226  IF( lsame( uplo, 'U' ) ) THEN
227 *
228 * Factorize the trailing columns of A using the upper triangle
229 * of A and working backwards, and compute the matrix W = U12*D
230 * for use in updating A11
231 *
232 * K is the main loop index, decreasing from N in steps of 1 or 2
233 *
234 * KW is the column of W which corresponds to column K of A
235 *
236  k = n
237  10 CONTINUE
238  kw = nb + k - n
239 *
240 * Exit from loop
241 *
242  IF( ( k.LE.n-nb+1 .AND. nb.LT.n ) .OR. k.LT.1 )
243  $ GO TO 30
244 *
245 * Copy column K of A to column KW of W and update it
246 *
247  CALL dcopy( k, a( 1, k ), 1, w( 1, kw ), 1 )
248  IF( k.LT.n )
249  $ CALL dgemv( 'No transpose', k, n-k, -one, a( 1, k+1 ), lda,
250  $ w( k, kw+1 ), ldw, one, w( 1, kw ), 1 )
251 *
252  kstep = 1
253 *
254 * Determine rows and columns to be interchanged and whether
255 * a 1-by-1 or 2-by-2 pivot block will be used
256 *
257  absakk = abs( w( k, kw ) )
258 *
259 * IMAX is the row-index of the largest off-diagonal element in
260 * column K, and COLMAX is its absolute value.
261 * Determine both COLMAX and IMAX.
262 *
263  IF( k.GT.1 ) THEN
264  imax = idamax( k-1, w( 1, kw ), 1 )
265  colmax = abs( w( imax, kw ) )
266  ELSE
267  colmax = zero
268  END IF
269 *
270  IF( max( absakk, colmax ).EQ.zero ) THEN
271 *
272 * Column K is zero or underflow: set INFO and continue
273 *
274  IF( info.EQ.0 )
275  $ info = k
276  kp = k
277  ELSE
278  IF( absakk.GE.alpha*colmax ) THEN
279 *
280 * no interchange, use 1-by-1 pivot block
281 *
282  kp = k
283  ELSE
284 *
285 * Copy column IMAX to column KW-1 of W and update it
286 *
287  CALL dcopy( imax, a( 1, imax ), 1, w( 1, kw-1 ), 1 )
288  CALL dcopy( k-imax, a( imax, imax+1 ), lda,
289  $ w( imax+1, kw-1 ), 1 )
290  IF( k.LT.n )
291  $ CALL dgemv( 'No transpose', k, n-k, -one, a( 1, k+1 ),
292  $ lda, w( imax, kw+1 ), ldw, one,
293  $ w( 1, kw-1 ), 1 )
294 *
295 * JMAX is the column-index of the largest off-diagonal
296 * element in row IMAX, and ROWMAX is its absolute value
297 *
298  jmax = imax + idamax( k-imax, w( imax+1, kw-1 ), 1 )
299  rowmax = abs( w( jmax, kw-1 ) )
300  IF( imax.GT.1 ) THEN
301  jmax = idamax( imax-1, w( 1, kw-1 ), 1 )
302  rowmax = max( rowmax, abs( w( jmax, kw-1 ) ) )
303  END IF
304 *
305  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
306 *
307 * no interchange, use 1-by-1 pivot block
308 *
309  kp = k
310  ELSE IF( abs( w( imax, kw-1 ) ).GE.alpha*rowmax ) THEN
311 *
312 * interchange rows and columns K and IMAX, use 1-by-1
313 * pivot block
314 *
315  kp = imax
316 *
317 * copy column KW-1 of W to column KW of W
318 *
319  CALL dcopy( k, w( 1, kw-1 ), 1, w( 1, kw ), 1 )
320  ELSE
321 *
322 * interchange rows and columns K-1 and IMAX, use 2-by-2
323 * pivot block
324 *
325  kp = imax
326  kstep = 2
327  END IF
328  END IF
329 *
330 * ============================================================
331 *
332 * KK is the column of A where pivoting step stopped
333 *
334  kk = k - kstep + 1
335 *
336 * KKW is the column of W which corresponds to column KK of A
337 *
338  kkw = nb + kk - n
339 *
340 * Interchange rows and columns KP and KK.
341 * Updated column KP is already stored in column KKW of W.
342 *
343  IF( kp.NE.kk ) THEN
344 *
345 * Copy non-updated column KK to column KP of submatrix A
346 * at step K. No need to copy element into column K
347 * (or K and K-1 for 2-by-2 pivot) of A, since these columns
348 * will be later overwritten.
349 *
350  a( kp, kp ) = a( kk, kk )
351  CALL dcopy( kk-1-kp, a( kp+1, kk ), 1, a( kp, kp+1 ),
352  $ lda )
353  IF( kp.GT.1 )
354  $ CALL dcopy( kp-1, a( 1, kk ), 1, a( 1, kp ), 1 )
355 *
356 * Interchange rows KK and KP in last K+1 to N columns of A
357 * (columns K (or K and K-1 for 2-by-2 pivot) of A will be
358 * later overwritten). Interchange rows KK and KP
359 * in last KKW to NB columns of W.
360 *
361  IF( k.LT.n )
362  $ CALL dswap( n-k, a( kk, k+1 ), lda, a( kp, k+1 ),
363  $ lda )
364  CALL dswap( n-kk+1, w( kk, kkw ), ldw, w( kp, kkw ),
365  $ ldw )
366  END IF
367 *
368  IF( kstep.EQ.1 ) THEN
369 *
370 * 1-by-1 pivot block D(k): column kw of W now holds
371 *
372 * W(kw) = U(k)*D(k),
373 *
374 * where U(k) is the k-th column of U
375 *
376 * Store subdiag. elements of column U(k)
377 * and 1-by-1 block D(k) in column k of A.
378 * NOTE: Diagonal element U(k,k) is a UNIT element
379 * and not stored.
380 * A(k,k) := D(k,k) = W(k,kw)
381 * A(1:k-1,k) := U(1:k-1,k) = W(1:k-1,kw)/D(k,k)
382 *
383  CALL dcopy( k, w( 1, kw ), 1, a( 1, k ), 1 )
384  r1 = one / a( k, k )
385  CALL dscal( k-1, r1, a( 1, k ), 1 )
386 *
387  ELSE
388 *
389 * 2-by-2 pivot block D(k): columns kw and kw-1 of W now hold
390 *
391 * ( W(kw-1) W(kw) ) = ( U(k-1) U(k) )*D(k)
392 *
393 * where U(k) and U(k-1) are the k-th and (k-1)-th columns
394 * of U
395 *
396 * Store U(1:k-2,k-1) and U(1:k-2,k) and 2-by-2
397 * block D(k-1:k,k-1:k) in columns k-1 and k of A.
398 * NOTE: 2-by-2 diagonal block U(k-1:k,k-1:k) is a UNIT
399 * block and not stored.
400 * A(k-1:k,k-1:k) := D(k-1:k,k-1:k) = W(k-1:k,kw-1:kw)
401 * A(1:k-2,k-1:k) := U(1:k-2,k:k-1:k) =
402 * = W(1:k-2,kw-1:kw) * ( D(k-1:k,k-1:k)**(-1) )
403 *
404  IF( k.GT.2 ) THEN
405 *
406 * Compose the columns of the inverse of 2-by-2 pivot
407 * block D in the following way to reduce the number
408 * of FLOPS when we myltiply panel ( W(kw-1) W(kw) ) by
409 * this inverse
410 *
411 * D**(-1) = ( d11 d21 )**(-1) =
412 * ( d21 d22 )
413 *
414 * = 1/(d11*d22-d21**2) * ( ( d22 ) (-d21 ) ) =
415 * ( (-d21 ) ( d11 ) )
416 *
417 * = 1/d21 * 1/((d11/d21)*(d22/d21)-1) *
418 *
419 * * ( ( d22/d21 ) ( -1 ) ) =
420 * ( ( -1 ) ( d11/d21 ) )
421 *
422 * = 1/d21 * 1/(D22*D11-1) * ( ( D11 ) ( -1 ) ) =
423 * ( ( -1 ) ( D22 ) )
424 *
425 * = 1/d21 * T * ( ( D11 ) ( -1 ) )
426 * ( ( -1 ) ( D22 ) )
427 *
428 * = D21 * ( ( D11 ) ( -1 ) )
429 * ( ( -1 ) ( D22 ) )
430 *
431  d21 = w( k-1, kw )
432  d11 = w( k, kw ) / d21
433  d22 = w( k-1, kw-1 ) / d21
434  t = one / ( d11*d22-one )
435  d21 = t / d21
436 *
437 * Update elements in columns A(k-1) and A(k) as
438 * dot products of rows of ( W(kw-1) W(kw) ) and columns
439 * of D**(-1)
440 *
441  DO 20 j = 1, k - 2
442  a( j, k-1 ) = d21*( d11*w( j, kw-1 )-w( j, kw ) )
443  a( j, k ) = d21*( d22*w( j, kw )-w( j, kw-1 ) )
444  20 CONTINUE
445  END IF
446 *
447 * Copy D(k) to A
448 *
449  a( k-1, k-1 ) = w( k-1, kw-1 )
450  a( k-1, k ) = w( k-1, kw )
451  a( k, k ) = w( k, kw )
452 *
453  END IF
454 *
455  END IF
456 *
457 * Store details of the interchanges in IPIV
458 *
459  IF( kstep.EQ.1 ) THEN
460  ipiv( k ) = kp
461  ELSE
462  ipiv( k ) = -kp
463  ipiv( k-1 ) = -kp
464  END IF
465 *
466 * Decrease K and return to the start of the main loop
467 *
468  k = k - kstep
469  GO TO 10
470 *
471  30 CONTINUE
472 *
473 * Update the upper triangle of A11 (= A(1:k,1:k)) as
474 *
475 * A11 := A11 - U12*D*U12**T = A11 - U12*W**T
476 *
477 * computing blocks of NB columns at a time
478 *
479  DO 50 j = ( ( k-1 ) / nb )*nb + 1, 1, -nb
480  jb = min( nb, k-j+1 )
481 *
482 * Update the upper triangle of the diagonal block
483 *
484  DO 40 jj = j, j + jb - 1
485  CALL dgemv( 'No transpose', jj-j+1, n-k, -one,
486  $ a( j, k+1 ), lda, w( jj, kw+1 ), ldw, one,
487  $ a( j, jj ), 1 )
488  40 CONTINUE
489 *
490 * Update the rectangular superdiagonal block
491 *
492  CALL dgemm( 'No transpose', 'Transpose', j-1, jb, n-k, -one,
493  $ a( 1, k+1 ), lda, w( j, kw+1 ), ldw, one,
494  $ a( 1, j ), lda )
495  50 CONTINUE
496 *
497 * Put U12 in standard form by partially undoing the interchanges
498 * in columns k+1:n looping backwards from k+1 to n
499 *
500  j = k + 1
501  60 CONTINUE
502 *
503 * Undo the interchanges (if any) of rows JJ and JP at each
504 * step J
505 *
506 * (Here, J is a diagonal index)
507  jj = j
508  jp = ipiv( j )
509  IF( jp.LT.0 ) THEN
510  jp = -jp
511 * (Here, J is a diagonal index)
512  j = j + 1
513  END IF
514 * (NOTE: Here, J is used to determine row length. Length N-J+1
515 * of the rows to swap back doesn't include diagonal element)
516  j = j + 1
517  IF( jp.NE.jj .AND. j.LE.n )
518  $ CALL dswap( n-j+1, a( jp, j ), lda, a( jj, j ), lda )
519  IF( j.LT.n )
520  $ GO TO 60
521 *
522 * Set KB to the number of columns factorized
523 *
524  kb = n - k
525 *
526  ELSE
527 *
528 * Factorize the leading columns of A using the lower triangle
529 * of A and working forwards, and compute the matrix W = L21*D
530 * for use in updating A22
531 *
532 * K is the main loop index, increasing from 1 in steps of 1 or 2
533 *
534  k = 1
535  70 CONTINUE
536 *
537 * Exit from loop
538 *
539  IF( ( k.GE.nb .AND. nb.LT.n ) .OR. k.GT.n )
540  $ GO TO 90
541 *
542 * Copy column K of A to column K of W and update it
543 *
544  CALL dcopy( n-k+1, a( k, k ), 1, w( k, k ), 1 )
545  CALL dgemv( 'No transpose', n-k+1, k-1, -one, a( k, 1 ), lda,
546  $ w( k, 1 ), ldw, one, w( k, k ), 1 )
547 *
548  kstep = 1
549 *
550 * Determine rows and columns to be interchanged and whether
551 * a 1-by-1 or 2-by-2 pivot block will be used
552 *
553  absakk = abs( w( k, k ) )
554 *
555 * IMAX is the row-index of the largest off-diagonal element in
556 * column K, and COLMAX is its absolute value.
557 * Determine both COLMAX and IMAX.
558 *
559  IF( k.LT.n ) THEN
560  imax = k + idamax( n-k, w( k+1, k ), 1 )
561  colmax = abs( w( imax, k ) )
562  ELSE
563  colmax = zero
564  END IF
565 *
566  IF( max( absakk, colmax ).EQ.zero ) THEN
567 *
568 * Column K is zero or underflow: set INFO and continue
569 *
570  IF( info.EQ.0 )
571  $ info = k
572  kp = k
573  ELSE
574  IF( absakk.GE.alpha*colmax ) THEN
575 *
576 * no interchange, use 1-by-1 pivot block
577 *
578  kp = k
579  ELSE
580 *
581 * Copy column IMAX to column K+1 of W and update it
582 *
583  CALL dcopy( imax-k, a( imax, k ), lda, w( k, k+1 ), 1 )
584  CALL dcopy( n-imax+1, a( imax, imax ), 1, w( imax, k+1 ),
585  $ 1 )
586  CALL dgemv( 'No transpose', n-k+1, k-1, -one, a( k, 1 ),
587  $ lda, w( imax, 1 ), ldw, one, w( k, k+1 ), 1 )
588 *
589 * JMAX is the column-index of the largest off-diagonal
590 * element in row IMAX, and ROWMAX is its absolute value
591 *
592  jmax = k - 1 + idamax( imax-k, w( k, k+1 ), 1 )
593  rowmax = abs( w( jmax, k+1 ) )
594  IF( imax.LT.n ) THEN
595  jmax = imax + idamax( n-imax, w( imax+1, k+1 ), 1 )
596  rowmax = max( rowmax, abs( w( jmax, k+1 ) ) )
597  END IF
598 *
599  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
600 *
601 * no interchange, use 1-by-1 pivot block
602 *
603  kp = k
604  ELSE IF( abs( w( imax, k+1 ) ).GE.alpha*rowmax ) THEN
605 *
606 * interchange rows and columns K and IMAX, use 1-by-1
607 * pivot block
608 *
609  kp = imax
610 *
611 * copy column K+1 of W to column K of W
612 *
613  CALL dcopy( n-k+1, w( k, k+1 ), 1, w( k, k ), 1 )
614  ELSE
615 *
616 * interchange rows and columns K+1 and IMAX, use 2-by-2
617 * pivot block
618 *
619  kp = imax
620  kstep = 2
621  END IF
622  END IF
623 *
624 * ============================================================
625 *
626 * KK is the column of A where pivoting step stopped
627 *
628  kk = k + kstep - 1
629 *
630 * Interchange rows and columns KP and KK.
631 * Updated column KP is already stored in column KK of W.
632 *
633  IF( kp.NE.kk ) THEN
634 *
635 * Copy non-updated column KK to column KP of submatrix A
636 * at step K. No need to copy element into column K
637 * (or K and K+1 for 2-by-2 pivot) of A, since these columns
638 * will be later overwritten.
639 *
640  a( kp, kp ) = a( kk, kk )
641  CALL dcopy( kp-kk-1, a( kk+1, kk ), 1, a( kp, kk+1 ),
642  $ lda )
643  IF( kp.LT.n )
644  $ CALL dcopy( n-kp, a( kp+1, kk ), 1, a( kp+1, kp ), 1 )
645 *
646 * Interchange rows KK and KP in first K-1 columns of A
647 * (columns K (or K and K+1 for 2-by-2 pivot) of A will be
648 * later overwritten). Interchange rows KK and KP
649 * in first KK columns of W.
650 *
651  IF( k.GT.1 )
652  $ CALL dswap( k-1, a( kk, 1 ), lda, a( kp, 1 ), lda )
653  CALL dswap( kk, w( kk, 1 ), ldw, w( kp, 1 ), ldw )
654  END IF
655 *
656  IF( kstep.EQ.1 ) THEN
657 *
658 * 1-by-1 pivot block D(k): column k of W now holds
659 *
660 * W(k) = L(k)*D(k),
661 *
662 * where L(k) is the k-th column of L
663 *
664 * Store subdiag. elements of column L(k)
665 * and 1-by-1 block D(k) in column k of A.
666 * (NOTE: Diagonal element L(k,k) is a UNIT element
667 * and not stored)
668 * A(k,k) := D(k,k) = W(k,k)
669 * A(k+1:N,k) := L(k+1:N,k) = W(k+1:N,k)/D(k,k)
670 *
671  CALL dcopy( n-k+1, w( k, k ), 1, a( k, k ), 1 )
672  IF( k.LT.n ) THEN
673  r1 = one / a( k, k )
674  CALL dscal( n-k, r1, a( k+1, k ), 1 )
675  END IF
676 *
677  ELSE
678 *
679 * 2-by-2 pivot block D(k): columns k and k+1 of W now hold
680 *
681 * ( W(k) W(k+1) ) = ( L(k) L(k+1) )*D(k)
682 *
683 * where L(k) and L(k+1) are the k-th and (k+1)-th columns
684 * of L
685 *
686 * Store L(k+2:N,k) and L(k+2:N,k+1) and 2-by-2
687 * block D(k:k+1,k:k+1) in columns k and k+1 of A.
688 * (NOTE: 2-by-2 diagonal block L(k:k+1,k:k+1) is a UNIT
689 * block and not stored)
690 * A(k:k+1,k:k+1) := D(k:k+1,k:k+1) = W(k:k+1,k:k+1)
691 * A(k+2:N,k:k+1) := L(k+2:N,k:k+1) =
692 * = W(k+2:N,k:k+1) * ( D(k:k+1,k:k+1)**(-1) )
693 *
694  IF( k.LT.n-1 ) THEN
695 *
696 * Compose the columns of the inverse of 2-by-2 pivot
697 * block D in the following way to reduce the number
698 * of FLOPS when we myltiply panel ( W(k) W(k+1) ) by
699 * this inverse
700 *
701 * D**(-1) = ( d11 d21 )**(-1) =
702 * ( d21 d22 )
703 *
704 * = 1/(d11*d22-d21**2) * ( ( d22 ) (-d21 ) ) =
705 * ( (-d21 ) ( d11 ) )
706 *
707 * = 1/d21 * 1/((d11/d21)*(d22/d21)-1) *
708 *
709 * * ( ( d22/d21 ) ( -1 ) ) =
710 * ( ( -1 ) ( d11/d21 ) )
711 *
712 * = 1/d21 * 1/(D22*D11-1) * ( ( D11 ) ( -1 ) ) =
713 * ( ( -1 ) ( D22 ) )
714 *
715 * = 1/d21 * T * ( ( D11 ) ( -1 ) )
716 * ( ( -1 ) ( D22 ) )
717 *
718 * = D21 * ( ( D11 ) ( -1 ) )
719 * ( ( -1 ) ( D22 ) )
720 *
721  d21 = w( k+1, k )
722  d11 = w( k+1, k+1 ) / d21
723  d22 = w( k, k ) / d21
724  t = one / ( d11*d22-one )
725  d21 = t / d21
726 *
727 * Update elements in columns A(k) and A(k+1) as
728 * dot products of rows of ( W(k) W(k+1) ) and columns
729 * of D**(-1)
730 *
731  DO 80 j = k + 2, n
732  a( j, k ) = d21*( d11*w( j, k )-w( j, k+1 ) )
733  a( j, k+1 ) = d21*( d22*w( j, k+1 )-w( j, k ) )
734  80 CONTINUE
735  END IF
736 *
737 * Copy D(k) to A
738 *
739  a( k, k ) = w( k, k )
740  a( k+1, k ) = w( k+1, k )
741  a( k+1, k+1 ) = w( k+1, k+1 )
742 *
743  END IF
744 *
745  END IF
746 *
747 * Store details of the interchanges in IPIV
748 *
749  IF( kstep.EQ.1 ) THEN
750  ipiv( k ) = kp
751  ELSE
752  ipiv( k ) = -kp
753  ipiv( k+1 ) = -kp
754  END IF
755 *
756 * Increase K and return to the start of the main loop
757 *
758  k = k + kstep
759  GO TO 70
760 *
761  90 CONTINUE
762 *
763 * Update the lower triangle of A22 (= A(k:n,k:n)) as
764 *
765 * A22 := A22 - L21*D*L21**T = A22 - L21*W**T
766 *
767 * computing blocks of NB columns at a time
768 *
769  DO 110 j = k, n, nb
770  jb = min( nb, n-j+1 )
771 *
772 * Update the lower triangle of the diagonal block
773 *
774  DO 100 jj = j, j + jb - 1
775  CALL dgemv( 'No transpose', j+jb-jj, k-1, -one,
776  $ a( jj, 1 ), lda, w( jj, 1 ), ldw, one,
777  $ a( jj, jj ), 1 )
778  100 CONTINUE
779 *
780 * Update the rectangular subdiagonal block
781 *
782  IF( j+jb.LE.n )
783  $ CALL dgemm( 'No transpose', 'Transpose', n-j-jb+1, jb,
784  $ k-1, -one, a( j+jb, 1 ), lda, w( j, 1 ), ldw,
785  $ one, a( j+jb, j ), lda )
786  110 CONTINUE
787 *
788 * Put L21 in standard form by partially undoing the interchanges
789 * of rows in columns 1:k-1 looping backwards from k-1 to 1
790 *
791  j = k - 1
792  120 CONTINUE
793 *
794 * Undo the interchanges (if any) of rows JJ and JP at each
795 * step J
796 *
797 * (Here, J is a diagonal index)
798  jj = j
799  jp = ipiv( j )
800  IF( jp.LT.0 ) THEN
801  jp = -jp
802 * (Here, J is a diagonal index)
803  j = j - 1
804  END IF
805 * (NOTE: Here, J is used to determine row length. Length J
806 * of the rows to swap back doesn't include diagonal element)
807  j = j - 1
808  IF( jp.NE.jj .AND. j.GE.1 )
809  $ CALL dswap( j, a( jp, 1 ), lda, a( jj, 1 ), lda )
810  IF( j.GT.1 )
811  $ GO TO 120
812 *
813 * Set KB to the number of columns factorized
814 *
815  kb = k - 1
816 *
817  END IF
818  RETURN
819 *
820 * End of DLASYF
821 *
822  END
subroutine dcopy(N, DX, INCX, DY, INCY)
DCOPY
Definition: dcopy.f:53
subroutine dlasyf(UPLO, N, NB, KB, A, LDA, IPIV, W, LDW, INFO)
DLASYF computes a partial factorization of a real symmetric matrix using the Bunch-Kaufman diagonal p...
Definition: dlasyf.f:178
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158
subroutine dgemm(TRANSA, TRANSB, M, N, K, ALPHA, A, LDA, B, LDB, BETA, C, LDC)
DGEMM
Definition: dgemm.f:189
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55