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