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