LAPACK 3.12.0
LAPACK: Linear Algebra PACKage
Loading...
Searching...
No Matches
csyr2k.f
Go to the documentation of this file.
1*> \brief \b CSYR2K
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 CSYR2K(UPLO,TRANS,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
12*
13* .. Scalar Arguments ..
14* COMPLEX ALPHA,BETA
15* INTEGER K,LDA,LDB,LDC,N
16* CHARACTER TRANS,UPLO
17* ..
18* .. Array Arguments ..
19* COMPLEX A(LDA,*),B(LDB,*),C(LDC,*)
20* ..
21*
22*
23*> \par Purpose:
24* =============
25*>
26*> \verbatim
27*>
28*> CSYR2K performs one of the symmetric rank 2k operations
29*>
30*> C := alpha*A*B**T + alpha*B*A**T + beta*C,
31*>
32*> or
33*>
34*> C := alpha*A**T*B + alpha*B**T*A + beta*C,
35*>
36*> where alpha and beta are scalars, C is an n by n symmetric matrix
37*> and A and B are n by k matrices in the first case and k by n
38*> matrices in the second case.
39*> \endverbatim
40*
41* Arguments:
42* ==========
43*
44*> \param[in] UPLO
45*> \verbatim
46*> UPLO is CHARACTER*1
47*> On entry, UPLO specifies whether the upper or lower
48*> triangular part of the array C is to be referenced as
49*> follows:
50*>
51*> UPLO = 'U' or 'u' Only the upper triangular part of C
52*> is to be referenced.
53*>
54*> UPLO = 'L' or 'l' Only the lower triangular part of C
55*> is to be referenced.
56*> \endverbatim
57*>
58*> \param[in] TRANS
59*> \verbatim
60*> TRANS is CHARACTER*1
61*> On entry, TRANS specifies the operation to be performed as
62*> follows:
63*>
64*> TRANS = 'N' or 'n' C := alpha*A*B**T + alpha*B*A**T +
65*> beta*C.
66*>
67*> TRANS = 'T' or 't' C := alpha*A**T*B + alpha*B**T*A +
68*> beta*C.
69*> \endverbatim
70*>
71*> \param[in] N
72*> \verbatim
73*> N is INTEGER
74*> On entry, N specifies the order of the matrix C. N must be
75*> at least zero.
76*> \endverbatim
77*>
78*> \param[in] K
79*> \verbatim
80*> K is INTEGER
81*> On entry with TRANS = 'N' or 'n', K specifies the number
82*> of columns of the matrices A and B, and on entry with
83*> TRANS = 'T' or 't', K specifies the number of rows of the
84*> matrices A and B. K must be at least zero.
85*> \endverbatim
86*>
87*> \param[in] ALPHA
88*> \verbatim
89*> ALPHA is COMPLEX
90*> On entry, ALPHA specifies the scalar alpha.
91*> \endverbatim
92*>
93*> \param[in] A
94*> \verbatim
95*> A is COMPLEX array, dimension ( LDA, ka ), where ka is
96*> k when TRANS = 'N' or 'n', and is n otherwise.
97*> Before entry with TRANS = 'N' or 'n', the leading n by k
98*> part of the array A must contain the matrix A, otherwise
99*> the leading k by n part of the array A must contain the
100*> matrix A.
101*> \endverbatim
102*>
103*> \param[in] LDA
104*> \verbatim
105*> LDA is INTEGER
106*> On entry, LDA specifies the first dimension of A as declared
107*> in the calling (sub) program. When TRANS = 'N' or 'n'
108*> then LDA must be at least max( 1, n ), otherwise LDA must
109*> be at least max( 1, k ).
110*> \endverbatim
111*>
112*> \param[in] B
113*> \verbatim
114*> B is COMPLEX array, dimension ( LDB, kb ), where kb is
115*> k when TRANS = 'N' or 'n', and is n otherwise.
116*> Before entry with TRANS = 'N' or 'n', the leading n by k
117*> part of the array B must contain the matrix B, otherwise
118*> the leading k by n part of the array B must contain the
119*> matrix B.
120*> \endverbatim
121*>
122*> \param[in] LDB
123*> \verbatim
124*> LDB is INTEGER
125*> On entry, LDB specifies the first dimension of B as declared
126*> in the calling (sub) program. When TRANS = 'N' or 'n'
127*> then LDB must be at least max( 1, n ), otherwise LDB must
128*> be at least max( 1, k ).
129*> \endverbatim
130*>
131*> \param[in] BETA
132*> \verbatim
133*> BETA is COMPLEX
134*> On entry, BETA specifies the scalar beta.
135*> \endverbatim
136*>
137*> \param[in,out] C
138*> \verbatim
139*> C is COMPLEX array, dimension ( LDC, N )
140*> Before entry with UPLO = 'U' or 'u', the leading n by n
141*> upper triangular part of the array C must contain the upper
142*> triangular part of the symmetric matrix and the strictly
143*> lower triangular part of C is not referenced. On exit, the
144*> upper triangular part of the array C is overwritten by the
145*> upper triangular part of the updated matrix.
146*> Before entry with UPLO = 'L' or 'l', the leading n by n
147*> lower triangular part of the array C must contain the lower
148*> triangular part of the symmetric matrix and the strictly
149*> upper triangular part of C is not referenced. On exit, the
150*> lower triangular part of the array C is overwritten by the
151*> lower triangular part of the updated matrix.
152*> \endverbatim
153*>
154*> \param[in] LDC
155*> \verbatim
156*> LDC is INTEGER
157*> On entry, LDC specifies the first dimension of C as declared
158*> in the calling (sub) program. LDC must be at least
159*> max( 1, n ).
160*> \endverbatim
161*
162* Authors:
163* ========
164*
165*> \author Univ. of Tennessee
166*> \author Univ. of California Berkeley
167*> \author Univ. of Colorado Denver
168*> \author NAG Ltd.
169*
170*> \ingroup her2k
171*
172*> \par Further Details:
173* =====================
174*>
175*> \verbatim
176*>
177*> Level 3 Blas routine.
178*>
179*> -- Written on 8-February-1989.
180*> Jack Dongarra, Argonne National Laboratory.
181*> Iain Duff, AERE Harwell.
182*> Jeremy Du Croz, Numerical Algorithms Group Ltd.
183*> Sven Hammarling, Numerical Algorithms Group Ltd.
184*> \endverbatim
185*>
186* =====================================================================
187 SUBROUTINE csyr2k(UPLO,TRANS,N,K,ALPHA,A,LDA,B,LDB,BETA,C,LDC)
188*
189* -- Reference BLAS level3 routine --
190* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
191* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
192*
193* .. Scalar Arguments ..
194 COMPLEX ALPHA,BETA
195 INTEGER K,LDA,LDB,LDC,N
196 CHARACTER TRANS,UPLO
197* ..
198* .. Array Arguments ..
199 COMPLEX A(LDA,*),B(LDB,*),C(LDC,*)
200* ..
201*
202* =====================================================================
203*
204* .. External Functions ..
205 LOGICAL LSAME
206 EXTERNAL lsame
207* ..
208* .. External Subroutines ..
209 EXTERNAL xerbla
210* ..
211* .. Intrinsic Functions ..
212 INTRINSIC max
213* ..
214* .. Local Scalars ..
215 COMPLEX TEMP1,TEMP2
216 INTEGER I,INFO,J,L,NROWA
217 LOGICAL UPPER
218* ..
219* .. Parameters ..
220 COMPLEX ONE
221 parameter(one= (1.0e+0,0.0e+0))
222 COMPLEX ZERO
223 parameter(zero= (0.0e+0,0.0e+0))
224* ..
225*
226* Test the input parameters.
227*
228 IF (lsame(trans,'N')) THEN
229 nrowa = n
230 ELSE
231 nrowa = k
232 END IF
233 upper = lsame(uplo,'U')
234*
235 info = 0
236 IF ((.NOT.upper) .AND. (.NOT.lsame(uplo,'L'))) THEN
237 info = 1
238 ELSE IF ((.NOT.lsame(trans,'N')) .AND.
239 + (.NOT.lsame(trans,'T'))) THEN
240 info = 2
241 ELSE IF (n.LT.0) THEN
242 info = 3
243 ELSE IF (k.LT.0) THEN
244 info = 4
245 ELSE IF (lda.LT.max(1,nrowa)) THEN
246 info = 7
247 ELSE IF (ldb.LT.max(1,nrowa)) THEN
248 info = 9
249 ELSE IF (ldc.LT.max(1,n)) THEN
250 info = 12
251 END IF
252 IF (info.NE.0) THEN
253 CALL xerbla('CSYR2K',info)
254 RETURN
255 END IF
256*
257* Quick return if possible.
258*
259 IF ((n.EQ.0) .OR. (((alpha.EQ.zero).OR.
260 + (k.EQ.0)).AND. (beta.EQ.one))) RETURN
261*
262* And when alpha.eq.zero.
263*
264 IF (alpha.EQ.zero) THEN
265 IF (upper) THEN
266 IF (beta.EQ.zero) THEN
267 DO 20 j = 1,n
268 DO 10 i = 1,j
269 c(i,j) = zero
270 10 CONTINUE
271 20 CONTINUE
272 ELSE
273 DO 40 j = 1,n
274 DO 30 i = 1,j
275 c(i,j) = beta*c(i,j)
276 30 CONTINUE
277 40 CONTINUE
278 END IF
279 ELSE
280 IF (beta.EQ.zero) THEN
281 DO 60 j = 1,n
282 DO 50 i = j,n
283 c(i,j) = zero
284 50 CONTINUE
285 60 CONTINUE
286 ELSE
287 DO 80 j = 1,n
288 DO 70 i = j,n
289 c(i,j) = beta*c(i,j)
290 70 CONTINUE
291 80 CONTINUE
292 END IF
293 END IF
294 RETURN
295 END IF
296*
297* Start the operations.
298*
299 IF (lsame(trans,'N')) THEN
300*
301* Form C := alpha*A*B**T + alpha*B*A**T + C.
302*
303 IF (upper) THEN
304 DO 130 j = 1,n
305 IF (beta.EQ.zero) THEN
306 DO 90 i = 1,j
307 c(i,j) = zero
308 90 CONTINUE
309 ELSE IF (beta.NE.one) THEN
310 DO 100 i = 1,j
311 c(i,j) = beta*c(i,j)
312 100 CONTINUE
313 END IF
314 DO 120 l = 1,k
315 IF ((a(j,l).NE.zero) .OR. (b(j,l).NE.zero)) THEN
316 temp1 = alpha*b(j,l)
317 temp2 = alpha*a(j,l)
318 DO 110 i = 1,j
319 c(i,j) = c(i,j) + a(i,l)*temp1 +
320 + b(i,l)*temp2
321 110 CONTINUE
322 END IF
323 120 CONTINUE
324 130 CONTINUE
325 ELSE
326 DO 180 j = 1,n
327 IF (beta.EQ.zero) THEN
328 DO 140 i = j,n
329 c(i,j) = zero
330 140 CONTINUE
331 ELSE IF (beta.NE.one) THEN
332 DO 150 i = j,n
333 c(i,j) = beta*c(i,j)
334 150 CONTINUE
335 END IF
336 DO 170 l = 1,k
337 IF ((a(j,l).NE.zero) .OR. (b(j,l).NE.zero)) THEN
338 temp1 = alpha*b(j,l)
339 temp2 = alpha*a(j,l)
340 DO 160 i = j,n
341 c(i,j) = c(i,j) + a(i,l)*temp1 +
342 + b(i,l)*temp2
343 160 CONTINUE
344 END IF
345 170 CONTINUE
346 180 CONTINUE
347 END IF
348 ELSE
349*
350* Form C := alpha*A**T*B + alpha*B**T*A + C.
351*
352 IF (upper) THEN
353 DO 210 j = 1,n
354 DO 200 i = 1,j
355 temp1 = zero
356 temp2 = zero
357 DO 190 l = 1,k
358 temp1 = temp1 + a(l,i)*b(l,j)
359 temp2 = temp2 + b(l,i)*a(l,j)
360 190 CONTINUE
361 IF (beta.EQ.zero) THEN
362 c(i,j) = alpha*temp1 + alpha*temp2
363 ELSE
364 c(i,j) = beta*c(i,j) + alpha*temp1 +
365 + alpha*temp2
366 END IF
367 200 CONTINUE
368 210 CONTINUE
369 ELSE
370 DO 240 j = 1,n
371 DO 230 i = j,n
372 temp1 = zero
373 temp2 = zero
374 DO 220 l = 1,k
375 temp1 = temp1 + a(l,i)*b(l,j)
376 temp2 = temp2 + b(l,i)*a(l,j)
377 220 CONTINUE
378 IF (beta.EQ.zero) THEN
379 c(i,j) = alpha*temp1 + alpha*temp2
380 ELSE
381 c(i,j) = beta*c(i,j) + alpha*temp1 +
382 + alpha*temp2
383 END IF
384 230 CONTINUE
385 240 CONTINUE
386 END IF
387 END IF
388*
389 RETURN
390*
391* End of CSYR2K
392*
393 END
subroutine xerbla(srname, info)
Definition cblat2.f:3285
subroutine csyr2k(uplo, trans, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
CSYR2K
Definition csyr2k.f:188