ctrsv Subroutine

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

NAME

ctrsv(3f) - [BLAS:COMPLEX_BLAS_LEVEL2]

CX := INVERSE(A)*CX, where A is a triangular matrix.

SYNOPSIS

 subroutine ctrsv(uplo,trans,diag,n,a,lda,x,incx)

   .. Scalar Arguments ..
   integer,intent(in)    :: incx,lda,n
   character,intent(in)  :: diag,trans,uplo
   ..
   .. Array Arguments ..
   complex,intent(in)    :: a(lda,*)
   complex,intent(inout) :: x(*)
   ..

DEFINITION

CTRSV solves one of the systems of equations

 A*x = b,   or   A**T*x = b,   or   A**H*x = 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.

OPTIONS

UPLO

       UPLO is CHARACTER*1
        On entry, UPLO specifies whether the matrix is an upper or
        lower triangular matrix as follows:

           UPLO = 'U' or 'u'   A is an upper triangular matrix.

           UPLO = 'L' or 'l'   A is a lower triangular matrix.

TRANS

       TRANS is CHARACTER*1
        On entry, TRANS specifies the equations to be solved as
        follows:

           TRANS = 'N' or 'n'   A*x = b.

           TRANS = 'T' or 't'   A**T*x = b.

           TRANS = 'C' or 'c'   A**H*x = b.

DIAG

       DIAG is CHARACTER*1
        On entry, DIAG specifies whether or not A is unit
        triangular as follows:

           DIAG = 'U' or 'u'   A is assumed to be unit triangular.

           DIAG = 'N' or 'n'   A is not assumed to be unit
                               triangular.

N

       N is INTEGER
        On entry, N specifies the order of the matrix A.
        N must be at least zero.

A

       A is COMPLEX array, dimension ( LDA, N )
        Before entry with UPLO = 'U' or 'u', the leading n by n
        upper triangular part of the array A must contain the upper
        triangular matrix and the strictly lower triangular part of
        A is not referenced.
        Before entry with UPLO = 'L' or 'l', the leading n by n
        lower triangular part of the array A must contain the lower
        triangular matrix and the strictly upper triangular part of
        A is not referenced.
        Note that when DIAG = 'U' or 'u', the diagonal elements of
        A are not referenced either, but are assumed to be unity.

LDA

       LDA is INTEGER
        On entry, LDA specifies the first dimension of A as declared
        in the calling (sub) program. LDA must be at least
        max( 1, n ).

X

       X is COMPLEX array, dimension at least
        ( 1 + ( n - 1 )*abs( INCX ) ).
        Before entry, the incremented array X must contain the n
        element right-hand side vector b. On exit, X is overwritten
        with the solution vector x.

INCX

       INCX is INTEGER
        On entry, INCX specifies the increment for the elements of
        X. INCX must not be zero.

AUTHORS

  • Univ. of Tennessee
  • Univ. of California Berkeley
  • Univ. of Colorado Denver
  • NAG Ltd.

date:December 2016

FURTHER DETAILS

Level 2 Blas routine.

– Written on 22-October-1986. Jack Dongarra, Argonne National Lab. Jeremy Du Croz, Nag Central Office. Sven Hammarling, Nag Central Office. Richard Hanson, Sandia National Labs.

SEE ALSO

Online html documentation available at
http://www.netlib.org/lapack/explore-html/

Arguments

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

Contents

Source Code


Variables

Type Visibility Attributes Name Initial
integer, public :: i
integer, public :: info
integer, public :: ix
integer, public :: j
integer, public :: jx
integer, public :: kx
logical, public :: noconj
logical, public :: nounit
complex, public :: temp
complex, public, parameter :: zero = (0.0e+0,0.0e+0)

Source Code

       subroutine ctrsv(uplo,trans,diag,n,a,lda,x,incx)
      implicit none
!
!  -- Reference BLAS level2 routine (version 3.7.0) --
!  -- Reference BLAS is a software package provided by Univ. of Tennessee,    --
!  -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
!     December 2016
!
!     .. Scalar Arguments ..
      integer,intent(in)    :: incx,lda,n
      character,intent(in)  :: diag,trans,uplo
!     ..
!     .. Array Arguments ..
      complex,intent(in)    :: a(lda,*)
      complex,intent(inout) :: x(*)
!     ..
!
!  =====================================================================
!
!     .. Parameters ..
      complex zero
      parameter (zero= (0.0e+0,0.0e+0))
!     ..
!     .. Local Scalars ..
      complex temp
      integer i,info,ix,j,jx,kx
      logical noconj,nounit
!     ..
!     .. External Functions ..
!     ..
!     .. External Subroutines ..
!     ..
!     .. Intrinsic Functions ..
      intrinsic conjg,max
!     ..
!
!     Test the input parameters.
!
      info = 0
      if (.not.lsame(uplo,'U') .and. .not.lsame(uplo,'L')) then
          info = 1
      elseif (.not.lsame(trans,'N') .and. .not.lsame(trans,'T') .and.  .not.lsame(trans,'C')) then
          info = 2
      elseif (.not.lsame(diag,'U') .and. .not.lsame(diag,'N')) then
          info = 3
      elseif (n.lt.0) then
          info = 4
      elseif (lda.lt.max(1,n)) then
          info = 6
      elseif (incx.eq.0) then
          info = 8
      endif
      if (info.ne.0) then
          call xerbla('CTRSV ',info)
          return
      endif
!
!     Quick return if possible.
!
      if (n.eq.0) return
!
      noconj = lsame(trans,'T')
      nounit = lsame(diag,'N')
!
!     Set up the start point in X if the increment is not unity. This
!     will be  ( N - 1 )*INCX  too small for descending loops.
!
      if (incx.le.0) then
          kx = 1 - (n-1)*incx
      elseif (incx.ne.1) then
          kx = 1
      endif
!
!     Start the operations. In this version the elements of A are
!     accessed sequentially with one pass through A.
!
      if (lsame(trans,'N')) then
!
!        Form  x := inv( A )*x.
!
          if (lsame(uplo,'U')) then
              if (incx.eq.1) then
                  do j = n,1,-1
                      if (x(j).ne.zero) then
                          if (nounit) x(j) = x(j)/a(j,j)
                          temp = x(j)
                          do i = j - 1,1,-1
                              x(i) = x(i) - temp*a(i,j)
                          enddo
                      endif
                  enddo
              else
                  jx = kx + (n-1)*incx
                  do j = n,1,-1
                      if (x(jx).ne.zero) then
                          if (nounit) x(jx) = x(jx)/a(j,j)
                          temp = x(jx)
                          ix = jx
                          do i = j - 1,1,-1
                              ix = ix - incx
                              x(ix) = x(ix) - temp*a(i,j)
                          enddo
                      endif
                      jx = jx - incx
                  enddo
              endif
          else
              if (incx.eq.1) then
                  do j = 1,n
                      if (x(j).ne.zero) then
                          if (nounit) x(j) = x(j)/a(j,j)
                          temp = x(j)
                          do i = j + 1,n
                              x(i) = x(i) - temp*a(i,j)
                          enddo
                      endif
                  enddo
              else
                  jx = kx
                  do j = 1,n
                      if (x(jx).ne.zero) then
                          if (nounit) x(jx) = x(jx)/a(j,j)
                          temp = x(jx)
                          ix = jx
                          do i = j + 1,n
                              ix = ix + incx
                              x(ix) = x(ix) - temp*a(i,j)
                          enddo
                      endif
                      jx = jx + incx
                  enddo
              endif
          endif
      else
!
!        Form  x := inv( A**T )*x  or  x := inv( A**H )*x.
!
          if (lsame(uplo,'U')) then
              if (incx.eq.1) then
                  do j = 1,n
                      temp = x(j)
                      if (noconj) then
                          do i = 1,j - 1
                              temp = temp - a(i,j)*x(i)
                          enddo
                          if (nounit) temp = temp/a(j,j)
                      else
                          do i = 1,j - 1
                              temp = temp - conjg(a(i,j))*x(i)
                          enddo
                          if (nounit) temp = temp/conjg(a(j,j))
                      endif
                      x(j) = temp
                  enddo
              else
                  jx = kx
                  do j = 1,n
                      ix = kx
                      temp = x(jx)
                      if (noconj) then
                          do i = 1,j - 1
                              temp = temp - a(i,j)*x(ix)
                              ix = ix + incx
                          enddo
                          if (nounit) temp = temp/a(j,j)
                      else
                          do i = 1,j - 1
                              temp = temp - conjg(a(i,j))*x(ix)
                              ix = ix + incx
                          enddo
                          if (nounit) temp = temp/conjg(a(j,j))
                      endif
                      x(jx) = temp
                      jx = jx + incx
                  enddo
              endif
          else
              if (incx.eq.1) then
                  do j = n,1,-1
                      temp = x(j)
                      if (noconj) then
                          do i = n,j + 1,-1
                              temp = temp - a(i,j)*x(i)
                          enddo
                          if (nounit) temp = temp/a(j,j)
                      else
                          do i = n,j + 1,-1
                              temp = temp - conjg(a(i,j))*x(i)
                          enddo
                          if (nounit) temp = temp/conjg(a(j,j))
                      endif
                      x(j) = temp
                  enddo
              else
                  kx = kx + (n-1)*incx
                  jx = kx
                  do j = n,1,-1
                      ix = kx
                      temp = x(jx)
                      if (noconj) then
                          do i = n,j + 1,-1
                              temp = temp - a(i,j)*x(ix)
                              ix = ix - incx
                          enddo
                          if (nounit) temp = temp/a(j,j)
                      else
                          do i = n,j + 1,-1
                              temp = temp - conjg(a(i,j))*x(ix)
                              ix = ix - incx
                          enddo
                          if (nounit) temp = temp/conjg(a(j,j))
                      endif
                      x(jx) = temp
                      jx = jx - incx
                  enddo
              endif
          endif
      endif
!
!     End of CTRSV .
!
      end subroutine ctrsv