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

◆ csytrf_rk()

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

CSYTRF_RK computes the factorization of a complex symmetric indefinite matrix using the bounded Bunch-Kaufman (rook) diagonal pivoting method (BLAS3 blocked algorithm).

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

Purpose:
!> CSYTRF_RK computes the factorization of a complex symmetric matrix A
!> using the bounded Bunch-Kaufman (rook) diagonal pivoting method:
!>
!>    A = P*U*D*(U**T)*(P**T) or A = P*L*D*(L**T)*(P**T),
!>
!> where U (or L) is unit upper (or lower) triangular matrix,
!> U**T (or L**T) is the transpose of U (or L), P is a permutation
!> matrix, P**T is the transpose of P, and D is symmetric 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.
!> For more information see Further Details section.
!> 
Parameters
[in]UPLO
!>          UPLO is CHARACTER*1
!>          Specifies whether the upper or lower triangular part of the
!>          symmetric matrix A is stored:
!>          = 'U':  Upper triangular
!>          = 'L':  Lower triangular
!> 
[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 symmetric 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, contains:
!>            a) ONLY diagonal elements of the symmetric block diagonal
!>               matrix D on the diagonal of A, i.e. D(k,k) = A(k,k);
!>               (superdiagonal (or subdiagonal) elements of D
!>                are stored on exit in array E), and
!>            b) If UPLO = 'U': factor U in the superdiagonal part of A.
!>               If UPLO = 'L': factor L in the subdiagonal part of A.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[out]E
!>          E is COMPLEX array, dimension (N)
!>          On exit, contains the superdiagonal (or subdiagonal)
!>          elements of the symmetric block diagonal matrix D
!>          with 1-by-1 or 2-by-2 diagonal blocks, where
!>          If UPLO = 'U': E(i) = D(i-1,i), i=2:N, E(1) is set to 0;
!>          If UPLO = 'L': E(i) = D(i+1,i), i=1:N-1, E(N) is set to 0.
!>
!>          NOTE: For 1-by-1 diagonal block D(k), where
!>          1 <= k <= N, the element E(k) is set to 0 in both
!>          UPLO = 'U' or UPLO = 'L' cases.
!> 
[out]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          IPIV describes the permutation matrix P in the factorization
!>          of matrix A as follows. The absolute value of IPIV(k)
!>          represents the index of row and column that were
!>          interchanged with the k-th row and column. The value of UPLO
!>          describes the order in which the interchanges were applied.
!>          Also, the sign of IPIV represents the block structure of
!>          the symmetric block diagonal matrix D with 1-by-1 or 2-by-2
!>          diagonal blocks which correspond to 1 or 2 interchanges
!>          at each factorization step. For more info see Further
!>          Details section.
!>
!>          If UPLO = 'U',
!>          ( in factorization order, k decreases from N to 1 ):
!>            a) A single positive entry IPIV(k) > 0 means:
!>               D(k,k) is a 1-by-1 diagonal block.
!>               If IPIV(k) != k, rows and columns k and IPIV(k) were
!>               interchanged in the matrix A(1:N,1:N);
!>               If IPIV(k) = k, no interchange occurred.
!>
!>            b) A pair of consecutive negative entries
!>               IPIV(k) < 0 and IPIV(k-1) < 0 means:
!>               D(k-1:k,k-1:k) is a 2-by-2 diagonal block.
!>               (NOTE: negative entries in IPIV appear ONLY in pairs).
!>               1) If -IPIV(k) != k, rows and columns
!>                  k and -IPIV(k) were interchanged
!>                  in the matrix A(1:N,1:N).
!>                  If -IPIV(k) = k, no interchange occurred.
!>               2) If -IPIV(k-1) != k-1, rows and columns
!>                  k-1 and -IPIV(k-1) were interchanged
!>                  in the matrix A(1:N,1:N).
!>                  If -IPIV(k-1) = k-1, no interchange occurred.
!>
!>            c) In both cases a) and b), always ABS( IPIV(k) ) <= k.
!>
!>            d) NOTE: Any entry IPIV(k) is always NONZERO on output.
!>
!>          If UPLO = 'L',
!>          ( in factorization order, k increases from 1 to N ):
!>            a) A single positive entry IPIV(k) > 0 means:
!>               D(k,k) is a 1-by-1 diagonal block.
!>               If IPIV(k) != k, rows and columns k and IPIV(k) were
!>               interchanged in the matrix A(1:N,1:N).
!>               If IPIV(k) = k, no interchange occurred.
!>
!>            b) A pair of consecutive negative entries
!>               IPIV(k) < 0 and IPIV(k+1) < 0 means:
!>               D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
!>               (NOTE: negative entries in IPIV appear ONLY in pairs).
!>               1) If -IPIV(k) != k, rows and columns
!>                  k and -IPIV(k) were interchanged
!>                  in the matrix A(1:N,1:N).
!>                  If -IPIV(k) = k, no interchange occurred.
!>               2) If -IPIV(k+1) != k+1, rows and columns
!>                  k-1 and -IPIV(k-1) were interchanged
!>                  in the matrix A(1:N,1:N).
!>                  If -IPIV(k+1) = k+1, no interchange occurred.
!>
!>            c) In both cases a) and b), always ABS( IPIV(k) ) >= k.
!>
!>            d) NOTE: Any entry IPIV(k) is always NONZERO on output.
!> 
[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 = -k, the k-th argument had an illegal value
!>
!>          > 0: If INFO = k, the matrix A is singular, because:
!>                 If UPLO = 'U': column k in the upper
!>                 triangular part of A contains all zeros.
!>                 If UPLO = 'L': column k in the lower
!>                 triangular part of A contains all zeros.
!>
!>               Therefore D(k,k) is exactly zero, and superdiagonal
!>               elements of column k of U (or subdiagonal elements of
!>               column k of L ) are all zeros. 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.
!>
!>               NOTE: INFO only stores the first occurrence of
!>               a singularity, any subsequent occurrence of singularity
!>               is not stored in INFO even though the factorization
!>               always completes.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!> TODO: put correct description
!> 
Contributors:
!>
!>  December 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 255 of file csytrf_rk.f.

257*
258* -- LAPACK computational routine --
259* -- LAPACK is a software package provided by Univ. of Tennessee, --
260* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
261*
262* .. Scalar Arguments ..
263 CHARACTER UPLO
264 INTEGER INFO, LDA, LWORK, N
265* ..
266* .. Array Arguments ..
267 INTEGER IPIV( * )
268 COMPLEX A( LDA, * ), E( * ), WORK( * )
269* ..
270*
271* =====================================================================
272*
273* .. Local Scalars ..
274 LOGICAL LQUERY, UPPER
275 INTEGER I, IINFO, IP, IWS, K, KB, LDWORK, LWKOPT,
276 $ NB, NBMIN
277* ..
278* .. External Functions ..
279 LOGICAL LSAME
280 INTEGER ILAENV
281 REAL SROUNDUP_LWORK
282 EXTERNAL lsame, ilaenv, sroundup_lwork
283* ..
284* .. External Subroutines ..
285 EXTERNAL clasyf_rk, csytf2_rk, cswap, xerbla
286* ..
287* .. Intrinsic Functions ..
288 INTRINSIC abs, max
289* ..
290* .. Executable Statements ..
291*
292* Test the input parameters.
293*
294 info = 0
295 upper = lsame( uplo, 'U' )
296 lquery = ( lwork.EQ.-1 )
297 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
298 info = -1
299 ELSE IF( n.LT.0 ) THEN
300 info = -2
301 ELSE IF( lda.LT.max( 1, n ) ) THEN
302 info = -4
303 ELSE IF( lwork.LT.1 .AND. .NOT.lquery ) THEN
304 info = -8
305 END IF
306*
307 IF( info.EQ.0 ) THEN
308*
309* Determine the block size
310*
311 nb = ilaenv( 1, 'CSYTRF_RK', uplo, n, -1, -1, -1 )
312 lwkopt = max( 1, n*nb )
313 work( 1 ) = sroundup_lwork(lwkopt)
314 END IF
315*
316 IF( info.NE.0 ) THEN
317 CALL xerbla( 'CSYTRF_RK', -info )
318 RETURN
319 ELSE IF( lquery ) THEN
320 RETURN
321 END IF
322*
323 nbmin = 2
324 ldwork = n
325 IF( nb.GT.1 .AND. nb.LT.n ) THEN
326 iws = ldwork*nb
327 IF( lwork.LT.iws ) THEN
328 nb = max( lwork / ldwork, 1 )
329 nbmin = max( 2, ilaenv( 2, 'CSYTRF_RK',
330 $ uplo, n, -1, -1, -1 ) )
331 END IF
332 ELSE
333 iws = 1
334 END IF
335 IF( nb.LT.nbmin )
336 $ nb = n
337*
338 IF( upper ) THEN
339*
340* Factorize A as U*D*U**T using the upper triangle of A
341*
342* K is the main loop index, decreasing from N to 1 in steps of
343* KB, where KB is the number of columns factorized by CLASYF_RK;
344* KB is either NB or NB-1, or K for the last block
345*
346 k = n
347 10 CONTINUE
348*
349* If K < 1, exit from loop
350*
351 IF( k.LT.1 )
352 $ GO TO 15
353*
354 IF( k.GT.nb ) THEN
355*
356* Factorize columns k-kb+1:k of A and use blocked code to
357* update columns 1:k-kb
358*
359 CALL clasyf_rk( uplo, k, nb, kb, a, lda, e,
360 $ ipiv, work, ldwork, iinfo )
361 ELSE
362*
363* Use unblocked code to factorize columns 1:k of A
364*
365 CALL csytf2_rk( uplo, k, a, lda, e, ipiv, iinfo )
366 kb = k
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
373*
374* No need to adjust IPIV
375*
376*
377* Apply permutations to the leading panel 1:k-1
378*
379* Read IPIV from the last block factored, i.e.
380* indices k-kb+1:k and apply row permutations to the
381* last k+1 colunms k+1:N after that block
382* (We can do the simple loop over IPIV with decrement -1,
383* since the ABS value of IPIV( I ) represents the row index
384* of the interchange with row i in both 1x1 and 2x2 pivot cases)
385*
386 IF( k.LT.n ) THEN
387 DO i = k, ( k - kb + 1 ), -1
388 ip = abs( ipiv( i ) )
389 IF( ip.NE.i ) THEN
390 CALL cswap( n-k, a( i, k+1 ), lda,
391 $ a( ip, k+1 ), lda )
392 END IF
393 END DO
394 END IF
395*
396* Decrease K and return to the start of the main loop
397*
398 k = k - kb
399 GO TO 10
400*
401* This label is the exit from main loop over K decreasing
402* from N to 1 in steps of KB
403*
404 15 CONTINUE
405*
406 ELSE
407*
408* Factorize A as L*D*L**T using the lower triangle of A
409*
410* K is the main loop index, increasing from 1 to N in steps of
411* KB, where KB is the number of columns factorized by CLASYF_RK;
412* KB is either NB or NB-1, or N-K+1 for the last block
413*
414 k = 1
415 20 CONTINUE
416*
417* If K > N, exit from loop
418*
419 IF( k.GT.n )
420 $ GO TO 35
421*
422 IF( k.LE.n-nb ) THEN
423*
424* Factorize columns k:k+kb-1 of A and use blocked code to
425* update columns k+kb:n
426*
427 CALL clasyf_rk( uplo, n-k+1, nb, kb, a( k, k ), lda,
428 $ e( k ),
429 $ ipiv( k ), work, ldwork, iinfo )
430
431
432 ELSE
433*
434* Use unblocked code to factorize columns k:n of A
435*
436 CALL csytf2_rk( uplo, n-k+1, a( k, k ), lda, e( k ),
437 $ ipiv( k ), iinfo )
438 kb = n - k + 1
439*
440 END IF
441*
442* Set INFO on the first occurrence of a zero pivot
443*
444 IF( info.EQ.0 .AND. iinfo.GT.0 )
445 $ info = iinfo + k - 1
446*
447* Adjust IPIV
448*
449 DO i = k, k + kb - 1
450 IF( ipiv( i ).GT.0 ) THEN
451 ipiv( i ) = ipiv( i ) + k - 1
452 ELSE
453 ipiv( i ) = ipiv( i ) - k + 1
454 END IF
455 END DO
456*
457* Apply permutations to the leading panel 1:k-1
458*
459* Read IPIV from the last block factored, i.e.
460* indices k:k+kb-1 and apply row permutations to the
461* first k-1 colunms 1:k-1 before that block
462* (We can do the simple loop over IPIV with increment 1,
463* since the ABS value of IPIV( I ) represents the row index
464* of the interchange with row i in both 1x1 and 2x2 pivot cases)
465*
466 IF( k.GT.1 ) THEN
467 DO i = k, ( k + kb - 1 ), 1
468 ip = abs( ipiv( i ) )
469 IF( ip.NE.i ) THEN
470 CALL cswap( k-1, a( i, 1 ), lda,
471 $ a( ip, 1 ), lda )
472 END IF
473 END DO
474 END IF
475*
476* Increase K and return to the start of the main loop
477*
478 k = k + kb
479 GO TO 20
480*
481* This label is the exit from main loop over K increasing
482* from 1 to N in steps of KB
483*
484 35 CONTINUE
485*
486* End Lower
487*
488 END IF
489*
490 work( 1 ) = sroundup_lwork(lwkopt)
491 RETURN
492*
493* End of CSYTRF_RK
494*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine csytf2_rk(uplo, n, a, lda, e, ipiv, info)
CSYTF2_RK computes the factorization of a complex symmetric indefinite matrix using the bounded Bunch...
Definition csytf2_rk.f:239
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
subroutine clasyf_rk(uplo, n, nb, kb, a, lda, e, ipiv, w, ldw, info)
CLASYF_RK computes a partial factorization of a complex symmetric indefinite matrix using bounded Bun...
Definition clasyf_rk.f:260
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
real function sroundup_lwork(lwork)
SROUNDUP_LWORK
subroutine cswap(n, cx, incx, cy, incy)
CSWAP
Definition cswap.f:81
Here is the call graph for this function:
Here is the caller graph for this function: