LAPACK  3.4.2
LAPACK: Linear Algebra PACKage
 All Files Functions Groups
dtbsv.f
Go to the documentation of this file.
1 *> \brief \b DTBSV
2 *
3 * =========== DOCUMENTATION ===========
4 *
5 * Online html documentation available at
6 * http://www.netlib.org/lapack/explore-html/
7 *
8 * Definition:
9 * ===========
10 *
11 * SUBROUTINE DTBSV(UPLO,TRANS,DIAG,N,K,A,LDA,X,INCX)
12 *
13 * .. Scalar Arguments ..
14 * INTEGER INCX,K,LDA,N
15 * CHARACTER DIAG,TRANS,UPLO
16 * ..
17 * .. Array Arguments ..
18 * DOUBLE PRECISION A(LDA,*),X(*)
19 * ..
20 *
21 *
22 *> \par Purpose:
23 * =============
24 *>
25 *> \verbatim
26 *>
27 *> DTBSV solves one of the systems of equations
28 *>
29 *> A*x = b, or A**T*x = b,
30 *>
31 *> where b and x are n element vectors and A is an n by n unit, or
32 *> non-unit, upper or lower triangular band matrix, with ( k + 1 )
33 *> diagonals.
34 *>
35 *> No test for singularity or near-singularity is included in this
36 *> routine. Such tests must be performed before calling this routine.
37 *> \endverbatim
38 *
39 * Arguments:
40 * ==========
41 *
42 *> \param[in] UPLO
43 *> \verbatim
44 *> UPLO is CHARACTER*1
45 *> On entry, UPLO specifies whether the matrix is an upper or
46 *> lower triangular matrix as follows:
47 *>
48 *> UPLO = 'U' or 'u' A is an upper triangular matrix.
49 *>
50 *> UPLO = 'L' or 'l' A is a lower triangular matrix.
51 *> \endverbatim
52 *>
53 *> \param[in] TRANS
54 *> \verbatim
55 *> TRANS is CHARACTER*1
56 *> On entry, TRANS specifies the equations to be solved as
57 *> follows:
58 *>
59 *> TRANS = 'N' or 'n' A*x = b.
60 *>
61 *> TRANS = 'T' or 't' A**T*x = b.
62 *>
63 *> TRANS = 'C' or 'c' A**T*x = b.
64 *> \endverbatim
65 *>
66 *> \param[in] DIAG
67 *> \verbatim
68 *> DIAG is CHARACTER*1
69 *> On entry, DIAG specifies whether or not A is unit
70 *> triangular as follows:
71 *>
72 *> DIAG = 'U' or 'u' A is assumed to be unit triangular.
73 *>
74 *> DIAG = 'N' or 'n' A is not assumed to be unit
75 *> triangular.
76 *> \endverbatim
77 *>
78 *> \param[in] N
79 *> \verbatim
80 *> N is INTEGER
81 *> On entry, N specifies the order of the matrix A.
82 *> N must be at least zero.
83 *> \endverbatim
84 *>
85 *> \param[in] K
86 *> \verbatim
87 *> K is INTEGER
88 *> On entry with UPLO = 'U' or 'u', K specifies the number of
89 *> super-diagonals of the matrix A.
90 *> On entry with UPLO = 'L' or 'l', K specifies the number of
91 *> sub-diagonals of the matrix A.
92 *> K must satisfy 0 .le. K.
93 *> \endverbatim
94 *>
95 *> \param[in] A
96 *> \verbatim
97 *> A is DOUBLE PRECISION array of DIMENSION ( LDA, n ).
98 *> Before entry with UPLO = 'U' or 'u', the leading ( k + 1 )
99 *> by n part of the array A must contain the upper triangular
100 *> band part of the matrix of coefficients, supplied column by
101 *> column, with the leading diagonal of the matrix in row
102 *> ( k + 1 ) of the array, the first super-diagonal starting at
103 *> position 2 in row k, and so on. The top left k by k triangle
104 *> of the array A is not referenced.
105 *> The following program segment will transfer an upper
106 *> triangular band matrix from conventional full matrix storage
107 *> to band storage:
108 *>
109 *> DO 20, J = 1, N
110 *> M = K + 1 - J
111 *> DO 10, I = MAX( 1, J - K ), J
112 *> A( M + I, J ) = matrix( I, J )
113 *> 10 CONTINUE
114 *> 20 CONTINUE
115 *>
116 *> Before entry with UPLO = 'L' or 'l', the leading ( k + 1 )
117 *> by n part of the array A must contain the lower triangular
118 *> band part of the matrix of coefficients, supplied column by
119 *> column, with the leading diagonal of the matrix in row 1 of
120 *> the array, the first sub-diagonal starting at position 1 in
121 *> row 2, and so on. The bottom right k by k triangle of the
122 *> array A is not referenced.
123 *> The following program segment will transfer a lower
124 *> triangular band matrix from conventional full matrix storage
125 *> to band storage:
126 *>
127 *> DO 20, J = 1, N
128 *> M = 1 - J
129 *> DO 10, I = J, MIN( N, J + K )
130 *> A( M + I, J ) = matrix( I, J )
131 *> 10 CONTINUE
132 *> 20 CONTINUE
133 *>
134 *> Note that when DIAG = 'U' or 'u' the elements of the array A
135 *> corresponding to the diagonal elements of the matrix are not
136 *> referenced, but are assumed to be unity.
137 *> \endverbatim
138 *>
139 *> \param[in] LDA
140 *> \verbatim
141 *> LDA is INTEGER
142 *> On entry, LDA specifies the first dimension of A as declared
143 *> in the calling (sub) program. LDA must be at least
144 *> ( k + 1 ).
145 *> \endverbatim
146 *>
147 *> \param[in,out] X
148 *> \verbatim
149 *> X is DOUBLE PRECISION array of dimension at least
150 *> ( 1 + ( n - 1 )*abs( INCX ) ).
151 *> Before entry, the incremented array X must contain the n
152 *> element right-hand side vector b. On exit, X is overwritten
153 *> with the solution vector x.
154 *> \endverbatim
155 *>
156 *> \param[in] INCX
157 *> \verbatim
158 *> INCX is INTEGER
159 *> On entry, INCX specifies the increment for the elements of
160 *> X. INCX must not be zero.
161 *> \endverbatim
162 *
163 * Authors:
164 * ========
165 *
166 *> \author Univ. of Tennessee
167 *> \author Univ. of California Berkeley
168 *> \author Univ. of Colorado Denver
169 *> \author NAG Ltd.
170 *
171 *> \date November 2011
172 *
173 *> \ingroup double_blas_level2
174 *
175 *> \par Further Details:
176 * =====================
177 *>
178 *> \verbatim
179 *>
180 *> Level 2 Blas routine.
181 *>
182 *> -- Written on 22-October-1986.
183 *> Jack Dongarra, Argonne National Lab.
184 *> Jeremy Du Croz, Nag Central Office.
185 *> Sven Hammarling, Nag Central Office.
186 *> Richard Hanson, Sandia National Labs.
187 *> \endverbatim
188 *>
189 * =====================================================================
190  SUBROUTINE dtbsv(UPLO,TRANS,DIAG,N,K,A,LDA,X,INCX)
191 *
192 * -- Reference BLAS level2 routine (version 3.4.0) --
193 * -- Reference BLAS is a software package provided by Univ. of Tennessee, --
194 * -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
195 * November 2011
196 *
197 * .. Scalar Arguments ..
198  INTEGER incx,k,lda,n
199  CHARACTER diag,trans,uplo
200 * ..
201 * .. Array Arguments ..
202  DOUBLE PRECISION a(lda,*),x(*)
203 * ..
204 *
205 * =====================================================================
206 *
207 * .. Parameters ..
208  DOUBLE PRECISION zero
209  parameter(zero=0.0d+0)
210 * ..
211 * .. Local Scalars ..
212  DOUBLE PRECISION temp
213  INTEGER i,info,ix,j,jx,kplus1,kx,l
214  LOGICAL nounit
215 * ..
216 * .. External Functions ..
217  LOGICAL lsame
218  EXTERNAL lsame
219 * ..
220 * .. External Subroutines ..
221  EXTERNAL xerbla
222 * ..
223 * .. Intrinsic Functions ..
224  INTRINSIC max,min
225 * ..
226 *
227 * Test the input parameters.
228 *
229  info = 0
230  IF (.NOT.lsame(uplo,'U') .AND. .NOT.lsame(uplo,'L')) THEN
231  info = 1
232  ELSE IF (.NOT.lsame(trans,'N') .AND. .NOT.lsame(trans,'T') .AND.
233  + .NOT.lsame(trans,'C')) THEN
234  info = 2
235  ELSE IF (.NOT.lsame(diag,'U') .AND. .NOT.lsame(diag,'N')) THEN
236  info = 3
237  ELSE IF (n.LT.0) THEN
238  info = 4
239  ELSE IF (k.LT.0) THEN
240  info = 5
241  ELSE IF (lda.LT. (k+1)) THEN
242  info = 7
243  ELSE IF (incx.EQ.0) THEN
244  info = 9
245  END IF
246  IF (info.NE.0) THEN
247  CALL xerbla('DTBSV ',info)
248  return
249  END IF
250 *
251 * Quick return if possible.
252 *
253  IF (n.EQ.0) return
254 *
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.
337 *
338  IF (lsame(uplo,'U')) THEN
339  kplus1 = k + 1
340  IF (incx.EQ.1) THEN
341  DO 100 j = 1,n
342  temp = x(j)
343  l = kplus1 - j
344  DO 90 i = max(1,j-k),j - 1
345  temp = temp - a(l+i,j)*x(i)
346  90 continue
347  IF (nounit) temp = temp/a(kplus1,j)
348  x(j) = temp
349  100 continue
350  ELSE
351  jx = kx
352  DO 120 j = 1,n
353  temp = x(jx)
354  ix = kx
355  l = kplus1 - j
356  DO 110 i = max(1,j-k),j - 1
357  temp = temp - a(l+i,j)*x(ix)
358  ix = ix + incx
359  110 continue
360  IF (nounit) temp = temp/a(kplus1,j)
361  x(jx) = temp
362  jx = jx + incx
363  IF (j.GT.k) kx = kx + incx
364  120 continue
365  END IF
366  ELSE
367  IF (incx.EQ.1) THEN
368  DO 140 j = n,1,-1
369  temp = x(j)
370  l = 1 - j
371  DO 130 i = min(n,j+k),j + 1,-1
372  temp = temp - a(l+i,j)*x(i)
373  130 continue
374  IF (nounit) temp = temp/a(1,j)
375  x(j) = temp
376  140 continue
377  ELSE
378  kx = kx + (n-1)*incx
379  jx = kx
380  DO 160 j = n,1,-1
381  temp = x(jx)
382  ix = kx
383  l = 1 - j
384  DO 150 i = min(n,j+k),j + 1,-1
385  temp = temp - a(l+i,j)*x(ix)
386  ix = ix - incx
387  150 continue
388  IF (nounit) temp = temp/a(1,j)
389  x(jx) = temp
390  jx = jx - incx
391  IF ((n-j).GE.k) kx = kx - incx
392  160 continue
393  END IF
394  END IF
395  END IF
396 *
397  return
398 *
399 * End of DTBSV .
400 *
401  END