LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
subroutine ztrsv ( character  UPLO,
character  TRANS,
character  DIAG,
integer  N,
complex*16, dimension(lda,*)  A,
integer  LDA,
complex*16, dimension(*)  X,
integer  INCX 
)

ZTRSV

Purpose:
 ZTRSV  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.
Parameters
[in]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.
[in]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.
[in]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.
[in]N
          N is INTEGER
           On entry, N specifies the order of the matrix A.
           N must be at least zero.
[in]A
          A is COMPLEX*16 array of 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.
[in]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 ).
[in,out]X
          X is COMPLEX*16 array of 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.
[in]INCX
          INCX is INTEGER
           On entry, INCX specifies the increment for the elements of
           X. INCX must not be zero.
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Date
November 2011
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.

Definition at line 151 of file ztrsv.f.

151 *
152 * -- Reference BLAS level2 routine (version 3.4.0) --
153 * -- Reference BLAS is a software package provided by Univ. of Tennessee, --
154 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
155 * November 2011
156 *
157 * .. Scalar Arguments ..
158  INTEGER incx,lda,n
159  CHARACTER diag,trans,uplo
160 * ..
161 * .. Array Arguments ..
162  COMPLEX*16 a(lda,*),x(*)
163 * ..
164 *
165 * =====================================================================
166 *
167 * .. Parameters ..
168  COMPLEX*16 zero
169  parameter(zero= (0.0d+0,0.0d+0))
170 * ..
171 * .. Local Scalars ..
172  COMPLEX*16 temp
173  INTEGER i,info,ix,j,jx,kx
174  LOGICAL noconj,nounit
175 * ..
176 * .. External Functions ..
177  LOGICAL lsame
178  EXTERNAL lsame
179 * ..
180 * .. External Subroutines ..
181  EXTERNAL xerbla
182 * ..
183 * .. Intrinsic Functions ..
184  INTRINSIC dconjg,max
185 * ..
186 *
187 * Test the input parameters.
188 *
189  info = 0
190  IF (.NOT.lsame(uplo,'U') .AND. .NOT.lsame(uplo,'L')) THEN
191  info = 1
192  ELSE IF (.NOT.lsame(trans,'N') .AND. .NOT.lsame(trans,'T') .AND.
193  + .NOT.lsame(trans,'C')) THEN
194  info = 2
195  ELSE IF (.NOT.lsame(diag,'U') .AND. .NOT.lsame(diag,'N')) THEN
196  info = 3
197  ELSE IF (n.LT.0) THEN
198  info = 4
199  ELSE IF (lda.LT.max(1,n)) THEN
200  info = 6
201  ELSE IF (incx.EQ.0) THEN
202  info = 8
203  END IF
204  IF (info.NE.0) THEN
205  CALL xerbla('ZTRSV ',info)
206  RETURN
207  END IF
208 *
209 * Quick return if possible.
210 *
211  IF (n.EQ.0) RETURN
212 *
213  noconj = lsame(trans,'T')
214  nounit = lsame(diag,'N')
215 *
216 * Set up the start point in X if the increment is not unity. This
217 * will be ( N - 1 )*INCX too small for descending loops.
218 *
219  IF (incx.LE.0) THEN
220  kx = 1 - (n-1)*incx
221  ELSE IF (incx.NE.1) THEN
222  kx = 1
223  END IF
224 *
225 * Start the operations. In this version the elements of A are
226 * accessed sequentially with one pass through A.
227 *
228  IF (lsame(trans,'N')) THEN
229 *
230 * Form x := inv( A )*x.
231 *
232  IF (lsame(uplo,'U')) THEN
233  IF (incx.EQ.1) THEN
234  DO 20 j = n,1,-1
235  IF (x(j).NE.zero) THEN
236  IF (nounit) x(j) = x(j)/a(j,j)
237  temp = x(j)
238  DO 10 i = j - 1,1,-1
239  x(i) = x(i) - temp*a(i,j)
240  10 CONTINUE
241  END IF
242  20 CONTINUE
243  ELSE
244  jx = kx + (n-1)*incx
245  DO 40 j = n,1,-1
246  IF (x(jx).NE.zero) THEN
247  IF (nounit) x(jx) = x(jx)/a(j,j)
248  temp = x(jx)
249  ix = jx
250  DO 30 i = j - 1,1,-1
251  ix = ix - incx
252  x(ix) = x(ix) - temp*a(i,j)
253  30 CONTINUE
254  END IF
255  jx = jx - incx
256  40 CONTINUE
257  END IF
258  ELSE
259  IF (incx.EQ.1) THEN
260  DO 60 j = 1,n
261  IF (x(j).NE.zero) THEN
262  IF (nounit) x(j) = x(j)/a(j,j)
263  temp = x(j)
264  DO 50 i = j + 1,n
265  x(i) = x(i) - temp*a(i,j)
266  50 CONTINUE
267  END IF
268  60 CONTINUE
269  ELSE
270  jx = kx
271  DO 80 j = 1,n
272  IF (x(jx).NE.zero) THEN
273  IF (nounit) x(jx) = x(jx)/a(j,j)
274  temp = x(jx)
275  ix = jx
276  DO 70 i = j + 1,n
277  ix = ix + incx
278  x(ix) = x(ix) - temp*a(i,j)
279  70 CONTINUE
280  END IF
281  jx = jx + incx
282  80 CONTINUE
283  END IF
284  END IF
285  ELSE
286 *
287 * Form x := inv( A**T )*x or x := inv( A**H )*x.
288 *
289  IF (lsame(uplo,'U')) THEN
290  IF (incx.EQ.1) THEN
291  DO 110 j = 1,n
292  temp = x(j)
293  IF (noconj) THEN
294  DO 90 i = 1,j - 1
295  temp = temp - a(i,j)*x(i)
296  90 CONTINUE
297  IF (nounit) temp = temp/a(j,j)
298  ELSE
299  DO 100 i = 1,j - 1
300  temp = temp - dconjg(a(i,j))*x(i)
301  100 CONTINUE
302  IF (nounit) temp = temp/dconjg(a(j,j))
303  END IF
304  x(j) = temp
305  110 CONTINUE
306  ELSE
307  jx = kx
308  DO 140 j = 1,n
309  ix = kx
310  temp = x(jx)
311  IF (noconj) THEN
312  DO 120 i = 1,j - 1
313  temp = temp - a(i,j)*x(ix)
314  ix = ix + incx
315  120 CONTINUE
316  IF (nounit) temp = temp/a(j,j)
317  ELSE
318  DO 130 i = 1,j - 1
319  temp = temp - dconjg(a(i,j))*x(ix)
320  ix = ix + incx
321  130 CONTINUE
322  IF (nounit) temp = temp/dconjg(a(j,j))
323  END IF
324  x(jx) = temp
325  jx = jx + incx
326  140 CONTINUE
327  END IF
328  ELSE
329  IF (incx.EQ.1) THEN
330  DO 170 j = n,1,-1
331  temp = x(j)
332  IF (noconj) THEN
333  DO 150 i = n,j + 1,-1
334  temp = temp - a(i,j)*x(i)
335  150 CONTINUE
336  IF (nounit) temp = temp/a(j,j)
337  ELSE
338  DO 160 i = n,j + 1,-1
339  temp = temp - dconjg(a(i,j))*x(i)
340  160 CONTINUE
341  IF (nounit) temp = temp/dconjg(a(j,j))
342  END IF
343  x(j) = temp
344  170 CONTINUE
345  ELSE
346  kx = kx + (n-1)*incx
347  jx = kx
348  DO 200 j = n,1,-1
349  ix = kx
350  temp = x(jx)
351  IF (noconj) THEN
352  DO 180 i = n,j + 1,-1
353  temp = temp - a(i,j)*x(ix)
354  ix = ix - incx
355  180 CONTINUE
356  IF (nounit) temp = temp/a(j,j)
357  ELSE
358  DO 190 i = n,j + 1,-1
359  temp = temp - dconjg(a(i,j))*x(ix)
360  ix = ix - incx
361  190 CONTINUE
362  IF (nounit) temp = temp/dconjg(a(j,j))
363  END IF
364  x(jx) = temp
365  jx = jx - incx
366  200 CONTINUE
367  END IF
368  END IF
369  END IF
370 *
371  RETURN
372 *
373 * End of ZTRSV .
374 *
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
logical function lsame(CA, CB)
LSAME
Definition: lsame.f:55

Here is the call graph for this function:

Here is the caller graph for this function: