LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
zhetf2.f
Go to the documentation of this file.
1*> \brief \b ZHETF2 computes the factorization of a complex Hermitian matrix, using the diagonal pivoting method (unblocked algorithm, calling Level 2 BLAS).
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> \htmlonly
9*> Download ZHETF2 + dependencies
10*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zhetf2.f">
11*> [TGZ]</a>
12*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zhetf2.f">
13*> [ZIP]</a>
14*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zhetf2.f">
15*> [TXT]</a>
16*> \endhtmlonly
17*
18* Definition:
19* ===========
20*
21* SUBROUTINE ZHETF2( UPLO, N, A, LDA, IPIV, INFO )
22*
23* .. Scalar Arguments ..
24* CHARACTER UPLO
25* INTEGER INFO, LDA, N
26* ..
27* .. Array Arguments ..
28* INTEGER IPIV( * )
29* COMPLEX*16 A( LDA, * )
30* ..
31*
32*
33*> \par Purpose:
34* =============
35*>
36*> \verbatim
37*>
38*> ZHETF2 computes the factorization of a complex Hermitian matrix A
39*> using the Bunch-Kaufman diagonal pivoting method:
40*>
41*> A = U*D*U**H or A = L*D*L**H
42*>
43*> where U (or L) is a product of permutation and unit upper (lower)
44*> triangular matrices, U**H is the conjugate transpose of U, and D is
45*> Hermitian and block diagonal with 1-by-1 and 2-by-2 diagonal blocks.
46*>
47*> This is the unblocked version of the algorithm, calling Level 2 BLAS.
48*> \endverbatim
49*
50* Arguments:
51* ==========
52*
53*> \param[in] UPLO
54*> \verbatim
55*> UPLO is CHARACTER*1
56*> Specifies whether the upper or lower triangular part of the
57*> Hermitian matrix A is stored:
58*> = 'U': Upper triangular
59*> = 'L': Lower triangular
60*> \endverbatim
61*>
62*> \param[in] N
63*> \verbatim
64*> N is INTEGER
65*> The order of the matrix A. N >= 0.
66*> \endverbatim
67*>
68*> \param[in,out] A
69*> \verbatim
70*> A is COMPLEX*16 array, dimension (LDA,N)
71*> On entry, the Hermitian matrix A. If UPLO = 'U', the leading
72*> n-by-n upper triangular part of A contains the upper
73*> triangular part of the matrix A, and the strictly lower
74*> triangular part of A is not referenced. If UPLO = 'L', the
75*> leading n-by-n lower triangular part of A contains the lower
76*> triangular part of the matrix A, and the strictly upper
77*> triangular part of A is not referenced.
78*>
79*> On exit, the block diagonal matrix D and the multipliers used
80*> to obtain the factor U or L (see below for further details).
81*> \endverbatim
82*>
83*> \param[in] LDA
84*> \verbatim
85*> LDA is INTEGER
86*> The leading dimension of the array A. LDA >= max(1,N).
87*> \endverbatim
88*>
89*> \param[out] IPIV
90*> \verbatim
91*> IPIV is INTEGER array, dimension (N)
92*> Details of the interchanges and the block structure of D.
93*>
94*> If UPLO = 'U':
95*> If IPIV(k) > 0, then rows and columns k and IPIV(k) were
96*> interchanged and D(k,k) is a 1-by-1 diagonal block.
97*>
98*> If IPIV(k) = IPIV(k-1) < 0, then rows and columns
99*> k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k)
100*> is a 2-by-2 diagonal block.
101*>
102*> If UPLO = 'L':
103*> If IPIV(k) > 0, then rows and columns k and IPIV(k) were
104*> interchanged and D(k,k) is a 1-by-1 diagonal block.
105*>
106*> If IPIV(k) = IPIV(k+1) < 0, then rows and columns
107*> k+1 and -IPIV(k) were interchanged and D(k:k+1,k:k+1)
108*> is a 2-by-2 diagonal block.
109*> \endverbatim
110*>
111*> \param[out] INFO
112*> \verbatim
113*> INFO is INTEGER
114*> = 0: successful exit
115*> < 0: if INFO = -k, the k-th argument had an illegal value
116*> > 0: if INFO = k, D(k,k) is exactly zero. The factorization
117*> has been completed, but the block diagonal matrix D is
118*> exactly singular, and division by zero will occur if it
119*> is used to solve a system of equations.
120*> \endverbatim
121*
122* Authors:
123* ========
124*
125*> \author Univ. of Tennessee
126*> \author Univ. of California Berkeley
127*> \author Univ. of Colorado Denver
128*> \author NAG Ltd.
129*
130*> \ingroup hetf2
131*
132*> \par Further Details:
133* =====================
134*>
135*> \verbatim
136*>
137*> If UPLO = 'U', then A = U*D*U**H, where
138*> U = P(n)*U(n)* ... *P(k)U(k)* ...,
139*> i.e., U is a product of terms P(k)*U(k), where k decreases from n to
140*> 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1
141*> and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as
142*> defined by IPIV(k), and U(k) is a unit upper triangular matrix, such
143*> that if the diagonal block D(k) is of order s (s = 1 or 2), then
144*>
145*> ( I v 0 ) k-s
146*> U(k) = ( 0 I 0 ) s
147*> ( 0 0 I ) n-k
148*> k-s s n-k
149*>
150*> If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k).
151*> If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k),
152*> and A(k,k), and v overwrites A(1:k-2,k-1:k).
153*>
154*> If UPLO = 'L', then A = L*D*L**H, where
155*> L = P(1)*L(1)* ... *P(k)*L(k)* ...,
156*> i.e., L is a product of terms P(k)*L(k), where k increases from 1 to
157*> n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1
158*> and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as
159*> defined by IPIV(k), and L(k) is a unit lower triangular matrix, such
160*> that if the diagonal block D(k) is of order s (s = 1 or 2), then
161*>
162*> ( I 0 0 ) k-1
163*> L(k) = ( 0 I 0 ) s
164*> ( 0 v I ) n-k-s+1
165*> k-1 s n-k-s+1
166*>
167*> If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k).
168*> If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k),
169*> and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
170*> \endverbatim
171*
172*> \par Contributors:
173* ==================
174*>
175*> \verbatim
176*> 09-29-06 - patch from
177*> Bobby Cheng, MathWorks
178*>
179*> Replace l.210 and l.393
180*> IF( MAX( ABSAKK, COLMAX ).EQ.ZERO ) THEN
181*> by
182*> IF( (MAX( ABSAKK, COLMAX ).EQ.ZERO) .OR. DISNAN(ABSAKK) ) THEN
183*>
184*> 01-01-96 - Based on modifications by
185*> J. Lewis, Boeing Computer Services Company
186*> A. Petitet, Computer Science Dept., Univ. of Tenn., Knoxville, USA
187*> \endverbatim
188*
189* =====================================================================
190 SUBROUTINE zhetf2( UPLO, N, A, LDA, IPIV, INFO )
191*
192* -- LAPACK computational routine --
193* -- LAPACK is a software package provided by Univ. of Tennessee, --
194* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
195*
196* .. Scalar Arguments ..
197 CHARACTER UPLO
198 INTEGER INFO, LDA, N
199* ..
200* .. Array Arguments ..
201 INTEGER IPIV( * )
202 COMPLEX*16 A( LDA, * )
203* ..
204*
205* =====================================================================
206*
207* .. Parameters ..
208 DOUBLE PRECISION ZERO, ONE
209 parameter( zero = 0.0d+0, one = 1.0d+0 )
210 DOUBLE PRECISION EIGHT, SEVTEN
211 parameter( eight = 8.0d+0, sevten = 17.0d+0 )
212* ..
213* .. Local Scalars ..
214 LOGICAL UPPER
215 INTEGER I, IMAX, J, JMAX, K, KK, KP, KSTEP
216 DOUBLE PRECISION ABSAKK, ALPHA, COLMAX, D, D11, D22, R1, ROWMAX,
217 $ TT
218 COMPLEX*16 D12, D21, T, WK, WKM1, WKP1, ZDUM
219* ..
220* .. External Functions ..
221 LOGICAL LSAME, DISNAN
222 INTEGER IZAMAX
223 DOUBLE PRECISION DLAPY2
224 EXTERNAL lsame, izamax, dlapy2, disnan
225* ..
226* .. External Subroutines ..
227 EXTERNAL xerbla, zdscal, zher, zswap
228* ..
229* .. Intrinsic Functions ..
230 INTRINSIC abs, dble, dcmplx, dconjg, dimag, max, sqrt
231* ..
232* .. Statement Functions ..
233 DOUBLE PRECISION CABS1
234* ..
235* .. Statement Function definitions ..
236 cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
237* ..
238* .. Executable Statements ..
239*
240* Test the input parameters.
241*
242 info = 0
243 upper = lsame( uplo, 'U' )
244 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
245 info = -1
246 ELSE IF( n.LT.0 ) THEN
247 info = -2
248 ELSE IF( lda.LT.max( 1, n ) ) THEN
249 info = -4
250 END IF
251 IF( info.NE.0 ) THEN
252 CALL xerbla( 'ZHETF2', -info )
253 RETURN
254 END IF
255*
256* Initialize ALPHA for use in choosing pivot block size.
257*
258 alpha = ( one+sqrt( sevten ) ) / eight
259*
260 IF( upper ) THEN
261*
262* Factorize A as U*D*U**H using the upper triangle of A
263*
264* K is the main loop index, decreasing from N to 1 in steps of
265* 1 or 2
266*
267 k = n
268 10 CONTINUE
269*
270* If K < 1, exit from loop
271*
272 IF( k.LT.1 )
273 $ GO TO 90
274 kstep = 1
275*
276* Determine rows and columns to be interchanged and whether
277* a 1-by-1 or 2-by-2 pivot block will be used
278*
279 absakk = abs( dble( a( k, k ) ) )
280*
281* IMAX is the row-index of the largest off-diagonal element in
282* column K, and COLMAX is its absolute value.
283* Determine both COLMAX and IMAX.
284*
285 IF( k.GT.1 ) THEN
286 imax = izamax( k-1, a( 1, k ), 1 )
287 colmax = cabs1( a( imax, k ) )
288 ELSE
289 colmax = zero
290 END IF
291*
292 IF( (max( absakk, colmax ).EQ.zero) .OR. disnan(absakk) ) THEN
293*
294* Column K is zero or underflow, or contains a NaN:
295* set INFO and continue
296*
297 IF( info.EQ.0 )
298 $ info = k
299 kp = k
300 a( k, k ) = dble( a( k, k ) )
301 ELSE
302*
303* ============================================================
304*
305* Test for interchange
306*
307 IF( absakk.GE.alpha*colmax ) THEN
308*
309* no interchange, use 1-by-1 pivot block
310*
311 kp = k
312 ELSE
313*
314* JMAX is the column-index of the largest off-diagonal
315* element in row IMAX, and ROWMAX is its absolute value.
316* Determine only ROWMAX.
317*
318 jmax = imax + izamax( k-imax, a( imax, imax+1 ), lda )
319 rowmax = cabs1( a( imax, jmax ) )
320 IF( imax.GT.1 ) THEN
321 jmax = izamax( imax-1, a( 1, imax ), 1 )
322 rowmax = max( rowmax, cabs1( a( jmax, imax ) ) )
323 END IF
324*
325 IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
326*
327* no interchange, use 1-by-1 pivot block
328*
329 kp = k
330*
331 ELSE IF( abs( dble( a( imax, imax ) ) ).GE.alpha*rowmax )
332 $ THEN
333*
334* interchange rows and columns K and IMAX, use 1-by-1
335* pivot block
336*
337 kp = imax
338 ELSE
339*
340* interchange rows and columns K-1 and IMAX, use 2-by-2
341* pivot block
342*
343 kp = imax
344 kstep = 2
345 END IF
346*
347 END IF
348*
349* ============================================================
350*
351 kk = k - kstep + 1
352 IF( kp.NE.kk ) THEN
353*
354* Interchange rows and columns KK and KP in the leading
355* submatrix A(1:k,1:k)
356*
357 CALL zswap( kp-1, a( 1, kk ), 1, a( 1, kp ), 1 )
358 DO 20 j = kp + 1, kk - 1
359 t = dconjg( a( j, kk ) )
360 a( j, kk ) = dconjg( a( kp, j ) )
361 a( kp, j ) = t
362 20 CONTINUE
363 a( kp, kk ) = dconjg( a( kp, kk ) )
364 r1 = dble( a( kk, kk ) )
365 a( kk, kk ) = dble( a( kp, kp ) )
366 a( kp, kp ) = r1
367 IF( kstep.EQ.2 ) THEN
368 a( k, k ) = dble( a( k, k ) )
369 t = a( k-1, k )
370 a( k-1, k ) = a( kp, k )
371 a( kp, k ) = t
372 END IF
373 ELSE
374 a( k, k ) = dble( a( k, k ) )
375 IF( kstep.EQ.2 )
376 $ a( k-1, k-1 ) = dble( a( k-1, k-1 ) )
377 END IF
378*
379* Update the leading submatrix
380*
381 IF( kstep.EQ.1 ) THEN
382*
383* 1-by-1 pivot block D(k): column k now holds
384*
385* W(k) = U(k)*D(k)
386*
387* where U(k) is the k-th column of U
388*
389* Perform a rank-1 update of A(1:k-1,1:k-1) as
390*
391* A := A - U(k)*D(k)*U(k)**H = A - W(k)*1/D(k)*W(k)**H
392*
393 r1 = one / dble( a( k, k ) )
394 CALL zher( uplo, k-1, -r1, a( 1, k ), 1, a, lda )
395*
396* Store U(k) in column k
397*
398 CALL zdscal( k-1, r1, a( 1, k ), 1 )
399 ELSE
400*
401* 2-by-2 pivot block D(k): columns k and k-1 now hold
402*
403* ( W(k-1) W(k) ) = ( U(k-1) U(k) )*D(k)
404*
405* where U(k) and U(k-1) are the k-th and (k-1)-th columns
406* of U
407*
408* Perform a rank-2 update of A(1:k-2,1:k-2) as
409*
410* A := A - ( U(k-1) U(k) )*D(k)*( U(k-1) U(k) )**H
411* = A - ( W(k-1) W(k) )*inv(D(k))*( W(k-1) W(k) )**H
412*
413 IF( k.GT.2 ) THEN
414*
415 d = dlapy2( dble( a( k-1, k ) ),
416 $ dimag( a( k-1, k ) ) )
417 d22 = dble( a( k-1, k-1 ) ) / d
418 d11 = dble( a( k, k ) ) / d
419 tt = one / ( d11*d22-one )
420 d12 = a( k-1, k ) / d
421 d = tt / d
422*
423 DO 40 j = k - 2, 1, -1
424 wkm1 = d*( d11*a( j, k-1 )-dconjg( d12 )*
425 $ a( j, k ) )
426 wk = d*( d22*a( j, k )-d12*a( j, k-1 ) )
427 DO 30 i = j, 1, -1
428 a( i, j ) = a( i, j ) - a( i, k )*dconjg( wk ) -
429 $ a( i, k-1 )*dconjg( wkm1 )
430 30 CONTINUE
431 a( j, k ) = wk
432 a( j, k-1 ) = wkm1
433 a( j, j ) = dcmplx( dble( a( j, j ) ), 0.0d+0 )
434 40 CONTINUE
435*
436 END IF
437*
438 END IF
439 END IF
440*
441* Store details of the interchanges in IPIV
442*
443 IF( kstep.EQ.1 ) THEN
444 ipiv( k ) = kp
445 ELSE
446 ipiv( k ) = -kp
447 ipiv( k-1 ) = -kp
448 END IF
449*
450* Decrease K and return to the start of the main loop
451*
452 k = k - kstep
453 GO TO 10
454*
455 ELSE
456*
457* Factorize A as L*D*L**H using the lower triangle of A
458*
459* K is the main loop index, increasing from 1 to N in steps of
460* 1 or 2
461*
462 k = 1
463 50 CONTINUE
464*
465* If K > N, exit from loop
466*
467 IF( k.GT.n )
468 $ GO TO 90
469 kstep = 1
470*
471* Determine rows and columns to be interchanged and whether
472* a 1-by-1 or 2-by-2 pivot block will be used
473*
474 absakk = abs( dble( a( k, k ) ) )
475*
476* IMAX is the row-index of the largest off-diagonal element in
477* column K, and COLMAX is its absolute value.
478* Determine both COLMAX and IMAX.
479*
480 IF( k.LT.n ) THEN
481 imax = k + izamax( n-k, a( k+1, k ), 1 )
482 colmax = cabs1( a( imax, k ) )
483 ELSE
484 colmax = zero
485 END IF
486*
487 IF( (max( absakk, colmax ).EQ.zero) .OR. disnan(absakk) ) THEN
488*
489* Column K is zero or underflow, or contains a NaN:
490* set INFO and continue
491*
492 IF( info.EQ.0 )
493 $ info = k
494 kp = k
495 a( k, k ) = dble( a( k, k ) )
496 ELSE
497*
498* ============================================================
499*
500* Test for interchange
501*
502 IF( absakk.GE.alpha*colmax ) THEN
503*
504* no interchange, use 1-by-1 pivot block
505*
506 kp = k
507 ELSE
508*
509* JMAX is the column-index of the largest off-diagonal
510* element in row IMAX, and ROWMAX is its absolute value.
511* Determine only ROWMAX.
512*
513 jmax = k - 1 + izamax( imax-k, a( imax, k ), lda )
514 rowmax = cabs1( a( imax, jmax ) )
515 IF( imax.LT.n ) THEN
516 jmax = imax + izamax( n-imax, a( imax+1, imax ), 1 )
517 rowmax = max( rowmax, cabs1( a( jmax, imax ) ) )
518 END IF
519*
520 IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
521*
522* no interchange, use 1-by-1 pivot block
523*
524 kp = k
525*
526 ELSE IF( abs( dble( a( imax, imax ) ) ).GE.alpha*rowmax )
527 $ THEN
528*
529* interchange rows and columns K and IMAX, use 1-by-1
530* pivot block
531*
532 kp = imax
533 ELSE
534*
535* interchange rows and columns K+1 and IMAX, use 2-by-2
536* pivot block
537*
538 kp = imax
539 kstep = 2
540 END IF
541*
542 END IF
543*
544* ============================================================
545*
546 kk = k + kstep - 1
547 IF( kp.NE.kk ) THEN
548*
549* Interchange rows and columns KK and KP in the trailing
550* submatrix A(k:n,k:n)
551*
552 IF( kp.LT.n )
553 $ CALL zswap( n-kp, a( kp+1, kk ), 1, a( kp+1, kp ), 1 )
554 DO 60 j = kk + 1, kp - 1
555 t = dconjg( a( j, kk ) )
556 a( j, kk ) = dconjg( a( kp, j ) )
557 a( kp, j ) = t
558 60 CONTINUE
559 a( kp, kk ) = dconjg( a( kp, kk ) )
560 r1 = dble( a( kk, kk ) )
561 a( kk, kk ) = dble( a( kp, kp ) )
562 a( kp, kp ) = r1
563 IF( kstep.EQ.2 ) THEN
564 a( k, k ) = dble( a( k, k ) )
565 t = a( k+1, k )
566 a( k+1, k ) = a( kp, k )
567 a( kp, k ) = t
568 END IF
569 ELSE
570 a( k, k ) = dble( a( k, k ) )
571 IF( kstep.EQ.2 )
572 $ a( k+1, k+1 ) = dble( a( k+1, k+1 ) )
573 END IF
574*
575* Update the trailing submatrix
576*
577 IF( kstep.EQ.1 ) THEN
578*
579* 1-by-1 pivot block D(k): column k now holds
580*
581* W(k) = L(k)*D(k)
582*
583* where L(k) is the k-th column of L
584*
585 IF( k.LT.n ) THEN
586*
587* Perform a rank-1 update of A(k+1:n,k+1:n) as
588*
589* A := A - L(k)*D(k)*L(k)**H = A - W(k)*(1/D(k))*W(k)**H
590*
591 r1 = one / dble( a( k, k ) )
592 CALL zher( uplo, n-k, -r1, a( k+1, k ), 1,
593 $ a( k+1, k+1 ), lda )
594*
595* Store L(k) in column K
596*
597 CALL zdscal( n-k, r1, a( k+1, k ), 1 )
598 END IF
599 ELSE
600*
601* 2-by-2 pivot block D(k)
602*
603 IF( k.LT.n-1 ) THEN
604*
605* Perform a rank-2 update of A(k+2:n,k+2:n) as
606*
607* A := A - ( L(k) L(k+1) )*D(k)*( L(k) L(k+1) )**H
608* = A - ( W(k) W(k+1) )*inv(D(k))*( W(k) W(k+1) )**H
609*
610* where L(k) and L(k+1) are the k-th and (k+1)-th
611* columns of L
612*
613 d = dlapy2( dble( a( k+1, k ) ),
614 $ dimag( a( k+1, k ) ) )
615 d11 = dble( a( k+1, k+1 ) ) / d
616 d22 = dble( a( k, k ) ) / d
617 tt = one / ( d11*d22-one )
618 d21 = a( k+1, k ) / d
619 d = tt / d
620*
621 DO 80 j = k + 2, n
622 wk = d*( d11*a( j, k )-d21*a( j, k+1 ) )
623 wkp1 = d*( d22*a( j, k+1 )-dconjg( d21 )*
624 $ a( j, k ) )
625 DO 70 i = j, n
626 a( i, j ) = a( i, j ) - a( i, k )*dconjg( wk ) -
627 $ a( i, k+1 )*dconjg( wkp1 )
628 70 CONTINUE
629 a( j, k ) = wk
630 a( j, k+1 ) = wkp1
631 a( j, j ) = dcmplx( dble( a( j, j ) ), 0.0d+0 )
632 80 CONTINUE
633 END IF
634 END IF
635 END IF
636*
637* Store details of the interchanges in IPIV
638*
639 IF( kstep.EQ.1 ) THEN
640 ipiv( k ) = kp
641 ELSE
642 ipiv( k ) = -kp
643 ipiv( k+1 ) = -kp
644 END IF
645*
646* Increase K and return to the start of the main loop
647*
648 k = k + kstep
649 GO TO 50
650*
651 END IF
652*
653 90 CONTINUE
654 RETURN
655*
656* End of ZHETF2
657*
658 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine zher(uplo, n, alpha, x, incx, a, lda)
ZHER
Definition zher.f:135
subroutine zhetf2(uplo, n, a, lda, ipiv, info)
ZHETF2 computes the factorization of a complex Hermitian matrix, using the diagonal pivoting method (...
Definition zhetf2.f:191
subroutine zdscal(n, da, zx, incx)
ZDSCAL
Definition zdscal.f:78
subroutine zswap(n, zx, incx, zy, incy)
ZSWAP
Definition zswap.f:81