LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
ctrsv.f
Go to the documentation of this file.
1*> \brief \b CTRSV
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 CTRSV(UPLO,TRANS,DIAG,N,A,LDA,X,INCX)
12*
13* .. Scalar Arguments ..
14* INTEGER INCX,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*> CTRSV 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 matrix.
33*>
34*> No test for singularity or near-singularity is included in this
35*> routine. Such tests must be performed before calling this routine.
36*> \endverbatim
37*
38* Arguments:
39* ==========
40*
41*> \param[in] UPLO
42*> \verbatim
43*> UPLO is CHARACTER*1
44*> On entry, UPLO specifies whether the matrix is an upper or
45*> lower triangular matrix as follows:
46*>
47*> UPLO = 'U' or 'u' A is an upper triangular matrix.
48*>
49*> UPLO = 'L' or 'l' A is a lower triangular matrix.
50*> \endverbatim
51*>
52*> \param[in] TRANS
53*> \verbatim
54*> TRANS is CHARACTER*1
55*> On entry, TRANS specifies the equations to be solved as
56*> follows:
57*>
58*> TRANS = 'N' or 'n' A*x = b.
59*>
60*> TRANS = 'T' or 't' A**T*x = b.
61*>
62*> TRANS = 'C' or 'c' A**H*x = b.
63*> \endverbatim
64*>
65*> \param[in] DIAG
66*> \verbatim
67*> DIAG is CHARACTER*1
68*> On entry, DIAG specifies whether or not A is unit
69*> triangular as follows:
70*>
71*> DIAG = 'U' or 'u' A is assumed to be unit triangular.
72*>
73*> DIAG = 'N' or 'n' A is not assumed to be unit
74*> triangular.
75*> \endverbatim
76*>
77*> \param[in] N
78*> \verbatim
79*> N is INTEGER
80*> On entry, N specifies the order of the matrix A.
81*> N must be at least zero.
82*> \endverbatim
83*>
84*> \param[in] A
85*> \verbatim
86*> A is COMPLEX array, dimension ( LDA, N )
87*> Before entry with UPLO = 'U' or 'u', the leading n by n
88*> upper triangular part of the array A must contain the upper
89*> triangular matrix and the strictly lower triangular part of
90*> A is not referenced.
91*> Before entry with UPLO = 'L' or 'l', the leading n by n
92*> lower triangular part of the array A must contain the lower
93*> triangular matrix and the strictly upper triangular part of
94*> A is not referenced.
95*> Note that when DIAG = 'U' or 'u', the diagonal elements of
96*> A are not referenced either, but are assumed to be unity.
97*> \endverbatim
98*>
99*> \param[in] LDA
100*> \verbatim
101*> LDA is INTEGER
102*> On entry, LDA specifies the first dimension of A as declared
103*> in the calling (sub) program. LDA must be at least
104*> max( 1, n ).
105*> \endverbatim
106*>
107*> \param[in,out] X
108*> \verbatim
109*> X is COMPLEX array, dimension at least
110*> ( 1 + ( n - 1 )*abs( INCX ) ).
111*> Before entry, the incremented array X must contain the n
112*> element right-hand side vector b. On exit, X is overwritten
113*> with the solution vector x.
114*> \endverbatim
115*>
116*> \param[in] INCX
117*> \verbatim
118*> INCX is INTEGER
119*> On entry, INCX specifies the increment for the elements of
120*> X. INCX must not be zero.
121*> \endverbatim
122*
123* Authors:
124* ========
125*
126*> \author Univ. of Tennessee
127*> \author Univ. of California Berkeley
128*> \author Univ. of Colorado Denver
129*> \author NAG Ltd.
130*
131*> \ingroup trsv
132*
133*> \par Further Details:
134* =====================
135*>
136*> \verbatim
137*>
138*> Level 2 Blas routine.
139*>
140*> -- Written on 22-October-1986.
141*> Jack Dongarra, Argonne National Lab.
142*> Jeremy Du Croz, Nag Central Office.
143*> Sven Hammarling, Nag Central Office.
144*> Richard Hanson, Sandia National Labs.
145*> \endverbatim
146*>
147* =====================================================================
148 SUBROUTINE ctrsv(UPLO,TRANS,DIAG,N,A,LDA,X,INCX)
149*
150* -- Reference BLAS level2 routine --
151* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
152* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
153*
154* .. Scalar Arguments ..
155 INTEGER INCX,LDA,N
156 CHARACTER DIAG,TRANS,UPLO
157* ..
158* .. Array Arguments ..
159 COMPLEX A(LDA,*),X(*)
160* ..
161*
162* =====================================================================
163*
164* .. Parameters ..
165 COMPLEX ZERO
166 parameter(zero= (0.0e+0,0.0e+0))
167* ..
168* .. Local Scalars ..
169 COMPLEX TEMP
170 INTEGER I,INFO,IX,J,JX,KX
171 LOGICAL NOCONJ,NOUNIT
172* ..
173* .. External Functions ..
174 LOGICAL LSAME
175 EXTERNAL lsame
176* ..
177* .. External Subroutines ..
178 EXTERNAL xerbla
179* ..
180* .. Intrinsic Functions ..
181 INTRINSIC conjg,max
182* ..
183*
184* Test the input parameters.
185*
186 info = 0
187 IF (.NOT.lsame(uplo,'U') .AND. .NOT.lsame(uplo,'L')) THEN
188 info = 1
189 ELSE IF (.NOT.lsame(trans,'N') .AND.
190 + .NOT.lsame(trans,'T') .AND.
191 + .NOT.lsame(trans,'C')) THEN
192 info = 2
193 ELSE IF (.NOT.lsame(diag,'U') .AND.
194 + .NOT.lsame(diag,'N')) THEN
195 info = 3
196 ELSE IF (n.LT.0) THEN
197 info = 4
198 ELSE IF (lda.LT.max(1,n)) THEN
199 info = 6
200 ELSE IF (incx.EQ.0) THEN
201 info = 8
202 END IF
203 IF (info.NE.0) THEN
204 CALL xerbla('CTRSV ',info)
205 RETURN
206 END IF
207*
208* Quick return if possible.
209*
210 IF (n.EQ.0) RETURN
211*
212 noconj = lsame(trans,'T')
213 nounit = lsame(diag,'N')
214*
215* Set up the start point in X if the increment is not unity. This
216* will be ( N - 1 )*INCX too small for descending loops.
217*
218 IF (incx.LE.0) THEN
219 kx = 1 - (n-1)*incx
220 ELSE IF (incx.NE.1) THEN
221 kx = 1
222 END IF
223*
224* Start the operations. In this version the elements of A are
225* accessed sequentially with one pass through A.
226*
227 IF (lsame(trans,'N')) THEN
228*
229* Form x := inv( A )*x.
230*
231 IF (lsame(uplo,'U')) THEN
232 IF (incx.EQ.1) THEN
233 DO 20 j = n,1,-1
234 IF (x(j).NE.zero) THEN
235 IF (nounit) x(j) = x(j)/a(j,j)
236 temp = x(j)
237 DO 10 i = j - 1,1,-1
238 x(i) = x(i) - temp*a(i,j)
239 10 CONTINUE
240 END IF
241 20 CONTINUE
242 ELSE
243 jx = kx + (n-1)*incx
244 DO 40 j = n,1,-1
245 IF (x(jx).NE.zero) THEN
246 IF (nounit) x(jx) = x(jx)/a(j,j)
247 temp = x(jx)
248 ix = jx
249 DO 30 i = j - 1,1,-1
250 ix = ix - incx
251 x(ix) = x(ix) - temp*a(i,j)
252 30 CONTINUE
253 END IF
254 jx = jx - incx
255 40 CONTINUE
256 END IF
257 ELSE
258 IF (incx.EQ.1) THEN
259 DO 60 j = 1,n
260 IF (x(j).NE.zero) THEN
261 IF (nounit) x(j) = x(j)/a(j,j)
262 temp = x(j)
263 DO 50 i = j + 1,n
264 x(i) = x(i) - temp*a(i,j)
265 50 CONTINUE
266 END IF
267 60 CONTINUE
268 ELSE
269 jx = kx
270 DO 80 j = 1,n
271 IF (x(jx).NE.zero) THEN
272 IF (nounit) x(jx) = x(jx)/a(j,j)
273 temp = x(jx)
274 ix = jx
275 DO 70 i = j + 1,n
276 ix = ix + incx
277 x(ix) = x(ix) - temp*a(i,j)
278 70 CONTINUE
279 END IF
280 jx = jx + incx
281 80 CONTINUE
282 END IF
283 END IF
284 ELSE
285*
286* Form x := inv( A**T )*x or x := inv( A**H )*x.
287*
288 IF (lsame(uplo,'U')) THEN
289 IF (incx.EQ.1) THEN
290 DO 110 j = 1,n
291 temp = x(j)
292 IF (noconj) THEN
293 DO 90 i = 1,j - 1
294 temp = temp - a(i,j)*x(i)
295 90 CONTINUE
296 IF (nounit) temp = temp/a(j,j)
297 ELSE
298 DO 100 i = 1,j - 1
299 temp = temp - conjg(a(i,j))*x(i)
300 100 CONTINUE
301 IF (nounit) temp = temp/conjg(a(j,j))
302 END IF
303 x(j) = temp
304 110 CONTINUE
305 ELSE
306 jx = kx
307 DO 140 j = 1,n
308 ix = kx
309 temp = x(jx)
310 IF (noconj) THEN
311 DO 120 i = 1,j - 1
312 temp = temp - a(i,j)*x(ix)
313 ix = ix + incx
314 120 CONTINUE
315 IF (nounit) temp = temp/a(j,j)
316 ELSE
317 DO 130 i = 1,j - 1
318 temp = temp - conjg(a(i,j))*x(ix)
319 ix = ix + incx
320 130 CONTINUE
321 IF (nounit) temp = temp/conjg(a(j,j))
322 END IF
323 x(jx) = temp
324 jx = jx + incx
325 140 CONTINUE
326 END IF
327 ELSE
328 IF (incx.EQ.1) THEN
329 DO 170 j = n,1,-1
330 temp = x(j)
331 IF (noconj) THEN
332 DO 150 i = n,j + 1,-1
333 temp = temp - a(i,j)*x(i)
334 150 CONTINUE
335 IF (nounit) temp = temp/a(j,j)
336 ELSE
337 DO 160 i = n,j + 1,-1
338 temp = temp - conjg(a(i,j))*x(i)
339 160 CONTINUE
340 IF (nounit) temp = temp/conjg(a(j,j))
341 END IF
342 x(j) = temp
343 170 CONTINUE
344 ELSE
345 kx = kx + (n-1)*incx
346 jx = kx
347 DO 200 j = n,1,-1
348 ix = kx
349 temp = x(jx)
350 IF (noconj) THEN
351 DO 180 i = n,j + 1,-1
352 temp = temp - a(i,j)*x(ix)
353 ix = ix - incx
354 180 CONTINUE
355 IF (nounit) temp = temp/a(j,j)
356 ELSE
357 DO 190 i = n,j + 1,-1
358 temp = temp - conjg(a(i,j))*x(ix)
359 ix = ix - incx
360 190 CONTINUE
361 IF (nounit) temp = temp/conjg(a(j,j))
362 END IF
363 x(jx) = temp
364 jx = jx - incx
365 200 CONTINUE
366 END IF
367 END IF
368 END IF
369*
370 RETURN
371*
372* End of CTRSV
373*
374 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine ctrsv(uplo, trans, diag, n, a, lda, x, incx)
CTRSV
Definition ctrsv.f:149