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

Functions

subroutine zgtsv (n, nrhs, dl, d, du, b, ldb, info)
  ZGTSV computes the solution to system of linear equations A * X = B for GT matrices
subroutine zgtsvx (fact, trans, n, nrhs, dl, d, du, dlf, df, duf, du2, ipiv, b, ldb, x, ldx, rcond, ferr, berr, work, rwork, info)
  ZGTSVX computes the solution to system of linear equations A * X = B for GT matrices

Detailed Description

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

Function Documentation

◆ zgtsv()

subroutine zgtsv ( integer n,
integer nrhs,
complex*16, dimension( * ) dl,
complex*16, dimension( * ) d,
complex*16, dimension( * ) du,
complex*16, dimension( ldb, * ) b,
integer ldb,
integer info )

ZGTSV computes the solution to system of linear equations A * X = B for GT matrices

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

Purpose:
!>
!> ZGTSV  solves the equation
!>
!>    A*X = B,
!>
!> where A is an N-by-N tridiagonal matrix, by Gaussian elimination with
!> partial pivoting.
!>
!> Note that the equation  A**T *X = B  may be solved by interchanging the
!> order of the arguments DU and DL.
!> 
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]DL
!>          DL is COMPLEX*16 array, dimension (N-1)
!>          On entry, DL must contain the (n-1) subdiagonal elements of
!>          A.
!>          On exit, DL is overwritten by the (n-2) elements of the
!>          second superdiagonal of the upper triangular matrix U from
!>          the LU factorization of A, in DL(1), ..., DL(n-2).
!> 
[in,out]D
!>          D is COMPLEX*16 array, dimension (N)
!>          On entry, D must contain the diagonal elements of A.
!>          On exit, D is overwritten by the n diagonal elements of U.
!> 
[in,out]DU
!>          DU is COMPLEX*16 array, dimension (N-1)
!>          On entry, DU must contain the (n-1) superdiagonal elements
!>          of A.
!>          On exit, DU is overwritten by the (n-1) elements of the first
!>          superdiagonal of U.
!> 
[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, U(i,i) is exactly zero, 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 123 of file zgtsv.f.

124*
125* -- LAPACK driver routine --
126* -- LAPACK is a software package provided by Univ. of Tennessee, --
127* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
128*
129* .. Scalar Arguments ..
130 INTEGER INFO, LDB, N, NRHS
131* ..
132* .. Array Arguments ..
133 COMPLEX*16 B( LDB, * ), D( * ), DL( * ), DU( * )
134* ..
135*
136* =====================================================================
137*
138* .. Parameters ..
139 COMPLEX*16 ZERO
140 parameter( zero = ( 0.0d+0, 0.0d+0 ) )
141* ..
142* .. Local Scalars ..
143 INTEGER J, K
144 COMPLEX*16 MULT, TEMP, ZDUM
145* ..
146* .. Intrinsic Functions ..
147 INTRINSIC abs, dble, dimag, max
148* ..
149* .. External Subroutines ..
150 EXTERNAL xerbla
151* ..
152* .. Statement Functions ..
153 DOUBLE PRECISION CABS1
154* ..
155* .. Statement Function definitions ..
156 cabs1( zdum ) = abs( dble( zdum ) ) + abs( dimag( zdum ) )
157* ..
158* .. Executable Statements ..
159*
160 info = 0
161 IF( n.LT.0 ) THEN
162 info = -1
163 ELSE IF( nrhs.LT.0 ) THEN
164 info = -2
165 ELSE IF( ldb.LT.max( 1, n ) ) THEN
166 info = -7
167 END IF
168 IF( info.NE.0 ) THEN
169 CALL xerbla( 'ZGTSV ', -info )
170 RETURN
171 END IF
172*
173 IF( n.EQ.0 )
174 $ RETURN
175*
176 DO 30 k = 1, n - 1
177 IF( dl( k ).EQ.zero ) THEN
178*
179* Subdiagonal is zero, no elimination is required.
180*
181 IF( d( k ).EQ.zero ) THEN
182*
183* Diagonal is zero: set INFO = K and return; a unique
184* solution can not be found.
185*
186 info = k
187 RETURN
188 END IF
189 ELSE IF( cabs1( d( k ) ).GE.cabs1( dl( k ) ) ) THEN
190*
191* No row interchange required
192*
193 mult = dl( k ) / d( k )
194 d( k+1 ) = d( k+1 ) - mult*du( k )
195 DO 10 j = 1, nrhs
196 b( k+1, j ) = b( k+1, j ) - mult*b( k, j )
197 10 CONTINUE
198 IF( k.LT.( n-1 ) )
199 $ dl( k ) = zero
200 ELSE
201*
202* Interchange rows K and K+1
203*
204 mult = d( k ) / dl( k )
205 d( k ) = dl( k )
206 temp = d( k+1 )
207 d( k+1 ) = du( k ) - mult*temp
208 IF( k.LT.( n-1 ) ) THEN
209 dl( k ) = du( k+1 )
210 du( k+1 ) = -mult*dl( k )
211 END IF
212 du( k ) = temp
213 DO 20 j = 1, nrhs
214 temp = b( k, j )
215 b( k, j ) = b( k+1, j )
216 b( k+1, j ) = temp - mult*b( k+1, j )
217 20 CONTINUE
218 END IF
219 30 CONTINUE
220 IF( d( n ).EQ.zero ) THEN
221 info = n
222 RETURN
223 END IF
224*
225* Back solve with the matrix U from the factorization.
226*
227 DO 50 j = 1, nrhs
228 b( n, j ) = b( n, j ) / d( n )
229 IF( n.GT.1 )
230 $ b( n-1, j ) = ( b( n-1, j )-du( n-1 )*b( n, j ) ) / d( n-1 )
231 DO 40 k = n - 2, 1, -1
232 b( k, j ) = ( b( k, j )-du( k )*b( k+1, j )-dl( k )*
233 $ b( k+2, j ) ) / d( k )
234 40 CONTINUE
235 50 CONTINUE
236*
237 RETURN
238*
239* End of ZGTSV
240*
subroutine xerbla(srname, info)
XERBLA
Definition xerbla.f:60
#define max(a, b)
Definition macros.h:21

◆ zgtsvx()

subroutine zgtsvx ( character fact,
character trans,
integer n,
integer nrhs,
complex*16, dimension( * ) dl,
complex*16, dimension( * ) d,
complex*16, dimension( * ) du,
complex*16, dimension( * ) dlf,
complex*16, dimension( * ) df,
complex*16, dimension( * ) duf,
complex*16, dimension( * ) du2,
integer, dimension( * ) ipiv,
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 )

ZGTSVX computes the solution to system of linear equations A * X = B for GT matrices

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

Purpose:
!>
!> ZGTSVX uses the LU factorization to compute the solution to a complex
!> system of linear equations A * X = B, A**T * X = B, or A**H * X = B,
!> where A is a tridiagonal matrix of order N 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 LU decomposition is used to factor the matrix A
!>    as A = L * U, where L is a product of permutation and unit lower
!>    bidiagonal matrices and U is upper triangular with nonzeros in
!>    only the main diagonal and first two superdiagonals.
!>
!> 2. If some U(i,i)=0, so that U is exactly singular, 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 A has been
!>          supplied on entry.
!>          = 'F':  DLF, DF, DUF, DU2, and IPIV contain the factored form
!>                  of A; DL, D, DU, DLF, DF, DUF, DU2 and IPIV will not
!>                  be modified.
!>          = 'N':  The matrix will be copied to DLF, DF, and DUF
!>                  and factored.
!> 
[in]TRANS
!>          TRANS is CHARACTER*1
!>          Specifies the form of the system of equations:
!>          = 'N':  A * X = B     (No transpose)
!>          = 'T':  A**T * X = B  (Transpose)
!>          = 'C':  A**H * X = B  (Conjugate transpose)
!> 
[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]DL
!>          DL is COMPLEX*16 array, dimension (N-1)
!>          The (n-1) subdiagonal elements of A.
!> 
[in]D
!>          D is COMPLEX*16 array, dimension (N)
!>          The n diagonal elements of A.
!> 
[in]DU
!>          DU is COMPLEX*16 array, dimension (N-1)
!>          The (n-1) superdiagonal elements of A.
!> 
[in,out]DLF
!>          DLF is COMPLEX*16 array, dimension (N-1)
!>          If FACT = 'F', then DLF is an input argument and on entry
!>          contains the (n-1) multipliers that define the matrix L from
!>          the LU factorization of A as computed by ZGTTRF.
!>
!>          If FACT = 'N', then DLF is an output argument and on exit
!>          contains the (n-1) multipliers that define the matrix L from
!>          the LU factorization of A.
!> 
[in,out]DF
!>          DF is COMPLEX*16 array, dimension (N)
!>          If FACT = 'F', then DF is an input argument and on entry
!>          contains the n diagonal elements of the upper triangular
!>          matrix U from the LU factorization of A.
!>
!>          If FACT = 'N', then DF is an output argument and on exit
!>          contains the n diagonal elements of the upper triangular
!>          matrix U from the LU factorization of A.
!> 
[in,out]DUF
!>          DUF is COMPLEX*16 array, dimension (N-1)
!>          If FACT = 'F', then DUF is an input argument and on entry
!>          contains the (n-1) elements of the first superdiagonal of U.
!>
!>          If FACT = 'N', then DUF is an output argument and on exit
!>          contains the (n-1) elements of the first superdiagonal of U.
!> 
[in,out]DU2
!>          DU2 is COMPLEX*16 array, dimension (N-2)
!>          If FACT = 'F', then DU2 is an input argument and on entry
!>          contains the (n-2) elements of the second superdiagonal of
!>          U.
!>
!>          If FACT = 'N', then DU2 is an output argument and on exit
!>          contains the (n-2) elements of the second superdiagonal of
!>          U.
!> 
[in,out]IPIV
!>          IPIV is INTEGER array, dimension (N)
!>          If FACT = 'F', then IPIV is an input argument and on entry
!>          contains the pivot indices from the LU factorization of A as
!>          computed by ZGTTRF.
!>
!>          If FACT = 'N', then IPIV is an output argument and on exit
!>          contains the pivot indices from the LU factorization of A;
!>          row i of the matrix was interchanged with row IPIV(i).
!>          IPIV(i) will always be either i or i+1; IPIV(i) = i indicates
!>          a row interchange was not required.
!> 
[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 estimate of 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 estimated 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).  The estimate is as reliable as
!>          the estimate for RCOND, and is almost always a slight
!>          overestimate of the true error.
!> 
[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 (2*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:  U(i,i) is exactly zero.  The factorization
!>                       has not been completed unless i = N, but the
!>                       factor U is exactly singular, so the solution
!>                       and error bounds could not be 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 291 of file zgtsvx.f.

294*
295* -- LAPACK driver routine --
296* -- LAPACK is a software package provided by Univ. of Tennessee, --
297* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
298*
299* .. Scalar Arguments ..
300 CHARACTER FACT, TRANS
301 INTEGER INFO, LDB, LDX, N, NRHS
302 DOUBLE PRECISION RCOND
303* ..
304* .. Array Arguments ..
305 INTEGER IPIV( * )
306 DOUBLE PRECISION BERR( * ), FERR( * ), RWORK( * )
307 COMPLEX*16 B( LDB, * ), D( * ), DF( * ), DL( * ),
308 $ DLF( * ), DU( * ), DU2( * ), DUF( * ),
309 $ WORK( * ), X( LDX, * )
310* ..
311*
312* =====================================================================
313*
314* .. Parameters ..
315 DOUBLE PRECISION ZERO
316 parameter( zero = 0.0d+0 )
317* ..
318* .. Local Scalars ..
319 LOGICAL NOFACT, NOTRAN
320 CHARACTER NORM
321 DOUBLE PRECISION ANORM
322* ..
323* .. External Functions ..
324 LOGICAL LSAME
325 DOUBLE PRECISION DLAMCH, ZLANGT
326 EXTERNAL lsame, dlamch, zlangt
327* ..
328* .. External Subroutines ..
329 EXTERNAL xerbla, zcopy, zgtcon, zgtrfs, zgttrf, zgttrs,
330 $ zlacpy
331* ..
332* .. Intrinsic Functions ..
333 INTRINSIC max
334* ..
335* .. Executable Statements ..
336*
337 info = 0
338 nofact = lsame( fact, 'N' )
339 notran = lsame( trans, 'N' )
340 IF( .NOT.nofact .AND. .NOT.lsame( fact, 'F' ) ) THEN
341 info = -1
342 ELSE IF( .NOT.notran .AND. .NOT.lsame( trans, 'T' ) .AND. .NOT.
343 $ lsame( trans, 'C' ) ) THEN
344 info = -2
345 ELSE IF( n.LT.0 ) THEN
346 info = -3
347 ELSE IF( nrhs.LT.0 ) THEN
348 info = -4
349 ELSE IF( ldb.LT.max( 1, n ) ) THEN
350 info = -14
351 ELSE IF( ldx.LT.max( 1, n ) ) THEN
352 info = -16
353 END IF
354 IF( info.NE.0 ) THEN
355 CALL xerbla( 'ZGTSVX', -info )
356 RETURN
357 END IF
358*
359 IF( nofact ) THEN
360*
361* Compute the LU factorization of A.
362*
363 CALL zcopy( n, d, 1, df, 1 )
364 IF( n.GT.1 ) THEN
365 CALL zcopy( n-1, dl, 1, dlf, 1 )
366 CALL zcopy( n-1, du, 1, duf, 1 )
367 END IF
368 CALL zgttrf( n, dlf, df, duf, du2, ipiv, info )
369*
370* Return if INFO is non-zero.
371*
372 IF( info.GT.0 )THEN
373 rcond = zero
374 RETURN
375 END IF
376 END IF
377*
378* Compute the norm of the matrix A.
379*
380 IF( notran ) THEN
381 norm = '1'
382 ELSE
383 norm = 'I'
384 END IF
385 anorm = zlangt( norm, n, dl, d, du )
386*
387* Compute the reciprocal of the condition number of A.
388*
389 CALL zgtcon( norm, n, dlf, df, duf, du2, ipiv, anorm, rcond, work,
390 $ info )
391*
392* Compute the solution vectors X.
393*
394 CALL zlacpy( 'Full', n, nrhs, b, ldb, x, ldx )
395 CALL zgttrs( trans, n, nrhs, dlf, df, duf, du2, ipiv, x, ldx,
396 $ info )
397*
398* Use iterative refinement to improve the computed solutions and
399* compute error bounds and backward error estimates for them.
400*
401 CALL zgtrfs( trans, n, nrhs, dl, d, du, dlf, df, duf, du2, ipiv,
402 $ b, ldb, x, ldx, ferr, berr, work, rwork, info )
403*
404* Set INFO = N+1 if the matrix is singular to working precision.
405*
406 IF( rcond.LT.dlamch( 'Epsilon' ) )
407 $ info = n + 1
408*
409 RETURN
410*
411* End of ZGTSVX
412*
norm(diag(diag(diag(inv(mat))) -id.SOL), 2) % destroy mumps instance id.JOB
logical function lsame(ca, cb)
LSAME
Definition lsame.f:53
subroutine zgtrfs(trans, n, nrhs, dl, d, du, dlf, df, duf, du2, ipiv, b, ldb, x, ldx, ferr, berr, work, rwork, info)
ZGTRFS
Definition zgtrfs.f:210
subroutine zgtcon(norm, n, dl, d, du, du2, ipiv, anorm, rcond, work, info)
ZGTCON
Definition zgtcon.f:141
subroutine zgttrs(trans, n, nrhs, dl, d, du, du2, ipiv, b, ldb, info)
ZGTTRS
Definition zgttrs.f:138
subroutine zgttrf(n, dl, d, du, du2, ipiv, info)
ZGTTRF
Definition zgttrf.f:124
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 zlangt(norm, n, dl, d, du)
ZLANGT returns the value of the 1-norm, Frobenius norm, infinity-norm, or the largest absolute value ...
Definition zlangt.f:106
subroutine zcopy(n, zx, incx, zy, incy)
ZCOPY
Definition zcopy.f:81
double precision function dlamch(cmach)
DLAMCH
Definition dlamch.f:69