stdlib_linalg_blas_d Module



Functions

public pure function stdlib_dasum(n, dx, incx)

DASUM takes the sum of the absolute values.

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: dx(*)
integer(kind=ilp), intent(in) :: incx

Return Value real(kind=dp)

public pure function stdlib_ddot(n, dx, incx, dy, incy)

DDOT forms the dot product of two vectors. uses unrolled loops for increments equal to one.

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: dx(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(in) :: dy(*)
integer(kind=ilp), intent(in) :: incy

Return Value real(kind=dp)

public pure function stdlib_dnrm2(n, x, incx)

DNRM2 returns the euclidean norm of a vector via the function name, so that DNRM2 := sqrt( x'*x )

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: x(*)
integer(kind=ilp), intent(in) :: incx

Return Value real(kind=dp)

public pure function stdlib_dsdot(n, sx, incx, sy, incy)

Compute the inner product of two vectors with extended precision accumulation and result. Returns D.P. dot product accumulated in D.P., for S.P. SX and SY DSDOT = sum for I = 0 to N-1 of SX(LX+IINCX) * SY(LY+IINCY), where LX = 1 if INCX >= 0, else LX = 1+(1-N)*INCX, and LY is defined in a similar way using INCY.

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: n
real(kind=sp), intent(in) :: sx(*)
integer(kind=ilp), intent(in) :: incx
real(kind=sp), intent(in) :: sy(*)
integer(kind=ilp), intent(in) :: incy

Return Value real(kind=dp)

public pure function stdlib_dzasum(n, zx, incx)

DZASUM takes the sum of the (|Re(.)| + |Im(.)|)'s of a complex vector and returns a double precision result.

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: n
complex(kind=dp), intent(in) :: zx(*)
integer(kind=ilp), intent(in) :: incx

Return Value real(kind=dp)

public pure function stdlib_dznrm2(n, x, incx)

DZNRM2 returns the euclidean norm of a vector via the function name, so that DZNRM2 := sqrt( x*Hx )

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: n
complex(kind=dp), intent(in) :: x(*)
integer(kind=ilp), intent(in) :: incx

Return Value real(kind=dp)


Subroutines

public pure subroutine stdlib_daxpy(n, da, dx, incx, dy, incy)

DAXPY constant times a vector plus a vector. uses unrolled loops for increments equal to one.

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: da
real(kind=dp), intent(in) :: dx(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(inout) :: dy(*)
integer(kind=ilp), intent(in) :: incy

public pure subroutine stdlib_dcopy(n, dx, incx, dy, incy)

DCOPY copies a vector, x, to a vector, y. uses unrolled loops for increments equal to 1.

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: dx(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(out) :: dy(*)
integer(kind=ilp), intent(in) :: incy

public pure subroutine stdlib_dgbmv(trans, m, n, kl, ku, alpha, a, lda, x, incx, beta, y, incy)

DGBMV performs one of the matrix-vector operations y := alphaAx + betay, or y := alphaATx + betay, where alpha and beta are scalars, x and y are vectors and A is an m by n band matrix, with kl sub-diagonals and ku super-diagonals.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: trans
integer(kind=ilp), intent(in) :: m
integer(kind=ilp), intent(in) :: n
integer(kind=ilp), intent(in) :: kl
integer(kind=ilp), intent(in) :: ku
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(in) :: x(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(in) :: beta
real(kind=dp), intent(inout) :: y(*)
integer(kind=ilp), intent(in) :: incy

public pure subroutine stdlib_dgemm(transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)

DGEMM performs one of the matrix-matrix operations C := alphaop( A )op( B ) + betaC, where op( X ) is one of op( X ) = X or op( X ) = X*T, alpha and beta are scalars, and A, B and C are matrices, with op( A ) an m by k matrix, op( B ) a k by n matrix and C an m by n matrix.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: transa
character(len=1), intent(in) :: transb
integer(kind=ilp), intent(in) :: m
integer(kind=ilp), intent(in) :: n
integer(kind=ilp), intent(in) :: k
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(in) :: b(ldb,*)
integer(kind=ilp), intent(in) :: ldb
real(kind=dp), intent(in) :: beta
real(kind=dp), intent(inout) :: c(ldc,*)
integer(kind=ilp), intent(in) :: ldc

public pure subroutine stdlib_dgemv(trans, m, n, alpha, a, lda, x, incx, beta, y, incy)

DGEMV performs one of the matrix-vector operations y := alphaAx + betay, or y := alphaATx + betay, where alpha and beta are scalars, x and y are vectors and A is an m by n matrix.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: trans
integer(kind=ilp), intent(in) :: m
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(in) :: x(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(in) :: beta
real(kind=dp), intent(inout) :: y(*)
integer(kind=ilp), intent(in) :: incy

public pure subroutine stdlib_dger(m, n, alpha, x, incx, y, incy, a, lda)

DGER performs the rank 1 operation A := alphaxy**T + A, where alpha is a scalar, x is an m element vector, y is an n element vector and A is an m by n matrix.

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: m
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: x(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(in) :: y(*)
integer(kind=ilp), intent(in) :: incy
real(kind=dp), intent(inout) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda

public pure subroutine stdlib_drot(n, dx, incx, dy, incy, c, s)

DROT applies a plane rotation.

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(inout) :: dx(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(inout) :: dy(*)
integer(kind=ilp), intent(in) :: incy
real(kind=dp), intent(in) :: c
real(kind=dp), intent(in) :: s

public pure subroutine stdlib_drotg(a, b, c, s)

The computation uses the formulas sigma = sgn(a) if |a| > |b| = sgn(b) if |b| >= |a| r = sigmasqrt( a2 + b2 ) c = 1; s = 0 if r = 0 c = a/r; s = b/r if r != 0 The subroutine also computes z = s if |a| > |b|, = 1/c if |b| >= |a| and c != 0 = 1 if c = 0 This allows c and s to be reconstructed from z as follows: If z = 1, set c = 0, s = 1. If |z| < 1, set c = sqrt(1 - z2) and s = z. If |z| > 1, set c = 1/z and s = sqrt( 1 - c*2).

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(inout) :: a
real(kind=dp), intent(inout) :: b
real(kind=dp), intent(out) :: c
real(kind=dp), intent(out) :: s

public pure subroutine stdlib_drotm(n, dx, incx, dy, incy, dparam)

DROTM applies the modified Givens transformation, , to the 2-by-N matrix where indicates transpose. The elements of are in DX(LX+IINCX), I = 0:N-1, where LX = 1 if INCX >= 0, else LX = (-INCX)N, and similarly for DY using LY and INCY. With DPARAM(1)=DFLAG, has one of the following forms:
See DROTMG for a description of data storage in DPARAM.

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(inout) :: dx(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(inout) :: dy(*)
integer(kind=ilp), intent(in) :: incy
real(kind=dp), intent(in) :: dparam(5)

public pure subroutine stdlib_drotmg(dd1, dd2, dx1, dy1, dparam)

DROTMG Constructs the modified Givens transformation matrix which zeros the second component of the 2-vector With DPARAM(1)=DFLAG, has one of the following forms:
Locations 2-4 of DPARAM contain DH11, DH21, DH12 and DH22 respectively. (Values of 1.0, -1.0, or 0.0 implied by the value of DPARAM(1) are not stored in DPARAM.) The values of parameters GAMSQ and RGAMSQ may be inexact. This is OK as they are only used for testing the size of DD1 and DD2. All actual scaling of data is done using GAM.

Arguments

Type IntentOptional Attributes Name
real(kind=dp), intent(inout) :: dd1
real(kind=dp), intent(inout) :: dd2
real(kind=dp), intent(inout) :: dx1
real(kind=dp), intent(in) :: dy1
real(kind=dp), intent(out) :: dparam(5)

public pure subroutine stdlib_dsbmv(uplo, n, k, alpha, a, lda, x, incx, beta, y, incy)

DSBMV performs the matrix-vector operation y := alphaAx + beta*y, where alpha and beta are scalars, x and y are n element vectors and A is an n by n symmetric band matrix, with k super-diagonals.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
integer(kind=ilp), intent(in) :: n
integer(kind=ilp), intent(in) :: k
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(in) :: x(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(in) :: beta
real(kind=dp), intent(inout) :: y(*)
integer(kind=ilp), intent(in) :: incy

public pure subroutine stdlib_dscal(n, da, dx, incx)

DSCAL scales a vector by a constant. uses unrolled loops for increment equal to 1.

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: da
real(kind=dp), intent(inout) :: dx(*)
integer(kind=ilp), intent(in) :: incx

public pure subroutine stdlib_dspmv(uplo, n, alpha, ap, x, incx, beta, y, incy)

DSPMV performs the matrix-vector operation y := alphaAx + beta*y, where alpha and beta are scalars, x and y are n element vectors and A is an n by n symmetric matrix, supplied in packed form.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: ap(*)
real(kind=dp), intent(in) :: x(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(in) :: beta
real(kind=dp), intent(inout) :: y(*)
integer(kind=ilp), intent(in) :: incy

public pure subroutine stdlib_dspr(uplo, n, alpha, x, incx, ap)

DSPR performs the symmetric rank 1 operation A := alphaxx**T + A, where alpha is a real scalar, x is an n element vector and A is an n by n symmetric matrix, supplied in packed form.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: x(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(inout) :: ap(*)

public pure subroutine stdlib_dspr2(uplo, n, alpha, x, incx, y, incy, ap)

DSPR2 performs the symmetric rank 2 operation A := alphaxyT + alphayxT + A, where alpha is a scalar, x and y are n element vectors and A is an n by n symmetric matrix, supplied in packed form.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: x(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(in) :: y(*)
integer(kind=ilp), intent(in) :: incy
real(kind=dp), intent(inout) :: ap(*)

public pure subroutine stdlib_dswap(n, dx, incx, dy, incy)

DSWAP interchanges two vectors. uses unrolled loops for increments equal to 1.

Arguments

Type IntentOptional Attributes Name
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(inout) :: dx(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(inout) :: dy(*)
integer(kind=ilp), intent(in) :: incy

public pure subroutine stdlib_dsymm(side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)

DSYMM performs one of the matrix-matrix operations C := alphaAB + betaC, or C := alphaBA + betaC, where alpha and beta are scalars, A is a symmetric matrix and B and C are m by n matrices.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: side
character(len=1), intent(in) :: uplo
integer(kind=ilp), intent(in) :: m
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(in) :: b(ldb,*)
integer(kind=ilp), intent(in) :: ldb
real(kind=dp), intent(in) :: beta
real(kind=dp), intent(inout) :: c(ldc,*)
integer(kind=ilp), intent(in) :: ldc

public pure subroutine stdlib_dsymv(uplo, n, alpha, a, lda, x, incx, beta, y, incy)

DSYMV performs the matrix-vector operation y := alphaAx + beta*y, where alpha and beta are scalars, x and y are n element vectors and A is an n by n symmetric matrix.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(in) :: x(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(in) :: beta
real(kind=dp), intent(inout) :: y(*)
integer(kind=ilp), intent(in) :: incy

public pure subroutine stdlib_dsyr(uplo, n, alpha, x, incx, a, lda)

DSYR performs the symmetric rank 1 operation A := alphaxx**T + A, where alpha is a real scalar, x is an n element vector and A is an n by n symmetric matrix.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: x(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(inout) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda

public pure subroutine stdlib_dsyr2(uplo, n, alpha, x, incx, y, incy, a, lda)

DSYR2 performs the symmetric rank 2 operation A := alphaxyT + alphayxT + A, where alpha is a scalar, x and y are n element vectors and A is an n by n symmetric matrix.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: x(*)
integer(kind=ilp), intent(in) :: incx
real(kind=dp), intent(in) :: y(*)
integer(kind=ilp), intent(in) :: incy
real(kind=dp), intent(inout) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda

public pure subroutine stdlib_dsyr2k(uplo, trans, n, k, alpha, a, lda, b, ldb, beta, c, ldc)

DSYR2K performs one of the symmetric rank 2k operations C := alphaABT + alphaBAT + betaC, or C := alphaATB + alphaBTA + betaC, where alpha and beta are scalars, C is an n by n symmetric matrix and A and B are n by k matrices in the first case and k by n matrices in the second case.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
character(len=1), intent(in) :: trans
integer(kind=ilp), intent(in) :: n
integer(kind=ilp), intent(in) :: k
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(in) :: b(ldb,*)
integer(kind=ilp), intent(in) :: ldb
real(kind=dp), intent(in) :: beta
real(kind=dp), intent(inout) :: c(ldc,*)
integer(kind=ilp), intent(in) :: ldc

public pure subroutine stdlib_dsyrk(uplo, trans, n, k, alpha, a, lda, beta, c, ldc)

DSYRK performs one of the symmetric rank k operations C := alphaAAT + betaC, or C := alphaATA + betaC, where alpha and beta are scalars, C is an n by n symmetric matrix and A is an n by k matrix in the first case and a k by n matrix in the second case.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
character(len=1), intent(in) :: trans
integer(kind=ilp), intent(in) :: n
integer(kind=ilp), intent(in) :: k
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(in) :: beta
real(kind=dp), intent(inout) :: c(ldc,*)
integer(kind=ilp), intent(in) :: ldc

public pure subroutine stdlib_dtbmv(uplo, trans, diag, n, k, a, lda, x, incx)

DTBMV performs one of the matrix-vector operations x := Ax, or x := ATx, where x is an n element vector and A is an n by n unit, or non-unit, upper or lower triangular band matrix, with ( k + 1 ) diagonals.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
character(len=1), intent(in) :: trans
character(len=1), intent(in) :: diag
integer(kind=ilp), intent(in) :: n
integer(kind=ilp), intent(in) :: k
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(inout) :: x(*)
integer(kind=ilp), intent(in) :: incx

public pure subroutine stdlib_dtbsv(uplo, trans, diag, n, k, a, lda, x, incx)

DTBSV solves one of the systems of equations Ax = b, or ATx = b, where b and x are n element vectors and A is an n by n unit, or non-unit, upper or lower triangular band matrix, with ( k + 1 ) diagonals. No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
character(len=1), intent(in) :: trans
character(len=1), intent(in) :: diag
integer(kind=ilp), intent(in) :: n
integer(kind=ilp), intent(in) :: k
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(inout) :: x(*)
integer(kind=ilp), intent(in) :: incx

public pure subroutine stdlib_dtpmv(uplo, trans, diag, n, ap, x, incx)

DTPMV performs one of the matrix-vector operations x := Ax, or x := ATx, where x is an n element vector and A is an n by n unit, or non-unit, upper or lower triangular matrix, supplied in packed form.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
character(len=1), intent(in) :: trans
character(len=1), intent(in) :: diag
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: ap(*)
real(kind=dp), intent(inout) :: x(*)
integer(kind=ilp), intent(in) :: incx

public pure subroutine stdlib_dtpsv(uplo, trans, diag, n, ap, x, incx)

DTPSV solves one of the systems of equations Ax = b, or ATx = b, where b and x are n element vectors and A is an n by n unit, or non-unit, upper or lower triangular matrix, supplied in packed form. No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
character(len=1), intent(in) :: trans
character(len=1), intent(in) :: diag
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: ap(*)
real(kind=dp), intent(inout) :: x(*)
integer(kind=ilp), intent(in) :: incx

public pure subroutine stdlib_dtrmm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)

DTRMM performs one of the matrix-matrix operations B := alphaop( A )B, or B := alphaBop( A ), where alpha is a scalar, B is an m by n matrix, A is a unit, or non-unit, upper or lower triangular matrix and op( A ) is one of op( A ) = A or op( A ) = A**T.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: side
character(len=1), intent(in) :: uplo
character(len=1), intent(in) :: transa
character(len=1), intent(in) :: diag
integer(kind=ilp), intent(in) :: m
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(inout) :: b(ldb,*)
integer(kind=ilp), intent(in) :: ldb

public pure subroutine stdlib_dtrmv(uplo, trans, diag, n, a, lda, x, incx)

DTRMV performs one of the matrix-vector operations x := Ax, or x := ATx, where x is an n element vector and A is an n by n unit, or non-unit, upper or lower triangular matrix.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
character(len=1), intent(in) :: trans
character(len=1), intent(in) :: diag
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(inout) :: x(*)
integer(kind=ilp), intent(in) :: incx

public pure subroutine stdlib_dtrsm(side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)

DTRSM solves one of the matrix equations op( A )X = alphaB, or Xop( A ) = alphaB, where alpha is a scalar, X and B are m by n matrices, A is a unit, or non-unit, upper or lower triangular matrix and op( A ) is one of op( A ) = A or op( A ) = A**T. The matrix X is overwritten on B.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: side
character(len=1), intent(in) :: uplo
character(len=1), intent(in) :: transa
character(len=1), intent(in) :: diag
integer(kind=ilp), intent(in) :: m
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: alpha
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(inout) :: b(ldb,*)
integer(kind=ilp), intent(in) :: ldb

public pure subroutine stdlib_dtrsv(uplo, trans, diag, n, a, lda, x, incx)

DTRSV solves one of the systems of equations Ax = b, or ATx = b, where b and x are n element vectors and A is an n by n unit, or non-unit, upper or lower triangular matrix. No test for singularity or near-singularity is included in this routine. Such tests must be performed before calling this routine.

Arguments

Type IntentOptional Attributes Name
character(len=1), intent(in) :: uplo
character(len=1), intent(in) :: trans
character(len=1), intent(in) :: diag
integer(kind=ilp), intent(in) :: n
real(kind=dp), intent(in) :: a(lda,*)
integer(kind=ilp), intent(in) :: lda
real(kind=dp), intent(inout) :: x(*)
integer(kind=ilp), intent(in) :: incx