3 * =========== DOCUMENTATION ===========
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
9 *> Download ZSPTRF + dependencies
10 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.tgz?format=tgz&filename=/lapack/lapack_routine/zsptrf.f">
12 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.zip?format=zip&filename=/lapack/lapack_routine/zsptrf.f">
14 *> <a href="http://www.netlib.org/cgi-bin/netlibfiles.txt?format=txt&filename=/lapack/lapack_routine/zsptrf.f">
21 * SUBROUTINE ZSPTRF( UPLO, N, AP, IPIV, INFO )
23 * .. Scalar Arguments ..
27 * .. Array Arguments ..
38 *> ZSPTRF computes the factorization of a complex symmetric matrix A
39 *> stored in packed format using the Bunch-Kaufman diagonal pivoting
42 *> A = U*D*U**T or A = L*D*L**T
44 *> where U (or L) is a product of permutation and unit upper (lower)
45 *> triangular matrices, and D is symmetric and block diagonal with
46 *> 1-by-1 and 2-by-2 diagonal blocks.
54 *> UPLO is CHARACTER*1
55 *> = 'U': Upper triangle of A is stored;
56 *> = 'L': Lower triangle of A is stored.
62 *> The order of the matrix A. N >= 0.
67 *> AP is COMPLEX*16 array, dimension (N*(N+1)/2)
68 *> On entry, the upper or lower triangle of the symmetric matrix
69 *> A, packed columnwise in a linear array. The j-th column of A
70 *> is stored in the array AP as follows:
71 *> if UPLO = 'U', AP(i + (j-1)*j/2) = A(i,j) for 1<=i<=j;
72 *> if UPLO = 'L', AP(i + (j-1)*(2n-j)/2) = A(i,j) for j<=i<=n.
74 *> On exit, the block diagonal matrix D and the multipliers used
75 *> to obtain the factor U or L, stored as a packed triangular
76 *> matrix overwriting A (see below for further details).
81 *> IPIV is INTEGER array, dimension (N)
82 *> Details of the interchanges and the block structure of D.
83 *> If IPIV(k) > 0, then rows and columns k and IPIV(k) were
84 *> interchanged and D(k,k) is a 1-by-1 diagonal block.
85 *> If UPLO = 'U' and IPIV(k) = IPIV(k-1) < 0, then rows and
86 *> columns k-1 and -IPIV(k) were interchanged and D(k-1:k,k-1:k)
87 *> is a 2-by-2 diagonal block. If UPLO = 'L' and IPIV(k) =
88 *> IPIV(k+1) < 0, then rows and columns k+1 and -IPIV(k) were
89 *> interchanged and D(k:k+1,k:k+1) is a 2-by-2 diagonal block.
95 *> = 0: successful exit
96 *> < 0: if INFO = -i, the i-th argument had an illegal value
97 *> > 0: if INFO = i, D(i,i) is exactly zero. The factorization
98 *> has been completed, but the block diagonal matrix D is
99 *> exactly singular, and division by zero will occur if it
100 *> is used to solve a system of equations.
106 *> \author Univ. of Tennessee
107 *> \author Univ. of California Berkeley
108 *> \author Univ. of Colorado Denver
111 *> \date November 2011
113 *> \ingroup complex16OTHERcomputational
115 *> \par Further Details:
116 * =====================
120 *> 5-96 - Based on modifications by J. Lewis, Boeing Computer Services
123 *> If UPLO = 'U', then A = U*D*U**T, where
124 *> U = P(n)*U(n)* ... *P(k)U(k)* ...,
125 *> i.e., U is a product of terms P(k)*U(k), where k decreases from n to
126 *> 1 in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1
127 *> and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as
128 *> defined by IPIV(k), and U(k) is a unit upper triangular matrix, such
129 *> that if the diagonal block D(k) is of order s (s = 1 or 2), then
132 *> U(k) = ( 0 I 0 ) s
136 *> If s = 1, D(k) overwrites A(k,k), and v overwrites A(1:k-1,k).
137 *> If s = 2, the upper triangle of D(k) overwrites A(k-1,k-1), A(k-1,k),
138 *> and A(k,k), and v overwrites A(1:k-2,k-1:k).
140 *> If UPLO = 'L', then A = L*D*L**T, where
141 *> L = P(1)*L(1)* ... *P(k)*L(k)* ...,
142 *> i.e., L is a product of terms P(k)*L(k), where k increases from 1 to
143 *> n in steps of 1 or 2, and D is a block diagonal matrix with 1-by-1
144 *> and 2-by-2 diagonal blocks D(k). P(k) is a permutation matrix as
145 *> defined by IPIV(k), and L(k) is a unit lower triangular matrix, such
146 *> that if the diagonal block D(k) is of order s (s = 1 or 2), then
149 *> L(k) = ( 0 I 0 ) s
153 *> If s = 1, D(k) overwrites A(k,k), and v overwrites A(k+1:n,k).
154 *> If s = 2, the lower triangle of D(k) overwrites A(k,k), A(k+1,k),
155 *> and A(k+1,k+1), and v overwrites A(k+2:n,k:k+1).
158 * =====================================================================
159 SUBROUTINE ZSPTRF( UPLO, N, AP, IPIV, INFO )
161 * -- LAPACK computational routine (version 3.4.0) --
162 * -- LAPACK is a software package provided by Univ. of Tennessee, --
163 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
166 * .. Scalar Arguments ..
170 * .. Array Arguments ..
175 * =====================================================================
178 DOUBLE PRECISION ZERO, ONE
179 PARAMETER ( ZERO = 0.0D+0, ONE = 1.0D+0 )
180 DOUBLE PRECISION EIGHT, SEVTEN
181 PARAMETER ( EIGHT = 8.0D+0, SEVTEN = 17.0D+0 )
183 PARAMETER ( CONE = ( 1.0D+0, 0.0D+0 ) )
185 * .. Local Scalars ..
187 INTEGER I, IMAX, J, JMAX, K, KC, KK, KNC, KP, KPC,
189 DOUBLE PRECISION ABSAKK, ALPHA, COLMAX, ROWMAX
190 COMPLEX*16 D11, D12, D21, D22, R1, T, WK, WKM1, WKP1, ZDUM
192 * .. External Functions ..
195 EXTERNAL LSAME, IZAMAX
197 * .. External Subroutines ..
198 EXTERNAL XERBLA, ZSCAL, ZSPR, ZSWAP
200 * .. Intrinsic Functions ..
201 INTRINSIC ABS, DBLE, DIMAG, MAX, SQRT
203 * .. Statement Functions ..
204 DOUBLE PRECISION CABS1
206 * .. Statement Function definitions ..
207 CABS1( ZDUM ) = ABS( DBLE( ZDUM ) ) + ABS( DIMAG( ZDUM ) )
209 * .. Executable Statements ..
211 * Test the input parameters.
214 UPPER = LSAME( UPLO, 'U' )
215 IF( .NOT.UPPER .AND. .NOT.LSAME( UPLO, 'L' ) ) THEN
217 ELSE IF( N.LT.0 ) THEN
221 CALL XERBLA( 'ZSPTRF', -INFO )
225 * Initialize ALPHA for use in choosing pivot block size.
227 ALPHA = ( ONE+SQRT( SEVTEN ) ) / EIGHT
231 * Factorize A as U*D*U**T using the upper triangle of A
233 * K is the main loop index, decreasing from N to 1 in steps of
237 KC = ( N-1 )*N / 2 + 1
241 * If K < 1, exit from loop
247 * Determine rows and columns to be interchanged and whether
248 * a 1-by-1 or 2-by-2 pivot block will be used
250 ABSAKK = CABS1( AP( KC+K-1 ) )
252 * IMAX is the row-index of the largest off-diagonal element in
253 * column K, and COLMAX is its absolute value
256 IMAX = IZAMAX( K-1, AP( KC ), 1 )
257 COLMAX = CABS1( AP( KC+IMAX-1 ) )
262 IF( MAX( ABSAKK, COLMAX ).EQ.ZERO ) THEN
264 * Column K is zero: set INFO and continue
270 IF( ABSAKK.GE.ALPHA*COLMAX ) THEN
272 * no interchange, use 1-by-1 pivot block
279 KX = IMAX*( IMAX+1 ) / 2 + IMAX
280 DO 20 J = IMAX + 1, K
281 IF( CABS1( AP( KX ) ).GT.ROWMAX ) THEN
282 ROWMAX = CABS1( AP( KX ) )
287 KPC = ( IMAX-1 )*IMAX / 2 + 1
289 JMAX = IZAMAX( IMAX-1, AP( KPC ), 1 )
290 ROWMAX = MAX( ROWMAX, CABS1( AP( KPC+JMAX-1 ) ) )
293 IF( ABSAKK.GE.ALPHA*COLMAX*( COLMAX / ROWMAX ) ) THEN
295 * no interchange, use 1-by-1 pivot block
298 ELSE IF( CABS1( AP( KPC+IMAX-1 ) ).GE.ALPHA*ROWMAX ) THEN
300 * interchange rows and columns K and IMAX, use 1-by-1
306 * interchange rows and columns K-1 and IMAX, use 2-by-2
319 * Interchange rows and columns KK and KP in the leading
320 * submatrix A(1:k,1:k)
322 CALL ZSWAP( KP-1, AP( KNC ), 1, AP( KPC ), 1 )
324 DO 30 J = KP + 1, KK - 1
327 AP( KNC+J-1 ) = AP( KX )
331 AP( KNC+KK-1 ) = AP( KPC+KP-1 )
333 IF( KSTEP.EQ.2 ) THEN
335 AP( KC+K-2 ) = AP( KC+KP-1 )
340 * Update the leading submatrix
342 IF( KSTEP.EQ.1 ) THEN
344 * 1-by-1 pivot block D(k): column k now holds
348 * where U(k) is the k-th column of U
350 * Perform a rank-1 update of A(1:k-1,1:k-1) as
352 * A := A - U(k)*D(k)*U(k)**T = A - W(k)*1/D(k)*W(k)**T
354 R1 = CONE / AP( KC+K-1 )
355 CALL ZSPR( UPLO, K-1, -R1, AP( KC ), 1, AP )
357 * Store U(k) in column k
359 CALL ZSCAL( K-1, R1, AP( KC ), 1 )
362 * 2-by-2 pivot block D(k): columns k and k-1 now hold
364 * ( W(k-1) W(k) ) = ( U(k-1) U(k) )*D(k)
366 * where U(k) and U(k-1) are the k-th and (k-1)-th columns
369 * Perform a rank-2 update of A(1:k-2,1:k-2) as
371 * A := A - ( U(k-1) U(k) )*D(k)*( U(k-1) U(k) )**T
372 * = A - ( W(k-1) W(k) )*inv(D(k))*( W(k-1) W(k) )**T
376 D12 = AP( K-1+( K-1 )*K / 2 )
377 D22 = AP( K-1+( K-2 )*( K-1 ) / 2 ) / D12
378 D11 = AP( K+( K-1 )*K / 2 ) / D12
379 T = CONE / ( D11*D22-CONE )
382 DO 50 J = K - 2, 1, -1
383 WKM1 = D12*( D11*AP( J+( K-2 )*( K-1 ) / 2 )-
384 $ AP( J+( K-1 )*K / 2 ) )
385 WK = D12*( D22*AP( J+( K-1 )*K / 2 )-
386 $ AP( J+( K-2 )*( K-1 ) / 2 ) )
388 AP( I+( J-1 )*J / 2 ) = AP( I+( J-1 )*J / 2 ) -
389 $ AP( I+( K-1 )*K / 2 )*WK -
390 $ AP( I+( K-2 )*( K-1 ) / 2 )*WKM1
392 AP( J+( K-1 )*K / 2 ) = WK
393 AP( J+( K-2 )*( K-1 ) / 2 ) = WKM1
400 * Store details of the interchanges in IPIV
402 IF( KSTEP.EQ.1 ) THEN
409 * Decrease K and return to the start of the main loop
417 * Factorize A as L*D*L**T using the lower triangle of A
419 * K is the main loop index, increasing from 1 to N in steps of
428 * If K > N, exit from loop
434 * Determine rows and columns to be interchanged and whether
435 * a 1-by-1 or 2-by-2 pivot block will be used
437 ABSAKK = CABS1( AP( KC ) )
439 * IMAX is the row-index of the largest off-diagonal element in
440 * column K, and COLMAX is its absolute value
443 IMAX = K + IZAMAX( N-K, AP( KC+1 ), 1 )
444 COLMAX = CABS1( AP( KC+IMAX-K ) )
449 IF( MAX( ABSAKK, COLMAX ).EQ.ZERO ) THEN
451 * Column K is zero: set INFO and continue
457 IF( ABSAKK.GE.ALPHA*COLMAX ) THEN
459 * no interchange, use 1-by-1 pivot block
464 * JMAX is the column-index of the largest off-diagonal
465 * element in row IMAX, and ROWMAX is its absolute value
469 DO 70 J = K, IMAX - 1
470 IF( CABS1( AP( KX ) ).GT.ROWMAX ) THEN
471 ROWMAX = CABS1( AP( KX ) )
476 KPC = NPP - ( N-IMAX+1 )*( N-IMAX+2 ) / 2 + 1
478 JMAX = IMAX + IZAMAX( N-IMAX, AP( KPC+1 ), 1 )
479 ROWMAX = MAX( ROWMAX, CABS1( AP( KPC+JMAX-IMAX ) ) )
482 IF( ABSAKK.GE.ALPHA*COLMAX*( COLMAX / ROWMAX ) ) THEN
484 * no interchange, use 1-by-1 pivot block
487 ELSE IF( CABS1( AP( KPC ) ).GE.ALPHA*ROWMAX ) THEN
489 * interchange rows and columns K and IMAX, use 1-by-1
495 * interchange rows and columns K+1 and IMAX, use 2-by-2
505 $ KNC = KNC + N - K + 1
508 * Interchange rows and columns KK and KP in the trailing
509 * submatrix A(k:n,k:n)
512 $ CALL ZSWAP( N-KP, AP( KNC+KP-KK+1 ), 1, AP( KPC+1 ),
515 DO 80 J = KK + 1, KP - 1
518 AP( KNC+J-KK ) = AP( KX )
522 AP( KNC ) = AP( KPC )
524 IF( KSTEP.EQ.2 ) THEN
526 AP( KC+1 ) = AP( KC+KP-K )
531 * Update the trailing submatrix
533 IF( KSTEP.EQ.1 ) THEN
535 * 1-by-1 pivot block D(k): column k now holds
539 * where L(k) is the k-th column of L
543 * Perform a rank-1 update of A(k+1:n,k+1:n) as
545 * A := A - L(k)*D(k)*L(k)**T = A - W(k)*(1/D(k))*W(k)**T
548 CALL ZSPR( UPLO, N-K, -R1, AP( KC+1 ), 1,
551 * Store L(k) in column K
553 CALL ZSCAL( N-K, R1, AP( KC+1 ), 1 )
557 * 2-by-2 pivot block D(k): columns K and K+1 now hold
559 * ( W(k) W(k+1) ) = ( L(k) L(k+1) )*D(k)
561 * where L(k) and L(k+1) are the k-th and (k+1)-th columns
566 * Perform a rank-2 update of A(k+2:n,k+2:n) as
568 * A := A - ( L(k) L(k+1) )*D(k)*( L(k) L(k+1) )**T
569 * = A - ( W(k) W(k+1) )*inv(D(k))*( W(k) W(k+1) )**T
571 * where L(k) and L(k+1) are the k-th and (k+1)-th
574 D21 = AP( K+1+( K-1 )*( 2*N-K ) / 2 )
575 D11 = AP( K+1+K*( 2*N-K-1 ) / 2 ) / D21
576 D22 = AP( K+( K-1 )*( 2*N-K ) / 2 ) / D21
577 T = CONE / ( D11*D22-CONE )
581 WK = D21*( D11*AP( J+( K-1 )*( 2*N-K ) / 2 )-
582 $ AP( J+K*( 2*N-K-1 ) / 2 ) )
583 WKP1 = D21*( D22*AP( J+K*( 2*N-K-1 ) / 2 )-
584 $ AP( J+( K-1 )*( 2*N-K ) / 2 ) )
586 AP( I+( J-1 )*( 2*N-J ) / 2 ) = AP( I+( J-1 )*
587 $ ( 2*N-J ) / 2 ) - AP( I+( K-1 )*( 2*N-K ) /
588 $ 2 )*WK - AP( I+K*( 2*N-K-1 ) / 2 )*WKP1
590 AP( J+( K-1 )*( 2*N-K ) / 2 ) = WK
591 AP( J+K*( 2*N-K-1 ) / 2 ) = WKP1
597 * Store details of the interchanges in IPIV
599 IF( KSTEP.EQ.1 ) THEN
606 * Increase K and return to the start of the main loop