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

◆ dsytrf()

subroutine dsytrf ( character  uplo,
integer  n,
double precision, dimension( lda, * )  a,
integer  lda,
integer, dimension( * )  ipiv,
double precision, dimension( * )  work,
integer  lwork,
integer  info 
)

DSYTRF

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

Purpose:
 DSYTRF computes the factorization of a real symmetric matrix A using
 the Bunch-Kaufman diagonal pivoting method.  The form of the
 factorization is

    A = U**T*D*U  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 DOUBLE PRECISION 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 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]WORK
          WORK is DOUBLE PRECISION 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**T*D*U, 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).

Definition at line 181 of file dsytrf.f.

182*
183* -- LAPACK computational routine --
184* -- LAPACK is a software package provided by Univ. of Tennessee, --
185* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
186*
187* .. Scalar Arguments ..
188 CHARACTER UPLO
189 INTEGER INFO, LDA, LWORK, N
190* ..
191* .. Array Arguments ..
192 INTEGER IPIV( * )
193 DOUBLE PRECISION A( LDA, * ), WORK( * )
194* ..
195*
196* =====================================================================
197*
198* .. Local Scalars ..
199 LOGICAL LQUERY, UPPER
200 INTEGER IINFO, IWS, J, K, KB, LDWORK, LWKOPT, NB, NBMIN
201* ..
202* .. External Functions ..
203 LOGICAL LSAME
204 INTEGER ILAENV
205 EXTERNAL lsame, ilaenv
206* ..
207* .. External Subroutines ..
208 EXTERNAL dlasyf, dsytf2, xerbla
209* ..
210* .. Intrinsic Functions ..
211 INTRINSIC max
212* ..
213* .. Executable Statements ..
214*
215* Test the input parameters.
216*
217 info = 0
218 upper = lsame( uplo, 'U' )
219 lquery = ( lwork.EQ.-1 )
220 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
221 info = -1
222 ELSE IF( n.LT.0 ) THEN
223 info = -2
224 ELSE IF( lda.LT.max( 1, n ) ) THEN
225 info = -4
226 ELSE IF( lwork.LT.1 .AND. .NOT.lquery ) THEN
227 info = -7
228 END IF
229*
230 IF( info.EQ.0 ) THEN
231*
232* Determine the block size
233*
234 nb = ilaenv( 1, 'DSYTRF', uplo, n, -1, -1, -1 )
235 lwkopt = max( 1, n*nb )
236 work( 1 ) = lwkopt
237 END IF
238*
239 IF( info.NE.0 ) THEN
240 CALL xerbla( 'DSYTRF', -info )
241 RETURN
242 ELSE IF( lquery ) THEN
243 RETURN
244 END IF
245*
246 nbmin = 2
247 ldwork = n
248 IF( nb.GT.1 .AND. nb.LT.n ) THEN
249 iws = ldwork*nb
250 IF( lwork.LT.iws ) THEN
251 nb = max( lwork / ldwork, 1 )
252 nbmin = max( 2, ilaenv( 2, 'DSYTRF', uplo, n, -1, -1, -1 ) )
253 END IF
254 ELSE
255 iws = 1
256 END IF
257 IF( nb.LT.nbmin )
258 $ nb = n
259*
260 IF( upper ) THEN
261*
262* Factorize A as U**T*D*U using the upper triangle of A
263*
264* K is the main loop index, decreasing from N to 1 in steps of
265* KB, where KB is the number of columns factorized by DLASYF;
266* KB is either NB or NB-1, or K for the last block
267*
268 k = n
269 10 CONTINUE
270*
271* If K < 1, exit from loop
272*
273 IF( k.LT.1 )
274 $ GO TO 40
275*
276 IF( k.GT.nb ) THEN
277*
278* Factorize columns k-kb+1:k of A and use blocked code to
279* update columns 1:k-kb
280*
281 CALL dlasyf( uplo, k, nb, kb, a, lda, ipiv, work, ldwork,
282 $ iinfo )
283 ELSE
284*
285* Use unblocked code to factorize columns 1:k of A
286*
287 CALL dsytf2( uplo, k, a, lda, ipiv, iinfo )
288 kb = k
289 END IF
290*
291* Set INFO on the first occurrence of a zero pivot
292*
293 IF( info.EQ.0 .AND. iinfo.GT.0 )
294 $ info = iinfo
295*
296* Decrease K and return to the start of the main loop
297*
298 k = k - kb
299 GO TO 10
300*
301 ELSE
302*
303* Factorize A as L*D*L**T using the lower triangle of A
304*
305* K is the main loop index, increasing from 1 to N in steps of
306* KB, where KB is the number of columns factorized by DLASYF;
307* KB is either NB or NB-1, or N-K+1 for the last block
308*
309 k = 1
310 20 CONTINUE
311*
312* If K > N, exit from loop
313*
314 IF( k.GT.n )
315 $ GO TO 40
316*
317 IF( k.LE.n-nb ) THEN
318*
319* Factorize columns k:k+kb-1 of A and use blocked code to
320* update columns k+kb:n
321*
322 CALL dlasyf( uplo, n-k+1, nb, kb, a( k, k ), lda, ipiv( k ),
323 $ work, ldwork, iinfo )
324 ELSE
325*
326* Use unblocked code to factorize columns k:n of A
327*
328 CALL dsytf2( uplo, n-k+1, a( k, k ), lda, ipiv( k ), iinfo )
329 kb = n - k + 1
330 END IF
331*
332* Set INFO on the first occurrence of a zero pivot
333*
334 IF( info.EQ.0 .AND. iinfo.GT.0 )
335 $ info = iinfo + k - 1
336*
337* Adjust IPIV
338*
339 DO 30 j = k, k + kb - 1
340 IF( ipiv( j ).GT.0 ) THEN
341 ipiv( j ) = ipiv( j ) + k - 1
342 ELSE
343 ipiv( j ) = ipiv( j ) - k + 1
344 END IF
345 30 CONTINUE
346*
347* Increase K and return to the start of the main loop
348*
349 k = k + kb
350 GO TO 20
351*
352 END IF
353*
354 40 CONTINUE
355 work( 1 ) = lwkopt
356 RETURN
357*
358* End of DSYTRF
359*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dsytf2(uplo, n, a, lda, ipiv, info)
DSYTF2 computes the factorization of a real symmetric indefinite matrix, using the diagonal pivoting ...
Definition dsytf2.f:194
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine dlasyf(uplo, n, nb, kb, a, lda, ipiv, w, ldw, info)
DLASYF computes a partial factorization of a real symmetric matrix using the Bunch-Kaufman diagonal p...
Definition dlasyf.f:176
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: