LAPACK 3.11.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ctbsv.f
Go to the documentation of this file.
1*> \brief \b CTBSV
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 CTBSV(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* COMPLEX A(LDA,*),X(*)
19* ..
20*
21*
22*> \par Purpose:
23* =============
24*>
25*> \verbatim
26*>
27*> CTBSV solves one of the systems of equations
28*>
29*> A*x = b, or A**T*x = b, or A**H*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**H*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 COMPLEX array, 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 COMPLEX array, 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*> \ingroup complex_blas_level2
172*
173*> \par Further Details:
174* =====================
175*>
176*> \verbatim
177*>
178*> Level 2 Blas routine.
179*>
180*> -- Written on 22-October-1986.
181*> Jack Dongarra, Argonne National Lab.
182*> Jeremy Du Croz, Nag Central Office.
183*> Sven Hammarling, Nag Central Office.
184*> Richard Hanson, Sandia National Labs.
185*> \endverbatim
186*>
187* =====================================================================
188 SUBROUTINE ctbsv(UPLO,TRANS,DIAG,N,K,A,LDA,X,INCX)
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 A(LDA,*),X(*)
200* ..
201*
202* =====================================================================
203*
204* .. Parameters ..
205 COMPLEX ZERO
206 parameter(zero= (0.0e+0,0.0e+0))
207* ..
208* .. Local Scalars ..
209 COMPLEX 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 conjg,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. .NOT.lsame(trans,'T') .AND.
230 + .NOT.lsame(trans,'C')) THEN
231 info = 2
232 ELSE IF (.NOT.lsame(diag,'U') .AND. .NOT.lsame(diag,'N')) THEN
233 info = 3
234 ELSE IF (n.LT.0) THEN
235 info = 4
236 ELSE IF (k.LT.0) THEN
237 info = 5
238 ELSE IF (lda.LT. (k+1)) THEN
239 info = 7
240 ELSE IF (incx.EQ.0) THEN
241 info = 9
242 END IF
243 IF (info.NE.0) THEN
244 CALL xerbla('CTBSV ',info)
245 RETURN
246 END IF
247*
248* Quick return if possible.
249*
250 IF (n.EQ.0) RETURN
251*
252 noconj = lsame(trans,'T')
253 nounit = lsame(diag,'N')
254*
255* Set up the start point in X if the increment is not unity. This
256* will be ( N - 1 )*INCX too small for descending loops.
257*
258 IF (incx.LE.0) THEN
259 kx = 1 - (n-1)*incx
260 ELSE IF (incx.NE.1) THEN
261 kx = 1
262 END IF
263*
264* Start the operations. In this version the elements of A are
265* accessed by sequentially with one pass through A.
266*
267 IF (lsame(trans,'N')) THEN
268*
269* Form x := inv( A )*x.
270*
271 IF (lsame(uplo,'U')) THEN
272 kplus1 = k + 1
273 IF (incx.EQ.1) THEN
274 DO 20 j = n,1,-1
275 IF (x(j).NE.zero) THEN
276 l = kplus1 - j
277 IF (nounit) x(j) = x(j)/a(kplus1,j)
278 temp = x(j)
279 DO 10 i = j - 1,max(1,j-k),-1
280 x(i) = x(i) - temp*a(l+i,j)
281 10 CONTINUE
282 END IF
283 20 CONTINUE
284 ELSE
285 kx = kx + (n-1)*incx
286 jx = kx
287 DO 40 j = n,1,-1
288 kx = kx - incx
289 IF (x(jx).NE.zero) THEN
290 ix = kx
291 l = kplus1 - j
292 IF (nounit) x(jx) = x(jx)/a(kplus1,j)
293 temp = x(jx)
294 DO 30 i = j - 1,max(1,j-k),-1
295 x(ix) = x(ix) - temp*a(l+i,j)
296 ix = ix - incx
297 30 CONTINUE
298 END IF
299 jx = jx - incx
300 40 CONTINUE
301 END IF
302 ELSE
303 IF (incx.EQ.1) THEN
304 DO 60 j = 1,n
305 IF (x(j).NE.zero) THEN
306 l = 1 - j
307 IF (nounit) x(j) = x(j)/a(1,j)
308 temp = x(j)
309 DO 50 i = j + 1,min(n,j+k)
310 x(i) = x(i) - temp*a(l+i,j)
311 50 CONTINUE
312 END IF
313 60 CONTINUE
314 ELSE
315 jx = kx
316 DO 80 j = 1,n
317 kx = kx + incx
318 IF (x(jx).NE.zero) THEN
319 ix = kx
320 l = 1 - j
321 IF (nounit) x(jx) = x(jx)/a(1,j)
322 temp = x(jx)
323 DO 70 i = j + 1,min(n,j+k)
324 x(ix) = x(ix) - temp*a(l+i,j)
325 ix = ix + incx
326 70 CONTINUE
327 END IF
328 jx = jx + incx
329 80 CONTINUE
330 END IF
331 END IF
332 ELSE
333*
334* Form x := inv( A**T )*x or x := inv( A**H )*x.
335*
336 IF (lsame(uplo,'U')) THEN
337 kplus1 = k + 1
338 IF (incx.EQ.1) THEN
339 DO 110 j = 1,n
340 temp = x(j)
341 l = kplus1 - j
342 IF (noconj) THEN
343 DO 90 i = max(1,j-k),j - 1
344 temp = temp - a(l+i,j)*x(i)
345 90 CONTINUE
346 IF (nounit) temp = temp/a(kplus1,j)
347 ELSE
348 DO 100 i = max(1,j-k),j - 1
349 temp = temp - conjg(a(l+i,j))*x(i)
350 100 CONTINUE
351 IF (nounit) temp = temp/conjg(a(kplus1,j))
352 END IF
353 x(j) = temp
354 110 CONTINUE
355 ELSE
356 jx = kx
357 DO 140 j = 1,n
358 temp = x(jx)
359 ix = kx
360 l = kplus1 - j
361 IF (noconj) THEN
362 DO 120 i = max(1,j-k),j - 1
363 temp = temp - a(l+i,j)*x(ix)
364 ix = ix + incx
365 120 CONTINUE
366 IF (nounit) temp = temp/a(kplus1,j)
367 ELSE
368 DO 130 i = max(1,j-k),j - 1
369 temp = temp - conjg(a(l+i,j))*x(ix)
370 ix = ix + incx
371 130 CONTINUE
372 IF (nounit) temp = temp/conjg(a(kplus1,j))
373 END IF
374 x(jx) = temp
375 jx = jx + incx
376 IF (j.GT.k) kx = kx + incx
377 140 CONTINUE
378 END IF
379 ELSE
380 IF (incx.EQ.1) THEN
381 DO 170 j = n,1,-1
382 temp = x(j)
383 l = 1 - j
384 IF (noconj) THEN
385 DO 150 i = min(n,j+k),j + 1,-1
386 temp = temp - a(l+i,j)*x(i)
387 150 CONTINUE
388 IF (nounit) temp = temp/a(1,j)
389 ELSE
390 DO 160 i = min(n,j+k),j + 1,-1
391 temp = temp - conjg(a(l+i,j))*x(i)
392 160 CONTINUE
393 IF (nounit) temp = temp/conjg(a(1,j))
394 END IF
395 x(j) = temp
396 170 CONTINUE
397 ELSE
398 kx = kx + (n-1)*incx
399 jx = kx
400 DO 200 j = n,1,-1
401 temp = x(jx)
402 ix = kx
403 l = 1 - j
404 IF (noconj) THEN
405 DO 180 i = min(n,j+k),j + 1,-1
406 temp = temp - a(l+i,j)*x(ix)
407 ix = ix - incx
408 180 CONTINUE
409 IF (nounit) temp = temp/a(1,j)
410 ELSE
411 DO 190 i = min(n,j+k),j + 1,-1
412 temp = temp - conjg(a(l+i,j))*x(ix)
413 ix = ix - incx
414 190 CONTINUE
415 IF (nounit) temp = temp/conjg(a(1,j))
416 END IF
417 x(jx) = temp
418 jx = jx - incx
419 IF ((n-j).GE.k) kx = kx - incx
420 200 CONTINUE
421 END IF
422 END IF
423 END IF
424*
425 RETURN
426*
427* End of CTBSV
428*
429 END
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine ctbsv(UPLO, TRANS, DIAG, N, K, A, LDA, X, INCX)
CTBSV
Definition: ctbsv.f:189