LAPACK  3.6.0
LAPACK: Linear Algebra PACKage
Collaboration diagram for complex16:

Functions

subroutine zptcon (N, D, E, ANORM, RCOND, RWORK, INFO)
 ZPTCON More...
 
subroutine zpteqr (COMPZ, N, D, E, Z, LDZ, WORK, INFO)
 ZPTEQR More...
 
subroutine zptrfs (UPLO, N, NRHS, D, E, DF, EF, B, LDB, X, LDX, FERR, BERR, WORK, RWORK, INFO)
 ZPTRFS More...
 
subroutine zpttrf (N, D, E, INFO)
 ZPTTRF More...
 
subroutine zpttrs (UPLO, N, NRHS, D, E, B, LDB, INFO)
 ZPTTRS More...
 
subroutine zptts2 (IUPLO, N, NRHS, D, E, B, LDB)
 ZPTTS2 solves a tridiagonal system of the form AX=B using the L D LH factorization computed by spttrf. More...
 

Detailed Description

This is the group of complex16 computational functions for PT matrices

Function Documentation

subroutine zptcon ( integer  N,
double precision, dimension( * )  D,
complex*16, dimension( * )  E,
double precision  ANORM,
double precision  RCOND,
double precision, dimension( * )  RWORK,
integer  INFO 
)

ZPTCON

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

Purpose:
 ZPTCON computes the reciprocal of the condition number (in the
 1-norm) of a complex Hermitian positive definite tridiagonal matrix
 using the factorization A = L*D*L**H or A = U**H*D*U computed by
 ZPTTRF.

 Norm(inv(A)) is computed by a direct method, and the reciprocal of
 the condition number is computed as
                  RCOND = 1 / (ANORM * norm(inv(A))).
Parameters
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]D
          D is DOUBLE PRECISION array, dimension (N)
          The n diagonal elements of the diagonal matrix D from the
          factorization of A, as computed by ZPTTRF.
[in]E
          E is COMPLEX*16 array, dimension (N-1)
          The (n-1) off-diagonal elements of the unit bidiagonal factor
          U or L from the factorization of A, as computed by ZPTTRF.
[in]ANORM
          ANORM is DOUBLE PRECISION
          The 1-norm of the original matrix A.
[out]RCOND
          RCOND is DOUBLE PRECISION
          The reciprocal of the condition number of the matrix A,
          computed as RCOND = 1/(ANORM * AINVNM), where AINVNM is the
          1-norm of inv(A) computed in this routine.
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012
Further Details:
  The method used is described in Nicholas J. Higham, "Efficient
  Algorithms for Computing the Condition Number of a Tridiagonal
  Matrix", SIAM J. Sci. Stat. Comput., Vol. 7, No. 1, January 1986.

Definition at line 121 of file zptcon.f.

121 *
122 * -- LAPACK computational routine (version 3.4.2) --
123 * -- LAPACK is a software package provided by Univ. of Tennessee, --
124 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
125 * September 2012
126 *
127 * .. Scalar Arguments ..
128  INTEGER info, n
129  DOUBLE PRECISION anorm, rcond
130 * ..
131 * .. Array Arguments ..
132  DOUBLE PRECISION d( * ), rwork( * )
133  COMPLEX*16 e( * )
134 * ..
135 *
136 * =====================================================================
137 *
138 * .. Parameters ..
139  DOUBLE PRECISION one, zero
140  parameter( one = 1.0d+0, zero = 0.0d+0 )
141 * ..
142 * .. Local Scalars ..
143  INTEGER i, ix
144  DOUBLE PRECISION ainvnm
145 * ..
146 * .. External Functions ..
147  INTEGER idamax
148  EXTERNAL idamax
149 * ..
150 * .. External Subroutines ..
151  EXTERNAL xerbla
152 * ..
153 * .. Intrinsic Functions ..
154  INTRINSIC abs
155 * ..
156 * .. Executable Statements ..
157 *
158 * Test the input arguments.
159 *
160  info = 0
161  IF( n.LT.0 ) THEN
162  info = -1
163  ELSE IF( anorm.LT.zero ) THEN
164  info = -4
165  END IF
166  IF( info.NE.0 ) THEN
167  CALL xerbla( 'ZPTCON', -info )
168  RETURN
169  END IF
170 *
171 * Quick return if possible
172 *
173  rcond = zero
174  IF( n.EQ.0 ) THEN
175  rcond = one
176  RETURN
177  ELSE IF( anorm.EQ.zero ) THEN
178  RETURN
179  END IF
180 *
181 * Check that D(1:N) is positive.
182 *
183  DO 10 i = 1, n
184  IF( d( i ).LE.zero )
185  $ RETURN
186  10 CONTINUE
187 *
188 * Solve M(A) * x = e, where M(A) = (m(i,j)) is given by
189 *
190 * m(i,j) = abs(A(i,j)), i = j,
191 * m(i,j) = -abs(A(i,j)), i .ne. j,
192 *
193 * and e = [ 1, 1, ..., 1 ]**T. Note M(A) = M(L)*D*M(L)**H.
194 *
195 * Solve M(L) * x = e.
196 *
197  rwork( 1 ) = one
198  DO 20 i = 2, n
199  rwork( i ) = one + rwork( i-1 )*abs( e( i-1 ) )
200  20 CONTINUE
201 *
202 * Solve D * M(L)**H * x = b.
203 *
204  rwork( n ) = rwork( n ) / d( n )
205  DO 30 i = n - 1, 1, -1
206  rwork( i ) = rwork( i ) / d( i ) + rwork( i+1 )*abs( e( i ) )
207  30 CONTINUE
208 *
209 * Compute AINVNM = max(x(i)), 1<=i<=n.
210 *
211  ix = idamax( n, rwork, 1 )
212  ainvnm = abs( rwork( ix ) )
213 *
214 * Compute the reciprocal condition number.
215 *
216  IF( ainvnm.NE.zero )
217  $ rcond = ( one / ainvnm ) / anorm
218 *
219  RETURN
220 *
221 * End of ZPTCON
222 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
integer function idamax(N, DX, INCX)
IDAMAX
Definition: idamax.f:53

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zpteqr ( character  COMPZ,
integer  N,
double precision, dimension( * )  D,
double precision, dimension( * )  E,
complex*16, dimension( ldz, * )  Z,
integer  LDZ,
double precision, dimension( * )  WORK,
integer  INFO 
)

ZPTEQR

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

Purpose:
 ZPTEQR computes all eigenvalues and, optionally, eigenvectors of a
 symmetric positive definite tridiagonal matrix by first factoring the
 matrix using DPTTRF and then calling ZBDSQR to compute the singular
 values of the bidiagonal factor.

 This routine computes the eigenvalues of the positive definite
 tridiagonal matrix to high relative accuracy.  This means that if the
 eigenvalues range over many orders of magnitude in size, then the
 small eigenvalues and corresponding eigenvectors will be computed
 more accurately than, for example, with the standard QR method.

 The eigenvectors of a full or band positive definite Hermitian matrix
 can also be found if ZHETRD, ZHPTRD, or ZHBTRD has been used to
 reduce this matrix to tridiagonal form.  (The reduction to
 tridiagonal form, however, may preclude the possibility of obtaining
 high relative accuracy in the small eigenvalues of the original
 matrix, if these eigenvalues range over many orders of magnitude.)
Parameters
[in]COMPZ
          COMPZ is CHARACTER*1
          = 'N':  Compute eigenvalues only.
          = 'V':  Compute eigenvectors of original Hermitian
                  matrix also.  Array Z contains the unitary matrix
                  used to reduce the original matrix to tridiagonal
                  form.
          = 'I':  Compute eigenvectors of tridiagonal matrix also.
[in]N
          N is INTEGER
          The order of the matrix.  N >= 0.
[in,out]D
          D is DOUBLE PRECISION array, dimension (N)
          On entry, the n diagonal elements of the tridiagonal matrix.
          On normal exit, D contains the eigenvalues, in descending
          order.
[in,out]E
          E is DOUBLE PRECISION array, dimension (N-1)
          On entry, the (n-1) subdiagonal elements of the tridiagonal
          matrix.
          On exit, E has been destroyed.
[in,out]Z
          Z is COMPLEX*16 array, dimension (LDZ, N)
          On entry, if COMPZ = 'V', the unitary matrix used in the
          reduction to tridiagonal form.
          On exit, if COMPZ = 'V', the orthonormal eigenvectors of the
          original Hermitian matrix;
          if COMPZ = 'I', the orthonormal eigenvectors of the
          tridiagonal matrix.
          If INFO > 0 on exit, Z contains the eigenvectors associated
          with only the stored eigenvalues.
          If  COMPZ = 'N', then Z is not referenced.
[in]LDZ
          LDZ is INTEGER
          The leading dimension of the array Z.  LDZ >= 1, and if
          COMPZ = 'V' or 'I', LDZ >= max(1,N).
[out]WORK
          WORK is DOUBLE PRECISION array, dimension (4*N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit.
          < 0:  if INFO = -i, the i-th argument had an illegal value.
          > 0:  if INFO = i, and i is:
                <= N  the Cholesky factorization of the matrix could
                      not be performed because the i-th principal minor
                      was not positive definite.
                > N   the SVD algorithm failed to converge;
                      if INFO = N+i, i off-diagonal elements of the
                      bidiagonal factor did not converge to zero.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 147 of file zpteqr.f.

147 *
148 * -- LAPACK computational routine (version 3.4.2) --
149 * -- LAPACK is a software package provided by Univ. of Tennessee, --
150 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
151 * September 2012
152 *
153 * .. Scalar Arguments ..
154  CHARACTER compz
155  INTEGER info, ldz, n
156 * ..
157 * .. Array Arguments ..
158  DOUBLE PRECISION d( * ), e( * ), work( * )
159  COMPLEX*16 z( ldz, * )
160 * ..
161 *
162 * ====================================================================
163 *
164 * .. Parameters ..
165  COMPLEX*16 czero, cone
166  parameter( czero = ( 0.0d+0, 0.0d+0 ),
167  $ cone = ( 1.0d+0, 0.0d+0 ) )
168 * ..
169 * .. External Functions ..
170  LOGICAL lsame
171  EXTERNAL lsame
172 * ..
173 * .. External Subroutines ..
174  EXTERNAL dpttrf, xerbla, zbdsqr, zlaset
175 * ..
176 * .. Local Arrays ..
177  COMPLEX*16 c( 1, 1 ), vt( 1, 1 )
178 * ..
179 * .. Local Scalars ..
180  INTEGER i, icompz, nru
181 * ..
182 * .. Intrinsic Functions ..
183  INTRINSIC max, sqrt
184 * ..
185 * .. Executable Statements ..
186 *
187 * Test the input parameters.
188 *
189  info = 0
190 *
191  IF( lsame( compz, 'N' ) ) THEN
192  icompz = 0
193  ELSE IF( lsame( compz, 'V' ) ) THEN
194  icompz = 1
195  ELSE IF( lsame( compz, 'I' ) ) THEN
196  icompz = 2
197  ELSE
198  icompz = -1
199  END IF
200  IF( icompz.LT.0 ) THEN
201  info = -1
202  ELSE IF( n.LT.0 ) THEN
203  info = -2
204  ELSE IF( ( ldz.LT.1 ) .OR. ( icompz.GT.0 .AND. ldz.LT.max( 1,
205  $ n ) ) ) THEN
206  info = -6
207  END IF
208  IF( info.NE.0 ) THEN
209  CALL xerbla( 'ZPTEQR', -info )
210  RETURN
211  END IF
212 *
213 * Quick return if possible
214 *
215  IF( n.EQ.0 )
216  $ RETURN
217 *
218  IF( n.EQ.1 ) THEN
219  IF( icompz.GT.0 )
220  $ z( 1, 1 ) = cone
221  RETURN
222  END IF
223  IF( icompz.EQ.2 )
224  $ CALL zlaset( 'Full', n, n, czero, cone, z, ldz )
225 *
226 * Call DPTTRF to factor the matrix.
227 *
228  CALL dpttrf( n, d, e, info )
229  IF( info.NE.0 )
230  $ RETURN
231  DO 10 i = 1, n
232  d( i ) = sqrt( d( i ) )
233  10 CONTINUE
234  DO 20 i = 1, n - 1
235  e( i ) = e( i )*d( i )
236  20 CONTINUE
237 *
238 * Call ZBDSQR to compute the singular values/vectors of the
239 * bidiagonal factor.
240 *
241  IF( icompz.GT.0 ) THEN
242  nru = n
243  ELSE
244  nru = 0
245  END IF
246  CALL zbdsqr( 'Lower', n, 0, nru, 0, d, e, vt, 1, z, ldz, c, 1,
247  $ work, info )
248 *
249 * Square the singular values.
250 *
251  IF( info.EQ.0 ) THEN
252  DO 30 i = 1, n
253  d( i ) = d( i )*d( i )
254  30 CONTINUE
255  ELSE
256  info = n + info
257  END IF
258 *
259  RETURN
260 *
261 * End of ZPTEQR
262 *
subroutine dpttrf(N, D, E, INFO)
DPTTRF
Definition: dpttrf.f:93
subroutine zbdsqr(UPLO, N, NCVT, NRU, NCC, D, E, VT, LDVT, U, LDU, C, LDC, RWORK, INFO)
ZBDSQR
Definition: zbdsqr.f:224
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zlaset(UPLO, M, N, ALPHA, BETA, A, LDA)
ZLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values...
Definition: zlaset.f:108
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:

subroutine zptrfs ( character  UPLO,
integer  N,
integer  NRHS,
double precision, dimension( * )  D,
complex*16, dimension( * )  E,
double precision, dimension( * )  DF,
complex*16, dimension( * )  EF,
complex*16, dimension( ldb, * )  B,
integer  LDB,
complex*16, dimension( ldx, * )  X,
integer  LDX,
double precision, dimension( * )  FERR,
double precision, dimension( * )  BERR,
complex*16, dimension( * )  WORK,
double precision, dimension( * )  RWORK,
integer  INFO 
)

ZPTRFS

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

Purpose:
 ZPTRFS improves the computed solution to a system of linear
 equations when the coefficient matrix is Hermitian positive definite
 and tridiagonal, and provides error bounds and backward error
 estimates for the solution.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          Specifies whether the superdiagonal or the subdiagonal of the
          tridiagonal matrix A is stored and the form of the
          factorization:
          = 'U':  E is the superdiagonal of A, and A = U**H*D*U;
          = 'L':  E is the subdiagonal of A, and A = L*D*L**H.
          (The two forms are equivalent if A is real.)
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in]NRHS
          NRHS is INTEGER
          The number of right hand sides, i.e., the number of columns
          of the matrix B.  NRHS >= 0.
[in]D
          D is DOUBLE PRECISION array, dimension (N)
          The n real diagonal elements of the tridiagonal matrix A.
[in]E
          E is COMPLEX*16 array, dimension (N-1)
          The (n-1) off-diagonal elements of the tridiagonal matrix A
          (see UPLO).
[in]DF
          DF is DOUBLE PRECISION array, dimension (N)
          The n diagonal elements of the diagonal matrix D from
          the factorization computed by ZPTTRF.
[in]EF
          EF is COMPLEX*16 array, dimension (N-1)
          The (n-1) off-diagonal elements of the unit bidiagonal
          factor U or L from the factorization computed by ZPTTRF
          (see UPLO).
[in]B
          B is COMPLEX*16 array, dimension (LDB,NRHS)
          The right hand side matrix B.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[in,out]X
          X is COMPLEX*16 array, dimension (LDX,NRHS)
          On entry, the solution matrix X, as computed by ZPTTRS.
          On exit, the improved solution matrix X.
[in]LDX
          LDX is INTEGER
          The leading dimension of the array X.  LDX >= max(1,N).
[out]FERR
          FERR is DOUBLE PRECISION array, dimension (NRHS)
          The forward error bound for each solution vector
          X(j) (the j-th column of the solution matrix X).
          If XTRUE is the true solution corresponding to X(j), FERR(j)
          is an estimated upper bound for the magnitude of the largest
          element in (X(j) - XTRUE) divided by the magnitude of the
          largest element in X(j).
[out]BERR
          BERR is DOUBLE PRECISION array, dimension (NRHS)
          The componentwise relative backward error of each solution
          vector X(j) (i.e., the smallest relative change in
          any element of A or B that makes X(j) an exact solution).
[out]WORK
          WORK is COMPLEX*16 array, dimension (N)
[out]RWORK
          RWORK is DOUBLE PRECISION array, dimension (N)
[out]INFO
          INFO is INTEGER
          = 0:  successful exit
          < 0:  if INFO = -i, the i-th argument had an illegal value
Internal Parameters:
  ITMAX is the maximum number of steps of iterative refinement.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 185 of file zptrfs.f.

185 *
186 * -- LAPACK computational routine (version 3.4.2) --
187 * -- LAPACK is a software package provided by Univ. of Tennessee, --
188 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
189 * September 2012
190 *
191 * .. Scalar Arguments ..
192  CHARACTER uplo
193  INTEGER info, ldb, ldx, n, nrhs
194 * ..
195 * .. Array Arguments ..
196  DOUBLE PRECISION berr( * ), d( * ), df( * ), ferr( * ),
197  $ rwork( * )
198  COMPLEX*16 b( ldb, * ), e( * ), ef( * ), work( * ),
199  $ x( ldx, * )
200 * ..
201 *
202 * =====================================================================
203 *
204 * .. Parameters ..
205  INTEGER itmax
206  parameter( itmax = 5 )
207  DOUBLE PRECISION zero
208  parameter( zero = 0.0d+0 )
209  DOUBLE PRECISION one
210  parameter( one = 1.0d+0 )
211  DOUBLE PRECISION two
212  parameter( two = 2.0d+0 )
213  DOUBLE PRECISION three
214  parameter( three = 3.0d+0 )
215 * ..
216 * .. Local Scalars ..
217  LOGICAL upper
218  INTEGER count, i, ix, j, nz
219  DOUBLE PRECISION eps, lstres, s, safe1, safe2, safmin
220  COMPLEX*16 bi, cx, dx, ex, zdum
221 * ..
222 * .. External Functions ..
223  LOGICAL lsame
224  INTEGER idamax
225  DOUBLE PRECISION dlamch
226  EXTERNAL lsame, idamax, dlamch
227 * ..
228 * .. External Subroutines ..
229  EXTERNAL xerbla, zaxpy, zpttrs
230 * ..
231 * .. Intrinsic Functions ..
232  INTRINSIC abs, dble, dcmplx, dconjg, dimag, max
233 * ..
234 * .. Statement Functions ..
235  DOUBLE PRECISION cabs1
236 * ..
237 * .. Statement Function definitions ..
238  cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
239 * ..
240 * .. Executable Statements ..
241 *
242 * Test the input parameters.
243 *
244  info = 0
245  upper = lsame( uplo, 'U' )
246  IF( .NOT.upper .AND. .NOT.lsame( uplo, 'L' ) ) THEN
247  info = -1
248  ELSE IF( n.LT.0 ) THEN
249  info = -2
250  ELSE IF( nrhs.LT.0 ) THEN
251  info = -3
252  ELSE IF( ldb.LT.max( 1, n ) ) THEN
253  info = -9
254  ELSE IF( ldx.LT.max( 1, n ) ) THEN
255  info = -11
256  END IF
257  IF( info.NE.0 ) THEN
258  CALL xerbla( 'ZPTRFS', -info )
259  RETURN
260  END IF
261 *
262 * Quick return if possible
263 *
264  IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
265  DO 10 j = 1, nrhs
266  ferr( j ) = zero
267  berr( j ) = zero
268  10 CONTINUE
269  RETURN
270  END IF
271 *
272 * NZ = maximum number of nonzero elements in each row of A, plus 1
273 *
274  nz = 4
275  eps = dlamch( 'Epsilon' )
276  safmin = dlamch( 'Safe minimum' )
277  safe1 = nz*safmin
278  safe2 = safe1 / eps
279 *
280 * Do for each right hand side
281 *
282  DO 100 j = 1, nrhs
283 *
284  count = 1
285  lstres = three
286  20 CONTINUE
287 *
288 * Loop until stopping criterion is satisfied.
289 *
290 * Compute residual R = B - A * X. Also compute
291 * abs(A)*abs(x) + abs(b) for use in the backward error bound.
292 *
293  IF( upper ) THEN
294  IF( n.EQ.1 ) THEN
295  bi = b( 1, j )
296  dx = d( 1 )*x( 1, j )
297  work( 1 ) = bi - dx
298  rwork( 1 ) = cabs1( bi ) + cabs1( dx )
299  ELSE
300  bi = b( 1, j )
301  dx = d( 1 )*x( 1, j )
302  ex = e( 1 )*x( 2, j )
303  work( 1 ) = bi - dx - ex
304  rwork( 1 ) = cabs1( bi ) + cabs1( dx ) +
305  $ cabs1( e( 1 ) )*cabs1( x( 2, j ) )
306  DO 30 i = 2, n - 1
307  bi = b( i, j )
308  cx = dconjg( e( i-1 ) )*x( i-1, j )
309  dx = d( i )*x( i, j )
310  ex = e( i )*x( i+1, j )
311  work( i ) = bi - cx - dx - ex
312  rwork( i ) = cabs1( bi ) +
313  $ cabs1( e( i-1 ) )*cabs1( x( i-1, j ) ) +
314  $ cabs1( dx ) + cabs1( e( i ) )*
315  $ cabs1( x( i+1, j ) )
316  30 CONTINUE
317  bi = b( n, j )
318  cx = dconjg( e( n-1 ) )*x( n-1, j )
319  dx = d( n )*x( n, j )
320  work( n ) = bi - cx - dx
321  rwork( n ) = cabs1( bi ) + cabs1( e( n-1 ) )*
322  $ cabs1( x( n-1, j ) ) + cabs1( dx )
323  END IF
324  ELSE
325  IF( n.EQ.1 ) THEN
326  bi = b( 1, j )
327  dx = d( 1 )*x( 1, j )
328  work( 1 ) = bi - dx
329  rwork( 1 ) = cabs1( bi ) + cabs1( dx )
330  ELSE
331  bi = b( 1, j )
332  dx = d( 1 )*x( 1, j )
333  ex = dconjg( e( 1 ) )*x( 2, j )
334  work( 1 ) = bi - dx - ex
335  rwork( 1 ) = cabs1( bi ) + cabs1( dx ) +
336  $ cabs1( e( 1 ) )*cabs1( x( 2, j ) )
337  DO 40 i = 2, n - 1
338  bi = b( i, j )
339  cx = e( i-1 )*x( i-1, j )
340  dx = d( i )*x( i, j )
341  ex = dconjg( e( i ) )*x( i+1, j )
342  work( i ) = bi - cx - dx - ex
343  rwork( i ) = cabs1( bi ) +
344  $ cabs1( e( i-1 ) )*cabs1( x( i-1, j ) ) +
345  $ cabs1( dx ) + cabs1( e( i ) )*
346  $ cabs1( x( i+1, j ) )
347  40 CONTINUE
348  bi = b( n, j )
349  cx = e( n-1 )*x( n-1, j )
350  dx = d( n )*x( n, j )
351  work( n ) = bi - cx - dx
352  rwork( n ) = cabs1( bi ) + cabs1( e( n-1 ) )*
353  $ cabs1( x( n-1, j ) ) + cabs1( dx )
354  END IF
355  END IF
356 *
357 * Compute componentwise relative backward error from formula
358 *
359 * max(i) ( abs(R(i)) / ( abs(A)*abs(X) + abs(B) )(i) )
360 *
361 * where abs(Z) is the componentwise absolute value of the matrix
362 * or vector Z. If the i-th component of the denominator is less
363 * than SAFE2, then SAFE1 is added to the i-th components of the
364 * numerator and denominator before dividing.
365 *
366  s = zero
367  DO 50 i = 1, n
368  IF( rwork( i ).GT.safe2 ) THEN
369  s = max( s, cabs1( work( i ) ) / rwork( i ) )
370  ELSE
371  s = max( s, ( cabs1( work( i ) )+safe1 ) /
372  $ ( rwork( i )+safe1 ) )
373  END IF
374  50 CONTINUE
375  berr( j ) = s
376 *
377 * Test stopping criterion. Continue iterating if
378 * 1) The residual BERR(J) is larger than machine epsilon, and
379 * 2) BERR(J) decreased by at least a factor of 2 during the
380 * last iteration, and
381 * 3) At most ITMAX iterations tried.
382 *
383  IF( berr( j ).GT.eps .AND. two*berr( j ).LE.lstres .AND.
384  $ count.LE.itmax ) THEN
385 *
386 * Update solution and try again.
387 *
388  CALL zpttrs( uplo, n, 1, df, ef, work, n, info )
389  CALL zaxpy( n, dcmplx( one ), work, 1, x( 1, j ), 1 )
390  lstres = berr( j )
391  count = count + 1
392  GO TO 20
393  END IF
394 *
395 * Bound error from formula
396 *
397 * norm(X - XTRUE) / norm(X) .le. FERR =
398 * norm( abs(inv(A))*
399 * ( abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) ))) / norm(X)
400 *
401 * where
402 * norm(Z) is the magnitude of the largest component of Z
403 * inv(A) is the inverse of A
404 * abs(Z) is the componentwise absolute value of the matrix or
405 * vector Z
406 * NZ is the maximum number of nonzeros in any row of A, plus 1
407 * EPS is machine epsilon
408 *
409 * The i-th component of abs(R)+NZ*EPS*(abs(A)*abs(X)+abs(B))
410 * is incremented by SAFE1 if the i-th component of
411 * abs(A)*abs(X) + abs(B) is less than SAFE2.
412 *
413  DO 60 i = 1, n
414  IF( rwork( i ).GT.safe2 ) THEN
415  rwork( i ) = cabs1( work( i ) ) + nz*eps*rwork( i )
416  ELSE
417  rwork( i ) = cabs1( work( i ) ) + nz*eps*rwork( i ) +
418  $ safe1
419  END IF
420  60 CONTINUE
421  ix = idamax( n, rwork, 1 )
422  ferr( j ) = rwork( ix )
423 *
424 * Estimate the norm of inv(A).
425 *
426 * Solve M(A) * x = e, where M(A) = (m(i,j)) is given by
427 *
428 * m(i,j) = abs(A(i,j)), i = j,
429 * m(i,j) = -abs(A(i,j)), i .ne. j,
430 *
431 * and e = [ 1, 1, ..., 1 ]**T. Note M(A) = M(L)*D*M(L)**H.
432 *
433 * Solve M(L) * x = e.
434 *
435  rwork( 1 ) = one
436  DO 70 i = 2, n
437  rwork( i ) = one + rwork( i-1 )*abs( ef( i-1 ) )
438  70 CONTINUE
439 *
440 * Solve D * M(L)**H * x = b.
441 *
442  rwork( n ) = rwork( n ) / df( n )
443  DO 80 i = n - 1, 1, -1
444  rwork( i ) = rwork( i ) / df( i ) +
445  $ rwork( i+1 )*abs( ef( i ) )
446  80 CONTINUE
447 *
448 * Compute norm(inv(A)) = max(x(i)), 1<=i<=n.
449 *
450  ix = idamax( n, rwork, 1 )
451  ferr( j ) = ferr( j )*abs( rwork( ix ) )
452 *
453 * Normalize error.
454 *
455  lstres = zero
456  DO 90 i = 1, n
457  lstres = max( lstres, abs( x( i, j ) ) )
458  90 CONTINUE
459  IF( lstres.NE.zero )
460  $ ferr( j ) = ferr( j ) / lstres
461 *
462  100 CONTINUE
463 *
464  RETURN
465 *
466 * End of ZPTRFS
467 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zpttrs(UPLO, N, NRHS, D, E, B, LDB, INFO)
ZPTTRS
Definition: zpttrs.f:123
double precision function dlamch(CMACH)
DLAMCH
Definition: dlamch.f:65
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55
subroutine zaxpy(N, ZA, ZX, INCX, ZY, INCY)
ZAXPY
Definition: zaxpy.f:53
integer function idamax(N, DX, INCX)
IDAMAX
Definition: idamax.f:53

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zpttrf ( integer  N,
double precision, dimension( * )  D,
complex*16, dimension( * )  E,
integer  INFO 
)

ZPTTRF

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

Purpose:
 ZPTTRF computes the L*D*L**H factorization of a complex Hermitian
 positive definite tridiagonal matrix A.  The factorization may also
 be regarded as having the form A = U**H *D*U.
Parameters
[in]N
          N is INTEGER
          The order of the matrix A.  N >= 0.
[in,out]D
          D is DOUBLE PRECISION array, dimension (N)
          On entry, the n diagonal elements of the tridiagonal matrix
          A.  On exit, the n diagonal elements of the diagonal matrix
          D from the L*D*L**H factorization of A.
[in,out]E
          E is COMPLEX*16 array, dimension (N-1)
          On entry, the (n-1) subdiagonal elements of the tridiagonal
          matrix A.  On exit, the (n-1) subdiagonal elements of the
          unit bidiagonal factor L from the L*D*L**H factorization of A.
          E can also be regarded as the superdiagonal of the unit
          bidiagonal factor U from the U**H *D*U factorization of A.
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -k, the k-th argument had an illegal value
          > 0: if INFO = k, the leading minor of order k is not
               positive definite; if k < N, the factorization could not
               be completed, while if k = N, the factorization was
               completed, but D(N) <= 0.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 94 of file zpttrf.f.

94 *
95 * -- LAPACK computational routine (version 3.4.2) --
96 * -- LAPACK is a software package provided by Univ. of Tennessee, --
97 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
98 * September 2012
99 *
100 * .. Scalar Arguments ..
101  INTEGER info, n
102 * ..
103 * .. Array Arguments ..
104  DOUBLE PRECISION d( * )
105  COMPLEX*16 e( * )
106 * ..
107 *
108 * =====================================================================
109 *
110 * .. Parameters ..
111  DOUBLE PRECISION zero
112  parameter( zero = 0.0d+0 )
113 * ..
114 * .. Local Scalars ..
115  INTEGER i, i4
116  DOUBLE PRECISION eii, eir, f, g
117 * ..
118 * .. External Subroutines ..
119  EXTERNAL xerbla
120 * ..
121 * .. Intrinsic Functions ..
122  INTRINSIC dble, dcmplx, dimag, mod
123 * ..
124 * .. Executable Statements ..
125 *
126 * Test the input parameters.
127 *
128  info = 0
129  IF( n.LT.0 ) THEN
130  info = -1
131  CALL xerbla( 'ZPTTRF', -info )
132  RETURN
133  END IF
134 *
135 * Quick return if possible
136 *
137  IF( n.EQ.0 )
138  $ RETURN
139 *
140 * Compute the L*D*L**H (or U**H *D*U) factorization of A.
141 *
142  i4 = mod( n-1, 4 )
143  DO 10 i = 1, i4
144  IF( d( i ).LE.zero ) THEN
145  info = i
146  GO TO 30
147  END IF
148  eir = dble( e( i ) )
149  eii = dimag( e( i ) )
150  f = eir / d( i )
151  g = eii / d( i )
152  e( i ) = dcmplx( f, g )
153  d( i+1 ) = d( i+1 ) - f*eir - g*eii
154  10 CONTINUE
155 *
156  DO 20 i = i4 + 1, n - 4, 4
157 *
158 * Drop out of the loop if d(i) <= 0: the matrix is not positive
159 * definite.
160 *
161  IF( d( i ).LE.zero ) THEN
162  info = i
163  GO TO 30
164  END IF
165 *
166 * Solve for e(i) and d(i+1).
167 *
168  eir = dble( e( i ) )
169  eii = dimag( e( i ) )
170  f = eir / d( i )
171  g = eii / d( i )
172  e( i ) = dcmplx( f, g )
173  d( i+1 ) = d( i+1 ) - f*eir - g*eii
174 *
175  IF( d( i+1 ).LE.zero ) THEN
176  info = i + 1
177  GO TO 30
178  END IF
179 *
180 * Solve for e(i+1) and d(i+2).
181 *
182  eir = dble( e( i+1 ) )
183  eii = dimag( e( i+1 ) )
184  f = eir / d( i+1 )
185  g = eii / d( i+1 )
186  e( i+1 ) = dcmplx( f, g )
187  d( i+2 ) = d( i+2 ) - f*eir - g*eii
188 *
189  IF( d( i+2 ).LE.zero ) THEN
190  info = i + 2
191  GO TO 30
192  END IF
193 *
194 * Solve for e(i+2) and d(i+3).
195 *
196  eir = dble( e( i+2 ) )
197  eii = dimag( e( i+2 ) )
198  f = eir / d( i+2 )
199  g = eii / d( i+2 )
200  e( i+2 ) = dcmplx( f, g )
201  d( i+3 ) = d( i+3 ) - f*eir - g*eii
202 *
203  IF( d( i+3 ).LE.zero ) THEN
204  info = i + 3
205  GO TO 30
206  END IF
207 *
208 * Solve for e(i+3) and d(i+4).
209 *
210  eir = dble( e( i+3 ) )
211  eii = dimag( e( i+3 ) )
212  f = eir / d( i+3 )
213  g = eii / d( i+3 )
214  e( i+3 ) = dcmplx( f, g )
215  d( i+4 ) = d( i+4 ) - f*eir - g*eii
216  20 CONTINUE
217 *
218 * Check d(n) for positive definiteness.
219 *
220  IF( d( n ).LE.zero )
221  $ info = n
222 *
223  30 CONTINUE
224  RETURN
225 *
226 * End of ZPTTRF
227 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zpttrs ( character  UPLO,
integer  N,
integer  NRHS,
double precision, dimension( * )  D,
complex*16, dimension( * )  E,
complex*16, dimension( ldb, * )  B,
integer  LDB,
integer  INFO 
)

ZPTTRS

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

Purpose:
 ZPTTRS solves a tridiagonal system of the form
    A * X = B
 using the factorization A = U**H *D* U or A = L*D*L**H computed by ZPTTRF.
 D is a diagonal matrix specified in the vector D, U (or L) is a unit
 bidiagonal matrix whose superdiagonal (subdiagonal) is specified in
 the vector E, and X and B are N by NRHS matrices.
Parameters
[in]UPLO
          UPLO is CHARACTER*1
          Specifies the form of the factorization and whether the
          vector E is the superdiagonal of the upper bidiagonal factor
          U or the subdiagonal of the lower bidiagonal factor L.
          = 'U':  A = U**H *D*U, E is the superdiagonal of U
          = 'L':  A = L*D*L**H, E is the subdiagonal of L
[in]N
          N is INTEGER
          The order of the tridiagonal matrix A.  N >= 0.
[in]NRHS
          NRHS is INTEGER
          The number of right hand sides, i.e., the number of columns
          of the matrix B.  NRHS >= 0.
[in]D
          D is DOUBLE PRECISION array, dimension (N)
          The n diagonal elements of the diagonal matrix D from the
          factorization A = U**H *D*U or A = L*D*L**H.
[in]E
          E is COMPLEX*16 array, dimension (N-1)
          If UPLO = 'U', the (n-1) superdiagonal elements of the unit
          bidiagonal factor U from the factorization A = U**H*D*U.
          If UPLO = 'L', the (n-1) subdiagonal elements of the unit
          bidiagonal factor L from the factorization A = L*D*L**H.
[in,out]B
          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
          On entry, the right hand side vectors B for the system of
          linear equations.
          On exit, the solution vectors, X.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
[out]INFO
          INFO is INTEGER
          = 0: successful exit
          < 0: if INFO = -k, the k-th argument had an illegal value
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 123 of file zpttrs.f.

123 *
124 * -- LAPACK computational routine (version 3.4.2) --
125 * -- LAPACK is a software package provided by Univ. of Tennessee, --
126 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
127 * September 2012
128 *
129 * .. Scalar Arguments ..
130  CHARACTER uplo
131  INTEGER info, ldb, n, nrhs
132 * ..
133 * .. Array Arguments ..
134  DOUBLE PRECISION d( * )
135  COMPLEX*16 b( ldb, * ), e( * )
136 * ..
137 *
138 * =====================================================================
139 *
140 * .. Local Scalars ..
141  LOGICAL upper
142  INTEGER iuplo, j, jb, nb
143 * ..
144 * .. External Functions ..
145  INTEGER ilaenv
146  EXTERNAL ilaenv
147 * ..
148 * .. External Subroutines ..
149  EXTERNAL xerbla, zptts2
150 * ..
151 * .. Intrinsic Functions ..
152  INTRINSIC max, min
153 * ..
154 * .. Executable Statements ..
155 *
156 * Test the input arguments.
157 *
158  info = 0
159  upper = ( uplo.EQ.'U' .OR. uplo.EQ.'u' )
160  IF( .NOT.upper .AND. .NOT.( uplo.EQ.'L' .OR. uplo.EQ.'l' ) ) THEN
161  info = -1
162  ELSE IF( n.LT.0 ) THEN
163  info = -2
164  ELSE IF( nrhs.LT.0 ) THEN
165  info = -3
166  ELSE IF( ldb.LT.max( 1, n ) ) THEN
167  info = -7
168  END IF
169  IF( info.NE.0 ) THEN
170  CALL xerbla( 'ZPTTRS', -info )
171  RETURN
172  END IF
173 *
174 * Quick return if possible
175 *
176  IF( n.EQ.0 .OR. nrhs.EQ.0 )
177  $ RETURN
178 *
179 * Determine the number of right-hand sides to solve at a time.
180 *
181  IF( nrhs.EQ.1 ) THEN
182  nb = 1
183  ELSE
184  nb = max( 1, ilaenv( 1, 'ZPTTRS', uplo, n, nrhs, -1, -1 ) )
185  END IF
186 *
187 * Decode UPLO
188 *
189  IF( upper ) THEN
190  iuplo = 1
191  ELSE
192  iuplo = 0
193  END IF
194 *
195  IF( nb.GE.nrhs ) THEN
196  CALL zptts2( iuplo, n, nrhs, d, e, b, ldb )
197  ELSE
198  DO 10 j = 1, nrhs, nb
199  jb = min( nrhs-j+1, nb )
200  CALL zptts2( iuplo, n, jb, d, e, b( 1, j ), ldb )
201  10 CONTINUE
202  END IF
203 *
204  RETURN
205 *
206 * End of ZPTTRS
207 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine zptts2(IUPLO, N, NRHS, D, E, B, LDB)
ZPTTS2 solves a tridiagonal system of the form AX=B using the L D LH factorization computed by spttrf...
Definition: zptts2.f:115
integer function ilaenv(ISPEC, NAME, OPTS, N1, N2, N3, N4)
Definition: tstiee.f:83

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine zptts2 ( integer  IUPLO,
integer  N,
integer  NRHS,
double precision, dimension( * )  D,
complex*16, dimension( * )  E,
complex*16, dimension( ldb, * )  B,
integer  LDB 
)

ZPTTS2 solves a tridiagonal system of the form AX=B using the L D LH factorization computed by spttrf.

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

Purpose:
 ZPTTS2 solves a tridiagonal system of the form
    A * X = B
 using the factorization A = U**H *D*U or A = L*D*L**H computed by ZPTTRF.
 D is a diagonal matrix specified in the vector D, U (or L) is a unit
 bidiagonal matrix whose superdiagonal (subdiagonal) is specified in
 the vector E, and X and B are N by NRHS matrices.
Parameters
[in]IUPLO
          IUPLO is INTEGER
          Specifies the form of the factorization and whether the
          vector E is the superdiagonal of the upper bidiagonal factor
          U or the subdiagonal of the lower bidiagonal factor L.
          = 1:  A = U**H *D*U, E is the superdiagonal of U
          = 0:  A = L*D*L**H, E is the subdiagonal of L
[in]N
          N is INTEGER
          The order of the tridiagonal matrix A.  N >= 0.
[in]NRHS
          NRHS is INTEGER
          The number of right hand sides, i.e., the number of columns
          of the matrix B.  NRHS >= 0.
[in]D
          D is DOUBLE PRECISION array, dimension (N)
          The n diagonal elements of the diagonal matrix D from the
          factorization A = U**H *D*U or A = L*D*L**H.
[in]E
          E is COMPLEX*16 array, dimension (N-1)
          If IUPLO = 1, the (n-1) superdiagonal elements of the unit
          bidiagonal factor U from the factorization A = U**H*D*U.
          If IUPLO = 0, the (n-1) subdiagonal elements of the unit
          bidiagonal factor L from the factorization A = L*D*L**H.
[in,out]B
          B is DOUBLE PRECISION array, dimension (LDB,NRHS)
          On entry, the right hand side vectors B for the system of
          linear equations.
          On exit, the solution vectors, X.
[in]LDB
          LDB is INTEGER
          The leading dimension of the array B.  LDB >= max(1,N).
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
September 2012

Definition at line 115 of file zptts2.f.

115 *
116 * -- LAPACK computational routine (version 3.4.2) --
117 * -- LAPACK is a software package provided by Univ. of Tennessee, --
118 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
119 * September 2012
120 *
121 * .. Scalar Arguments ..
122  INTEGER iuplo, ldb, n, nrhs
123 * ..
124 * .. Array Arguments ..
125  DOUBLE PRECISION d( * )
126  COMPLEX*16 b( ldb, * ), e( * )
127 * ..
128 *
129 * =====================================================================
130 *
131 * .. Local Scalars ..
132  INTEGER i, j
133 * ..
134 * .. External Subroutines ..
135  EXTERNAL zdscal
136 * ..
137 * .. Intrinsic Functions ..
138  INTRINSIC dconjg
139 * ..
140 * .. Executable Statements ..
141 *
142 * Quick return if possible
143 *
144  IF( n.LE.1 ) THEN
145  IF( n.EQ.1 )
146  $ CALL zdscal( nrhs, 1.d0 / d( 1 ), b, ldb )
147  RETURN
148  END IF
149 *
150  IF( iuplo.EQ.1 ) THEN
151 *
152 * Solve A * X = B using the factorization A = U**H *D*U,
153 * overwriting each right hand side vector with its solution.
154 *
155  IF( nrhs.LE.2 ) THEN
156  j = 1
157  10 CONTINUE
158 *
159 * Solve U**H * x = b.
160 *
161  DO 20 i = 2, n
162  b( i, j ) = b( i, j ) - b( i-1, j )*dconjg( e( i-1 ) )
163  20 CONTINUE
164 *
165 * Solve D * U * x = b.
166 *
167  DO 30 i = 1, n
168  b( i, j ) = b( i, j ) / d( i )
169  30 CONTINUE
170  DO 40 i = n - 1, 1, -1
171  b( i, j ) = b( i, j ) - b( i+1, j )*e( i )
172  40 CONTINUE
173  IF( j.LT.nrhs ) THEN
174  j = j + 1
175  GO TO 10
176  END IF
177  ELSE
178  DO 70 j = 1, nrhs
179 *
180 * Solve U**H * x = b.
181 *
182  DO 50 i = 2, n
183  b( i, j ) = b( i, j ) - b( i-1, j )*dconjg( e( i-1 ) )
184  50 CONTINUE
185 *
186 * Solve D * U * x = b.
187 *
188  b( n, j ) = b( n, j ) / d( n )
189  DO 60 i = n - 1, 1, -1
190  b( i, j ) = b( i, j ) / d( i ) - b( i+1, j )*e( i )
191  60 CONTINUE
192  70 CONTINUE
193  END IF
194  ELSE
195 *
196 * Solve A * X = B using the factorization A = L*D*L**H,
197 * overwriting each right hand side vector with its solution.
198 *
199  IF( nrhs.LE.2 ) THEN
200  j = 1
201  80 CONTINUE
202 *
203 * Solve L * x = b.
204 *
205  DO 90 i = 2, n
206  b( i, j ) = b( i, j ) - b( i-1, j )*e( i-1 )
207  90 CONTINUE
208 *
209 * Solve D * L**H * x = b.
210 *
211  DO 100 i = 1, n
212  b( i, j ) = b( i, j ) / d( i )
213  100 CONTINUE
214  DO 110 i = n - 1, 1, -1
215  b( i, j ) = b( i, j ) - b( i+1, j )*dconjg( e( i ) )
216  110 CONTINUE
217  IF( j.LT.nrhs ) THEN
218  j = j + 1
219  GO TO 80
220  END IF
221  ELSE
222  DO 140 j = 1, nrhs
223 *
224 * Solve L * x = b.
225 *
226  DO 120 i = 2, n
227  b( i, j ) = b( i, j ) - b( i-1, j )*e( i-1 )
228  120 CONTINUE
229 *
230 * Solve D * L**H * x = b.
231 *
232  b( n, j ) = b( n, j ) / d( n )
233  DO 130 i = n - 1, 1, -1
234  b( i, j ) = b( i, j ) / d( i ) -
235  $ b( i+1, j )*dconjg( e( i ) )
236  130 CONTINUE
237  140 CONTINUE
238  END IF
239  END IF
240 *
241  RETURN
242 *
243 * End of ZPTTS2
244 *
subroutine zdscal(N, DA, ZX, INCX)
ZDSCAL
Definition: zdscal.f:54

Here is the call graph for this function:

Here is the caller graph for this function: