LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
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).
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 *> 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 *> If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0, then rows and
96 *> columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k)
97 *> is a 2-by-2 diagonal block. If UPLO = 'L' and IPIV(k) =
98 *> IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were
99 *> interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
100 *> \endverbatim
101 *>
102 *> \param[out] INFO
103 *> \verbatim
104 *> INFO is INTEGER
105 *> = 0: successful exit
106 *> < 0: if INFO = -k, the k-th argument had an illegal value
107 *> > 0: if INFO = k, D(k,k) is exactly zero. The factorization
108 *> has been completed, but the block diagonal matrix D is
109 *> exactly singular, and division by zero will occur if it
110 *> is used to solve a system of equations.
111 *> \endverbatim
112 *
113 * Authors:
114 * ========
115 *
116 *> \author Univ. of Tennessee
117 *> \author Univ. of California Berkeley
118 *> \author Univ. of Colorado Denver
119 *> \author NAG Ltd.
120 *
121 *> \date September 2012
122 *
123 *> \ingroup complex16HEcomputational
124 *
125 *> \par Further Details:
126 * =====================
127 *>
128 *> \verbatim
129 *>
130 *> If UPLO = 'U', then A = U*D*U**H, where
131 *> U = P(n)*U(n)* ... *P(k)U(k)* ...,
132 *> i.e., U is a product of terms P(k)*U(k), where k decreases from n to
133 *> 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1
134 *> and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as
135 *> defined by IPIV(k), and U(k) is a unit upper triangular matrix, such
136 *> that if the diagonal block D(k) is of order s (s = 1 or 2), then
137 *>
138 *> ( I v 0 ) k-s
139 *> U(k) = ( 0 I 0 ) s
140 *> ( 0 0 I ) n-k
141 *> k-s s n-k
142 *>
143 *> If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k).
144 *> If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k),
145 *> and A(k,k), and v overwrites A(1:k-2,k-1:k).
146 *>
147 *> If UPLO = 'L', then A = L*D*L**H, where
148 *> L = P(1)*L(1)* ... *P(k)*L(k)* ...,
149 *> i.e., L is a product of terms P(k)*L(k), where k increases from 1 to
150 *> n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1
151 *> and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as
152 *> defined by IPIV(k), and L(k) is a unit lower triangular matrix, such
153 *> that if the diagonal block D(k) is of order s (s = 1 or 2), then
154 *>
155 *> ( I 0 0 ) k-1
156 *> L(k) = ( 0 I 0 ) s
157 *> ( 0 v I ) n-k-s+1
158 *> k-1 s n-k-s+1
159 *>
160 *> If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k).
161 *> If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k),
162 *> and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
163 *> \endverbatim
164 *
165 *> \par Contributors:
166 * ==================
167 *>
168 *> \verbatim
169 *> 09-29-06 - patch from
170 *> Bobby Cheng, MathWorks
171 *>
172 *> Replace l.210 and l.393
173 *> IF( MAX( ABSAKK, COLMAX ).EQ.ZERO ) THEN
174 *> by
175 *> IF( (MAX( ABSAKK, COLMAX ).EQ.ZERO) .OR. DISNAN(ABSAKK) ) THEN
176 *>
177 *> 01-01-96 - Based on modifications by
178 *> J. Lewis, Boeing Computer Services Company
179 *> A. Petitet, Computer Science Dept., Univ. of Tenn., Knoxville, USA
180 *> \endverbatim
181 *
182 * =====================================================================
183  SUBROUTINE zhetf2( UPLO, N, A, LDA, IPIV, INFO )
184 *
185 * -- LAPACK computational routine (version 3.4.2) --
186 * -- LAPACK is a software package provided by Univ. of Tennessee, --
187 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
188 * September 2012
189 *
190 * .. Scalar Arguments ..
191  CHARACTER uplo
192  INTEGER info, lda, n
193 * ..
194 * .. Array Arguments ..
195  INTEGER ipiv( * )
196  COMPLEX*16 a( lda, * )
197 * ..
198 *
199 * =====================================================================
200 *
201 * .. Parameters ..
202  DOUBLE PRECISION zero, one
203  parameter( zero = 0.0d+0, one = 1.0d+0 )
204  DOUBLE PRECISION eight, sevten
205  parameter( eight = 8.0d+0, sevten = 17.0d+0 )
206 * ..
207 * .. Local Scalars ..
208  LOGICAL upper
209  INTEGER i, imax, j, jmax, k, kk, kp, kstep
210  DOUBLE PRECISION absakk, alpha, colmax, d, d11, d22, r1, rowmax,
211  $ tt
212  COMPLEX*16 d12, d21, t, wk, wkm1, wkp1, zdum
213 * ..
214 * .. External Functions ..
215  LOGICAL lsame, disnan
216  INTEGER izamax
217  DOUBLE PRECISION dlapy2
218  EXTERNAL lsame, izamax, dlapy2, disnan
219 * ..
220 * .. External Subroutines ..
221  EXTERNAL xerbla, zdscal, zher, zswap
222 * ..
223 * .. Intrinsic Functions ..
224  INTRINSIC abs, dble, dcmplx, dconjg, dimag, max, sqrt
225 * ..
226 * .. Statement Functions ..
227  DOUBLE PRECISION cabs1
228 * ..
229 * .. Statement Function definitions ..
230  cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
231 * ..
232 * .. Executable Statements ..
233 *
234 * Test the input parameters.
235 *
236  info = 0
237  upper = lsame( uplo, 'U' )
238  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
239  info = -1
240  ELSE IF( n.LT.0 ) THEN
241  info = -2
242  ELSE IF( lda.LT.max( 1, n ) ) THEN
243  info = -4
244  END IF
245  IF( info.NE.0 ) THEN
246  CALL xerbla( 'ZHETF2', -info )
247  return
248  END IF
249 *
250 * Initialize ALPHA for use in choosing pivot block size.
251 *
252  alpha = ( one+sqrt( sevten ) ) / eight
253 *
254  IF( upper ) THEN
255 *
256 * Factorize A as U*D*U**H using the upper triangle of A
257 *
258 * K is the main loop index, decreasing from N to 1 in steps of
259 * 1 or 2
260 *
261  k = n
262  10 continue
263 *
264 * If K < 1, exit from loop
265 *
266  IF( k.LT.1 )
267  $ go to 90
268  kstep = 1
269 *
270 * Determine rows and columns to be interchanged and whether
271 * a 1-by-1 or 2-by-2 pivot block will be used
272 *
273  absakk = abs( dble( a( k, k ) ) )
274 *
275 * IMAX is the row-index of the largest off-diagonal element in
276 * column K, and COLMAX is its absolute value
277 *
278  IF( k.GT.1 ) THEN
279  imax = izamax( k-1, a( 1, k ), 1 )
280  colmax = cabs1( a( imax, k ) )
281  ELSE
282  colmax = zero
283  END IF
284 *
285  IF( (max( absakk, colmax ).EQ.zero) .OR. disnan(absakk) ) THEN
286 *
287 * Column K is zero or contains a NaN: set INFO and continue
288 *
289  IF( info.EQ.0 )
290  $ info = k
291  kp = k
292  a( k, k ) = dble( a( k, k ) )
293  ELSE
294  IF( absakk.GE.alpha*colmax ) THEN
295 *
296 * no interchange, use 1-by-1 pivot block
297 *
298  kp = k
299  ELSE
300 *
301 * JMAX is the column-index of the largest off-diagonal
302 * element in row IMAX, and ROWMAX is its absolute value
303 *
304  jmax = imax + izamax( k-imax, a( imax, imax+1 ), lda )
305  rowmax = cabs1( a( imax, jmax ) )
306  IF( imax.GT.1 ) THEN
307  jmax = izamax( imax-1, a( 1, imax ), 1 )
308  rowmax = max( rowmax, cabs1( a( jmax, imax ) ) )
309  END IF
310 *
311  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
312 *
313 * no interchange, use 1-by-1 pivot block
314 *
315  kp = k
316  ELSE IF( abs( dble( a( imax, imax ) ) ).GE.alpha*rowmax )
317  $ THEN
318 *
319 * interchange rows and columns K and IMAX, use 1-by-1
320 * pivot block
321 *
322  kp = imax
323  ELSE
324 *
325 * interchange rows and columns K-1 and IMAX, use 2-by-2
326 * pivot block
327 *
328  kp = imax
329  kstep = 2
330  END IF
331  END IF
332 *
333  kk = k - kstep + 1
334  IF( kp.NE.kk ) THEN
335 *
336 * Interchange rows and columns KK and KP in the leading
337 * submatrix A(1:k,1:k)
338 *
339  CALL zswap( kp-1, a( 1, kk ), 1, a( 1, kp ), 1 )
340  DO 20 j = kp + 1, kk - 1
341  t = dconjg( a( j, kk ) )
342  a( j, kk ) = dconjg( a( kp, j ) )
343  a( kp, j ) = t
344  20 continue
345  a( kp, kk ) = dconjg( a( kp, kk ) )
346  r1 = dble( a( kk, kk ) )
347  a( kk, kk ) = dble( a( kp, kp ) )
348  a( kp, kp ) = r1
349  IF( kstep.EQ.2 ) THEN
350  a( k, k ) = dble( a( k, k ) )
351  t = a( k-1, k )
352  a( k-1, k ) = a( kp, k )
353  a( kp, k ) = t
354  END IF
355  ELSE
356  a( k, k ) = dble( a( k, k ) )
357  IF( kstep.EQ.2 )
358  $ a( k-1, k-1 ) = dble( a( k-1, k-1 ) )
359  END IF
360 *
361 * Update the leading submatrix
362 *
363  IF( kstep.EQ.1 ) THEN
364 *
365 * 1-by-1 pivot block D(k): column k now holds
366 *
367 * W(k) = U(k)*D(k)
368 *
369 * where U(k) is the k-th column of U
370 *
371 * Perform a rank-1 update of A(1:k-1,1:k-1) as
372 *
373 * A := A - U(k)*D(k)*U(k)**H = A - W(k)*1/D(k)*W(k)**H
374 *
375  r1 = one / dble( a( k, k ) )
376  CALL zher( uplo, k-1, -r1, a( 1, k ), 1, a, lda )
377 *
378 * Store U(k) in column k
379 *
380  CALL zdscal( k-1, r1, a( 1, k ), 1 )
381  ELSE
382 *
383 * 2-by-2 pivot block D(k): columns k and k-1 now hold
384 *
385 * ( W(k-1) W(k) ) = ( U(k-1) U(k) )*D(k)
386 *
387 * where U(k) and U(k-1) are the k-th and (k-1)-th columns
388 * of U
389 *
390 * Perform a rank-2 update of A(1:k-2,1:k-2) as
391 *
392 * A := A - ( U(k-1) U(k) )*D(k)*( U(k-1) U(k) )**H
393 * = A - ( W(k-1) W(k) )*inv(D(k))*( W(k-1) W(k) )**H
394 *
395  IF( k.GT.2 ) THEN
396 *
397  d = dlapy2( dble( a( k-1, k ) ),
398  $ dimag( a( k-1, k ) ) )
399  d22 = dble( a( k-1, k-1 ) ) / d
400  d11 = dble( a( k, k ) ) / d
401  tt = one / ( d11*d22-one )
402  d12 = a( k-1, k ) / d
403  d = tt / d
404 *
405  DO 40 j = k - 2, 1, -1
406  wkm1 = d*( d11*a( j, k-1 )-dconjg( d12 )*
407  $ a( j, k ) )
408  wk = d*( d22*a( j, k )-d12*a( j, k-1 ) )
409  DO 30 i = j, 1, -1
410  a( i, j ) = a( i, j ) - a( i, k )*dconjg( wk ) -
411  $ a( i, k-1 )*dconjg( wkm1 )
412  30 continue
413  a( j, k ) = wk
414  a( j, k-1 ) = wkm1
415  a( j, j ) = dcmplx( dble( a( j, j ) ), 0.0d+0 )
416  40 continue
417 *
418  END IF
419 *
420  END IF
421  END IF
422 *
423 * Store details of the interchanges in IPIV
424 *
425  IF( kstep.EQ.1 ) THEN
426  ipiv( k ) = kp
427  ELSE
428  ipiv( k ) = -kp
429  ipiv( k-1 ) = -kp
430  END IF
431 *
432 * Decrease K and return to the start of the main loop
433 *
434  k = k - kstep
435  go to 10
436 *
437  ELSE
438 *
439 * Factorize A as L*D*L**H using the lower triangle of A
440 *
441 * K is the main loop index, increasing from 1 to N in steps of
442 * 1 or 2
443 *
444  k = 1
445  50 continue
446 *
447 * If K > N, exit from loop
448 *
449  IF( k.GT.n )
450  $ go to 90
451  kstep = 1
452 *
453 * Determine rows and columns to be interchanged and whether
454 * a 1-by-1 or 2-by-2 pivot block will be used
455 *
456  absakk = abs( dble( a( k, k ) ) )
457 *
458 * IMAX is the row-index of the largest off-diagonal element in
459 * column K, and COLMAX is its absolute value
460 *
461  IF( k.LT.n ) THEN
462  imax = k + izamax( n-k, a( k+1, k ), 1 )
463  colmax = cabs1( a( imax, k ) )
464  ELSE
465  colmax = zero
466  END IF
467 *
468  IF( (max( absakk, colmax ).EQ.zero) .OR. disnan(absakk) ) THEN
469 *
470 * Column K is zero or contains a NaN: set INFO and continue
471 *
472  IF( info.EQ.0 )
473  $ info = k
474  kp = k
475  a( k, k ) = dble( a( k, k ) )
476  ELSE
477  IF( absakk.GE.alpha*colmax ) THEN
478 *
479 * no interchange, use 1-by-1 pivot block
480 *
481  kp = k
482  ELSE
483 *
484 * JMAX is the column-index of the largest off-diagonal
485 * element in row IMAX, and ROWMAX is its absolute value
486 *
487  jmax = k - 1 + izamax( imax-k, a( imax, k ), lda )
488  rowmax = cabs1( a( imax, jmax ) )
489  IF( imax.LT.n ) THEN
490  jmax = imax + izamax( n-imax, a( imax+1, imax ), 1 )
491  rowmax = max( rowmax, cabs1( a( jmax, imax ) ) )
492  END IF
493 *
494  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
495 *
496 * no interchange, use 1-by-1 pivot block
497 *
498  kp = k
499  ELSE IF( abs( dble( a( imax, imax ) ) ).GE.alpha*rowmax )
500  $ THEN
501 *
502 * interchange rows and columns K and IMAX, use 1-by-1
503 * pivot block
504 *
505  kp = imax
506  ELSE
507 *
508 * interchange rows and columns K+1 and IMAX, use 2-by-2
509 * pivot block
510 *
511  kp = imax
512  kstep = 2
513  END IF
514  END IF
515 *
516  kk = k + kstep - 1
517  IF( kp.NE.kk ) THEN
518 *
519 * Interchange rows and columns KK and KP in the trailing
520 * submatrix A(k:n,k:n)
521 *
522  IF( kp.LT.n )
523  $ CALL zswap( n-kp, a( kp+1, kk ), 1, a( kp+1, kp ), 1 )
524  DO 60 j = kk + 1, kp - 1
525  t = dconjg( a( j, kk ) )
526  a( j, kk ) = dconjg( a( kp, j ) )
527  a( kp, j ) = t
528  60 continue
529  a( kp, kk ) = dconjg( a( kp, kk ) )
530  r1 = dble( a( kk, kk ) )
531  a( kk, kk ) = dble( a( kp, kp ) )
532  a( kp, kp ) = r1
533  IF( kstep.EQ.2 ) THEN
534  a( k, k ) = dble( a( k, k ) )
535  t = a( k+1, k )
536  a( k+1, k ) = a( kp, k )
537  a( kp, k ) = t
538  END IF
539  ELSE
540  a( k, k ) = dble( a( k, k ) )
541  IF( kstep.EQ.2 )
542  $ a( k+1, k+1 ) = dble( a( k+1, k+1 ) )
543  END IF
544 *
545 * Update the trailing submatrix
546 *
547  IF( kstep.EQ.1 ) THEN
548 *
549 * 1-by-1 pivot block D(k): column k now holds
550 *
551 * W(k) = L(k)*D(k)
552 *
553 * where L(k) is the k-th column of L
554 *
555  IF( k.LT.n ) THEN
556 *
557 * Perform a rank-1 update of A(k+1:n,k+1:n) as
558 *
559 * A := A - L(k)*D(k)*L(k)**H = A - W(k)*(1/D(k))*W(k)**H
560 *
561  r1 = one / dble( a( k, k ) )
562  CALL zher( uplo, n-k, -r1, a( k+1, k ), 1,
563  $ a( k+1, k+1 ), lda )
564 *
565 * Store L(k) in column K
566 *
567  CALL zdscal( n-k, r1, a( k+1, k ), 1 )
568  END IF
569  ELSE
570 *
571 * 2-by-2 pivot block D(k)
572 *
573  IF( k.LT.n-1 ) THEN
574 *
575 * Perform a rank-2 update of A(k+2:n,k+2:n) as
576 *
577 * A := A - ( L(k) L(k+1) )*D(k)*( L(k) L(k+1) )**H
578 * = A - ( W(k) W(k+1) )*inv(D(k))*( W(k) W(k+1) )**H
579 *
580 * where L(k) and L(k+1) are the k-th and (k+1)-th
581 * columns of L
582 *
583  d = dlapy2( dble( a( k+1, k ) ),
584  $ dimag( a( k+1, k ) ) )
585  d11 = dble( a( k+1, k+1 ) ) / d
586  d22 = dble( a( k, k ) ) / d
587  tt = one / ( d11*d22-one )
588  d21 = a( k+1, k ) / d
589  d = tt / d
590 *
591  DO 80 j = k + 2, n
592  wk = d*( d11*a( j, k )-d21*a( j, k+1 ) )
593  wkp1 = d*( d22*a( j, k+1 )-dconjg( d21 )*
594  $ a( j, k ) )
595  DO 70 i = j, n
596  a( i, j ) = a( i, j ) - a( i, k )*dconjg( wk ) -
597  $ a( i, k+1 )*dconjg( wkp1 )
598  70 continue
599  a( j, k ) = wk
600  a( j, k+1 ) = wkp1
601  a( j, j ) = dcmplx( dble( a( j, j ) ), 0.0d+0 )
602  80 continue
603  END IF
604  END IF
605  END IF
606 *
607 * Store details of the interchanges in IPIV
608 *
609  IF( kstep.EQ.1 ) THEN
610  ipiv( k ) = kp
611  ELSE
612  ipiv( k ) = -kp
613  ipiv( k+1 ) = -kp
614  END IF
615 *
616 * Increase K and return to the start of the main loop
617 *
618  k = k + kstep
619  go to 50
620 *
621  END IF
622 *
623  90 continue
624  return
625 *
626 * End of ZHETF2
627 *
628  END