LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ cpstf2()

subroutine cpstf2 ( character  uplo,
integer  n,
complex, dimension( lda, * )  a,
integer  lda,
integer, dimension( n )  piv,
integer  rank,
real  tol,
real, dimension( 2*n )  work,
integer  info 
)

CPSTF2 computes the Cholesky factorization with complete pivoting of complex Hermitian positive semidefinite matrix.

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

Purpose:
 CPSTF2 computes the Cholesky factorization with complete
 pivoting of a complex Hermitian positive semidefinite matrix A.

 The factorization has the form
    P**T * A * P = U**H * U ,  if UPLO = 'U',
    P**T * A * P = L  * L**H,  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 COMPLEX 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.

Definition at line 141 of file cpstf2.f.

142*
143* -- LAPACK computational routine --
144* -- LAPACK is a software package provided by Univ. of Tennessee, --
145* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
146*
147* .. Scalar Arguments ..
148 REAL TOL
149 INTEGER INFO, LDA, N, RANK
150 CHARACTER UPLO
151* ..
152* .. Array Arguments ..
153 COMPLEX A( LDA, * )
154 REAL WORK( 2*N )
155 INTEGER PIV( N )
156* ..
157*
158* =====================================================================
159*
160* .. Parameters ..
161 REAL ONE, ZERO
162 parameter( one = 1.0e+0, zero = 0.0e+0 )
163 COMPLEX CONE
164 parameter( cone = ( 1.0e+0, 0.0e+0 ) )
165* ..
166* .. Local Scalars ..
167 COMPLEX CTEMP
168 REAL AJJ, SSTOP, STEMP
169 INTEGER I, ITEMP, J, PVT
170 LOGICAL UPPER
171* ..
172* .. External Functions ..
173 REAL SLAMCH
174 LOGICAL LSAME, SISNAN
175 EXTERNAL slamch, lsame, sisnan
176* ..
177* .. External Subroutines ..
178 EXTERNAL cgemv, clacgv, csscal, cswap, 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 = 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 ), lda,
304 $ a( 1, j ), 1, cone, a( j, j+1 ), lda )
305 CALL clacgv( j-1, a( 1, j ), 1 )
306 CALL csscal( n-j, one / ajj, a( j, j+1 ), lda )
307 END IF
308*
309 150 CONTINUE
310*
311 ELSE
312*
313* Compute the Cholesky factorization P**T * A * P = L * L**H
314*
315 DO 180 j = 1, n
316*
317* Find pivot, test for exit, else swap rows and columns
318* Update dot products, compute possible pivots which are
319* stored in the second half of WORK
320*
321 DO 160 i = j, n
322*
323 IF( j.GT.1 ) THEN
324 work( i ) = work( i ) +
325 $ real( conjg( a( i, j-1 ) )*
326 $ a( i, j-1 ) )
327 END IF
328 work( n+i ) = real( a( i, i ) ) - work( i )
329*
330 160 CONTINUE
331*
332 IF( j.GT.1 ) THEN
333 itemp = maxloc( work( (n+j):(2*n) ), 1 )
334 pvt = itemp + j - 1
335 ajj = work( n+pvt )
336 IF( ajj.LE.sstop.OR.sisnan( ajj ) ) THEN
337 a( j, j ) = ajj
338 GO TO 190
339 END IF
340 END IF
341*
342 IF( j.NE.pvt ) THEN
343*
344* Pivot OK, so can now swap pivot rows and columns
345*
346 a( pvt, pvt ) = a( j, j )
347 CALL cswap( j-1, a( j, 1 ), lda, a( pvt, 1 ), lda )
348 IF( pvt.LT.n )
349 $ CALL cswap( n-pvt, a( pvt+1, j ), 1, a( pvt+1, pvt ),
350 $ 1 )
351 DO 170 i = j + 1, pvt - 1
352 ctemp = conjg( a( i, j ) )
353 a( i, j ) = conjg( a( pvt, i ) )
354 a( pvt, i ) = ctemp
355 170 CONTINUE
356 a( pvt, j ) = conjg( a( pvt, j ) )
357*
358* Swap dot products and PIV
359*
360 stemp = work( j )
361 work( j ) = work( pvt )
362 work( pvt ) = stemp
363 itemp = piv( pvt )
364 piv( pvt ) = piv( j )
365 piv( j ) = itemp
366 END IF
367*
368 ajj = sqrt( ajj )
369 a( j, j ) = ajj
370*
371* Compute elements J+1:N of column J
372*
373 IF( j.LT.n ) THEN
374 CALL clacgv( j-1, a( j, 1 ), lda )
375 CALL cgemv( 'No Trans', n-j, j-1, -cone, a( j+1, 1 ),
376 $ lda, a( j, 1 ), lda, cone, a( j+1, j ), 1 )
377 CALL clacgv( j-1, a( j, 1 ), lda )
378 CALL csscal( n-j, one / ajj, a( j+1, j ), 1 )
379 END IF
380*
381 180 CONTINUE
382*
383 END IF
384*
385* Ran to completion, A has full rank
386*
387 rank = n
388*
389 GO TO 200
390 190 CONTINUE
391*
392* Rank is number of steps completed. Set INFO = 1 to signal
393* that the factorization cannot be used to solve a system.
394*
395 rank = j - 1
396 info = 1
397*
398 200 CONTINUE
399 RETURN
400*
401* End of CPSTF2
402*
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
logical function sisnan(sin)
SISNAN tests input for NaN.
Definition sisnan.f:59
subroutine clacgv(n, x, incx)
CLACGV conjugates a complex vector.
Definition clacgv.f:74
real function slamch(cmach)
SLAMCH
Definition slamch.f:68
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine csscal(n, sa, cx, incx)
CSSCAL
Definition csscal.f:78
subroutine cswap(n, cx, incx, cy, incy)
CSWAP
Definition cswap.f:81
Here is the call graph for this function:
Here is the caller graph for this function: