LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zhptrf ( character  UPLO,
integer  N,
complex*16, dimension( * )  AP,
integer, dimension( * )  IPIV,
integer  INFO 
)

ZHPTRF

Download ZHPTRF + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 ZHPTRF computes the factorization of a complex Hermitian packed
 matrix A using the Bunch-Kaufman diagonal pivoting method:

    A = U*D*U**H  or  A = L*D*L**H

 where U (or L) is a product of permutation and unit upper (lower)
 triangular matrices, and D is Hermitian and block diagonal with
 1-by-1 and 2-by-2 diagonal blocks.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]AP
          AP is COMPLEX*16 array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the Hermitian matrix
          A, packed columnwise in a linear array.  The j-th column of A
          is stored in the array AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.

          On exit, the block diagonal matrix D and the multipliers used
          to obtain the factor U or L, stored as a packed triangular
          matrix overwriting A (see below for further details).
[out]IPIV
          IPIV is INTEGER array, dimension (N)
          Details of the interchanges and the block structure of D.
          If IPIV(k) > 0, then rows and columns k and IPIV(k) were
          interchanged and D(k,k) is a 1-by-1 diagonal block.
          If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0, then rows and
          columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k)
          is a 2-by-2 diagonal block.  If UPLO = 'L' and IPIV(k) =
          IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were
          interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -i, the i-th argument had an illegal value
          > 0: if INFO = i, D(i,i) is exactly zero.  The factorization
               has been completed, but the block diagonal matrix D is
               exactly singular, and division by zero will occur if it
               is used to solve a system of equations.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011
Further Details:
  If UPLO = 'U', then A = U*D*U**H, where
     U = P(n)*U(n)* ... *P(k)U(k)* ...,
  i.e., U is a product of terms P(k)*U(k), where k decreases from n to
  1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1
  and 2-by-2 diagonal blocks D(k).  P(k) is a permutation matrix as
  defined by IPIV(k), and U(k) is a unit upper triangular matrix, such
  that if the diagonal block D(k) is of order s (s = 1 or 2), then

             (   I    v    0   )   k-s
     U(k) =  (   0    I    0   )   s
             (   0    0    I   )   n-k
                k-s   s   n-k

  If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k).
  If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k),
  and A(k,k), and v overwrites A(1:k-2,k-1:k).

  If UPLO = 'L', then A = L*D*L**H, where
     L = P(1)*L(1)* ... *P(k)*L(k)* ...,
  i.e., L is a product of terms P(k)*L(k), where k increases from 1 to
  n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1
  and 2-by-2 diagonal blocks D(k).  P(k) is a permutation matrix as
  defined by IPIV(k), and L(k) is a unit lower triangular matrix, such
  that if the diagonal block D(k) is of order s (s = 1 or 2), then

             (   I    0     0   )  k-1
     L(k) =  (   0    I     0   )  s
             (   0    v     I   )  n-k-s+1
                k-1   s  n-k-s+1

  If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k).
  If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k),
  and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
Contributors:
J. Lewis, Boeing Computer Services Company

Definition at line 161 of file zhptrf.f.

161 *
162 * -- LAPACK computational routine (version 3.4.0) --
163 * -- LAPACK is a software package provided by Univ. of Tennessee, --
164 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
165 * November 2011
166 *
167 * .. Scalar Arguments ..
168  CHARACTER uplo
169  INTEGER info, n
170 * ..
171 * .. Array Arguments ..
172  INTEGER ipiv( * )
173  COMPLEX*16 ap( * )
174 * ..
175 *
176 * =====================================================================
177 *
178 * .. Parameters ..
179  DOUBLE PRECISION zero, one
180  parameter ( zero = 0.0d+0, one = 1.0d+0 )
181  DOUBLE PRECISION eight, sevten
182  parameter ( eight = 8.0d+0, sevten = 17.0d+0 )
183 * ..
184 * .. Local Scalars ..
185  LOGICAL upper
186  INTEGER i, imax, j, jmax, k, kc, kk, knc, kp, kpc,
187  $ kstep, kx, npp
188  DOUBLE PRECISION absakk, alpha, colmax, d, d11, d22, r1, rowmax,
189  $ tt
190  COMPLEX*16 d12, d21, t, wk, wkm1, wkp1, zdum
191 * ..
192 * .. External Functions ..
193  LOGICAL lsame
194  INTEGER izamax
195  DOUBLE PRECISION dlapy2
196  EXTERNAL lsame, izamax, dlapy2
197 * ..
198 * .. External Subroutines ..
199  EXTERNAL xerbla, zdscal, zhpr, zswap
200 * ..
201 * .. Intrinsic Functions ..
202  INTRINSIC abs, dble, dcmplx, dconjg, dimag, max, sqrt
203 * ..
204 * .. Statement Functions ..
205  DOUBLE PRECISION cabs1
206 * ..
207 * .. Statement Function definitions ..
208  cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
209 * ..
210 * .. Executable Statements ..
211 *
212 * Test the input parameters.
213 *
214  info = 0
215  upper = lsame( uplo, 'U' )
216  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
217  info = -1
218  ELSE IF( n.LT.0 ) THEN
219  info = -2
220  END IF
221  IF( info.NE.0 ) THEN
222  CALL xerbla( 'ZHPTRF', -info )
223  RETURN
224  END IF
225 *
226 * Initialize ALPHA for use in choosing pivot block size.
227 *
228  alpha = ( one+sqrt( sevten ) ) / eight
229 *
230  IF( upper ) THEN
231 *
232 * Factorize A as U*D*U**H using the upper triangle of A
233 *
234 * K is the main loop index, decreasing from N to 1 in steps of
235 * 1 or 2
236 *
237  k = n
238  kc = ( n-1 )*n / 2 + 1
239  10 CONTINUE
240  knc = kc
241 *
242 * If K < 1, exit from loop
243 *
244  IF( k.LT.1 )
245  $ GO TO 110
246  kstep = 1
247 *
248 * Determine rows and columns to be interchanged and whether
249 * a 1-by-1 or 2-by-2 pivot block will be used
250 *
251  absakk = abs( dble( ap( kc+k-1 ) ) )
252 *
253 * IMAX is the row-index of the largest off-diagonal element in
254 * column K, and COLMAX is its absolute value
255 *
256  IF( k.GT.1 ) THEN
257  imax = izamax( k-1, ap( kc ), 1 )
258  colmax = cabs1( ap( kc+imax-1 ) )
259  ELSE
260  colmax = zero
261  END IF
262 *
263  IF( max( absakk, colmax ).EQ.zero ) THEN
264 *
265 * Column K is zero: set INFO and continue
266 *
267  IF( info.EQ.0 )
268  $ info = k
269  kp = k
270  ap( kc+k-1 ) = dble( ap( kc+k-1 ) )
271  ELSE
272  IF( absakk.GE.alpha*colmax ) THEN
273 *
274 * no interchange, use 1-by-1 pivot block
275 *
276  kp = k
277  ELSE
278 *
279 * JMAX is the column-index of the largest off-diagonal
280 * element in row IMAX, and ROWMAX is its absolute value
281 *
282  rowmax = zero
283  jmax = imax
284  kx = imax*( imax+1 ) / 2 + imax
285  DO 20 j = imax + 1, k
286  IF( cabs1( ap( kx ) ).GT.rowmax ) THEN
287  rowmax = cabs1( ap( kx ) )
288  jmax = j
289  END IF
290  kx = kx + j
291  20 CONTINUE
292  kpc = ( imax-1 )*imax / 2 + 1
293  IF( imax.GT.1 ) THEN
294  jmax = izamax( imax-1, ap( kpc ), 1 )
295  rowmax = max( rowmax, cabs1( ap( kpc+jmax-1 ) ) )
296  END IF
297 *
298  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
299 *
300 * no interchange, use 1-by-1 pivot block
301 *
302  kp = k
303  ELSE IF( abs( dble( ap( kpc+imax-1 ) ) ).GE.alpha*
304  $ rowmax ) THEN
305 *
306 * interchange rows and columns K and IMAX, use 1-by-1
307 * pivot block
308 *
309  kp = imax
310  ELSE
311 *
312 * interchange rows and columns K-1 and IMAX, use 2-by-2
313 * pivot block
314 *
315  kp = imax
316  kstep = 2
317  END IF
318  END IF
319 *
320  kk = k - kstep + 1
321  IF( kstep.EQ.2 )
322  $ knc = knc - k + 1
323  IF( kp.NE.kk ) THEN
324 *
325 * Interchange rows and columns KK and KP in the leading
326 * submatrix A(1:k,1:k)
327 *
328  CALL zswap( kp-1, ap( knc ), 1, ap( kpc ), 1 )
329  kx = kpc + kp - 1
330  DO 30 j = kp + 1, kk - 1
331  kx = kx + j - 1
332  t = dconjg( ap( knc+j-1 ) )
333  ap( knc+j-1 ) = dconjg( ap( kx ) )
334  ap( kx ) = t
335  30 CONTINUE
336  ap( kx+kk-1 ) = dconjg( ap( kx+kk-1 ) )
337  r1 = dble( ap( knc+kk-1 ) )
338  ap( knc+kk-1 ) = dble( ap( kpc+kp-1 ) )
339  ap( kpc+kp-1 ) = r1
340  IF( kstep.EQ.2 ) THEN
341  ap( kc+k-1 ) = dble( ap( kc+k-1 ) )
342  t = ap( kc+k-2 )
343  ap( kc+k-2 ) = ap( kc+kp-1 )
344  ap( kc+kp-1 ) = t
345  END IF
346  ELSE
347  ap( kc+k-1 ) = dble( ap( kc+k-1 ) )
348  IF( kstep.EQ.2 )
349  $ ap( kc-1 ) = dble( ap( kc-1 ) )
350  END IF
351 *
352 * Update the leading submatrix
353 *
354  IF( kstep.EQ.1 ) THEN
355 *
356 * 1-by-1 pivot block D(k): column k now holds
357 *
358 * W(k) = U(k)*D(k)
359 *
360 * where U(k) is the k-th column of U
361 *
362 * Perform a rank-1 update of A(1:k-1,1:k-1) as
363 *
364 * A := A - U(k)*D(k)*U(k)**H = A - W(k)*1/D(k)*W(k)**H
365 *
366  r1 = one / dble( ap( kc+k-1 ) )
367  CALL zhpr( uplo, k-1, -r1, ap( kc ), 1, ap )
368 *
369 * Store U(k) in column k
370 *
371  CALL zdscal( k-1, r1, ap( kc ), 1 )
372  ELSE
373 *
374 * 2-by-2 pivot block D(k): columns k and k-1 now hold
375 *
376 * ( W(k-1) W(k) ) = ( U(k-1) U(k) )*D(k)
377 *
378 * where U(k) and U(k-1) are the k-th and (k-1)-th columns
379 * of U
380 *
381 * Perform a rank-2 update of A(1:k-2,1:k-2) as
382 *
383 * A := A - ( U(k-1) U(k) )*D(k)*( U(k-1) U(k) )**H
384 * = A - ( W(k-1) W(k) )*inv(D(k))*( W(k-1) W(k) )**H
385 *
386  IF( k.GT.2 ) THEN
387 *
388  d = dlapy2( dble( ap( k-1+( k-1 )*k / 2 ) ),
389  $ dimag( ap( k-1+( k-1 )*k / 2 ) ) )
390  d22 = dble( ap( k-1+( k-2 )*( k-1 ) / 2 ) ) / d
391  d11 = dble( ap( k+( k-1 )*k / 2 ) ) / d
392  tt = one / ( d11*d22-one )
393  d12 = ap( k-1+( k-1 )*k / 2 ) / d
394  d = tt / d
395 *
396  DO 50 j = k - 2, 1, -1
397  wkm1 = d*( d11*ap( j+( k-2 )*( k-1 ) / 2 )-
398  $ dconjg( d12 )*ap( j+( k-1 )*k / 2 ) )
399  wk = d*( d22*ap( j+( k-1 )*k / 2 )-d12*
400  $ ap( j+( k-2 )*( k-1 ) / 2 ) )
401  DO 40 i = j, 1, -1
402  ap( i+( j-1 )*j / 2 ) = ap( i+( j-1 )*j / 2 ) -
403  $ ap( i+( k-1 )*k / 2 )*dconjg( wk ) -
404  $ ap( i+( k-2 )*( k-1 ) / 2 )*dconjg( wkm1 )
405  40 CONTINUE
406  ap( j+( k-1 )*k / 2 ) = wk
407  ap( j+( k-2 )*( k-1 ) / 2 ) = wkm1
408  ap( j+( j-1 )*j / 2 ) = dcmplx( dble( ap( j+( j-
409  $ 1 )*j / 2 ) ), 0.0d+0 )
410  50 CONTINUE
411 *
412  END IF
413 *
414  END IF
415  END IF
416 *
417 * Store details of the interchanges in IPIV
418 *
419  IF( kstep.EQ.1 ) THEN
420  ipiv( k ) = kp
421  ELSE
422  ipiv( k ) = -kp
423  ipiv( k-1 ) = -kp
424  END IF
425 *
426 * Decrease K and return to the start of the main loop
427 *
428  k = k - kstep
429  kc = knc - k
430  GO TO 10
431 *
432  ELSE
433 *
434 * Factorize A as L*D*L**H using the lower triangle of A
435 *
436 * K is the main loop index, increasing from 1 to N in steps of
437 * 1 or 2
438 *
439  k = 1
440  kc = 1
441  npp = n*( n+1 ) / 2
442  60 CONTINUE
443  knc = kc
444 *
445 * If K > N, exit from loop
446 *
447  IF( k.GT.n )
448  $ GO TO 110
449  kstep = 1
450 *
451 * Determine rows and columns to be interchanged and whether
452 * a 1-by-1 or 2-by-2 pivot block will be used
453 *
454  absakk = abs( dble( ap( kc ) ) )
455 *
456 * IMAX is the row-index of the largest off-diagonal element in
457 * column K, and COLMAX is its absolute value
458 *
459  IF( k.LT.n ) THEN
460  imax = k + izamax( n-k, ap( kc+1 ), 1 )
461  colmax = cabs1( ap( kc+imax-k ) )
462  ELSE
463  colmax = zero
464  END IF
465 *
466  IF( max( absakk, colmax ).EQ.zero ) THEN
467 *
468 * Column K is zero: set INFO and continue
469 *
470  IF( info.EQ.0 )
471  $ info = k
472  kp = k
473  ap( kc ) = dble( ap( kc ) )
474  ELSE
475  IF( absakk.GE.alpha*colmax ) THEN
476 *
477 * no interchange, use 1-by-1 pivot block
478 *
479  kp = k
480  ELSE
481 *
482 * JMAX is the column-index of the largest off-diagonal
483 * element in row IMAX, and ROWMAX is its absolute value
484 *
485  rowmax = zero
486  kx = kc + imax - k
487  DO 70 j = k, imax - 1
488  IF( cabs1( ap( kx ) ).GT.rowmax ) THEN
489  rowmax = cabs1( ap( kx ) )
490  jmax = j
491  END IF
492  kx = kx + n - j
493  70 CONTINUE
494  kpc = npp - ( n-imax+1 )*( n-imax+2 ) / 2 + 1
495  IF( imax.LT.n ) THEN
496  jmax = imax + izamax( n-imax, ap( kpc+1 ), 1 )
497  rowmax = max( rowmax, cabs1( ap( kpc+jmax-imax ) ) )
498  END IF
499 *
500  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
501 *
502 * no interchange, use 1-by-1 pivot block
503 *
504  kp = k
505  ELSE IF( abs( dble( ap( kpc ) ) ).GE.alpha*rowmax ) THEN
506 *
507 * interchange rows and columns K and IMAX, use 1-by-1
508 * pivot block
509 *
510  kp = imax
511  ELSE
512 *
513 * interchange rows and columns K+1 and IMAX, use 2-by-2
514 * pivot block
515 *
516  kp = imax
517  kstep = 2
518  END IF
519  END IF
520 *
521  kk = k + kstep - 1
522  IF( kstep.EQ.2 )
523  $ knc = knc + n - k + 1
524  IF( kp.NE.kk ) THEN
525 *
526 * Interchange rows and columns KK and KP in the trailing
527 * submatrix A(k:n,k:n)
528 *
529  IF( kp.LT.n )
530  $ CALL zswap( n-kp, ap( knc+kp-kk+1 ), 1, ap( kpc+1 ),
531  $ 1 )
532  kx = knc + kp - kk
533  DO 80 j = kk + 1, kp - 1
534  kx = kx + n - j + 1
535  t = dconjg( ap( knc+j-kk ) )
536  ap( knc+j-kk ) = dconjg( ap( kx ) )
537  ap( kx ) = t
538  80 CONTINUE
539  ap( knc+kp-kk ) = dconjg( ap( knc+kp-kk ) )
540  r1 = dble( ap( knc ) )
541  ap( knc ) = dble( ap( kpc ) )
542  ap( kpc ) = r1
543  IF( kstep.EQ.2 ) THEN
544  ap( kc ) = dble( ap( kc ) )
545  t = ap( kc+1 )
546  ap( kc+1 ) = ap( kc+kp-k )
547  ap( kc+kp-k ) = t
548  END IF
549  ELSE
550  ap( kc ) = dble( ap( kc ) )
551  IF( kstep.EQ.2 )
552  $ ap( knc ) = dble( ap( knc ) )
553  END IF
554 *
555 * Update the trailing submatrix
556 *
557  IF( kstep.EQ.1 ) THEN
558 *
559 * 1-by-1 pivot block D(k): column k now holds
560 *
561 * W(k) = L(k)*D(k)
562 *
563 * where L(k) is the k-th column of L
564 *
565  IF( k.LT.n ) THEN
566 *
567 * Perform a rank-1 update of A(k+1:n,k+1:n) as
568 *
569 * A := A - L(k)*D(k)*L(k)**H = A - W(k)*(1/D(k))*W(k)**H
570 *
571  r1 = one / dble( ap( kc ) )
572  CALL zhpr( uplo, n-k, -r1, ap( kc+1 ), 1,
573  $ ap( kc+n-k+1 ) )
574 *
575 * Store L(k) in column K
576 *
577  CALL zdscal( n-k, r1, ap( kc+1 ), 1 )
578  END IF
579  ELSE
580 *
581 * 2-by-2 pivot block D(k): columns K and K+1 now hold
582 *
583 * ( W(k) W(k+1) ) = ( L(k) L(k+1) )*D(k)
584 *
585 * where L(k) and L(k+1) are the k-th and (k+1)-th columns
586 * of L
587 *
588  IF( k.LT.n-1 ) THEN
589 *
590 * Perform a rank-2 update of A(k+2:n,k+2:n) as
591 *
592 * A := A - ( L(k) L(k+1) )*D(k)*( L(k) L(k+1) )**H
593 * = A - ( W(k) W(k+1) )*inv(D(k))*( W(k) W(k+1) )**H
594 *
595 * where L(k) and L(k+1) are the k-th and (k+1)-th
596 * columns of L
597 *
598  d = dlapy2( dble( ap( k+1+( k-1 )*( 2*n-k ) / 2 ) ),
599  $ dimag( ap( k+1+( k-1 )*( 2*n-k ) / 2 ) ) )
600  d11 = dble( ap( k+1+k*( 2*n-k-1 ) / 2 ) ) / d
601  d22 = dble( ap( k+( k-1 )*( 2*n-k ) / 2 ) ) / d
602  tt = one / ( d11*d22-one )
603  d21 = ap( k+1+( k-1 )*( 2*n-k ) / 2 ) / d
604  d = tt / d
605 *
606  DO 100 j = k + 2, n
607  wk = d*( d11*ap( j+( k-1 )*( 2*n-k ) / 2 )-d21*
608  $ ap( j+k*( 2*n-k-1 ) / 2 ) )
609  wkp1 = d*( d22*ap( j+k*( 2*n-k-1 ) / 2 )-
610  $ dconjg( d21 )*ap( j+( k-1 )*( 2*n-k ) /
611  $ 2 ) )
612  DO 90 i = j, n
613  ap( i+( j-1 )*( 2*n-j ) / 2 ) = ap( i+( j-1 )*
614  $ ( 2*n-j ) / 2 ) - ap( i+( k-1 )*( 2*n-k ) /
615  $ 2 )*dconjg( wk ) - ap( i+k*( 2*n-k-1 ) / 2 )*
616  $ dconjg( wkp1 )
617  90 CONTINUE
618  ap( j+( k-1 )*( 2*n-k ) / 2 ) = wk
619  ap( j+k*( 2*n-k-1 ) / 2 ) = wkp1
620  ap( j+( j-1 )*( 2*n-j ) / 2 )
621  $ = dcmplx( dble( ap( j+( j-1 )*( 2*n-j ) / 2 ) ),
622  $ 0.0d+0 )
623  100 CONTINUE
624  END IF
625  END IF
626  END IF
627 *
628 * Store details of the interchanges in IPIV
629 *
630  IF( kstep.EQ.1 ) THEN
631  ipiv( k ) = kp
632  ELSE
633  ipiv( k ) = -kp
634  ipiv( k+1 ) = -kp
635  END IF
636 *
637 * Increase K and return to the start of the main loop
638 *
639  k = k + kstep
640  kc = knc + n - k + 2
641  GO TO 60
642 *
643  END IF
644 *
645  110 CONTINUE
646  RETURN
647 *
648 * End of ZHPTRF
649 *
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:52
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zhpr(UPLO, N, ALPHA, X, INCX, AP)
ZHPR
Definition: zhpr.f:132
integer function izamax(N, ZX, INCX)
IZAMAX
Definition: izamax.f:53
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:54
double precision function dlapy2(X, Y)
DLAPY2 returns sqrt(x2+y2).
Definition: dlapy2.f:65
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: