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