LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine zpstf2 ( character  UPLO,
integer  N,
complex*16, dimension( lda, * )  A,
integer  LDA,
integer, dimension( n )  PIV,
integer  RANK,
double precision  TOL,
double precision, dimension( 2*n )  WORK,
integer  INFO 
)

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

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

Purpose:
 ZPSTF2 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*16 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 DOUBLE PRECISION
          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 DOUBLE PRECISION 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 144 of file zpstf2.f.

144 *
145 * -- LAPACK computational routine (version 3.6.0) --
146 * -- LAPACK is a software package provided by Univ. of Tennessee, --
147 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
148 * November 2015
149 *
150 * .. Scalar Arguments ..
151  DOUBLE PRECISION tol
152  INTEGER info, lda, n, rank
153  CHARACTER uplo
154 * ..
155 * .. Array Arguments ..
156  COMPLEX*16 a( lda, * )
157  DOUBLE PRECISION work( 2*n )
158  INTEGER piv( n )
159 * ..
160 *
161 * =====================================================================
162 *
163 * .. Parameters ..
164  DOUBLE PRECISION one, zero
165  parameter ( one = 1.0d+0, zero = 0.0d+0 )
166  COMPLEX*16 cone
167  parameter ( cone = ( 1.0d+0, 0.0d+0 ) )
168 * ..
169 * .. Local Scalars ..
170  COMPLEX*16 ztemp
171  DOUBLE PRECISION ajj, dstop, dtemp
172  INTEGER i, itemp, j, pvt
173  LOGICAL upper
174 * ..
175 * .. External Functions ..
176  DOUBLE PRECISION dlamch
177  LOGICAL lsame, disnan
178  EXTERNAL dlamch, lsame, disnan
179 * ..
180 * .. External Subroutines ..
181  EXTERNAL zdscal, zgemv, zlacgv, zswap, xerbla
182 * ..
183 * .. Intrinsic Functions ..
184  INTRINSIC dble, dconjg, max, sqrt
185 * ..
186 * .. Executable Statements ..
187 *
188 * Test the input parameters
189 *
190  info = 0
191  upper = lsame( uplo, 'U' )
192  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
193  info = -1
194  ELSE IF( n.LT.0 ) THEN
195  info = -2
196  ELSE IF( lda.LT.max( 1, n ) ) THEN
197  info = -4
198  END IF
199  IF( info.NE.0 ) THEN
200  CALL xerbla( 'ZPSTF2', -info )
201  RETURN
202  END IF
203 *
204 * Quick return if possible
205 *
206  IF( n.EQ.0 )
207  $ RETURN
208 *
209 * Initialize PIV
210 *
211  DO 100 i = 1, n
212  piv( i ) = i
213  100 CONTINUE
214 *
215 * Compute stopping value
216 *
217  DO 110 i = 1, n
218  work( i ) = dble( a( i, i ) )
219  110 CONTINUE
220  pvt = maxloc( work( 1:n ), 1 )
221  ajj = dble( a( pvt, pvt ) )
222  IF( ajj.LE.zero.OR.disnan( ajj ) ) THEN
223  rank = 0
224  info = 1
225  GO TO 200
226  END IF
227 *
228 * Compute stopping value if not supplied
229 *
230  IF( tol.LT.zero ) THEN
231  dstop = n * dlamch( 'Epsilon' ) * ajj
232  ELSE
233  dstop = tol
234  END IF
235 *
236 * Set first half of WORK to zero, holds dot products
237 *
238  DO 120 i = 1, n
239  work( i ) = 0
240  120 CONTINUE
241 *
242  IF( upper ) THEN
243 *
244 * Compute the Cholesky factorization P**T * A * P = U**H* U
245 *
246  DO 150 j = 1, n
247 *
248 * Find pivot, test for exit, else swap rows and columns
249 * Update dot products, compute possible pivots which are
250 * stored in the second half of WORK
251 *
252  DO 130 i = j, n
253 *
254  IF( j.GT.1 ) THEN
255  work( i ) = work( i ) +
256  $ dble( dconjg( a( j-1, i ) )*
257  $ a( j-1, i ) )
258  END IF
259  work( n+i ) = dble( a( i, i ) ) - work( i )
260 *
261  130 CONTINUE
262 *
263  IF( j.GT.1 ) THEN
264  itemp = maxloc( work( (n+j):(2*n) ), 1 )
265  pvt = itemp + j - 1
266  ajj = work( n+pvt )
267  IF( ajj.LE.dstop.OR.disnan( ajj ) ) THEN
268  a( j, j ) = ajj
269  GO TO 190
270  END IF
271  END IF
272 *
273  IF( j.NE.pvt ) THEN
274 *
275 * Pivot OK, so can now swap pivot rows and columns
276 *
277  a( pvt, pvt ) = a( j, j )
278  CALL zswap( j-1, a( 1, j ), 1, a( 1, pvt ), 1 )
279  IF( pvt.LT.n )
280  $ CALL zswap( n-pvt, a( j, pvt+1 ), lda,
281  $ a( pvt, pvt+1 ), lda )
282  DO 140 i = j + 1, pvt - 1
283  ztemp = dconjg( a( j, i ) )
284  a( j, i ) = dconjg( a( i, pvt ) )
285  a( i, pvt ) = ztemp
286  140 CONTINUE
287  a( j, pvt ) = dconjg( a( j, pvt ) )
288 *
289 * Swap dot products and PIV
290 *
291  dtemp = work( j )
292  work( j ) = work( pvt )
293  work( pvt ) = dtemp
294  itemp = piv( pvt )
295  piv( pvt ) = piv( j )
296  piv( j ) = itemp
297  END IF
298 *
299  ajj = sqrt( ajj )
300  a( j, j ) = ajj
301 *
302 * Compute elements J+1:N of row J
303 *
304  IF( j.LT.n ) THEN
305  CALL zlacgv( j-1, a( 1, j ), 1 )
306  CALL zgemv( 'Trans', j-1, n-j, -cone, a( 1, j+1 ), lda,
307  $ a( 1, j ), 1, cone, a( j, j+1 ), lda )
308  CALL zlacgv( j-1, a( 1, j ), 1 )
309  CALL zdscal( n-j, one / ajj, a( j, j+1 ), lda )
310  END IF
311 *
312  150 CONTINUE
313 *
314  ELSE
315 *
316 * Compute the Cholesky factorization P**T * A * P = L * L**H
317 *
318  DO 180 j = 1, n
319 *
320 * Find pivot, test for exit, else swap rows and columns
321 * Update dot products, compute possible pivots which are
322 * stored in the second half of WORK
323 *
324  DO 160 i = j, n
325 *
326  IF( j.GT.1 ) THEN
327  work( i ) = work( i ) +
328  $ dble( dconjg( a( i, j-1 ) )*
329  $ a( i, j-1 ) )
330  END IF
331  work( n+i ) = dble( a( i, i ) ) - work( i )
332 *
333  160 CONTINUE
334 *
335  IF( j.GT.1 ) THEN
336  itemp = maxloc( work( (n+j):(2*n) ), 1 )
337  pvt = itemp + j - 1
338  ajj = work( n+pvt )
339  IF( ajj.LE.dstop.OR.disnan( ajj ) ) THEN
340  a( j, j ) = ajj
341  GO TO 190
342  END IF
343  END IF
344 *
345  IF( j.NE.pvt ) THEN
346 *
347 * Pivot OK, so can now swap pivot rows and columns
348 *
349  a( pvt, pvt ) = a( j, j )
350  CALL zswap( j-1, a( j, 1 ), lda, a( pvt, 1 ), lda )
351  IF( pvt.LT.n )
352  $ CALL zswap( n-pvt, a( pvt+1, j ), 1, a( pvt+1, pvt ),
353  $ 1 )
354  DO 170 i = j + 1, pvt - 1
355  ztemp = dconjg( a( i, j ) )
356  a( i, j ) = dconjg( a( pvt, i ) )
357  a( pvt, i ) = ztemp
358  170 CONTINUE
359  a( pvt, j ) = dconjg( a( pvt, j ) )
360 *
361 * Swap dot products and PIV
362 *
363  dtemp = work( j )
364  work( j ) = work( pvt )
365  work( pvt ) = dtemp
366  itemp = piv( pvt )
367  piv( pvt ) = piv( j )
368  piv( j ) = itemp
369  END IF
370 *
371  ajj = sqrt( ajj )
372  a( j, j ) = ajj
373 *
374 * Compute elements J+1:N of column J
375 *
376  IF( j.LT.n ) THEN
377  CALL zlacgv( j-1, a( j, 1 ), lda )
378  CALL zgemv( 'No Trans', n-j, j-1, -cone, a( j+1, 1 ),
379  $ lda, a( j, 1 ), lda, cone, a( j+1, j ), 1 )
380  CALL zlacgv( j-1, a( j, 1 ), lda )
381  CALL zdscal( n-j, one / ajj, a( j+1, j ), 1 )
382  END IF
383 *
384  180 CONTINUE
385 *
386  END IF
387 *
388 * Ran to completion, A has full rank
389 *
390  rank = n
391 *
392  GO TO 200
393  190 CONTINUE
394 *
395 * Rank is number of steps completed. Set INFO = 1 to signal
396 * that the factorization cannot be used to solve a system.
397 *
398  rank = j - 1
399  info = 1
400 *
401  200 CONTINUE
402  RETURN
403 *
404 * End of ZPSTF2
405 *
logical function disnan(DIN)
DISNAN tests input for NaN.
Definition: disnan.f:61
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
subroutine zgemv(TRANS, M, N, ALPHA, A, LDA, X, INCX, BETA, Y, INCY)
ZGEMV
Definition: zgemv.f:160
subroutine zswap(N, ZX, INCX, ZY, INCY)
ZSWAP
Definition: zswap.f:52
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:54
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zlacgv(N, X, INCX)
ZLACGV conjugates a complex vector.
Definition: zlacgv.f:76

Here is the call graph for this function:

Here is the caller graph for this function: