LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ chetrf_rook()

subroutine chetrf_rook ( character  uplo,
integer  n,
complex, dimension( lda, * )  a,
integer  lda,
integer, dimension( * )  ipiv,
complex, dimension( * )  work,
integer  lwork,
integer  info 
)

CHETRF_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 CHETRF_ROOK + dependencies [TGZ] [ZIP] [TXT]

Purpose:
 CHETRF_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 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 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.
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 211 of file chetrf_rook.f.

212*
213* -- LAPACK computational routine --
214* -- LAPACK is a software package provided by Univ. of Tennessee, --
215* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
216*
217* .. Scalar Arguments ..
218 CHARACTER UPLO
219 INTEGER INFO, LDA, LWORK, N
220* ..
221* .. Array Arguments ..
222 INTEGER IPIV( * )
223 COMPLEX A( LDA, * ), WORK( * )
224* ..
225*
226* =====================================================================
227*
228* .. Local Scalars ..
229 LOGICAL LQUERY, UPPER
230 INTEGER IINFO, IWS, J, K, KB, LDWORK, LWKOPT, NB, NBMIN
231* ..
232* .. External Functions ..
233 LOGICAL LSAME
234 INTEGER ILAENV
235 REAL SROUNDUP_LWORK
236 EXTERNAL lsame, ilaenv, sroundup_lwork
237* ..
238* .. External Subroutines ..
240* ..
241* .. Intrinsic Functions ..
242 INTRINSIC max
243* ..
244* .. Executable Statements ..
245*
246* Test the input parameters.
247*
248 info = 0
249 upper = lsame( uplo, 'U' )
250 lquery = ( lwork.EQ.-1 )
251 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
252 info = -1
253 ELSE IF( n.LT.0 ) THEN
254 info = -2
255 ELSE IF( lda.LT.max( 1, n ) ) THEN
256 info = -4
257 ELSE IF( lwork.LT.1 .AND. .NOT.lquery ) THEN
258 info = -7
259 END IF
260*
261 IF( info.EQ.0 ) THEN
262*
263* Determine the block size
264*
265 nb = ilaenv( 1, 'CHETRF_ROOK', uplo, n, -1, -1, -1 )
266 lwkopt = max( 1, n*nb )
267 work( 1 ) = sroundup_lwork(lwkopt)
268 END IF
269*
270 IF( info.NE.0 ) THEN
271 CALL xerbla( 'CHETRF_ROOK', -info )
272 RETURN
273 ELSE IF( lquery ) THEN
274 RETURN
275 END IF
276*
277 nbmin = 2
278 ldwork = n
279 IF( nb.GT.1 .AND. nb.LT.n ) THEN
280 iws = ldwork*nb
281 IF( lwork.LT.iws ) THEN
282 nb = max( lwork / ldwork, 1 )
283 nbmin = max( 2, ilaenv( 2, 'CHETRF_ROOK',
284 $ uplo, n, -1, -1, -1 ) )
285 END IF
286 ELSE
287 iws = 1
288 END IF
289 IF( nb.LT.nbmin )
290 $ nb = n
291*
292 IF( upper ) THEN
293*
294* Factorize A as U*D*U**T using the upper triangle of A
295*
296* K is the main loop index, decreasing from N to 1 in steps of
297* KB, where KB is the number of columns factorized by CLAHEF_ROOK;
298* KB is either NB or NB-1, or K for the last block
299*
300 k = n
301 10 CONTINUE
302*
303* If K < 1, exit from loop
304*
305 IF( k.LT.1 )
306 $ GO TO 40
307*
308 IF( k.GT.nb ) THEN
309*
310* Factorize columns k-kb+1:k of A and use blocked code to
311* update columns 1:k-kb
312*
313 CALL clahef_rook( uplo, k, nb, kb, a, lda,
314 $ ipiv, work, ldwork, iinfo )
315 ELSE
316*
317* Use unblocked code to factorize columns 1:k of A
318*
319 CALL chetf2_rook( uplo, k, a, lda, ipiv, iinfo )
320 kb = k
321 END IF
322*
323* Set INFO on the first occurrence of a zero pivot
324*
325 IF( info.EQ.0 .AND. iinfo.GT.0 )
326 $ info = iinfo
327*
328* No need to adjust IPIV
329*
330* Decrease K and return to the start of the main loop
331*
332 k = k - kb
333 GO TO 10
334*
335 ELSE
336*
337* Factorize A as L*D*L**T using the lower triangle of A
338*
339* K is the main loop index, increasing from 1 to N in steps of
340* KB, where KB is the number of columns factorized by CLAHEF_ROOK;
341* KB is either NB or NB-1, or N-K+1 for the last block
342*
343 k = 1
344 20 CONTINUE
345*
346* If K > N, exit from loop
347*
348 IF( k.GT.n )
349 $ GO TO 40
350*
351 IF( k.LE.n-nb ) THEN
352*
353* Factorize columns k:k+kb-1 of A and use blocked code to
354* update columns k+kb:n
355*
356 CALL clahef_rook( uplo, n-k+1, nb, kb, a( k, k ), lda,
357 $ ipiv( k ), work, ldwork, iinfo )
358 ELSE
359*
360* Use unblocked code to factorize columns k:n of A
361*
362 CALL chetf2_rook( uplo, n-k+1, a( k, k ), lda, ipiv( k ),
363 $ iinfo )
364 kb = n - k + 1
365 END IF
366*
367* Set INFO on the first occurrence of a zero pivot
368*
369 IF( info.EQ.0 .AND. iinfo.GT.0 )
370 $ info = iinfo + k - 1
371*
372* Adjust IPIV
373*
374 DO 30 j = k, k + kb - 1
375 IF( ipiv( j ).GT.0 ) THEN
376 ipiv( j ) = ipiv( j ) + k - 1
377 ELSE
378 ipiv( j ) = ipiv( j ) - k + 1
379 END IF
380 30 CONTINUE
381*
382* Increase K and return to the start of the main loop
383*
384 k = k + kb
385 GO TO 20
386*
387 END IF
388*
389 40 CONTINUE
390 work( 1 ) = sroundup_lwork(lwkopt)
391 RETURN
392*
393* End of CHETRF_ROOK
394*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine chetf2_rook(uplo, n, a, lda, ipiv, info)
CHETF2_ROOK computes the factorization of a complex Hermitian indefinite matrix using the bounded Bun...
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine clahef_rook(uplo, n, nb, kb, a, lda, ipiv, w, ldw, info)
Download CLAHEF_ROOK + dependencies <a href="http://www.netlib.org/cgi-bin/netlibfiles....
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
Here is the call graph for this function:
Here is the caller graph for this function: