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

◆ ssytrf_rook()

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

SSYTRF_ROOK

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

Purpose:
 SSYTRF_ROOK computes the factorization of a real symmetric 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 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.
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 REAL 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, 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':
               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':
               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 REAL 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 207 of file ssytrf_rook.f.

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