OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches

Functions

subroutine dptcon (n, d, e, anorm, rcond, work, info)
 DPTCON
subroutine dpteqr (compz, n, d, e, z, ldz, work, info)
 DPTEQR
subroutine dptrfs (n, nrhs, d, e, df, ef, b, ldb, x, ldx, ferr, berr, work, info)
 DPTRFS
subroutine dpttrf (n, d, e, info)
 DPTTRF
subroutine dpttrs (n, nrhs, d, e, b, ldb, info)
 DPTTRS
subroutine dptts2 (n, nrhs, d, e, b, ldb)
 DPTTS2 solves a tridiagonal system of the form AX=B using the L D LH factorization computed by spttrf.

Detailed Description

This is the group of double computational functions for PT matrices

Function Documentation

◆ dptcon()

subroutine dptcon ( integer n,
double precision, dimension( * ) d,
double precision, dimension( * ) e,
double precision anorm,
double precision rcond,
double precision, dimension( * ) work,
integer info )

DPTCON

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

Purpose:
!>
!> DPTCON computes the reciprocal of the condition number (in the
!> 1-norm) of a real symmetric positive definite tridiagonal matrix
!> using the factorization A = L*D*L**T or A = U**T*D*U computed by
!> DPTTRF.
!>
!> 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 DPTTRF.
!> 
[in]E
!>          E is DOUBLE PRECISION 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 DPTTRF.
!> 
[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]WORK
!>          WORK 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.
Further Details:
!>
!>  The method used is described in Nicholas J. Higham, , SIAM J. Sci. Stat. Comput., Vol. 7, No. 1, January 1986.
!> 

Definition at line 117 of file dptcon.f.

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

◆ dpteqr()

subroutine dpteqr ( character compz,
integer n,
double precision, dimension( * ) d,
double precision, dimension( * ) e,
double precision, dimension( ldz, * ) z,
integer ldz,
double precision, dimension( * ) work,
integer info )

DPTEQR

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

Purpose:
!>
!> DPTEQR computes all eigenvalues and, optionally, eigenvectors of a
!> symmetric positive definite tridiagonal matrix by first factoring the
!> matrix using DPTTRF, and then calling DBDSQR 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 symmetric positive definite matrix
!> can also be found if DSYTRD, DSPTRD, or DSBTRD 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 symmetric
!>                  matrix also.  Array Z contains the orthogonal
!>                  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 DOUBLE PRECISION array, dimension (LDZ, N)
!>          On entry, if COMPZ = 'V', the orthogonal matrix used in the
!>          reduction to tridiagonal form.
!>          On exit, if COMPZ = 'V', the orthonormal eigenvectors of the
!>          original symmetric 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.

Definition at line 144 of file dpteqr.f.

145*
146* -- LAPACK computational routine --
147* -- LAPACK is a software package provided by Univ. of Tennessee, --
148* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
149*
150* .. Scalar Arguments ..
151 CHARACTER COMPZ
152 INTEGER INFO, LDZ, N
153* ..
154* .. Array Arguments ..
155 DOUBLE PRECISION D( * ), E( * ), WORK( * ), Z( LDZ, * )
156* ..
157*
158* =====================================================================
159*
160* .. Parameters ..
161 DOUBLE PRECISION ZERO, ONE
162 parameter( zero = 0.0d0, one = 1.0d0 )
163* ..
164* .. External Functions ..
165 LOGICAL LSAME
166 EXTERNAL lsame
167* ..
168* .. External Subroutines ..
169 EXTERNAL dbdsqr, dlaset, dpttrf, xerbla
170* ..
171* .. Local Arrays ..
172 DOUBLE PRECISION C( 1, 1 ), VT( 1, 1 )
173* ..
174* .. Local Scalars ..
175 INTEGER I, ICOMPZ, NRU
176* ..
177* .. Intrinsic Functions ..
178 INTRINSIC max, sqrt
179* ..
180* .. Executable Statements ..
181*
182* Test the input parameters.
183*
184 info = 0
185*
186 IF( lsame( compz, 'N' ) ) THEN
187 icompz = 0
188 ELSE IF( lsame( compz, 'V' ) ) THEN
189 icompz = 1
190 ELSE IF( lsame( compz, 'I' ) ) THEN
191 icompz = 2
192 ELSE
193 icompz = -1
194 END IF
195 IF( icompz.LT.0 ) THEN
196 info = -1
197 ELSE IF( n.LT.0 ) THEN
198 info = -2
199 ELSE IF( ( ldz.LT.1 ) .OR. ( icompz.GT.0 .AND. ldz.LT.max( 1,
200 $ n ) ) ) THEN
201 info = -6
202 END IF
203 IF( info.NE.0 ) THEN
204 CALL xerbla( 'DPTEQR', -info )
205 RETURN
206 END IF
207*
208* Quick return if possible
209*
210 IF( n.EQ.0 )
211 $ RETURN
212*
213 IF( n.EQ.1 ) THEN
214 IF( icompz.GT.0 )
215 $ z( 1, 1 ) = one
216 RETURN
217 END IF
218 IF( icompz.EQ.2 )
219 $ CALL dlaset( 'Full', n, n, zero, one, z, ldz )
220*
221* Call DPTTRF to factor the matrix.
222*
223 CALL dpttrf( n, d, e, info )
224 IF( info.NE.0 )
225 $ RETURN
226 DO 10 i = 1, n
227 d( i ) = sqrt( d( i ) )
228 10 CONTINUE
229 DO 20 i = 1, n - 1
230 e( i ) = e( i )*d( i )
231 20 CONTINUE
232*
233* Call DBDSQR to compute the singular values/vectors of the
234* bidiagonal factor.
235*
236 IF( icompz.GT.0 ) THEN
237 nru = n
238 ELSE
239 nru = 0
240 END IF
241 CALL dbdsqr( 'Lower', n, 0, nru, 0, d, e, vt, 1, z, ldz, c, 1,
242 $ work, info )
243*
244* Square the singular values.
245*
246 IF( info.EQ.0 ) THEN
247 DO 30 i = 1, n
248 d( i ) = d( i )*d( i )
249 30 CONTINUE
250 ELSE
251 info = n + info
252 END IF
253*
254 RETURN
255*
256* End of DPTEQR
257*
subroutine dlaset(uplo, m, n, alpha, beta, a, lda)
DLASET initializes the off-diagonal elements and the diagonal elements of a matrix to given values.
Definition dlaset.f:110
subroutine dbdsqr(uplo, n, ncvt, nru, ncc, d, e, vt, ldvt, u, ldu, c, ldc, work, info)
DBDSQR
Definition dbdsqr.f:241
logical function lsame(ca, cb)
LSAME
Definition lsame.f:53
subroutine dpttrf(n, d, e, info)
DPTTRF
Definition dpttrf.f:91
#define max(a, b)
Definition macros.h:21

◆ dptrfs()

subroutine dptrfs ( integer n,
integer nrhs,
double precision, dimension( * ) d,
double precision, dimension( * ) e,
double precision, dimension( * ) df,
double precision, dimension( * ) ef,
double precision, dimension( ldb, * ) b,
integer ldb,
double precision, dimension( ldx, * ) x,
integer ldx,
double precision, dimension( * ) ferr,
double precision, dimension( * ) berr,
double precision, dimension( * ) work,
integer info )

DPTRFS

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

Purpose:
!>
!> DPTRFS improves the computed solution to a system of linear
!> equations when the coefficient matrix is symmetric positive definite
!> and tridiagonal, and provides error bounds and backward error
!> estimates for the solution.
!> 
Parameters
[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 diagonal elements of the tridiagonal matrix A.
!> 
[in]E
!>          E is DOUBLE PRECISION array, dimension (N-1)
!>          The (n-1) subdiagonal elements of the tridiagonal matrix A.
!> 
[in]DF
!>          DF is DOUBLE PRECISION array, dimension (N)
!>          The n diagonal elements of the diagonal matrix D from the
!>          factorization computed by DPTTRF.
!> 
[in]EF
!>          EF is DOUBLE PRECISION array, dimension (N-1)
!>          The (n-1) subdiagonal elements of the unit bidiagonal factor
!>          L from the factorization computed by DPTTRF.
!> 
[in]B
!>          B is DOUBLE PRECISION 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 DOUBLE PRECISION array, dimension (LDX,NRHS)
!>          On entry, the solution matrix X, as computed by DPTTRS.
!>          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 DOUBLE PRECISION array, dimension (2*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.

Definition at line 161 of file dptrfs.f.

163*
164* -- LAPACK computational routine --
165* -- LAPACK is a software package provided by Univ. of Tennessee, --
166* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
167*
168* .. Scalar Arguments ..
169 INTEGER INFO, LDB, LDX, N, NRHS
170* ..
171* .. Array Arguments ..
172 DOUBLE PRECISION B( LDB, * ), BERR( * ), D( * ), DF( * ),
173 $ E( * ), EF( * ), FERR( * ), WORK( * ),
174 $ X( LDX, * )
175* ..
176*
177* =====================================================================
178*
179* .. Parameters ..
180 INTEGER ITMAX
181 parameter( itmax = 5 )
182 DOUBLE PRECISION ZERO
183 parameter( zero = 0.0d+0 )
184 DOUBLE PRECISION ONE
185 parameter( one = 1.0d+0 )
186 DOUBLE PRECISION TWO
187 parameter( two = 2.0d+0 )
188 DOUBLE PRECISION THREE
189 parameter( three = 3.0d+0 )
190* ..
191* .. Local Scalars ..
192 INTEGER COUNT, I, IX, J, NZ
193 DOUBLE PRECISION BI, CX, DX, EPS, EX, LSTRES, S, SAFE1, SAFE2,
194 $ SAFMIN
195* ..
196* .. External Subroutines ..
197 EXTERNAL daxpy, dpttrs, xerbla
198* ..
199* .. Intrinsic Functions ..
200 INTRINSIC abs, max
201* ..
202* .. External Functions ..
203 INTEGER IDAMAX
204 DOUBLE PRECISION DLAMCH
205 EXTERNAL idamax, dlamch
206* ..
207* .. Executable Statements ..
208*
209* Test the input parameters.
210*
211 info = 0
212 IF( n.LT.0 ) THEN
213 info = -1
214 ELSE IF( nrhs.LT.0 ) THEN
215 info = -2
216 ELSE IF( ldb.LT.max( 1, n ) ) THEN
217 info = -8
218 ELSE IF( ldx.LT.max( 1, n ) ) THEN
219 info = -10
220 END IF
221 IF( info.NE.0 ) THEN
222 CALL xerbla( 'DPTRFS', -info )
223 RETURN
224 END IF
225*
226* Quick return if possible
227*
228 IF( n.EQ.0 .OR. nrhs.EQ.0 ) THEN
229 DO 10 j = 1, nrhs
230 ferr( j ) = zero
231 berr( j ) = zero
232 10 CONTINUE
233 RETURN
234 END IF
235*
236* NZ = maximum number of nonzero elements in each row of A, plus 1
237*
238 nz = 4
239 eps = dlamch( 'Epsilon' )
240 safmin = dlamch( 'Safe minimum' )
241 safe1 = nz*safmin
242 safe2 = safe1 / eps
243*
244* Do for each right hand side
245*
246 DO 90 j = 1, nrhs
247*
248 count = 1
249 lstres = three
250 20 CONTINUE
251*
252* Loop until stopping criterion is satisfied.
253*
254* Compute residual R = B - A * X. Also compute
255* abs(A)*abs(x) + abs(b) for use in the backward error bound.
256*
257 IF( n.EQ.1 ) THEN
258 bi = b( 1, j )
259 dx = d( 1 )*x( 1, j )
260 work( n+1 ) = bi - dx
261 work( 1 ) = abs( bi ) + abs( dx )
262 ELSE
263 bi = b( 1, j )
264 dx = d( 1 )*x( 1, j )
265 ex = e( 1 )*x( 2, j )
266 work( n+1 ) = bi - dx - ex
267 work( 1 ) = abs( bi ) + abs( dx ) + abs( ex )
268 DO 30 i = 2, n - 1
269 bi = b( i, j )
270 cx = e( i-1 )*x( i-1, j )
271 dx = d( i )*x( i, j )
272 ex = e( i )*x( i+1, j )
273 work( n+i ) = bi - cx - dx - ex
274 work( i ) = abs( bi ) + abs( cx ) + abs( dx ) + abs( ex )
275 30 CONTINUE
276 bi = b( n, j )
277 cx = e( n-1 )*x( n-1, j )
278 dx = d( n )*x( n, j )
279 work( n+n ) = bi - cx - dx
280 work( n ) = abs( bi ) + abs( cx ) + abs( dx )
281 END IF
282*
283* Compute componentwise relative backward error from formula
284*
285* max(i) ( abs(R(i)) / ( abs(A)*abs(X) + abs(B) )(i) )
286*
287* where abs(Z) is the componentwise absolute value of the matrix
288* or vector Z. If the i-th component of the denominator is less
289* than SAFE2, then SAFE1 is added to the i-th components of the
290* numerator and denominator before dividing.
291*
292 s = zero
293 DO 40 i = 1, n
294 IF( work( i ).GT.safe2 ) THEN
295 s = max( s, abs( work( n+i ) ) / work( i ) )
296 ELSE
297 s = max( s, ( abs( work( n+i ) )+safe1 ) /
298 $ ( work( i )+safe1 ) )
299 END IF
300 40 CONTINUE
301 berr( j ) = s
302*
303* Test stopping criterion. Continue iterating if
304* 1) The residual BERR(J) is larger than machine epsilon, and
305* 2) BERR(J) decreased by at least a factor of 2 during the
306* last iteration, and
307* 3) At most ITMAX iterations tried.
308*
309 IF( berr( j ).GT.eps .AND. two*berr( j ).LE.lstres .AND.
310 $ count.LE.itmax ) THEN
311*
312* Update solution and try again.
313*
314 CALL dpttrs( n, 1, df, ef, work( n+1 ), n, info )
315 CALL daxpy( n, one, work( n+1 ), 1, x( 1, j ), 1 )
316 lstres = berr( j )
317 count = count + 1
318 GO TO 20
319 END IF
320*
321* Bound error from formula
322*
323* norm(X - XTRUE) / norm(X) .le. FERR =
324* norm( abs(inv(A))*
325* ( abs(R) + NZ*EPS*( abs(A)*abs(X)+abs(B) ))) / norm(X)
326*
327* where
328* norm(Z) is the magnitude of the largest component of Z
329* inv(A) is the inverse of A
330* abs(Z) is the componentwise absolute value of the matrix or
331* vector Z
332* NZ is the maximum number of nonzeros in any row of A, plus 1
333* EPS is machine epsilon
334*
335* The i-th component of abs(R)+NZ*EPS*(abs(A)*abs(X)+abs(B))
336* is incremented by SAFE1 if the i-th component of
337* abs(A)*abs(X) + abs(B) is less than SAFE2.
338*
339 DO 50 i = 1, n
340 IF( work( i ).GT.safe2 ) THEN
341 work( i ) = abs( work( n+i ) ) + nz*eps*work( i )
342 ELSE
343 work( i ) = abs( work( n+i ) ) + nz*eps*work( i ) + safe1
344 END IF
345 50 CONTINUE
346 ix = idamax( n, work, 1 )
347 ferr( j ) = work( ix )
348*
349* Estimate the norm of inv(A).
350*
351* Solve M(A) * x = e, where M(A) = (m(i,j)) is given by
352*
353* m(i,j) = abs(A(i,j)), i = j,
354* m(i,j) = -abs(A(i,j)), i .ne. j,
355*
356* and e = [ 1, 1, ..., 1 ]**T. Note M(A) = M(L)*D*M(L)**T.
357*
358* Solve M(L) * x = e.
359*
360 work( 1 ) = one
361 DO 60 i = 2, n
362 work( i ) = one + work( i-1 )*abs( ef( i-1 ) )
363 60 CONTINUE
364*
365* Solve D * M(L)**T * x = b.
366*
367 work( n ) = work( n ) / df( n )
368 DO 70 i = n - 1, 1, -1
369 work( i ) = work( i ) / df( i ) + work( i+1 )*abs( ef( i ) )
370 70 CONTINUE
371*
372* Compute norm(inv(A)) = max(x(i)), 1<=i<=n.
373*
374 ix = idamax( n, work, 1 )
375 ferr( j ) = ferr( j )*abs( work( ix ) )
376*
377* Normalize error.
378*
379 lstres = zero
380 DO 80 i = 1, n
381 lstres = max( lstres, abs( x( i, j ) ) )
382 80 CONTINUE
383 IF( lstres.NE.zero )
384 $ ferr( j ) = ferr( j ) / lstres
385*
386 90 CONTINUE
387*
388 RETURN
389*
390* End of DPTRFS
391*
subroutine dpttrs(n, nrhs, d, e, b, ldb, info)
DPTTRS
Definition dpttrs.f:109
subroutine daxpy(n, da, dx, incx, dy, incy)
DAXPY
Definition daxpy.f:89
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69

◆ dpttrf()

subroutine dpttrf ( integer n,
double precision, dimension( * ) d,
double precision, dimension( * ) e,
integer info )

DPTTRF

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

Purpose:
!>
!> DPTTRF computes the L*D*L**T factorization of a real symmetric
!> positive definite tridiagonal matrix A.  The factorization may also
!> be regarded as having the form A = U**T*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**T factorization of A.
!> 
[in,out]E
!>          E is DOUBLE PRECISION 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**T factorization of A.
!>          E can also be regarded as the superdiagonal of the unit
!>          bidiagonal factor U from the U**T*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.

Definition at line 90 of file dpttrf.f.

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

◆ dpttrs()

subroutine dpttrs ( integer n,
integer nrhs,
double precision, dimension( * ) d,
double precision, dimension( * ) e,
double precision, dimension( ldb, * ) b,
integer ldb,
integer info )

DPTTRS

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

Purpose:
!>
!> DPTTRS solves a tridiagonal system of the form
!>    A * X = B
!> using the L*D*L**T factorization of A computed by DPTTRF.  D is a
!> diagonal matrix specified in the vector D, L is a unit bidiagonal
!> matrix whose subdiagonal is specified in the vector E, and X and B
!> are N by NRHS matrices.
!> 
Parameters
[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
!>          L*D*L**T factorization of A.
!> 
[in]E
!>          E is DOUBLE PRECISION array, dimension (N-1)
!>          The (n-1) subdiagonal elements of the unit bidiagonal factor
!>          L from the L*D*L**T factorization of A.  E can also be regarded
!>          as the superdiagonal of the unit bidiagonal factor U from the
!>          factorization A = U**T*D*U.
!> 
[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.

Definition at line 108 of file dpttrs.f.

109*
110* -- LAPACK computational routine --
111* -- LAPACK is a software package provided by Univ. of Tennessee, --
112* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
113*
114* .. Scalar Arguments ..
115 INTEGER INFO, LDB, N, NRHS
116* ..
117* .. Array Arguments ..
118 DOUBLE PRECISION B( LDB, * ), D( * ), E( * )
119* ..
120*
121* =====================================================================
122*
123* .. Local Scalars ..
124 INTEGER J, JB, NB
125* ..
126* .. External Functions ..
127 INTEGER ILAENV
128 EXTERNAL ilaenv
129* ..
130* .. External Subroutines ..
131 EXTERNAL dptts2, xerbla
132* ..
133* .. Intrinsic Functions ..
134 INTRINSIC max, min
135* ..
136* .. Executable Statements ..
137*
138* Test the input arguments.
139*
140 info = 0
141 IF( n.LT.0 ) THEN
142 info = -1
143 ELSE IF( nrhs.LT.0 ) THEN
144 info = -2
145 ELSE IF( ldb.LT.max( 1, n ) ) THEN
146 info = -6
147 END IF
148 IF( info.NE.0 ) THEN
149 CALL xerbla( 'DPTTRS', -info )
150 RETURN
151 END IF
152*
153* Quick return if possible
154*
155 IF( n.EQ.0 .OR. nrhs.EQ.0 )
156 $ RETURN
157*
158* Determine the number of right-hand sides to solve at a time.
159*
160 IF( nrhs.EQ.1 ) THEN
161 nb = 1
162 ELSE
163 nb = max( 1, ilaenv( 1, 'DPTTRS', ' ', n, nrhs, -1, -1 ) )
164 END IF
165*
166 IF( nb.GE.nrhs ) THEN
167 CALL dptts2( n, nrhs, d, e, b, ldb )
168 ELSE
169 DO 10 j = 1, nrhs, nb
170 jb = min( nrhs-j+1, nb )
171 CALL dptts2( n, jb, d, e, b( 1, j ), ldb )
172 10 CONTINUE
173 END IF
174*
175 RETURN
176*
177* End of DPTTRS
178*
integer function ilaenv(ispec, name, opts, n1, n2, n3, n4)
ILAENV
Definition ilaenv.f:162
subroutine dptts2(n, nrhs, d, e, b, ldb)
DPTTS2 solves a tridiagonal system of the form AX=B using the L D LH factorization computed by spttrf...
Definition dptts2.f:102
#define min(a, b)
Definition macros.h:20

◆ dptts2()

subroutine dptts2 ( integer n,
integer nrhs,
double precision, dimension( * ) d,
double precision, dimension( * ) e,
double precision, dimension( ldb, * ) b,
integer ldb )

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

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

Purpose:
!>
!> DPTTS2 solves a tridiagonal system of the form
!>    A * X = B
!> using the L*D*L**T factorization of A computed by DPTTRF.  D is a
!> diagonal matrix specified in the vector D, L is a unit bidiagonal
!> matrix whose subdiagonal is specified in the vector E, and X and B
!> are N by NRHS matrices.
!> 
Parameters
[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
!>          L*D*L**T factorization of A.
!> 
[in]E
!>          E is DOUBLE PRECISION array, dimension (N-1)
!>          The (n-1) subdiagonal elements of the unit bidiagonal factor
!>          L from the L*D*L**T factorization of A.  E can also be regarded
!>          as the superdiagonal of the unit bidiagonal factor U from the
!>          factorization A = U**T*D*U.
!> 
[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.

Definition at line 101 of file dptts2.f.

102*
103* -- LAPACK computational routine --
104* -- LAPACK is a software package provided by Univ. of Tennessee, --
105* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
106*
107* .. Scalar Arguments ..
108 INTEGER LDB, N, NRHS
109* ..
110* .. Array Arguments ..
111 DOUBLE PRECISION B( LDB, * ), D( * ), E( * )
112* ..
113*
114* =====================================================================
115*
116* .. Local Scalars ..
117 INTEGER I, J
118* ..
119* .. External Subroutines ..
120 EXTERNAL dscal
121* ..
122* .. Executable Statements ..
123*
124* Quick return if possible
125*
126 IF( n.LE.1 ) THEN
127 IF( n.EQ.1 )
128 $ CALL dscal( nrhs, 1.d0 / d( 1 ), b, ldb )
129 RETURN
130 END IF
131*
132* Solve A * X = B using the factorization A = L*D*L**T,
133* overwriting each right hand side vector with its solution.
134*
135 DO 30 j = 1, nrhs
136*
137* Solve L * x = b.
138*
139 DO 10 i = 2, n
140 b( i, j ) = b( i, j ) - b( i-1, j )*e( i-1 )
141 10 CONTINUE
142*
143* Solve D * L**T * x = b.
144*
145 b( n, j ) = b( n, j ) / d( n )
146 DO 20 i = n - 1, 1, -1
147 b( i, j ) = b( i, j ) / d( i ) - b( i+1, j )*e( i )
148 20 CONTINUE
149 30 CONTINUE
150*
151 RETURN
152*
153* End of DPTTS2
154*
subroutine dscal(n, da, dx, incx)
DSCAL
Definition dscal.f:79