LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
spstf2.f
Go to the documentation of this file.
1*> \brief \b SPSTF2 computes the Cholesky factorization with complete pivoting of a real symmetric positive semidefinite matrix.
2*
3* =========== DOCUMENTATION ===========
4*
5* Online html documentation available at
6* http://www.netlib.org/lapack/explore-html/
7*
8*> Download SPSTF2 + dependencies
9*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/spstf2.f">
10*> [TGZ]</a>
11*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/spstf2.f">
12*> [ZIP]</a>
13*> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/spstf2.f">
14*> [TXT]</a>
15*
16* Definition:
17* ===========
18*
19* SUBROUTINE SPSTF2( UPLO, N, A, LDA, PIV, RANK, TOL, WORK, INFO )
20*
21* .. Scalar Arguments ..
22* REAL TOL
23* INTEGER INFO, LDA, N, RANK
24* CHARACTER UPLO
25* ..
26* .. Array Arguments ..
27* REAL A( LDA, * ), WORK( 2*N )
28* INTEGER PIV( N )
29* ..
30*
31*
32*> \par Purpose:
33* =============
34*>
35*> \verbatim
36*>
37*> SPSTF2 computes the Cholesky factorization with complete
38*> pivoting of a real symmetric positive semidefinite matrix A.
39*>
40*> The factorization has the form
41*> P**T * A * P = U**T * U , if UPLO = 'U',
42*> P**T * A * P = L * L**T, if UPLO = 'L',
43*> where U is an upper triangular matrix and L is lower triangular, and
44*> P is stored as vector PIV.
45*>
46*> This algorithm does not attempt to check that A is positive
47*> semidefinite. This version of the algorithm calls level 2 BLAS.
48*> \endverbatim
49*
50* Arguments:
51* ==========
52*
53*> \param[in] UPLO
54*> \verbatim
55*> UPLO is CHARACTER*1
56*> Specifies whether the upper or lower triangular part of the
57*> symmetric matrix A is stored.
58*> = 'U': Upper triangular
59*> = 'L': Lower triangular
60*> \endverbatim
61*>
62*> \param[in] N
63*> \verbatim
64*> N is INTEGER
65*> The order of the matrix A. N >= 0.
66*> \endverbatim
67*>
68*> \param[in,out] A
69*> \verbatim
70*> A is REAL array, dimension (LDA,N)
71*> On entry, the symmetric matrix A. If UPLO = 'U', the leading
72*> n by n upper triangular part of A contains the upper
73*> triangular part of the matrix A, and the strictly lower
74*> triangular part of A is not referenced. If UPLO = 'L', the
75*> leading n by n lower triangular part of A contains the lower
76*> triangular part of the matrix A, and the strictly upper
77*> triangular part of A is not referenced.
78*>
79*> On exit, if INFO = 0, the factor U or L from the Cholesky
80*> factorization as above.
81*> \endverbatim
82*>
83*> \param[out] PIV
84*> \verbatim
85*> PIV is INTEGER array, dimension (N)
86*> PIV is such that the nonzero entries are P( PIV(K), K ) = 1.
87*> \endverbatim
88*>
89*> \param[out] RANK
90*> \verbatim
91*> RANK is INTEGER
92*> The rank of A given by the number of steps the algorithm
93*> completed.
94*> \endverbatim
95*>
96*> \param[in] TOL
97*> \verbatim
98*> TOL is REAL
99*> User defined tolerance. If TOL < 0, then N*U*MAX( A( K,K ) )
100*> will be used. The algorithm terminates at the (K-1)st step
101*> if the pivot <= TOL.
102*> \endverbatim
103*>
104*> \param[in] LDA
105*> \verbatim
106*> LDA is INTEGER
107*> The leading dimension of the array A. LDA >= max(1,N).
108*> \endverbatim
109*>
110*> \param[out] WORK
111*> \verbatim
112*> WORK is REAL array, dimension (2*N)
113*> Work space.
114*> \endverbatim
115*>
116*> \param[out] INFO
117*> \verbatim
118*> INFO is INTEGER
119*> < 0: If INFO = -K, the K-th argument had an illegal value,
120*> = 0: algorithm completed successfully, and
121*> > 0: the matrix A is either rank deficient with computed rank
122*> as returned in RANK, or is not positive semidefinite. See
123*> Section 7 of LAPACK Working Note #161 for further
124*> information.
125*> \endverbatim
126*
127* Authors:
128* ========
129*
130*> \author Univ. of Tennessee
131*> \author Univ. of California Berkeley
132*> \author Univ. of Colorado Denver
133*> \author NAG Ltd.
134*
135*> \ingroup pstf2
136*
137* =====================================================================
138 SUBROUTINE spstf2( UPLO, N, A, LDA, PIV, RANK, TOL, WORK,
139 $ INFO )
140*
141* -- LAPACK computational routine --
142* -- LAPACK is a software package provided by Univ. of Tennessee, --
143* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
144*
145* .. Scalar Arguments ..
146 REAL TOL
147 INTEGER INFO, LDA, N, RANK
148 CHARACTER UPLO
149* ..
150* .. Array Arguments ..
151 REAL A( LDA, * ), WORK( 2*N )
152 INTEGER PIV( N )
153* ..
154*
155* =====================================================================
156*
157* .. Parameters ..
158 REAL ONE, ZERO
159 parameter( one = 1.0e+0, zero = 0.0e+0 )
160* ..
161* .. Local Scalars ..
162 REAL AJJ, SSTOP, STEMP
163 INTEGER I, ITEMP, J, PVT
164 LOGICAL UPPER
165* ..
166* .. External Functions ..
167 REAL SLAMCH
168 LOGICAL LSAME, SISNAN
169 EXTERNAL slamch, lsame, sisnan
170* ..
171* .. External Subroutines ..
172 EXTERNAL sgemv, sscal, sswap, xerbla
173* ..
174* .. Intrinsic Functions ..
175 INTRINSIC max, sqrt, maxloc
176* ..
177* .. Executable Statements ..
178*
179* Test the input parameters
180*
181 info = 0
182 upper = lsame( uplo, 'U' )
183 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
184 info = -1
185 ELSE IF( n.LT.0 ) THEN
186 info = -2
187 ELSE IF( lda.LT.max( 1, n ) ) THEN
188 info = -4
189 END IF
190 IF( info.NE.0 ) THEN
191 CALL xerbla( 'SPSTF2', -info )
192 RETURN
193 END IF
194*
195* Quick return if possible
196*
197 IF( n.EQ.0 )
198 $ RETURN
199*
200* Initialize PIV
201*
202 DO 100 i = 1, n
203 piv( i ) = i
204 100 CONTINUE
205*
206* Compute stopping value
207*
208 pvt = 1
209 ajj = a( pvt, pvt )
210 DO i = 2, n
211 IF( a( i, i ).GT.ajj ) THEN
212 pvt = i
213 ajj = a( pvt, pvt )
214 END IF
215 END DO
216 IF( ajj.LE.zero.OR.sisnan( ajj ) ) THEN
217 rank = 0
218 info = 1
219 GO TO 170
220 END IF
221*
222* Compute stopping value if not supplied
223*
224 IF( tol.LT.zero ) THEN
225 sstop = real( n ) * slamch( 'Epsilon' ) * ajj
226 ELSE
227 sstop = tol
228 END IF
229*
230* Set first half of WORK to zero, holds dot products
231*
232 DO 110 i = 1, n
233 work( i ) = 0
234 110 CONTINUE
235*
236 IF( upper ) THEN
237*
238* Compute the Cholesky factorization P**T * A * P = U**T * U
239*
240 DO 130 j = 1, n
241*
242* Find pivot, test for exit, else swap rows and columns
243* Update dot products, compute possible pivots which are
244* stored in the second half of WORK
245*
246 DO 120 i = j, n
247*
248 IF( j.GT.1 ) THEN
249 work( i ) = work( i ) + a( j-1, i )**2
250 END IF
251 work( n+i ) = a( i, i ) - work( i )
252*
253 120 CONTINUE
254*
255 IF( j.GT.1 ) THEN
256 itemp = maxloc( work( (n+j):(2*n) ), 1 )
257 pvt = itemp + j - 1
258 ajj = work( n+pvt )
259 IF( ajj.LE.sstop.OR.sisnan( ajj ) ) THEN
260 a( j, j ) = ajj
261 GO TO 160
262 END IF
263 END IF
264*
265 IF( j.NE.pvt ) THEN
266*
267* Pivot OK, so can now swap pivot rows and columns
268*
269 a( pvt, pvt ) = a( j, j )
270 CALL sswap( j-1, a( 1, j ), 1, a( 1, pvt ), 1 )
271 IF( pvt.LT.n )
272 $ CALL sswap( n-pvt, a( j, pvt+1 ), lda,
273 $ a( pvt, pvt+1 ), lda )
274 CALL sswap( pvt-j-1, a( j, j+1 ), lda, a( j+1, pvt ),
275 $ 1 )
276*
277* Swap dot products and PIV
278*
279 stemp = work( j )
280 work( j ) = work( pvt )
281 work( pvt ) = stemp
282 itemp = piv( pvt )
283 piv( pvt ) = piv( j )
284 piv( j ) = itemp
285 END IF
286*
287 ajj = sqrt( ajj )
288 a( j, j ) = ajj
289*
290* Compute elements J+1:N of row J
291*
292 IF( j.LT.n ) THEN
293 CALL sgemv( 'Trans', j-1, n-j, -one, a( 1, j+1 ), lda,
294 $ a( 1, j ), 1, one, a( j, j+1 ), lda )
295 CALL sscal( n-j, one / ajj, a( j, j+1 ), lda )
296 END IF
297*
298 130 CONTINUE
299*
300 ELSE
301*
302* Compute the Cholesky factorization P**T * A * P = L * L**T
303*
304 DO 150 j = 1, n
305*
306* Find pivot, test for exit, else swap rows and columns
307* Update dot products, compute possible pivots which are
308* stored in the second half of WORK
309*
310 DO 140 i = j, n
311*
312 IF( j.GT.1 ) THEN
313 work( i ) = work( i ) + a( i, j-1 )**2
314 END IF
315 work( n+i ) = a( i, i ) - work( i )
316*
317 140 CONTINUE
318*
319 IF( j.GT.1 ) THEN
320 itemp = maxloc( work( (n+j):(2*n) ), 1 )
321 pvt = itemp + j - 1
322 ajj = work( n+pvt )
323 IF( ajj.LE.sstop.OR.sisnan( ajj ) ) THEN
324 a( j, j ) = ajj
325 GO TO 160
326 END IF
327 END IF
328*
329 IF( j.NE.pvt ) THEN
330*
331* Pivot OK, so can now swap pivot rows and columns
332*
333 a( pvt, pvt ) = a( j, j )
334 CALL sswap( j-1, a( j, 1 ), lda, a( pvt, 1 ), lda )
335 IF( pvt.LT.n )
336 $ CALL sswap( n-pvt, a( pvt+1, j ), 1, a( pvt+1,
337 $ pvt ),
338 $ 1 )
339 CALL sswap( pvt-j-1, a( j+1, j ), 1, a( pvt, j+1 ),
340 $ lda )
341*
342* Swap dot products and PIV
343*
344 stemp = work( j )
345 work( j ) = work( pvt )
346 work( pvt ) = stemp
347 itemp = piv( pvt )
348 piv( pvt ) = piv( j )
349 piv( j ) = itemp
350 END IF
351*
352 ajj = sqrt( ajj )
353 a( j, j ) = ajj
354*
355* Compute elements J+1:N of column J
356*
357 IF( j.LT.n ) THEN
358 CALL sgemv( 'No Trans', n-j, j-1, -one, a( j+1, 1 ),
359 $ 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*
386 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine sgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
SGEMV
Definition sgemv.f:158
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:140
subroutine sscal(n, sa, sx, incx)
SSCAL
Definition sscal.f:79
subroutine sswap(n, sx, incx, sy, incy)
SSWAP
Definition sswap.f:82