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

◆ dpstrf()

subroutine dpstrf ( character uplo,
integer n,
double precision, dimension( lda, * ) a,
integer lda,
integer, dimension( n ) piv,
integer rank,
double precision tol,
double precision, dimension( 2*n ) work,
integer info )

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

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

Purpose:
!>
!> DPSTRF 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 3 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 DOUBLE PRECISION 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.
!> 
[in]LDA
!>          LDA is INTEGER
!>          The leading dimension of the array A.  LDA >= max(1,N).
!> 
[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.
!> 
[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.

Definition at line 139 of file dpstrf.f.

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 DOUBLE PRECISION TOL
148 INTEGER INFO, LDA, N, RANK
149 CHARACTER UPLO
150* ..
151* .. Array Arguments ..
152 DOUBLE PRECISION A( LDA, * ), WORK( 2*N )
153 INTEGER PIV( N )
154* ..
155*
156* =====================================================================
157*
158* .. Parameters ..
159 DOUBLE PRECISION ONE, ZERO
160 parameter( one = 1.0d+0, zero = 0.0d+0 )
161* ..
162* .. Local Scalars ..
163 DOUBLE PRECISION AJJ, DSTOP, DTEMP
164 INTEGER I, ITEMP, J, JB, K, NB, PVT
165 LOGICAL UPPER
166* ..
167* .. External Functions ..
168 DOUBLE PRECISION DLAMCH
169 INTEGER ILAENV
170 LOGICAL LSAME, DISNAN
171 EXTERNAL dlamch, ilaenv, lsame, disnan
172* ..
173* .. External Subroutines ..
174 EXTERNAL dgemv, dpstf2, dscal, dswap, dsyrk,
175 $ xerbla
176* ..
177* .. Intrinsic Functions ..
178 INTRINSIC max, min, sqrt, maxloc
179* ..
180* .. Executable Statements ..
181*
182* Test the input parameters.
183*
184 info = 0
185 upper = lsame( uplo, 'U' )
186 IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
187 info = -1
188 ELSE IF( n.LT.0 ) THEN
189 info = -2
190 ELSE IF( lda.LT.max( 1, n ) ) THEN
191 info = -4
192 END IF
193 IF( info.NE.0 ) THEN
194 CALL xerbla( 'DPSTRF', -info )
195 RETURN
196 END IF
197*
198* Quick return if possible
199*
200 IF( n.EQ.0 )
201 $ RETURN
202*
203* Get block size
204*
205 nb = ilaenv( 1, 'DPOTRF', uplo, n, -1, -1, -1 )
206 IF( nb.LE.1 .OR. nb.GE.n ) THEN
207*
208* Use unblocked code
209*
210 CALL dpstf2( uplo, n, a( 1, 1 ), lda, piv, rank, tol, work,
211 $ info )
212 GO TO 200
213*
214 ELSE
215*
216* Initialize PIV
217*
218 DO 100 i = 1, n
219 piv( i ) = i
220 100 CONTINUE
221*
222* Compute stopping value
223*
224 pvt = 1
225 ajj = a( pvt, pvt )
226 DO i = 2, n
227 IF( a( i, i ).GT.ajj ) THEN
228 pvt = i
229 ajj = a( pvt, pvt )
230 END IF
231 END DO
232 IF( ajj.LE.zero.OR.disnan( ajj ) ) THEN
233 rank = 0
234 info = 1
235 GO TO 200
236 END IF
237*
238* Compute stopping value if not supplied
239*
240 IF( tol.LT.zero ) THEN
241 dstop = n * dlamch( 'Epsilon' ) * ajj
242 ELSE
243 dstop = tol
244 END IF
245*
246*
247 IF( upper ) THEN
248*
249* Compute the Cholesky factorization P**T * A * P = U**T * U
250*
251 DO 140 k = 1, n, nb
252*
253* Account for last block not being NB wide
254*
255 jb = min( nb, n-k+1 )
256*
257* Set relevant part of first half of WORK to zero,
258* holds dot products
259*
260 DO 110 i = k, n
261 work( i ) = 0
262 110 CONTINUE
263*
264 DO 130 j = k, k + jb - 1
265*
266* Find pivot, test for exit, else swap rows and columns
267* Update dot products, compute possible pivots which are
268* stored in the second half of WORK
269*
270 DO 120 i = j, n
271*
272 IF( j.GT.k ) THEN
273 work( i ) = work( i ) + a( j-1, i )**2
274 END IF
275 work( n+i ) = a( i, i ) - work( i )
276*
277 120 CONTINUE
278*
279 IF( j.GT.1 ) THEN
280 itemp = maxloc( work( (n+j):(2*n) ), 1 )
281 pvt = itemp + j - 1
282 ajj = work( n+pvt )
283 IF( ajj.LE.dstop.OR.disnan( ajj ) ) THEN
284 a( j, j ) = ajj
285 GO TO 190
286 END IF
287 END IF
288*
289 IF( j.NE.pvt ) THEN
290*
291* Pivot OK, so can now swap pivot rows and columns
292*
293 a( pvt, pvt ) = a( j, j )
294 CALL dswap( j-1, a( 1, j ), 1, a( 1, pvt ), 1 )
295 IF( pvt.LT.n )
296 $ CALL dswap( n-pvt, a( j, pvt+1 ), lda,
297 $ a( pvt, pvt+1 ), lda )
298 CALL dswap( pvt-j-1, a( j, j+1 ), lda,
299 $ a( j+1, pvt ), 1 )
300*
301* Swap dot products and PIV
302*
303 dtemp = work( j )
304 work( j ) = work( pvt )
305 work( pvt ) = dtemp
306 itemp = piv( pvt )
307 piv( pvt ) = piv( j )
308 piv( j ) = itemp
309 END IF
310*
311 ajj = sqrt( ajj )
312 a( j, j ) = ajj
313*
314* Compute elements J+1:N of row J.
315*
316 IF( j.LT.n ) THEN
317 CALL dgemv( 'Trans', j-k, n-j, -one, a( k,
318 $ j+1 ),
319 $ lda, a( k, j ), 1, one, a( j, j+1 ),
320 $ lda )
321 CALL dscal( n-j, one / ajj, a( j, j+1 ), lda )
322 END IF
323*
324 130 CONTINUE
325*
326* Update trailing matrix, J already incremented
327*
328 IF( k+jb.LE.n ) THEN
329 CALL dsyrk( 'Upper', 'Trans', n-j+1, jb, -one,
330 $ a( k, j ), lda, one, a( j, j ), lda )
331 END IF
332*
333 140 CONTINUE
334*
335 ELSE
336*
337* Compute the Cholesky factorization P**T * A * P = L * L**T
338*
339 DO 180 k = 1, n, nb
340*
341* Account for last block not being NB wide
342*
343 jb = min( nb, n-k+1 )
344*
345* Set relevant part of first half of WORK to zero,
346* holds dot products
347*
348 DO 150 i = k, n
349 work( i ) = 0
350 150 CONTINUE
351*
352 DO 170 j = k, k + jb - 1
353*
354* Find pivot, test for exit, else swap rows and columns
355* Update dot products, compute possible pivots which are
356* stored in the second half of WORK
357*
358 DO 160 i = j, n
359*
360 IF( j.GT.k ) THEN
361 work( i ) = work( i ) + a( i, j-1 )**2
362 END IF
363 work( n+i ) = a( i, i ) - work( i )
364*
365 160 CONTINUE
366*
367 IF( j.GT.1 ) THEN
368 itemp = maxloc( work( (n+j):(2*n) ), 1 )
369 pvt = itemp + j - 1
370 ajj = work( n+pvt )
371 IF( ajj.LE.dstop.OR.disnan( ajj ) ) THEN
372 a( j, j ) = ajj
373 GO TO 190
374 END IF
375 END IF
376*
377 IF( j.NE.pvt ) THEN
378*
379* Pivot OK, so can now swap pivot rows and columns
380*
381 a( pvt, pvt ) = a( j, j )
382 CALL dswap( j-1, a( j, 1 ), lda, a( pvt, 1 ),
383 $ lda )
384 IF( pvt.LT.n )
385 $ CALL dswap( n-pvt, a( pvt+1, j ), 1,
386 $ a( pvt+1, pvt ), 1 )
387 CALL dswap( pvt-j-1, a( j+1, j ), 1, a( pvt,
388 $ j+1 ),
389 $ lda )
390*
391* Swap dot products and PIV
392*
393 dtemp = work( j )
394 work( j ) = work( pvt )
395 work( pvt ) = dtemp
396 itemp = piv( pvt )
397 piv( pvt ) = piv( j )
398 piv( j ) = itemp
399 END IF
400*
401 ajj = sqrt( ajj )
402 a( j, j ) = ajj
403*
404* Compute elements J+1:N of column J.
405*
406 IF( j.LT.n ) THEN
407 CALL dgemv( 'No Trans', n-j, j-k, -one,
408 $ a( j+1, k ), lda, a( j, k ), lda, one,
409 $ a( j+1, j ), 1 )
410 CALL dscal( n-j, one / ajj, a( j+1, j ), 1 )
411 END IF
412*
413 170 CONTINUE
414*
415* Update trailing matrix, J already incremented
416*
417 IF( k+jb.LE.n ) THEN
418 CALL dsyrk( 'Lower', 'No Trans', n-j+1, jb, -one,
419 $ a( j, k ), lda, one, a( j, j ), lda )
420 END IF
421*
422 180 CONTINUE
423*
424 END IF
425 END IF
426*
427* Ran to completion, A has full rank
428*
429 rank = n
430*
431 GO TO 200
432 190 CONTINUE
433*
434* Rank is the number of steps completed. Set INFO = 1 to signal
435* that the factorization cannot be used to solve a system.
436*
437 rank = j - 1
438 info = 1
439*
440 200 CONTINUE
441 RETURN
442*
443* End of DPSTRF
444*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)
DGEMV
Definition dgemv.f:158
subroutine dsyrk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
DSYRK
Definition dsyrk.f:169
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:160
logical function disnan(din)
DISNAN tests input for NaN.
Definition disnan.f:57
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
subroutine dpstf2(uplo, n, a, lda, piv, rank, tol, work, info)
DPSTF2 computes the Cholesky factorization with complete pivoting of a real symmetric positive semide...
Definition dpstf2.f:140
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79
subroutine dswap(n, dx, incx, dy, incy)
DSWAP
Definition dswap.f:82
Here is the call graph for this function:
Here is the caller graph for this function: