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

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

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

Purpose:
 SPSTRF 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 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, 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 REAL
          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 REAL 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 143 of file spstrf.f.

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