LAPACK 3.12.1
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches

◆ ztbsv()

subroutine ztbsv ( character uplo,
character trans,
character diag,
integer n,
integer k,
complex*16, dimension(lda,*) a,
integer lda,
complex*16, dimension(*) x,
integer incx )

ZTBSV

Purpose:
!>
!> ZTBSV  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 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.
!> 
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]K
!>          K is INTEGER
!>           On entry with UPLO = 'U' or 'u', K specifies the number of
!>           super-diagonals of the matrix A.
!>           On entry with UPLO = 'L' or 'l', K specifies the number of
!>           sub-diagonals of the matrix A.
!>           K must satisfy  0 .le. K.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension ( LDA, N )
!>           Before entry with UPLO = 'U' or 'u', the leading ( k + 1 )
!>           by n part of the array A must contain the upper triangular
!>           band part of the matrix of coefficients, supplied column by
!>           column, with the leading diagonal of the matrix in row
!>           ( k + 1 ) of the array, the first super-diagonal starting at
!>           position 2 in row k, and so on. The top left k by k triangle
!>           of the array A is not referenced.
!>           The following program segment will transfer an upper
!>           triangular band matrix from conventional full matrix storage
!>           to band storage:
!>
!>                 DO 20, J = 1, N
!>                    M = K + 1 - J
!>                    DO 10, I = MAX( 1, J - K ), J
!>                       A( M + I, J ) = matrix( I, J )
!>              10    CONTINUE
!>              20 CONTINUE
!>
!>           Before entry with UPLO = 'L' or 'l', the leading ( k + 1 )
!>           by n part of the array A must contain the lower triangular
!>           band part of the matrix of coefficients, supplied column by
!>           column, with the leading diagonal of the matrix in row 1 of
!>           the array, the first sub-diagonal starting at position 1 in
!>           row 2, and so on. The bottom right k by k triangle of the
!>           array A is not referenced.
!>           The following program segment will transfer a lower
!>           triangular band matrix from conventional full matrix storage
!>           to band storage:
!>
!>                 DO 20, J = 1, N
!>                    M = 1 - J
!>                    DO 10, I = J, MIN( N, J + K )
!>                       A( M + I, J ) = matrix( I, J )
!>              10    CONTINUE
!>              20 CONTINUE
!>
!>           Note that when DIAG = 'U' or 'u' the elements of the array A
!>           corresponding to the diagonal elements of the matrix are not
!>           referenced, 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
!>           ( k + 1 ).
!> 
[in,out]X
!>          X is COMPLEX*16 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.
!> 
[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.
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 188 of file ztbsv.f.

189*
190* -- Reference BLAS level2 routine --
191* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
192* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
193*
194* .. Scalar Arguments ..
195 INTEGER INCX,K,LDA,N
196 CHARACTER DIAG,TRANS,UPLO
197* ..
198* .. Array Arguments ..
199 COMPLEX*16 A(LDA,*),X(*)
200* ..
201*
202* =====================================================================
203*
204* .. Parameters ..
205 COMPLEX*16 ZERO
206 parameter(zero= (0.0d+0,0.0d+0))
207* ..
208* .. Local Scalars ..
209 COMPLEX*16 TEMP
210 INTEGER I,INFO,IX,J,JX,KPLUS1,KX,L
211 LOGICAL NOCONJ,NOUNIT
212* ..
213* .. External Functions ..
214 LOGICAL LSAME
215 EXTERNAL lsame
216* ..
217* .. External Subroutines ..
218 EXTERNAL xerbla
219* ..
220* .. Intrinsic Functions ..
221 INTRINSIC dconjg,max,min
222* ..
223*
224* Test the input parameters.
225*
226 info = 0
227 IF (.NOT.lsame(uplo,'U') .AND. .NOT.lsame(uplo,'L')) THEN
228 info = 1
229 ELSE IF (.NOT.lsame(trans,'N') .AND.
230 + .NOT.lsame(trans,'T') .AND.
231 + .NOT.lsame(trans,'C')) THEN
232 info = 2
233 ELSE IF (.NOT.lsame(diag,'U') .AND.
234 + .NOT.lsame(diag,'N')) THEN
235 info = 3
236 ELSE IF (n.LT.0) THEN
237 info = 4
238 ELSE IF (k.LT.0) THEN
239 info = 5
240 ELSE IF (lda.LT. (k+1)) THEN
241 info = 7
242 ELSE IF (incx.EQ.0) THEN
243 info = 9
244 END IF
245 IF (info.NE.0) THEN
246 CALL xerbla('ZTBSV ',info)
247 RETURN
248 END IF
249*
250* Quick return if possible.
251*
252 IF (n.EQ.0) RETURN
253*
254 noconj = lsame(trans,'T')
255 nounit = lsame(diag,'N')
256*
257* Set up the start point in X if the increment is not unity. This
258* will be ( N - 1 )*INCX too small for descending loops.
259*
260 IF (incx.LE.0) THEN
261 kx = 1 - (n-1)*incx
262 ELSE IF (incx.NE.1) THEN
263 kx = 1
264 END IF
265*
266* Start the operations. In this version the elements of A are
267* accessed by sequentially with one pass through A.
268*
269 IF (lsame(trans,'N')) THEN
270*
271* Form x := inv( A )*x.
272*
273 IF (lsame(uplo,'U')) THEN
274 kplus1 = k + 1
275 IF (incx.EQ.1) THEN
276 DO 20 j = n,1,-1
277 IF (x(j).NE.zero) THEN
278 l = kplus1 - j
279 IF (nounit) x(j) = x(j)/a(kplus1,j)
280 temp = x(j)
281 DO 10 i = j - 1,max(1,j-k),-1
282 x(i) = x(i) - temp*a(l+i,j)
283 10 CONTINUE
284 END IF
285 20 CONTINUE
286 ELSE
287 kx = kx + (n-1)*incx
288 jx = kx
289 DO 40 j = n,1,-1
290 kx = kx - incx
291 IF (x(jx).NE.zero) THEN
292 ix = kx
293 l = kplus1 - j
294 IF (nounit) x(jx) = x(jx)/a(kplus1,j)
295 temp = x(jx)
296 DO 30 i = j - 1,max(1,j-k),-1
297 x(ix) = x(ix) - temp*a(l+i,j)
298 ix = ix - incx
299 30 CONTINUE
300 END IF
301 jx = jx - incx
302 40 CONTINUE
303 END IF
304 ELSE
305 IF (incx.EQ.1) THEN
306 DO 60 j = 1,n
307 IF (x(j).NE.zero) THEN
308 l = 1 - j
309 IF (nounit) x(j) = x(j)/a(1,j)
310 temp = x(j)
311 DO 50 i = j + 1,min(n,j+k)
312 x(i) = x(i) - temp*a(l+i,j)
313 50 CONTINUE
314 END IF
315 60 CONTINUE
316 ELSE
317 jx = kx
318 DO 80 j = 1,n
319 kx = kx + incx
320 IF (x(jx).NE.zero) THEN
321 ix = kx
322 l = 1 - j
323 IF (nounit) x(jx) = x(jx)/a(1,j)
324 temp = x(jx)
325 DO 70 i = j + 1,min(n,j+k)
326 x(ix) = x(ix) - temp*a(l+i,j)
327 ix = ix + incx
328 70 CONTINUE
329 END IF
330 jx = jx + incx
331 80 CONTINUE
332 END IF
333 END IF
334 ELSE
335*
336* Form x := inv( A**T )*x or x := inv( A**H )*x.
337*
338 IF (lsame(uplo,'U')) THEN
339 kplus1 = k + 1
340 IF (incx.EQ.1) THEN
341 DO 110 j = 1,n
342 temp = x(j)
343 l = kplus1 - j
344 IF (noconj) THEN
345 DO 90 i = max(1,j-k),j - 1
346 temp = temp - a(l+i,j)*x(i)
347 90 CONTINUE
348 IF (nounit) temp = temp/a(kplus1,j)
349 ELSE
350 DO 100 i = max(1,j-k),j - 1
351 temp = temp - dconjg(a(l+i,j))*x(i)
352 100 CONTINUE
353 IF (nounit) temp = temp/dconjg(a(kplus1,j))
354 END IF
355 x(j) = temp
356 110 CONTINUE
357 ELSE
358 jx = kx
359 DO 140 j = 1,n
360 temp = x(jx)
361 ix = kx
362 l = kplus1 - j
363 IF (noconj) THEN
364 DO 120 i = max(1,j-k),j - 1
365 temp = temp - a(l+i,j)*x(ix)
366 ix = ix + incx
367 120 CONTINUE
368 IF (nounit) temp = temp/a(kplus1,j)
369 ELSE
370 DO 130 i = max(1,j-k),j - 1
371 temp = temp - dconjg(a(l+i,j))*x(ix)
372 ix = ix + incx
373 130 CONTINUE
374 IF (nounit) temp = temp/dconjg(a(kplus1,j))
375 END IF
376 x(jx) = temp
377 jx = jx + incx
378 IF (j.GT.k) kx = kx + incx
379 140 CONTINUE
380 END IF
381 ELSE
382 IF (incx.EQ.1) THEN
383 DO 170 j = n,1,-1
384 temp = x(j)
385 l = 1 - j
386 IF (noconj) THEN
387 DO 150 i = min(n,j+k),j + 1,-1
388 temp = temp - a(l+i,j)*x(i)
389 150 CONTINUE
390 IF (nounit) temp = temp/a(1,j)
391 ELSE
392 DO 160 i = min(n,j+k),j + 1,-1
393 temp = temp - dconjg(a(l+i,j))*x(i)
394 160 CONTINUE
395 IF (nounit) temp = temp/dconjg(a(1,j))
396 END IF
397 x(j) = temp
398 170 CONTINUE
399 ELSE
400 kx = kx + (n-1)*incx
401 jx = kx
402 DO 200 j = n,1,-1
403 temp = x(jx)
404 ix = kx
405 l = 1 - j
406 IF (noconj) THEN
407 DO 180 i = min(n,j+k),j + 1,-1
408 temp = temp - a(l+i,j)*x(ix)
409 ix = ix - incx
410 180 CONTINUE
411 IF (nounit) temp = temp/a(1,j)
412 ELSE
413 DO 190 i = min(n,j+k),j + 1,-1
414 temp = temp - dconjg(a(l+i,j))*x(ix)
415 ix = ix - incx
416 190 CONTINUE
417 IF (nounit) temp = temp/dconjg(a(1,j))
418 END IF
419 x(jx) = temp
420 jx = jx - incx
421 IF ((n-j).GE.k) kx = kx - incx
422 200 CONTINUE
423 END IF
424 END IF
425 END IF
426*
427 RETURN
428*
429* End of ZTBSV
430*
subroutine xerbla(srname, info)
Definition cblat2.f:3285
logical function lsame(ca, cb)
LSAME
Definition lsame.f:48
Here is the call graph for this function:
Here is the caller graph for this function: