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