LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine dsptrf ( character  UPLO,
integer  N,
double precision, dimension( * )  AP,
integer, dimension( * )  IPIV,
integer  INFO 
)

DSPTRF

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

Purpose:
 DSPTRF computes the factorization of a real symmetric matrix A stored
 in packed format using the Bunch-Kaufman diagonal pivoting method:

    A = U*D*U**T  or  A = L*D*L**T

 where U (or L) is a product of permutation and unit upper (lower)
 triangular matrices, and D is symmetric and block diagonal with
 1-by-1 and 2-by-2 diagonal blocks.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          = 'U':  Upper triangle of A is stored;
          = 'L':  Lower triangle of A is stored.
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]AP
          AP is DOUBLE PRECISION array, dimension (N*(N+1)/2)
          On entry, the upper or lower triangle of the symmetric matrix
          A, packed columnwise in a linear array.  The j-th column of A
          is stored in the array AP as follows:
          if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
          if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.

          On exit, the block diagonal matrix D and the multipliers used
          to obtain the factor U or L, stored as a packed triangular
          matrix overwriting A (see below for further details).
[out]IPIV
          IPIV is INTEGER array, dimension (N)
          Details of the interchanges and the block structure of D.
          If IPIV(k) > 0, then rows and columns k and IPIV(k) were
          interchanged and D(k,k) is a 1-by-1 diagonal block.
          If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0, then rows and
          columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k)
          is a 2-by-2 diagonal block.  If UPLO = 'L' and IPIV(k) =
          IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were
          interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -i, the i-th argument had an illegal value
          > 0: if INFO = i, D(i,i) is exactly zero.  The factorization
               has been completed, but the block diagonal matrix D is
               exactly singular, and division by zero will occur if it
               is used to solve a system of equations.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011
Further Details:
  If UPLO = 'U', then A = U*D*U**T, where
     U = P(n)*U(n)* ... *P(k)U(k)* ...,
  i.e., U is a product of terms P(k)*U(k), where k decreases from n to
  1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1
  and 2-by-2 diagonal blocks D(k).  P(k) is a permutation matrix as
  defined by IPIV(k), and U(k) is a unit upper triangular matrix, such
  that if the diagonal block D(k) is of order s (s = 1 or 2), then

             (   I    v    0   )   k-s
     U(k) =  (   0    I    0   )   s
             (   0    0    I   )   n-k
                k-s   s   n-k

  If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k).
  If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k),
  and A(k,k), and v overwrites A(1:k-2,k-1:k).

  If UPLO = 'L', then A = L*D*L**T, where
     L = P(1)*L(1)* ... *P(k)*L(k)* ...,
  i.e., L is a product of terms P(k)*L(k), where k increases from 1 to
  n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1
  and 2-by-2 diagonal blocks D(k).  P(k) is a permutation matrix as
  defined by IPIV(k), and L(k) is a unit lower triangular matrix, such
  that if the diagonal block D(k) is of order s (s = 1 or 2), then

             (   I    0     0   )  k-1
     L(k) =  (   0    I     0   )  s
             (   0    v     I   )  n-k-s+1
                k-1   s  n-k-s+1

  If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k).
  If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k),
  and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
Contributors:
J. Lewis, Boeing Computer Services Company

Definition at line 161 of file dsptrf.f.

161 *
162 * -- LAPACK computational routine (version 3.4.0) --
163 * -- LAPACK is a software package provided by Univ. of Tennessee, --
164 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
165 * November 2011
166 *
167 * .. Scalar Arguments ..
168  CHARACTER uplo
169  INTEGER info, n
170 * ..
171 * .. Array Arguments ..
172  INTEGER ipiv( * )
173  DOUBLE PRECISION ap( * )
174 * ..
175 *
176 * =====================================================================
177 *
178 * .. Parameters ..
179  DOUBLE PRECISION zero, one
180  parameter ( zero = 0.0d+0, one = 1.0d+0 )
181  DOUBLE PRECISION eight, sevten
182  parameter ( eight = 8.0d+0, sevten = 17.0d+0 )
183 * ..
184 * .. Local Scalars ..
185  LOGICAL upper
186  INTEGER i, imax, j, jmax, k, kc, kk, knc, kp, kpc,
187  $ kstep, kx, npp
188  DOUBLE PRECISION absakk, alpha, colmax, d11, d12, d21, d22, r1,
189  $ rowmax, t, wk, wkm1, wkp1
190 * ..
191 * .. External Functions ..
192  LOGICAL lsame
193  INTEGER idamax
194  EXTERNAL lsame, idamax
195 * ..
196 * .. External Subroutines ..
197  EXTERNAL dscal, dspr, dswap, xerbla
198 * ..
199 * .. Intrinsic Functions ..
200  INTRINSIC abs, max, sqrt
201 * ..
202 * .. Executable Statements ..
203 *
204 * Test the input parameters.
205 *
206  info = 0
207  upper = lsame( uplo, 'U' )
208  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
209  info = -1
210  ELSE IF( n.LT.0 ) THEN
211  info = -2
212  END IF
213  IF( info.NE.0 ) THEN
214  CALL xerbla( 'DSPTRF', -info )
215  RETURN
216  END IF
217 *
218 * Initialize ALPHA for use in choosing pivot block size.
219 *
220  alpha = ( one+sqrt( sevten ) ) / eight
221 *
222  IF( upper ) THEN
223 *
224 * Factorize A as U*D*U**T using the upper triangle of A
225 *
226 * K is the main loop index, decreasing from N to 1 in steps of
227 * 1 or 2
228 *
229  k = n
230  kc = ( n-1 )*n / 2 + 1
231  10 CONTINUE
232  knc = kc
233 *
234 * If K < 1, exit from loop
235 *
236  IF( k.LT.1 )
237  $ GO TO 110
238  kstep = 1
239 *
240 * Determine rows and columns to be interchanged and whether
241 * a 1-by-1 or 2-by-2 pivot block will be used
242 *
243  absakk = abs( ap( kc+k-1 ) )
244 *
245 * IMAX is the row-index of the largest off-diagonal element in
246 * column K, and COLMAX is its absolute value
247 *
248  IF( k.GT.1 ) THEN
249  imax = idamax( k-1, ap( kc ), 1 )
250  colmax = abs( ap( kc+imax-1 ) )
251  ELSE
252  colmax = zero
253  END IF
254 *
255  IF( max( absakk, colmax ).EQ.zero ) THEN
256 *
257 * Column K is zero: set INFO and continue
258 *
259  IF( info.EQ.0 )
260  $ info = k
261  kp = k
262  ELSE
263  IF( absakk.GE.alpha*colmax ) THEN
264 *
265 * no interchange, use 1-by-1 pivot block
266 *
267  kp = k
268  ELSE
269 *
270  rowmax = zero
271  jmax = imax
272  kx = imax*( imax+1 ) / 2 + imax
273  DO 20 j = imax + 1, k
274  IF( abs( ap( kx ) ).GT.rowmax ) THEN
275  rowmax = abs( ap( kx ) )
276  jmax = j
277  END IF
278  kx = kx + j
279  20 CONTINUE
280  kpc = ( imax-1 )*imax / 2 + 1
281  IF( imax.GT.1 ) THEN
282  jmax = idamax( imax-1, ap( kpc ), 1 )
283  rowmax = max( rowmax, abs( ap( kpc+jmax-1 ) ) )
284  END IF
285 *
286  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
287 *
288 * no interchange, use 1-by-1 pivot block
289 *
290  kp = k
291  ELSE IF( abs( ap( kpc+imax-1 ) ).GE.alpha*rowmax ) THEN
292 *
293 * interchange rows and columns K and IMAX, use 1-by-1
294 * pivot block
295 *
296  kp = imax
297  ELSE
298 *
299 * interchange rows and columns K-1 and IMAX, use 2-by-2
300 * pivot block
301 *
302  kp = imax
303  kstep = 2
304  END IF
305  END IF
306 *
307  kk = k - kstep + 1
308  IF( kstep.EQ.2 )
309  $ knc = knc - k + 1
310  IF( kp.NE.kk ) THEN
311 *
312 * Interchange rows and columns KK and KP in the leading
313 * submatrix A(1:k,1:k)
314 *
315  CALL dswap( kp-1, ap( knc ), 1, ap( kpc ), 1 )
316  kx = kpc + kp - 1
317  DO 30 j = kp + 1, kk - 1
318  kx = kx + j - 1
319  t = ap( knc+j-1 )
320  ap( knc+j-1 ) = ap( kx )
321  ap( kx ) = t
322  30 CONTINUE
323  t = ap( knc+kk-1 )
324  ap( knc+kk-1 ) = ap( kpc+kp-1 )
325  ap( kpc+kp-1 ) = t
326  IF( kstep.EQ.2 ) THEN
327  t = ap( kc+k-2 )
328  ap( kc+k-2 ) = ap( kc+kp-1 )
329  ap( kc+kp-1 ) = t
330  END IF
331  END IF
332 *
333 * Update the leading submatrix
334 *
335  IF( kstep.EQ.1 ) THEN
336 *
337 * 1-by-1 pivot block D(k): column k now holds
338 *
339 * W(k) = U(k)*D(k)
340 *
341 * where U(k) is the k-th column of U
342 *
343 * Perform a rank-1 update of A(1:k-1,1:k-1) as
344 *
345 * A := A - U(k)*D(k)*U(k)**T = A - W(k)*1/D(k)*W(k)**T
346 *
347  r1 = one / ap( kc+k-1 )
348  CALL dspr( uplo, k-1, -r1, ap( kc ), 1, ap )
349 *
350 * Store U(k) in column k
351 *
352  CALL dscal( k-1, r1, ap( kc ), 1 )
353  ELSE
354 *
355 * 2-by-2 pivot block D(k): columns k and k-1 now hold
356 *
357 * ( W(k-1) W(k) ) = ( U(k-1) U(k) )*D(k)
358 *
359 * where U(k) and U(k-1) are the k-th and (k-1)-th columns
360 * of U
361 *
362 * Perform a rank-2 update of A(1:k-2,1:k-2) as
363 *
364 * A := A - ( U(k-1) U(k) )*D(k)*( U(k-1) U(k) )**T
365 * = A - ( W(k-1) W(k) )*inv(D(k))*( W(k-1) W(k) )**T
366 *
367  IF( k.GT.2 ) THEN
368 *
369  d12 = ap( k-1+( k-1 )*k / 2 )
370  d22 = ap( k-1+( k-2 )*( k-1 ) / 2 ) / d12
371  d11 = ap( k+( k-1 )*k / 2 ) / d12
372  t = one / ( d11*d22-one )
373  d12 = t / d12
374 *
375  DO 50 j = k - 2, 1, -1
376  wkm1 = d12*( d11*ap( j+( k-2 )*( k-1 ) / 2 )-
377  $ ap( j+( k-1 )*k / 2 ) )
378  wk = d12*( d22*ap( j+( k-1 )*k / 2 )-
379  $ ap( j+( k-2 )*( k-1 ) / 2 ) )
380  DO 40 i = j, 1, -1
381  ap( i+( j-1 )*j / 2 ) = ap( i+( j-1 )*j / 2 ) -
382  $ ap( i+( k-1 )*k / 2 )*wk -
383  $ ap( i+( k-2 )*( k-1 ) / 2 )*wkm1
384  40 CONTINUE
385  ap( j+( k-1 )*k / 2 ) = wk
386  ap( j+( k-2 )*( k-1 ) / 2 ) = wkm1
387  50 CONTINUE
388 *
389  END IF
390 *
391  END IF
392  END IF
393 *
394 * Store details of the interchanges in IPIV
395 *
396  IF( kstep.EQ.1 ) THEN
397  ipiv( k ) = kp
398  ELSE
399  ipiv( k ) = -kp
400  ipiv( k-1 ) = -kp
401  END IF
402 *
403 * Decrease K and return to the start of the main loop
404 *
405  k = k - kstep
406  kc = knc - k
407  GO TO 10
408 *
409  ELSE
410 *
411 * Factorize A as L*D*L**T using the lower triangle of A
412 *
413 * K is the main loop index, increasing from 1 to N in steps of
414 * 1 or 2
415 *
416  k = 1
417  kc = 1
418  npp = n*( n+1 ) / 2
419  60 CONTINUE
420  knc = kc
421 *
422 * If K > N, exit from loop
423 *
424  IF( k.GT.n )
425  $ GO TO 110
426  kstep = 1
427 *
428 * Determine rows and columns to be interchanged and whether
429 * a 1-by-1 or 2-by-2 pivot block will be used
430 *
431  absakk = abs( ap( kc ) )
432 *
433 * IMAX is the row-index of the largest off-diagonal element in
434 * column K, and COLMAX is its absolute value
435 *
436  IF( k.LT.n ) THEN
437  imax = k + idamax( n-k, ap( kc+1 ), 1 )
438  colmax = abs( ap( kc+imax-k ) )
439  ELSE
440  colmax = zero
441  END IF
442 *
443  IF( max( absakk, colmax ).EQ.zero ) THEN
444 *
445 * Column K is zero: set INFO and continue
446 *
447  IF( info.EQ.0 )
448  $ info = k
449  kp = k
450  ELSE
451  IF( absakk.GE.alpha*colmax ) THEN
452 *
453 * no interchange, use 1-by-1 pivot block
454 *
455  kp = k
456  ELSE
457 *
458 * JMAX is the column-index of the largest off-diagonal
459 * element in row IMAX, and ROWMAX is its absolute value
460 *
461  rowmax = zero
462  kx = kc + imax - k
463  DO 70 j = k, imax - 1
464  IF( abs( ap( kx ) ).GT.rowmax ) THEN
465  rowmax = abs( ap( kx ) )
466  jmax = j
467  END IF
468  kx = kx + n - j
469  70 CONTINUE
470  kpc = npp - ( n-imax+1 )*( n-imax+2 ) / 2 + 1
471  IF( imax.LT.n ) THEN
472  jmax = imax + idamax( n-imax, ap( kpc+1 ), 1 )
473  rowmax = max( rowmax, abs( ap( kpc+jmax-imax ) ) )
474  END IF
475 *
476  IF( absakk.GE.alpha*colmax*( colmax / rowmax ) ) THEN
477 *
478 * no interchange, use 1-by-1 pivot block
479 *
480  kp = k
481  ELSE IF( abs( ap( kpc ) ).GE.alpha*rowmax ) THEN
482 *
483 * interchange rows and columns K and IMAX, use 1-by-1
484 * pivot block
485 *
486  kp = imax
487  ELSE
488 *
489 * interchange rows and columns K+1 and IMAX, use 2-by-2
490 * pivot block
491 *
492  kp = imax
493  kstep = 2
494  END IF
495  END IF
496 *
497  kk = k + kstep - 1
498  IF( kstep.EQ.2 )
499  $ knc = knc + n - k + 1
500  IF( kp.NE.kk ) THEN
501 *
502 * Interchange rows and columns KK and KP in the trailing
503 * submatrix A(k:n,k:n)
504 *
505  IF( kp.LT.n )
506  $ CALL dswap( n-kp, ap( knc+kp-kk+1 ), 1, ap( kpc+1 ),
507  $ 1 )
508  kx = knc + kp - kk
509  DO 80 j = kk + 1, kp - 1
510  kx = kx + n - j + 1
511  t = ap( knc+j-kk )
512  ap( knc+j-kk ) = ap( kx )
513  ap( kx ) = t
514  80 CONTINUE
515  t = ap( knc )
516  ap( knc ) = ap( kpc )
517  ap( kpc ) = t
518  IF( kstep.EQ.2 ) THEN
519  t = ap( kc+1 )
520  ap( kc+1 ) = ap( kc+kp-k )
521  ap( kc+kp-k ) = t
522  END IF
523  END IF
524 *
525 * Update the trailing submatrix
526 *
527  IF( kstep.EQ.1 ) THEN
528 *
529 * 1-by-1 pivot block D(k): column k now holds
530 *
531 * W(k) = L(k)*D(k)
532 *
533 * where L(k) is the k-th column of L
534 *
535  IF( k.LT.n ) THEN
536 *
537 * Perform a rank-1 update of A(k+1:n,k+1:n) as
538 *
539 * A := A - L(k)*D(k)*L(k)**T = A - W(k)*(1/D(k))*W(k)**T
540 *
541  r1 = one / ap( kc )
542  CALL dspr( uplo, n-k, -r1, ap( kc+1 ), 1,
543  $ ap( kc+n-k+1 ) )
544 *
545 * Store L(k) in column K
546 *
547  CALL dscal( n-k, r1, ap( kc+1 ), 1 )
548  END IF
549  ELSE
550 *
551 * 2-by-2 pivot block D(k): columns K and K+1 now hold
552 *
553 * ( W(k) W(k+1) ) = ( L(k) L(k+1) )*D(k)
554 *
555 * where L(k) and L(k+1) are the k-th and (k+1)-th columns
556 * of L
557 *
558  IF( k.LT.n-1 ) THEN
559 *
560 * Perform a rank-2 update of A(k+2:n,k+2:n) as
561 *
562 * A := A - ( L(k) L(k+1) )*D(k)*( L(k) L(k+1) )**T
563 * = A - ( W(k) W(k+1) )*inv(D(k))*( W(k) W(k+1) )**T
564 *
565 * where L(k) and L(k+1) are the k-th and (k+1)-th
566 * columns of L
567 *
568  d21 = ap( k+1+( k-1 )*( 2*n-k ) / 2 )
569  d11 = ap( k+1+k*( 2*n-k-1 ) / 2 ) / d21
570  d22 = ap( k+( k-1 )*( 2*n-k ) / 2 ) / d21
571  t = one / ( d11*d22-one )
572  d21 = t / d21
573 *
574  DO 100 j = k + 2, n
575  wk = d21*( d11*ap( j+( k-1 )*( 2*n-k ) / 2 )-
576  $ ap( j+k*( 2*n-k-1 ) / 2 ) )
577  wkp1 = d21*( d22*ap( j+k*( 2*n-k-1 ) / 2 )-
578  $ ap( j+( k-1 )*( 2*n-k ) / 2 ) )
579 *
580  DO 90 i = j, n
581  ap( i+( j-1 )*( 2*n-j ) / 2 ) = ap( i+( j-1 )*
582  $ ( 2*n-j ) / 2 ) - ap( i+( k-1 )*( 2*n-k ) /
583  $ 2 )*wk - ap( i+k*( 2*n-k-1 ) / 2 )*wkp1
584  90 CONTINUE
585 *
586  ap( j+( k-1 )*( 2*n-k ) / 2 ) = wk
587  ap( j+k*( 2*n-k-1 ) / 2 ) = wkp1
588 *
589  100 CONTINUE
590  END IF
591  END IF
592  END IF
593 *
594 * Store details of the interchanges in IPIV
595 *
596  IF( kstep.EQ.1 ) THEN
597  ipiv( k ) = kp
598  ELSE
599  ipiv( k ) = -kp
600  ipiv( k+1 ) = -kp
601  END IF
602 *
603 * Increase K and return to the start of the main loop
604 *
605  k = k + kstep
606  kc = knc + n - k + 2
607  GO TO 60
608 *
609  END IF
610 *
611  110 CONTINUE
612  RETURN
613 *
614 * End of DSPTRF
615 *
integer function idamax(N, DX, INCX)
IDAMAX
Definition: idamax.f:53
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine dswap(N, DX, INCX, DY, INCY)
DSWAP
Definition: dswap.f:53
subroutine dscal(N, DA, DX, INCX)
DSCAL
Definition: dscal.f:55
subroutine dspr(UPLO, N, ALPHA, X, INCX, AP)
DSPR
Definition: dspr.f:129
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: