LAPACK  3.6.1
LAPACK: Linear Algebra PACKage
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 of 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 of 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 *> tranformed 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 *> \date November 2011
129 *
130 *> \ingroup complex_blas_level2
131 *
132 *> \par Further Details:
133 * =====================
134 *>
135 *> \verbatim
136 *>
137 *> Level 2 Blas routine.
138 *> The vector and matrix arguments are not referenced when N = 0, or M = 0
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 ctrmv(UPLO,TRANS,DIAG,N,A,LDA,X,INCX)
149 *
150 * -- Reference BLAS level2 routine (version 3.4.0) --
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 * November 2011
154 *
155 * .. Scalar Arguments ..
156  INTEGER INCX,LDA,N
157  CHARACTER DIAG,TRANS,UPLO
158 * ..
159 * .. Array Arguments ..
160  COMPLEX A(lda,*),X(*)
161 * ..
162 *
163 * =====================================================================
164 *
165 * .. Parameters ..
166  COMPLEX ZERO
167  parameter(zero= (0.0e+0,0.0e+0))
168 * ..
169 * .. Local Scalars ..
170  COMPLEX TEMP
171  INTEGER I,INFO,IX,J,JX,KX
172  LOGICAL NOCONJ,NOUNIT
173 * ..
174 * .. External Functions ..
175  LOGICAL LSAME
176  EXTERNAL lsame
177 * ..
178 * .. External Subroutines ..
179  EXTERNAL xerbla
180 * ..
181 * .. Intrinsic Functions ..
182  INTRINSIC conjg,max
183 * ..
184 *
185 * Test the input parameters.
186 *
187  info = 0
188  IF (.NOT.lsame(uplo,'U') .AND. .NOT.lsame(uplo,'L')) THEN
189  info = 1
190  ELSE IF (.NOT.lsame(trans,'N') .AND. .NOT.lsame(trans,'T') .AND.
191  + .NOT.lsame(trans,'C')) THEN
192  info = 2
193  ELSE IF (.NOT.lsame(diag,'U') .AND. .NOT.lsame(diag,'N')) THEN
194  info = 3
195  ELSE IF (n.LT.0) THEN
196  info = 4
197  ELSE IF (lda.LT.max(1,n)) THEN
198  info = 6
199  ELSE IF (incx.EQ.0) THEN
200  info = 8
201  END IF
202  IF (info.NE.0) THEN
203  CALL xerbla('CTRMV ',info)
204  RETURN
205  END IF
206 *
207 * Quick return if possible.
208 *
209  IF (n.EQ.0) RETURN
210 *
211  noconj = lsame(trans,'T')
212  nounit = lsame(diag,'N')
213 *
214 * Set up the start point in X if the increment is not unity. This
215 * will be ( N - 1 )*INCX too small for descending loops.
216 *
217  IF (incx.LE.0) THEN
218  kx = 1 - (n-1)*incx
219  ELSE IF (incx.NE.1) THEN
220  kx = 1
221  END IF
222 *
223 * Start the operations. In this version the elements of A are
224 * accessed sequentially with one pass through A.
225 *
226  IF (lsame(trans,'N')) THEN
227 *
228 * Form x := A*x.
229 *
230  IF (lsame(uplo,'U')) THEN
231  IF (incx.EQ.1) THEN
232  DO 20 j = 1,n
233  IF (x(j).NE.zero) THEN
234  temp = x(j)
235  DO 10 i = 1,j - 1
236  x(i) = x(i) + temp*a(i,j)
237  10 CONTINUE
238  IF (nounit) x(j) = x(j)*a(j,j)
239  END IF
240  20 CONTINUE
241  ELSE
242  jx = kx
243  DO 40 j = 1,n
244  IF (x(jx).NE.zero) THEN
245  temp = x(jx)
246  ix = kx
247  DO 30 i = 1,j - 1
248  x(ix) = x(ix) + temp*a(i,j)
249  ix = ix + incx
250  30 CONTINUE
251  IF (nounit) x(jx) = x(jx)*a(j,j)
252  END IF
253  jx = jx + incx
254  40 CONTINUE
255  END IF
256  ELSE
257  IF (incx.EQ.1) THEN
258  DO 60 j = n,1,-1
259  IF (x(j).NE.zero) THEN
260  temp = x(j)
261  DO 50 i = n,j + 1,-1
262  x(i) = x(i) + temp*a(i,j)
263  50 CONTINUE
264  IF (nounit) x(j) = x(j)*a(j,j)
265  END IF
266  60 CONTINUE
267  ELSE
268  kx = kx + (n-1)*incx
269  jx = kx
270  DO 80 j = n,1,-1
271  IF (x(jx).NE.zero) THEN
272  temp = x(jx)
273  ix = kx
274  DO 70 i = n,j + 1,-1
275  x(ix) = x(ix) + temp*a(i,j)
276  ix = ix - incx
277  70 CONTINUE
278  IF (nounit) x(jx) = x(jx)*a(j,j)
279  END IF
280  jx = jx - incx
281  80 CONTINUE
282  END IF
283  END IF
284  ELSE
285 *
286 * Form x := A**T*x or x := A**H*x.
287 *
288  IF (lsame(uplo,'U')) THEN
289  IF (incx.EQ.1) THEN
290  DO 110 j = n,1,-1
291  temp = x(j)
292  IF (noconj) THEN
293  IF (nounit) temp = temp*a(j,j)
294  DO 90 i = j - 1,1,-1
295  temp = temp + a(i,j)*x(i)
296  90 CONTINUE
297  ELSE
298  IF (nounit) temp = temp*conjg(a(j,j))
299  DO 100 i = j - 1,1,-1
300  temp = temp + conjg(a(i,j))*x(i)
301  100 CONTINUE
302  END IF
303  x(j) = temp
304  110 CONTINUE
305  ELSE
306  jx = kx + (n-1)*incx
307  DO 140 j = n,1,-1
308  temp = x(jx)
309  ix = jx
310  IF (noconj) THEN
311  IF (nounit) temp = temp*a(j,j)
312  DO 120 i = j - 1,1,-1
313  ix = ix - incx
314  temp = temp + a(i,j)*x(ix)
315  120 CONTINUE
316  ELSE
317  IF (nounit) temp = temp*conjg(a(j,j))
318  DO 130 i = j - 1,1,-1
319  ix = ix - incx
320  temp = temp + conjg(a(i,j))*x(ix)
321  130 CONTINUE
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 = 1,n
330  temp = x(j)
331  IF (noconj) THEN
332  IF (nounit) temp = temp*a(j,j)
333  DO 150 i = j + 1,n
334  temp = temp + a(i,j)*x(i)
335  150 CONTINUE
336  ELSE
337  IF (nounit) temp = temp*conjg(a(j,j))
338  DO 160 i = j + 1,n
339  temp = temp + conjg(a(i,j))*x(i)
340  160 CONTINUE
341  END IF
342  x(j) = temp
343  170 CONTINUE
344  ELSE
345  jx = kx
346  DO 200 j = 1,n
347  temp = x(jx)
348  ix = jx
349  IF (noconj) THEN
350  IF (nounit) temp = temp*a(j,j)
351  DO 180 i = j + 1,n
352  ix = ix + incx
353  temp = temp + a(i,j)*x(ix)
354  180 CONTINUE
355  ELSE
356  IF (nounit) temp = temp*conjg(a(j,j))
357  DO 190 i = j + 1,n
358  ix = ix + incx
359  temp = temp + conjg(a(i,j))*x(ix)
360  190 CONTINUE
361  END IF
362  x(jx) = temp
363  jx = jx + incx
364  200 CONTINUE
365  END IF
366  END IF
367  END IF
368 *
369  RETURN
370 *
371 * End of CTRMV .
372 *
373  END
subroutine xerbla(SRNAME, INFO)
XERBLA
Definition: xerbla.f:62
subroutine ctrmv(UPLO, TRANS, DIAG, N, A, LDA, X, INCX)
CTRMV
Definition: ctrmv.f:149