LAPACK 3.11.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 complex_blas_level2
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. .NOT.lsame(trans,'T') .AND.
190 + .NOT.lsame(trans,'C')) THEN
191 info = 2
192 ELSE IF (.NOT.lsame(diag,'U') .AND. .NOT.lsame(diag,'N')) THEN
193 info = 3
194 ELSE IF (n.LT.0) THEN
195 info = 4
196 ELSE IF (lda.LT.max(1,n)) THEN
197 info = 6
198 ELSE IF (incx.EQ.0) THEN
199 info = 8
200 END IF
201 IF (info.NE.0) THEN
202 CALL xerbla('CTRSV ',info)
203 RETURN
204 END IF
205*
206* Quick return if possible.
207*
208 IF (n.EQ.0) RETURN
209*
210 noconj = lsame(trans,'T')
211 nounit = lsame(diag,'N')
212*
213* Set up the start point in X if the increment is not unity. This
214* will be ( N - 1 )*INCX too small for descending loops.
215*
216 IF (incx.LE.0) THEN
217 kx = 1 - (n-1)*incx
218 ELSE IF (incx.NE.1) THEN
219 kx = 1
220 END IF
221*
222* Start the operations. In this version the elements of A are
223* accessed sequentially with one pass through A.
224*
225 IF (lsame(trans,'N')) THEN
226*
227* Form x := inv( A )*x.
228*
229 IF (lsame(uplo,'U')) THEN
230 IF (incx.EQ.1) THEN
231 DO 20 j = n,1,-1
232 IF (x(j).NE.zero) THEN
233 IF (nounit) x(j) = x(j)/a(j,j)
234 temp = x(j)
235 DO 10 i = j - 1,1,-1
236 x(i) = x(i) - temp*a(i,j)
237 10 CONTINUE
238 END IF
239 20 CONTINUE
240 ELSE
241 jx = kx + (n-1)*incx
242 DO 40 j = n,1,-1
243 IF (x(jx).NE.zero) THEN
244 IF (nounit) x(jx) = x(jx)/a(j,j)
245 temp = x(jx)
246 ix = jx
247 DO 30 i = j - 1,1,-1
248 ix = ix - incx
249 x(ix) = x(ix) - temp*a(i,j)
250 30 CONTINUE
251 END IF
252 jx = jx - incx
253 40 CONTINUE
254 END IF
255 ELSE
256 IF (incx.EQ.1) THEN
257 DO 60 j = 1,n
258 IF (x(j).NE.zero) THEN
259 IF (nounit) x(j) = x(j)/a(j,j)
260 temp = x(j)
261 DO 50 i = j + 1,n
262 x(i) = x(i) - temp*a(i,j)
263 50 CONTINUE
264 END IF
265 60 CONTINUE
266 ELSE
267 jx = kx
268 DO 80 j = 1,n
269 IF (x(jx).NE.zero) THEN
270 IF (nounit) x(jx) = x(jx)/a(j,j)
271 temp = x(jx)
272 ix = jx
273 DO 70 i = j + 1,n
274 ix = ix + incx
275 x(ix) = x(ix) - temp*a(i,j)
276 70 CONTINUE
277 END IF
278 jx = jx + incx
279 80 CONTINUE
280 END IF
281 END IF
282 ELSE
283*
284* Form x := inv( A**T )*x or x := inv( A**H )*x.
285*
286 IF (lsame(uplo,'U')) THEN
287 IF (incx.EQ.1) THEN
288 DO 110 j = 1,n
289 temp = x(j)
290 IF (noconj) THEN
291 DO 90 i = 1,j - 1
292 temp = temp - a(i,j)*x(i)
293 90 CONTINUE
294 IF (nounit) temp = temp/a(j,j)
295 ELSE
296 DO 100 i = 1,j - 1
297 temp = temp - conjg(a(i,j))*x(i)
298 100 CONTINUE
299 IF (nounit) temp = temp/conjg(a(j,j))
300 END IF
301 x(j) = temp
302 110 CONTINUE
303 ELSE
304 jx = kx
305 DO 140 j = 1,n
306 ix = kx
307 temp = x(jx)
308 IF (noconj) THEN
309 DO 120 i = 1,j - 1
310 temp = temp - a(i,j)*x(ix)
311 ix = ix + incx
312 120 CONTINUE
313 IF (nounit) temp = temp/a(j,j)
314 ELSE
315 DO 130 i = 1,j - 1
316 temp = temp - conjg(a(i,j))*x(ix)
317 ix = ix + incx
318 130 CONTINUE
319 IF (nounit) temp = temp/conjg(a(j,j))
320 END IF
321 x(jx) = temp
322 jx = jx + incx
323 140 CONTINUE
324 END IF
325 ELSE
326 IF (incx.EQ.1) THEN
327 DO 170 j = n,1,-1
328 temp = x(j)
329 IF (noconj) THEN
330 DO 150 i = n,j + 1,-1
331 temp = temp - a(i,j)*x(i)
332 150 CONTINUE
333 IF (nounit) temp = temp/a(j,j)
334 ELSE
335 DO 160 i = n,j + 1,-1
336 temp = temp - conjg(a(i,j))*x(i)
337 160 CONTINUE
338 IF (nounit) temp = temp/conjg(a(j,j))
339 END IF
340 x(j) = temp
341 170 CONTINUE
342 ELSE
343 kx = kx + (n-1)*incx
344 jx = kx
345 DO 200 j = n,1,-1
346 ix = kx
347 temp = x(jx)
348 IF (noconj) THEN
349 DO 180 i = n,j + 1,-1
350 temp = temp - a(i,j)*x(ix)
351 ix = ix - incx
352 180 CONTINUE
353 IF (nounit) temp = temp/a(j,j)
354 ELSE
355 DO 190 i = n,j + 1,-1
356 temp = temp - conjg(a(i,j))*x(ix)
357 ix = ix - incx
358 190 CONTINUE
359 IF (nounit) temp = temp/conjg(a(j,j))
360 END IF
361 x(jx) = temp
362 jx = jx - incx
363 200 CONTINUE
364 END IF
365 END IF
366 END IF
367*
368 RETURN
369*
370* End of CTRSV
371*
372 END
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:60
subroutine ctrsv(UPLO, TRANS, DIAG, N, A, LDA, X, INCX)
CTRSV
Definition: ctrsv.f:149