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

Functions

subroutine zptsv (n, nrhs, d, e, b, ldb, info)
  ZPTSV computes the solution to system of linear equations A * X = B for PT matrices
subroutine zptsvx (fact, n, nrhs, d, e, df, ef, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
  ZPTSVX computes the solution to system of linear equations A * X = B for PT matrices

Detailed Description

This is the group of complex16 solve driver functions for PT matrices

Function Documentation

◆ zptsv()

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

ZPTSV computes the solution to system of linear equations A * X = B for PT matrices

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

Purpose:
!>
!> ZPTSV computes the solution to a complex system of linear equations
!> A*X = B, where A is an N-by-N Hermitian positive definite tridiagonal
!> matrix, and X and B are N-by-NRHS matrices.
!>
!> A is factored as A = L*D*L**H, and the factored form of A is then
!> used to solve the system of equations.
!> 
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,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 factorization A = L*D*L**H.
!> 
[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.
!> 
[in,out]B
!>          B is COMPLEX*16 array, dimension (LDB,NRHS)
!>          On entry, the N-by-NRHS right hand side matrix B.
!>          On exit, if INFO = 0, the N-by-NRHS solution matrix 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 = -i, the i-th argument had an illegal value
!>          > 0:  if INFO = i, the leading minor of order i is not
!>                positive definite, and the solution has not been
!>                computed.  The factorization has not been completed
!>                unless i = N.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 114 of file zptsv.f.

115*
116* -- LAPACK driver routine --
117* -- LAPACK is a software package provided by Univ. of Tennessee, --
118* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
119*
120* .. Scalar Arguments ..
121 INTEGER INFO, LDB, N, NRHS
122* ..
123* .. Array Arguments ..
124 DOUBLE PRECISION D( * )
125 COMPLEX*16 B( LDB, * ), E( * )
126* ..
127*
128* =====================================================================
129*
130* .. External Subroutines ..
131 EXTERNAL xerbla, zpttrf, zpttrs
132* ..
133* .. Intrinsic Functions ..
134 INTRINSIC max
135* ..
136* .. Executable Statements ..
137*
138* Test the input parameters.
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( 'ZPTSV ', -info )
150 RETURN
151 END IF
152*
153* Compute the L*D*L**H (or U**H*D*U) factorization of A.
154*
155 CALL zpttrf( n, d, e, info )
156 IF( info.EQ.0 ) THEN
157*
158* Solve the system A*X = B, overwriting B with X.
159*
160 CALL zpttrs( 'Lower', n, nrhs, d, e, b, ldb, info )
161 END IF
162 RETURN
163*
164* End of ZPTSV
165*
subroutine xerbla(srname, info)
XERBLA
Definition xerbla.f:60
subroutine zpttrf(n, d, e, info)
ZPTTRF
Definition zpttrf.f:92
subroutine zpttrs(uplo, n, nrhs, d, e, b, ldb, info)
ZPTTRS
Definition zpttrs.f:121
#define max(a, b)
Definition macros.h:21

◆ zptsvx()

subroutine zptsvx ( character fact,
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 rcond,
double precision, dimension( * ) ferr,
double precision, dimension( * ) berr,
complex*16, dimension( * ) work,
double precision, dimension( * ) rwork,
integer info )

ZPTSVX computes the solution to system of linear equations A * X = B for PT matrices

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

Purpose:
!>
!> ZPTSVX uses the factorization A = L*D*L**H to compute the solution
!> to a complex system of linear equations A*X = B, where A is an
!> N-by-N Hermitian positive definite tridiagonal matrix and X and B
!> are N-by-NRHS matrices.
!>
!> Error bounds on the solution and a condition estimate are also
!> provided.
!> 
Description:
!>
!> The following steps are performed:
!>
!> 1. If FACT = 'N', the matrix A is factored as A = L*D*L**H, where L
!>    is a unit lower bidiagonal matrix and D is diagonal.  The
!>    factorization can also be regarded as having the form
!>    A = U**H*D*U.
!>
!> 2. If the leading i-by-i principal minor is not positive definite,
!>    then the routine returns with INFO = i. Otherwise, the factored
!>    form of A is used to estimate the condition number of the matrix
!>    A.  If the reciprocal of the condition number is less than machine
!>    precision, INFO = N+1 is returned as a warning, but the routine
!>    still goes on to solve for X and compute error bounds as
!>    described below.
!>
!> 3. The system of equations is solved for X using the factored form
!>    of A.
!>
!> 4. Iterative refinement is applied to improve the computed solution
!>    matrix and calculate error bounds and backward error estimates
!>    for it.
!> 
Parameters
[in]FACT
!>          FACT is CHARACTER*1
!>          Specifies whether or not the factored form of the matrix
!>          A is supplied on entry.
!>          = 'F':  On entry, DF and EF contain the factored form of A.
!>                  D, E, DF, and EF will not be modified.
!>          = 'N':  The matrix A will be copied to DF and EF and
!>                  factored.
!> 
[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 matrices B and X.  NRHS >= 0.
!> 
[in]D
!>          D is DOUBLE PRECISION array, dimension (N)
!>          The n diagonal elements of the tridiagonal matrix A.
!> 
[in]E
!>          E is COMPLEX*16 array, dimension (N-1)
!>          The (n-1) subdiagonal elements of the tridiagonal matrix A.
!> 
[in,out]DF
!>          DF is DOUBLE PRECISION array, dimension (N)
!>          If FACT = 'F', then DF is an input argument and on entry
!>          contains the n diagonal elements of the diagonal matrix D
!>          from the L*D*L**H factorization of A.
!>          If FACT = 'N', then DF is an output argument and on exit
!>          contains the n diagonal elements of the diagonal matrix D
!>          from the L*D*L**H factorization of A.
!> 
[in,out]EF
!>          EF is COMPLEX*16 array, dimension (N-1)
!>          If FACT = 'F', then EF is an input argument and on entry
!>          contains the (n-1) subdiagonal elements of the unit
!>          bidiagonal factor L from the L*D*L**H factorization of A.
!>          If FACT = 'N', then EF is an output argument and on exit
!>          contains the (n-1) subdiagonal elements of the unit
!>          bidiagonal factor L from the L*D*L**H factorization of A.
!> 
[in]B
!>          B is COMPLEX*16 array, dimension (LDB,NRHS)
!>          The N-by-NRHS right hand side matrix B.
!> 
[in]LDB
!>          LDB is INTEGER
!>          The leading dimension of the array B.  LDB >= max(1,N).
!> 
[out]X
!>          X is COMPLEX*16 array, dimension (LDX,NRHS)
!>          If INFO = 0 or INFO = N+1, the N-by-NRHS solution matrix X.
!> 
[in]LDX
!>          LDX is INTEGER
!>          The leading dimension of the array X.  LDX >= max(1,N).
!> 
[out]RCOND
!>          RCOND is DOUBLE PRECISION
!>          The reciprocal condition number of the matrix A.  If RCOND
!>          is less than the machine precision (in particular, if
!>          RCOND = 0), the matrix is singular to working precision.
!>          This condition is indicated by a return code of INFO > 0.
!> 
[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
!>          > 0:  if INFO = i, and i is
!>                <= N:  the leading minor of order i of A is
!>                       not positive definite, so the factorization
!>                       could not be completed, and the solution has not
!>                       been computed. RCOND = 0 is returned.
!>                = N+1: U is nonsingular, but RCOND is less than machine
!>                       precision, meaning that the matrix is singular
!>                       to working precision.  Nevertheless, the
!>                       solution and error bounds are computed because
!>                       there are a number of situations where the
!>                       computed solution can be more accurate than the
!>                       value of RCOND would suggest.
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.

Definition at line 232 of file zptsvx.f.

234*
235* -- LAPACK driver routine --
236* -- LAPACK is a software package provided by Univ. of Tennessee, --
237* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
238*
239* .. Scalar Arguments ..
240 CHARACTER FACT
241 INTEGER INFO, LDB, LDX, N, NRHS
242 DOUBLE PRECISION RCOND
243* ..
244* .. Array Arguments ..
245 DOUBLE PRECISION BERR( * ), D( * ), DF( * ), FERR( * ),
246 $ RWORK( * )
247 COMPLEX*16 B( LDB, * ), E( * ), EF( * ), WORK( * ),
248 $ X( LDX, * )
249* ..
250*
251* =====================================================================
252*
253* .. Parameters ..
254 DOUBLE PRECISION ZERO
255 parameter( zero = 0.0d+0 )
256* ..
257* .. Local Scalars ..
258 LOGICAL NOFACT
259 DOUBLE PRECISION ANORM
260* ..
261* .. External Functions ..
262 LOGICAL LSAME
263 DOUBLE PRECISION DLAMCH, ZLANHT
264 EXTERNAL lsame, dlamch, zlanht
265* ..
266* .. External Subroutines ..
267 EXTERNAL dcopy, xerbla, zcopy, zlacpy, zptcon, zptrfs,
268 $ zpttrf, zpttrs
269* ..
270* .. Intrinsic Functions ..
271 INTRINSIC max
272* ..
273* .. Executable Statements ..
274*
275* Test the input parameters.
276*
277 info = 0
278 nofact = lsame( fact, 'N' )
279 IF( .NOT.nofact .AND. .NOT.lsame( fact, 'F' ) ) THEN
280 info = -1
281 ELSE IF( n.LT.0 ) THEN
282 info = -2
283 ELSE IF( nrhs.LT.0 ) THEN
284 info = -3
285 ELSE IF( ldb.LT.max( 1, n ) ) THEN
286 info = -9
287 ELSE IF( ldx.LT.max( 1, n ) ) THEN
288 info = -11
289 END IF
290 IF( info.NE.0 ) THEN
291 CALL xerbla( 'ZPTSVX', -info )
292 RETURN
293 END IF
294*
295 IF( nofact ) THEN
296*
297* Compute the L*D*L**H (or U**H*D*U) factorization of A.
298*
299 CALL dcopy( n, d, 1, df, 1 )
300 IF( n.GT.1 )
301 $ CALL zcopy( n-1, e, 1, ef, 1 )
302 CALL zpttrf( n, df, ef, info )
303*
304* Return if INFO is non-zero.
305*
306 IF( info.GT.0 )THEN
307 rcond = zero
308 RETURN
309 END IF
310 END IF
311*
312* Compute the norm of the matrix A.
313*
314 anorm = zlanht( '1', n, d, e )
315*
316* Compute the reciprocal of the condition number of A.
317*
318 CALL zptcon( n, df, ef, anorm, rcond, rwork, info )
319*
320* Compute the solution vectors X.
321*
322 CALL zlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
323 CALL zpttrs( 'Lower', n, nrhs, df, ef, x, ldx, info )
324*
325* Use iterative refinement to improve the computed solutions and
326* compute error bounds and backward error estimates for them.
327*
328 CALL zptrfs( 'Lower', n, nrhs, d, e, df, ef, b, ldb, x, ldx, ferr,
329 $ berr, work, rwork, info )
330*
331* Set INFO = N+1 if the matrix is singular to working precision.
332*
333 IF( rcond.LT.dlamch( 'Epsilon' ) )
334 $ info = n + 1
335*
336 RETURN
337*
338* End of ZPTSVX
339*
logical function lsame(ca, cb)
LSAME
Definition lsame.f:53
subroutine zlacpy(uplo, m, n, a, lda, b, ldb)
ZLACPY copies all or part of one two-dimensional array to another.
Definition zlacpy.f:103
double precision function zlanht(norm, n, d, e)
ZLANHT returns the value of the 1-norm, or the Frobenius norm, or the infinity norm,...
Definition zlanht.f:101
subroutine zptrfs(uplo, n, nrhs, d, e, df, ef, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZPTRFS
Definition zptrfs.f:183
subroutine zptcon(n, d, e, anorm, rcond, rwork, info)
ZPTCON
Definition zptcon.f:119
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
subroutine dcopy(n, dx, incx, dy, incy)
DCOPY
Definition dcopy.f:82
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69