LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dpstrf ( character  UPLO,
integer  N,
double precision, dimension( lda, * )  A,
integer  LDA,
integer, dimension( n )  PIV,
integer  RANK,
double precision  TOL,
double precision, dimension( 2*n )  WORK,
integer  INFO 
)

DPSTRF computes the Cholesky factorization with complete pivoting of a real symmetric positive semidefinite matrix.

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

Purpose:
 DPSTRF computes the Cholesky factorization with complete
 pivoting of a real symmetric positive semidefinite matrix A.

 The factorization has the form
    P**T * A * P = U**T * U ,  if UPLO = 'U',
    P**T * A * P = L  * L**T,  if UPLO = 'L',
 where U is an upper triangular matrix and L is lower triangular, and
 P is stored as vector PIV.

 This algorithm does not attempt to check that A is positive
 semidefinite. This version of the algorithm calls level 3 BLAS.
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 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, if INFO = 0, the factor U or L from the Cholesky
          factorization as above.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[out]PIV
          PIV is INTEGER array, dimension (N)
          PIV is such that the nonzero entries are P( PIV(K), K ) = 1.
[out]RANK
          RANK is INTEGER
          The rank of A given by the number of steps the algorithm
          completed.
[in]TOL
          TOL is DOUBLE PRECISION
          User defined tolerance. If TOL < 0, then N*U*MAX( A(K,K) )
          will be used. The algorithm terminates at the (K-1)st step
          if the pivot <= TOL.
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (2*N)
          Work space.
[out]INFO
          INFO is INTEGER
          < 0: If INFO = -K, the K-th argument had an illegal value,
          = 0: algorithm completed successfully, and
          > 0: the matrix A is either rank deficient with computed rank
               as returned in RANK, or is not positive semidefinite. See
               Section 7 of LAPACK Working Note #161 for further
               information.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2015

Definition at line 144 of file dpstrf.f.

144 *
145 * -- LAPACK computational routine (version 3.6.0) --
146 * -- LAPACK is a software package provided by Univ. of Tennessee, --
147 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
148 * November 2015
149 *
150 * .. Scalar Arguments ..
151  DOUBLE PRECISION tol
152  INTEGER info, lda, n, rank
153  CHARACTER uplo
154 * ..
155 * .. Array Arguments ..
156  DOUBLE PRECISION a( lda, * ), work( 2*n )
157  INTEGER piv( n )
158 * ..
159 *
160 * =====================================================================
161 *
162 * .. Parameters ..
163  DOUBLE PRECISION one, zero
164  parameter ( one = 1.0d+0, zero = 0.0d+0 )
165 * ..
166 * .. Local Scalars ..
167  DOUBLE PRECISION ajj, dstop, dtemp
168  INTEGER i, itemp, j, jb, k, nb, pvt
169  LOGICAL upper
170 * ..
171 * .. External Functions ..
172  DOUBLE PRECISION dlamch
173  INTEGER ilaenv
174  LOGICAL lsame, disnan
175  EXTERNAL dlamch, ilaenv, lsame, disnan
176 * ..
177 * .. External Subroutines ..
178  EXTERNAL dgemv, dpstf2, dscal, dswap, dsyrk, xerbla
179 * ..
180 * .. Intrinsic Functions ..
181  INTRINSIC max, min, sqrt, maxloc
182 * ..
183 * .. Executable Statements ..
184 *
185 * Test the input parameters.
186 *
187  info = 0
188  upper = lsame( uplo, 'U' )
189  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
190  info = -1
191  ELSE IF( n.LT.0 ) THEN
192  info = -2
193  ELSE IF( lda.LT.max( 1, n ) ) THEN
194  info = -4
195  END IF
196  IF( info.NE.0 ) THEN
197  CALL xerbla( 'DPSTRF', -info )
198  RETURN
199  END IF
200 *
201 * Quick return if possible
202 *
203  IF( n.EQ.0 )
204  $ RETURN
205 *
206 * Get block size
207 *
208  nb = ilaenv( 1, 'DPOTRF', uplo, n, -1, -1, -1 )
209  IF( nb.LE.1 .OR. nb.GE.n ) THEN
210 *
211 * Use unblocked code
212 *
213  CALL dpstf2( uplo, n, a( 1, 1 ), lda, piv, rank, tol, work,
214  $ info )
215  GO TO 200
216 *
217  ELSE
218 *
219 * Initialize PIV
220 *
221  DO 100 i = 1, n
222  piv( i ) = i
223  100 CONTINUE
224 *
225 * Compute stopping value
226 *
227  pvt = 1
228  ajj = a( pvt, pvt )
229  DO i = 2, n
230  IF( a( i, i ).GT.ajj ) THEN
231  pvt = i
232  ajj = a( pvt, pvt )
233  END IF
234  END DO
235  IF( ajj.LE.zero.OR.disnan( ajj ) ) THEN
236  rank = 0
237  info = 1
238  GO TO 200
239  END IF
240 *
241 * Compute stopping value if not supplied
242 *
243  IF( tol.LT.zero ) THEN
244  dstop = n * dlamch( 'Epsilon' ) * ajj
245  ELSE
246  dstop = tol
247  END IF
248 *
249 *
250  IF( upper ) THEN
251 *
252 * Compute the Cholesky factorization P**T * A * P = U**T * U
253 *
254  DO 140 k = 1, n, nb
255 *
256 * Account for last block not being NB wide
257 *
258  jb = min( nb, n-k+1 )
259 *
260 * Set relevant part of first half of WORK to zero,
261 * holds dot products
262 *
263  DO 110 i = k, n
264  work( i ) = 0
265  110 CONTINUE
266 *
267  DO 130 j = k, k + jb - 1
268 *
269 * Find pivot, test for exit, else swap rows and columns
270 * Update dot products, compute possible pivots which are
271 * stored in the second half of WORK
272 *
273  DO 120 i = j, n
274 *
275  IF( j.GT.k ) THEN
276  work( i ) = work( i ) + a( j-1, i )**2
277  END IF
278  work( n+i ) = a( i, i ) - work( i )
279 *
280  120 CONTINUE
281 *
282  IF( j.GT.1 ) THEN
283  itemp = maxloc( work( (n+j):(2*n) ), 1 )
284  pvt = itemp + j - 1
285  ajj = work( n+pvt )
286  IF( ajj.LE.dstop.OR.disnan( ajj ) ) THEN
287  a( j, j ) = ajj
288  GO TO 190
289  END IF
290  END IF
291 *
292  IF( j.NE.pvt ) THEN
293 *
294 * Pivot OK, so can now swap pivot rows and columns
295 *
296  a( pvt, pvt ) = a( j, j )
297  CALL dswap( j-1, a( 1, j ), 1, a( 1, pvt ), 1 )
298  IF( pvt.LT.n )
299  $ CALL dswap( n-pvt, a( j, pvt+1 ), lda,
300  $ a( pvt, pvt+1 ), lda )
301  CALL dswap( pvt-j-1, a( j, j+1 ), lda,
302  $ a( j+1, pvt ), 1 )
303 *
304 * Swap dot products and PIV
305 *
306  dtemp = work( j )
307  work( j ) = work( pvt )
308  work( pvt ) = dtemp
309  itemp = piv( pvt )
310  piv( pvt ) = piv( j )
311  piv( j ) = itemp
312  END IF
313 *
314  ajj = sqrt( ajj )
315  a( j, j ) = ajj
316 *
317 * Compute elements J+1:N of row J.
318 *
319  IF( j.LT.n ) THEN
320  CALL dgemv( 'Trans', j-k, n-j, -one, a( k, j+1 ),
321  $ lda, a( k, j ), 1, one, a( j, j+1 ),
322  $ lda )
323  CALL dscal( n-j, one / ajj, a( j, j+1 ), lda )
324  END IF
325 *
326  130 CONTINUE
327 *
328 * Update trailing matrix, J already incremented
329 *
330  IF( k+jb.LE.n ) THEN
331  CALL dsyrk( 'Upper', 'Trans', n-j+1, jb, -one,
332  $ a( k, j ), lda, one, a( j, j ), lda )
333  END IF
334 *
335  140 CONTINUE
336 *
337  ELSE
338 *
339 * Compute the Cholesky factorization P**T * A * P = L * L**T
340 *
341  DO 180 k = 1, n, nb
342 *
343 * Account for last block not being NB wide
344 *
345  jb = min( nb, n-k+1 )
346 *
347 * Set relevant part of first half of WORK to zero,
348 * holds dot products
349 *
350  DO 150 i = k, n
351  work( i ) = 0
352  150 CONTINUE
353 *
354  DO 170 j = k, k + jb - 1
355 *
356 * Find pivot, test for exit, else swap rows and columns
357 * Update dot products, compute possible pivots which are
358 * stored in the second half of WORK
359 *
360  DO 160 i = j, n
361 *
362  IF( j.GT.k ) THEN
363  work( i ) = work( i ) + a( i, j-1 )**2
364  END IF
365  work( n+i ) = a( i, i ) - work( i )
366 *
367  160 CONTINUE
368 *
369  IF( j.GT.1 ) THEN
370  itemp = maxloc( work( (n+j):(2*n) ), 1 )
371  pvt = itemp + j - 1
372  ajj = work( n+pvt )
373  IF( ajj.LE.dstop.OR.disnan( ajj ) ) THEN
374  a( j, j ) = ajj
375  GO TO 190
376  END IF
377  END IF
378 *
379  IF( j.NE.pvt ) THEN
380 *
381 * Pivot OK, so can now swap pivot rows and columns
382 *
383  a( pvt, pvt ) = a( j, j )
384  CALL dswap( j-1, a( j, 1 ), lda, a( pvt, 1 ), lda )
385  IF( pvt.LT.n )
386  $ CALL dswap( n-pvt, a( pvt+1, j ), 1,
387  $ a( pvt+1, pvt ), 1 )
388  CALL dswap( pvt-j-1, a( j+1, j ), 1, a( pvt, j+1 ),
389  $ lda )
390 *
391 * Swap dot products and PIV
392 *
393  dtemp = work( j )
394  work( j ) = work( pvt )
395  work( pvt ) = dtemp
396  itemp = piv( pvt )
397  piv( pvt ) = piv( j )
398  piv( j ) = itemp
399  END IF
400 *
401  ajj = sqrt( ajj )
402  a( j, j ) = ajj
403 *
404 * Compute elements J+1:N of column J.
405 *
406  IF( j.LT.n ) THEN
407  CALL dgemv( 'No Trans', n-j, j-k, -one,
408  $ a( j+1, k ), lda, a( j, k ), lda, one,
409  $ a( j+1, j ), 1 )
410  CALL dscal( n-j, one / ajj, a( j+1, j ), 1 )
411  END IF
412 *
413  170 CONTINUE
414 *
415 * Update trailing matrix, J already incremented
416 *
417  IF( k+jb.LE.n ) THEN
418  CALL dsyrk( 'Lower', 'No Trans', n-j+1, jb, -one,
419  $ a( j, k ), lda, one, a( j, j ), lda )
420  END IF
421 *
422  180 CONTINUE
423 *
424  END IF
425  END IF
426 *
427 * Ran to completion, A has full rank
428 *
429  rank = n
430 *
431  GO TO 200
432  190 CONTINUE
433 *
434 * Rank is the number of steps completed. Set INFO = 1 to signal
435 * that the factorization cannot be used to solve a system.
436 *
437  rank = j - 1
438  info = 1
439 *
440  200 CONTINUE
441  RETURN
442 *
443 * End of DPSTRF
444 *
logical function disnan(DIN)
DISNAN tests input for NaN.
Definition: disnan.f:61
subroutine dpstf2(UPLO, N, A, LDA, PIV, RANK, TOL, WORK, INFO)
DPSTF2 computes the Cholesky factorization with complete pivoting of a real symmetric positive semide...
Definition: dpstf2.f:143
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine dgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
DGEMV
Definition: dgemv.f:158
subroutine dsyrk(UPLO, TRANS, N, K, ALPHA, A, LDA, BETA, C, LDC)
DSYRK
Definition: dsyrk.f:171
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
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: