LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine spstf2 ( 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 
)

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

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

Purpose:
 SPSTF2 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 2 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.
[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.
[in]LDA
          LDA is INTEGER
          The leading dimension of the array A.  LDA >= max(1,N).
[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 spstf2.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, pvt
168  LOGICAL upper
169 * ..
170 * .. External Functions ..
171  REAL slamch
172  LOGICAL lsame, sisnan
173  EXTERNAL slamch, lsame, sisnan
174 * ..
175 * .. External Subroutines ..
176  EXTERNAL sgemv, sscal, sswap, xerbla
177 * ..
178 * .. Intrinsic Functions ..
179  INTRINSIC max, sqrt, maxloc
180 * ..
181 * .. Executable Statements ..
182 *
183 * Test the input parameters
184 *
185  info = 0
186  upper = lsame( uplo, 'U' )
187  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
188  info = -1
189  ELSE IF( n.LT.0 ) THEN
190  info = -2
191  ELSE IF( lda.LT.max( 1, n ) ) THEN
192  info = -4
193  END IF
194  IF( info.NE.0 ) THEN
195  CALL xerbla( 'SPSTF2', -info )
196  RETURN
197  END IF
198 *
199 * Quick return if possible
200 *
201  IF( n.EQ.0 )
202  $ RETURN
203 *
204 * Initialize PIV
205 *
206  DO 100 i = 1, n
207  piv( i ) = i
208  100 CONTINUE
209 *
210 * Compute stopping value
211 *
212  pvt = 1
213  ajj = a( pvt, pvt )
214  DO i = 2, n
215  IF( a( i, i ).GT.ajj ) THEN
216  pvt = i
217  ajj = a( pvt, pvt )
218  END IF
219  END DO
220  IF( ajj.LE.zero.OR.sisnan( ajj ) ) THEN
221  rank = 0
222  info = 1
223  GO TO 170
224  END IF
225 *
226 * Compute stopping value if not supplied
227 *
228  IF( tol.LT.zero ) THEN
229  sstop = n * slamch( 'Epsilon' ) * ajj
230  ELSE
231  sstop = tol
232  END IF
233 *
234 * Set first half of WORK to zero, holds dot products
235 *
236  DO 110 i = 1, n
237  work( i ) = 0
238  110 CONTINUE
239 *
240  IF( upper ) THEN
241 *
242 * Compute the Cholesky factorization P**T * A * P = U**T * U
243 *
244  DO 130 j = 1, n
245 *
246 * Find pivot, test for exit, else swap rows and columns
247 * Update dot products, compute possible pivots which are
248 * stored in the second half of WORK
249 *
250  DO 120 i = j, n
251 *
252  IF( j.GT.1 ) THEN
253  work( i ) = work( i ) + a( j-1, i )**2
254  END IF
255  work( n+i ) = a( i, i ) - work( i )
256 *
257  120 CONTINUE
258 *
259  IF( j.GT.1 ) THEN
260  itemp = maxloc( work( (n+j):(2*n) ), 1 )
261  pvt = itemp + j - 1
262  ajj = work( n+pvt )
263  IF( ajj.LE.sstop.OR.sisnan( ajj ) ) THEN
264  a( j, j ) = ajj
265  GO TO 160
266  END IF
267  END IF
268 *
269  IF( j.NE.pvt ) THEN
270 *
271 * Pivot OK, so can now swap pivot rows and columns
272 *
273  a( pvt, pvt ) = a( j, j )
274  CALL sswap( j-1, a( 1, j ), 1, a( 1, pvt ), 1 )
275  IF( pvt.LT.n )
276  $ CALL sswap( n-pvt, a( j, pvt+1 ), lda,
277  $ a( pvt, pvt+1 ), lda )
278  CALL sswap( pvt-j-1, a( j, j+1 ), lda, a( j+1, pvt ), 1 )
279 *
280 * Swap dot products and PIV
281 *
282  stemp = work( j )
283  work( j ) = work( pvt )
284  work( pvt ) = stemp
285  itemp = piv( pvt )
286  piv( pvt ) = piv( j )
287  piv( j ) = itemp
288  END IF
289 *
290  ajj = sqrt( ajj )
291  a( j, j ) = ajj
292 *
293 * Compute elements J+1:N of row J
294 *
295  IF( j.LT.n ) THEN
296  CALL sgemv( 'Trans', j-1, n-j, -one, a( 1, j+1 ), lda,
297  $ a( 1, j ), 1, one, a( j, j+1 ), lda )
298  CALL sscal( n-j, one / ajj, a( j, j+1 ), lda )
299  END IF
300 *
301  130 CONTINUE
302 *
303  ELSE
304 *
305 * Compute the Cholesky factorization P**T * A * P = L * L**T
306 *
307  DO 150 j = 1, n
308 *
309 * Find pivot, test for exit, else swap rows and columns
310 * Update dot products, compute possible pivots which are
311 * stored in the second half of WORK
312 *
313  DO 140 i = j, n
314 *
315  IF( j.GT.1 ) THEN
316  work( i ) = work( i ) + a( i, j-1 )**2
317  END IF
318  work( n+i ) = a( i, i ) - work( i )
319 *
320  140 CONTINUE
321 *
322  IF( j.GT.1 ) THEN
323  itemp = maxloc( work( (n+j):(2*n) ), 1 )
324  pvt = itemp + j - 1
325  ajj = work( n+pvt )
326  IF( ajj.LE.sstop.OR.sisnan( ajj ) ) THEN
327  a( j, j ) = ajj
328  GO TO 160
329  END IF
330  END IF
331 *
332  IF( j.NE.pvt ) THEN
333 *
334 * Pivot OK, so can now swap pivot rows and columns
335 *
336  a( pvt, pvt ) = a( j, j )
337  CALL sswap( j-1, a( j, 1 ), lda, a( pvt, 1 ), lda )
338  IF( pvt.LT.n )
339  $ CALL sswap( n-pvt, a( pvt+1, j ), 1, a( pvt+1, pvt ),
340  $ 1 )
341  CALL sswap( pvt-j-1, a( j+1, j ), 1, a( pvt, j+1 ), lda )
342 *
343 * Swap dot products and PIV
344 *
345  stemp = work( j )
346  work( j ) = work( pvt )
347  work( pvt ) = stemp
348  itemp = piv( pvt )
349  piv( pvt ) = piv( j )
350  piv( j ) = itemp
351  END IF
352 *
353  ajj = sqrt( ajj )
354  a( j, j ) = ajj
355 *
356 * Compute elements J+1:N of column J
357 *
358  IF( j.LT.n ) THEN
359  CALL sgemv( 'No Trans', n-j, j-1, -one, a( j+1, 1 ), lda,
360  $ a( j, 1 ), lda, one, a( j+1, j ), 1 )
361  CALL sscal( n-j, one / ajj, a( j+1, j ), 1 )
362  END IF
363 *
364  150 CONTINUE
365 *
366  END IF
367 *
368 * Ran to completion, A has full rank
369 *
370  rank = n
371 *
372  GO TO 170
373  160 CONTINUE
374 *
375 * Rank is number of steps completed. Set INFO = 1 to signal
376 * that the factorization cannot be used to solve a system.
377 *
378  rank = j - 1
379  info = 1
380 *
381  170 CONTINUE
382  RETURN
383 *
384 * End of SPSTF2
385 *
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
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: