LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zhetrf_rook ( character  UPLO,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
integer, dimension( * )  IPIV,
complex*16, dimension( * )  WORK,
integer  LWORK,
integer  INFO 
)

ZHETRF_ROOK computes the factorization of a complex Hermitian indefinite matrix using the bounded Bunch-Kaufman ("rook") diagonal pivoting method (blocked algorithm, calling Level 3 BLAS).

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

Purpose:
 ZHETRF_ROOK computes the factorization of a complex Hermitian matrix A
 using the bounded Bunch-Kaufman ("rook") diagonal pivoting method.
 The form of the factorization is

    A = U*D*U**T  or  A = L*D*L**T

 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.

 This is the blocked version of the algorithm, calling Level 3 BLAS.
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]A
          A is COMPLEX*16 array, dimension (LDA,N)
          On entry, the Hermitian matrix A.  If UPLO = 'U', the leading
          N-by-N upper triangular part of A contains the upper
          triangular part of the matrix A, and the strictly lower
          triangular part of A is not referenced.  If UPLO = 'L', the
          leading N-by-N lower triangular part of A contains the lower
          triangular part of the matrix A, and the strictly upper
          triangular part of A is not referenced.

          On exit, the block diagonal matrix D and the multipliers used
          to obtain the factor U or L (see below for further details).
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]IPIV
          IPIV is INTEGER array, dimension (N)
          Details of the interchanges and the block structure of D.

          If UPLO = 'U':
             Only the last KB elements of IPIV are set.

             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 IPIV(k) < 0 and IPIV(k-1) < 0, then rows and
             columns k and -IPIV(k) were interchanged and rows and
             columns k-1 and -IPIV(k-1) were inerchaged,
             D(k-1:k,k-1:k) is a 2-by-2 diagonal block.

          If UPLO = 'L':
             Only the first KB elements of IPIV are set.

             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 IPIV(k) < 0 and IPIV(k+1) < 0, then rows and
             columns k and -IPIV(k) were interchanged and rows and
             columns k+1 and -IPIV(k+1) were inerchaged,
             D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
[out]WORK
          WORK is COMPLEX*16 array, dimension (MAX(1,LWORK)).
          On exit, if INFO = 0, WORK(1) returns the optimal LWORK.
[in]LWORK
          LWORK is INTEGER
          The length of WORK.  LWORK >=1.  For best performance
          LWORK >= N*NB, where NB is the block size returned by ILAENV.

          If LWORK = -1, then a workspace query is assumed; the routine
          only calculates the optimal size of the WORK array, returns
          this value as the first entry of the WORK array, and no error
          message related to LWORK is issued by XERBLA.
[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
June 2016
Further Details:
  If UPLO = 'U', then A = U*D*U**T, 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**T, 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:
  June 2016,  Igor Kozachenko,
                  Computer Science Division,
                  University of California, Berkeley

  September 2007, Sven Hammarling, Nicholas J. Higham, Craig Lucas,
                  School of Mathematics,
                  University of Manchester

Definition at line 214 of file zhetrf_rook.f.

214 *
215 * -- LAPACK computational routine (version 3.6.1) --
216 * -- LAPACK is a software package provided by Univ. of Tennessee, --
217 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
218 * June 2016
219 *
220 * .. Scalar Arguments ..
221  CHARACTER uplo
222  INTEGER info, lda, lwork, n
223 * ..
224 * .. Array Arguments ..
225  INTEGER ipiv( * )
226  COMPLEX*16 a( lda, * ), work( * )
227 * ..
228 *
229 * =====================================================================
230 *
231 * .. Local Scalars ..
232  LOGICAL lquery, upper
233  INTEGER iinfo, iws, j, k, kb, ldwork, lwkopt, nb, nbmin
234 * ..
235 * .. External Functions ..
236  LOGICAL lsame
237  INTEGER ilaenv
238  EXTERNAL lsame, ilaenv
239 * ..
240 * .. External Subroutines ..
241  EXTERNAL zlahef_rook, zhetf2_rook, xerbla
242 * ..
243 * .. Intrinsic Functions ..
244  INTRINSIC max
245 * ..
246 * .. Executable Statements ..
247 *
248 * Test the input parameters.
249 *
250  info = 0
251  upper = lsame( uplo, 'U' )
252  lquery = ( lwork.EQ.-1 )
253  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
254  info = -1
255  ELSE IF( n.LT.0 ) THEN
256  info = -2
257  ELSE IF( lda.LT.max( 1, n ) ) THEN
258  info = -4
259  ELSE IF( lwork.LT.1 .AND. .NOT.lquery ) THEN
260  info = -7
261  END IF
262 *
263  IF( info.EQ.0 ) THEN
264 *
265 * Determine the block size
266 *
267  nb = ilaenv( 1, 'ZHETRF_ROOK', uplo, n, -1, -1, -1 )
268  lwkopt = max( 1, n*nb )
269  work( 1 ) = lwkopt
270  END IF
271 *
272  IF( info.NE.0 ) THEN
273  CALL xerbla( 'ZHETRF_ROOK', -info )
274  RETURN
275  ELSE IF( lquery ) THEN
276  RETURN
277  END IF
278 *
279  nbmin = 2
280  ldwork = n
281  IF( nb.GT.1 .AND. nb.LT.n ) THEN
282  iws = ldwork*nb
283  IF( lwork.LT.iws ) THEN
284  nb = max( lwork / ldwork, 1 )
285  nbmin = max( 2, ilaenv( 2, 'ZHETRF_ROOK',
286  $ uplo, n, -1, -1, -1 ) )
287  END IF
288  ELSE
289  iws = 1
290  END IF
291  IF( nb.LT.nbmin )
292  $ nb = n
293 *
294  IF( upper ) THEN
295 *
296 * Factorize A as U*D*U**T using the upper triangle of A
297 *
298 * K is the main loop index, decreasing from N to 1 in steps of
299 * KB, where KB is the number of columns factorized by ZLAHEF_ROOK;
300 * KB is either NB or NB-1, or K for the last block
301 *
302  k = n
303  10 CONTINUE
304 *
305 * If K < 1, exit from loop
306 *
307  IF( k.LT.1 )
308  $ GO TO 40
309 *
310  IF( k.GT.nb ) THEN
311 *
312 * Factorize columns k-kb+1:k of A and use blocked code to
313 * update columns 1:k-kb
314 *
315  CALL zlahef_rook( uplo, k, nb, kb, a, lda,
316  $ ipiv, work, ldwork, iinfo )
317  ELSE
318 *
319 * Use unblocked code to factorize columns 1:k of A
320 *
321  CALL zhetf2_rook( uplo, k, a, lda, ipiv, iinfo )
322  kb = k
323  END IF
324 *
325 * Set INFO on the first occurrence of a zero pivot
326 *
327  IF( info.EQ.0 .AND. iinfo.GT.0 )
328  $ info = iinfo
329 *
330 * No need to adjust IPIV
331 *
332 * Decrease K and return to the start of the main loop
333 *
334  k = k - kb
335  GO TO 10
336 *
337  ELSE
338 *
339 * Factorize A as L*D*L**T using the lower triangle of A
340 *
341 * K is the main loop index, increasing from 1 to N in steps of
342 * KB, where KB is the number of columns factorized by ZLAHEF_ROOK;
343 * KB is either NB or NB-1, or N-K+1 for the last block
344 *
345  k = 1
346  20 CONTINUE
347 *
348 * If K > N, exit from loop
349 *
350  IF( k.GT.n )
351  $ GO TO 40
352 *
353  IF( k.LE.n-nb ) THEN
354 *
355 * Factorize columns k:k+kb-1 of A and use blocked code to
356 * update columns k+kb:n
357 *
358  CALL zlahef_rook( uplo, n-k+1, nb, kb, a( k, k ), lda,
359  $ ipiv( k ), work, ldwork, iinfo )
360  ELSE
361 *
362 * Use unblocked code to factorize columns k:n of A
363 *
364  CALL zhetf2_rook( uplo, n-k+1, a( k, k ), lda, ipiv( k ),
365  $ iinfo )
366  kb = n - k + 1
367  END IF
368 *
369 * Set INFO on the first occurrence of a zero pivot
370 *
371  IF( info.EQ.0 .AND. iinfo.GT.0 )
372  $ info = iinfo + k - 1
373 *
374 * Adjust IPIV
375 *
376  DO 30 j = k, k + kb - 1
377  IF( ipiv( j ).GT.0 ) THEN
378  ipiv( j ) = ipiv( j ) + k - 1
379  ELSE
380  ipiv( j ) = ipiv( j ) - k + 1
381  END IF
382  30 CONTINUE
383 *
384 * Increase K and return to the start of the main loop
385 *
386  k = k + kb
387  GO TO 20
388 *
389  END IF
390 *
391  40 CONTINUE
392  work( 1 ) = lwkopt
393  RETURN
394 *
395 * End of ZHETRF_ROOK
396 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zlahef_rook(UPLO, N, NB, KB, A, LDA, IPIV, W, LDW, INFO)
Download ZLAHEF_ROOK + dependencies [TGZ] [ZIP] [TXT]
Definition: zlahef_rook.f:186
subroutine zhetf2_rook(UPLO, N, A, LDA, IPIV, INFO)
ZHETF2_ROOK computes the factorization of a complex Hermitian indefinite matrix using the bounded Bun...
Definition: zhetf2_rook.f:196
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83
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: