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