OpenRadioss 2025.1.11
OpenRadioss project
Loading...
Searching...
No Matches
complex16

Functions

subroutine zgemm (transa, transb, m, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
 ZGEMM
subroutine zhemm (side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
 ZHEMM
subroutine zher2k (uplo, trans, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
 ZHER2K
subroutine zherk (uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
 ZHERK
subroutine zsymm (side, uplo, m, n, alpha, a, lda, b, ldb, beta, c, ldc)
 ZSYMM
subroutine zsyr2k (uplo, trans, n, k, alpha, a, lda, b, ldb, beta, c, ldc)
 ZSYR2K
subroutine zsyrk (uplo, trans, n, k, alpha, a, lda, beta, c, ldc)
 ZSYRK
subroutine ztrmm (side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
 ZTRMM
subroutine ztrsm (side, uplo, transa, diag, m, n, alpha, a, lda, b, ldb)
 ZTRSM

Detailed Description

This is the group of complex16 LEVEL 3 BLAS routines.

Function Documentation

◆ zgemm()

subroutine zgemm ( character transa,
character transb,
integer m,
integer n,
integer k,
complex*16 alpha,
complex*16, dimension(lda,*) a,
integer lda,
complex*16, dimension(ldb,*) b,
integer ldb,
complex*16 beta,
complex*16, dimension(ldc,*) c,
integer ldc )

ZGEMM

Purpose:
!>
!> ZGEMM  performs one of the matrix-matrix operations
!>
!>    C := alpha*op( A )*op( B ) + beta*C,
!>
!> where  op( X ) is one of
!>
!>    op( X ) = X   or   op( X ) = X**T   or   op( X ) = X**H,
!>
!> alpha and beta are scalars, and A, B and C are matrices, with op( A )
!> an m by k matrix,  op( B )  a  k by n matrix and  C an m by n matrix.
!> 
Parameters
[in]TRANSA
!>          TRANSA is CHARACTER*1
!>           On entry, TRANSA specifies the form of op( A ) to be used in
!>           the matrix multiplication as follows:
!>
!>              TRANSA = 'N' or 'n',  op( A ) = A.
!>
!>              TRANSA = 'T' or 't',  op( A ) = A**T.
!>
!>              TRANSA = 'C' or 'c',  op( A ) = A**H.
!> 
[in]TRANSB
!>          TRANSB is CHARACTER*1
!>           On entry, TRANSB specifies the form of op( B ) to be used in
!>           the matrix multiplication as follows:
!>
!>              TRANSB = 'N' or 'n',  op( B ) = B.
!>
!>              TRANSB = 'T' or 't',  op( B ) = B**T.
!>
!>              TRANSB = 'C' or 'c',  op( B ) = B**H.
!> 
[in]M
!>          M is INTEGER
!>           On entry,  M  specifies  the number  of rows  of the  matrix
!>           op( A )  and of the  matrix  C.  M  must  be at least  zero.
!> 
[in]N
!>          N is INTEGER
!>           On entry,  N  specifies the number  of columns of the matrix
!>           op( B ) and the number of columns of the matrix C. N must be
!>           at least zero.
!> 
[in]K
!>          K is INTEGER
!>           On entry,  K  specifies  the number of columns of the matrix
!>           op( A ) and the number of rows of the matrix op( B ). K must
!>           be at least  zero.
!> 
[in]ALPHA
!>          ALPHA is COMPLEX*16
!>           On entry, ALPHA specifies the scalar alpha.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension ( LDA, ka ), where ka is
!>           k  when  TRANSA = 'N' or 'n',  and is  m  otherwise.
!>           Before entry with  TRANSA = 'N' or 'n',  the leading  m by k
!>           part of the array  A  must contain the matrix  A,  otherwise
!>           the leading  k by m  part of the array  A  must contain  the
!>           matrix A.
!> 
[in]LDA
!>          LDA is INTEGER
!>           On entry, LDA specifies the first dimension of A as declared
!>           in the calling (sub) program. When  TRANSA = 'N' or 'n' then
!>           LDA must be at least  max( 1, m ), otherwise  LDA must be at
!>           least  max( 1, k ).
!> 
[in]B
!>          B is COMPLEX*16 array, dimension ( LDB, kb ), where kb is
!>           n  when  TRANSB = 'N' or 'n',  and is  k  otherwise.
!>           Before entry with  TRANSB = 'N' or 'n',  the leading  k by n
!>           part of the array  B  must contain the matrix  B,  otherwise
!>           the leading  n by k  part of the array  B  must contain  the
!>           matrix B.
!> 
[in]LDB
!>          LDB is INTEGER
!>           On entry, LDB specifies the first dimension of B as declared
!>           in the calling (sub) program. When  TRANSB = 'N' or 'n' then
!>           LDB must be at least  max( 1, k ), otherwise  LDB must be at
!>           least  max( 1, n ).
!> 
[in]BETA
!>          BETA is COMPLEX*16
!>           On entry,  BETA  specifies the scalar  beta.  When  BETA  is
!>           supplied as zero then C need not be set on input.
!> 
[in,out]C
!>          C is COMPLEX*16 array, dimension ( LDC, N )
!>           Before entry, the leading  m by n  part of the array  C must
!>           contain the matrix  C,  except when  beta  is zero, in which
!>           case C need not be set on entry.
!>           On exit, the array  C  is overwritten by the  m by n  matrix
!>           ( alpha*op( A )*op( B ) + beta*C ).
!> 
[in]LDC
!>          LDC is INTEGER
!>           On entry, LDC specifies the first dimension of C as declared
!>           in  the  calling  (sub)  program.   LDC  must  be  at  least
!>           max( 1, m ).
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  Level 3 Blas routine.
!>
!>  -- Written on 8-February-1989.
!>     Jack Dongarra, Argonne National Laboratory.
!>     Iain Duff, AERE Harwell.
!>     Jeremy Du Croz, Numerical Algorithms Group Ltd.
!>     Sven Hammarling, Numerical Algorithms Group Ltd.
!> 

Definition at line 186 of file zgemm.f.

187*
188* -- Reference BLAS level3 routine --
189* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
190* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
191*
192* .. Scalar Arguments ..
193 COMPLEX*16 ALPHA,BETA
194 INTEGER K,LDA,LDB,LDC,M,N
195 CHARACTER TRANSA,TRANSB
196* ..
197* .. Array Arguments ..
198 COMPLEX*16 A(LDA,*),B(LDB,*),C(LDC,*)
199* ..
200*
201* =====================================================================
202*
203* .. External Functions ..
204 LOGICAL LSAME
205 EXTERNAL lsame
206* ..
207* .. External Subroutines ..
208 EXTERNAL xerbla
209* ..
210* .. Intrinsic Functions ..
211 INTRINSIC dconjg,max
212* ..
213* .. Local Scalars ..
214 COMPLEX*16 TEMP
215 INTEGER I,INFO,J,L,NROWA,NROWB
216 LOGICAL CONJA,CONJB,NOTA,NOTB
217* ..
218* .. Parameters ..
219 COMPLEX*16 ONE
220 parameter(one= (1.0d+0,0.0d+0))
221 COMPLEX*16 ZERO
222 parameter(zero= (0.0d+0,0.0d+0))
223* ..
224*
225* Set NOTA and NOTB as true if A and B respectively are not
226* conjugated or transposed, set CONJA and CONJB as true if A and
227* B respectively are to be transposed but not conjugated and set
228* NROWA and NROWB as the number of rows of A and B respectively.
229*
230 nota = lsame(transa,'N')
231 notb = lsame(transb,'N')
232 conja = lsame(transa,'C')
233 conjb = lsame(transb,'C')
234 IF (nota) THEN
235 nrowa = m
236 ELSE
237 nrowa = k
238 END IF
239 IF (notb) THEN
240 nrowb = k
241 ELSE
242 nrowb = n
243 END IF
244*
245* Test the input parameters.
246*
247 info = 0
248 IF ((.NOT.nota) .AND. (.NOT.conja) .AND.
249 + (.NOT.lsame(transa,'T'))) THEN
250 info = 1
251 ELSE IF ((.NOT.notb) .AND. (.NOT.conjb) .AND.
252 + (.NOT.lsame(transb,'T'))) THEN
253 info = 2
254 ELSE IF (m.LT.0) THEN
255 info = 3
256 ELSE IF (n.LT.0) THEN
257 info = 4
258 ELSE IF (k.LT.0) THEN
259 info = 5
260 ELSE IF (lda.LT.max(1,nrowa)) THEN
261 info = 8
262 ELSE IF (ldb.LT.max(1,nrowb)) THEN
263 info = 10
264 ELSE IF (ldc.LT.max(1,m)) THEN
265 info = 13
266 END IF
267 IF (info.NE.0) THEN
268 CALL xerbla('ZGEMM ',info)
269 RETURN
270 END IF
271*
272* Quick return if possible.
273*
274 IF ((m.EQ.0) .OR. (n.EQ.0) .OR.
275 + (((alpha.EQ.zero).OR. (k.EQ.0)).AND. (beta.EQ.one))) RETURN
276*
277* And when alpha.eq.zero.
278*
279 IF (alpha.EQ.zero) THEN
280 IF (beta.EQ.zero) THEN
281 DO 20 j = 1,n
282 DO 10 i = 1,m
283 c(i,j) = zero
284 10 CONTINUE
285 20 CONTINUE
286 ELSE
287 DO 40 j = 1,n
288 DO 30 i = 1,m
289 c(i,j) = beta*c(i,j)
290 30 CONTINUE
291 40 CONTINUE
292 END IF
293 RETURN
294 END IF
295*
296* Start the operations.
297*
298 IF (notb) THEN
299 IF (nota) THEN
300*
301* Form C := alpha*A*B + beta*C.
302*
303 DO 90 j = 1,n
304 IF (beta.EQ.zero) THEN
305 DO 50 i = 1,m
306 c(i,j) = zero
307 50 CONTINUE
308 ELSE IF (beta.NE.one) THEN
309 DO 60 i = 1,m
310 c(i,j) = beta*c(i,j)
311 60 CONTINUE
312 END IF
313 DO 80 l = 1,k
314 temp = alpha*b(l,j)
315 DO 70 i = 1,m
316 c(i,j) = c(i,j) + temp*a(i,l)
317 70 CONTINUE
318 80 CONTINUE
319 90 CONTINUE
320 ELSE IF (conja) THEN
321*
322* Form C := alpha*A**H*B + beta*C.
323*
324 DO 120 j = 1,n
325 DO 110 i = 1,m
326 temp = zero
327 DO 100 l = 1,k
328 temp = temp + dconjg(a(l,i))*b(l,j)
329 100 CONTINUE
330 IF (beta.EQ.zero) THEN
331 c(i,j) = alpha*temp
332 ELSE
333 c(i,j) = alpha*temp + beta*c(i,j)
334 END IF
335 110 CONTINUE
336 120 CONTINUE
337 ELSE
338*
339* Form C := alpha*A**T*B + beta*C
340*
341 DO 150 j = 1,n
342 DO 140 i = 1,m
343 temp = zero
344 DO 130 l = 1,k
345 temp = temp + a(l,i)*b(l,j)
346 130 CONTINUE
347 IF (beta.EQ.zero) THEN
348 c(i,j) = alpha*temp
349 ELSE
350 c(i,j) = alpha*temp + beta*c(i,j)
351 END IF
352 140 CONTINUE
353 150 CONTINUE
354 END IF
355 ELSE IF (nota) THEN
356 IF (conjb) THEN
357*
358* Form C := alpha*A*B**H + beta*C.
359*
360 DO 200 j = 1,n
361 IF (beta.EQ.zero) THEN
362 DO 160 i = 1,m
363 c(i,j) = zero
364 160 CONTINUE
365 ELSE IF (beta.NE.one) THEN
366 DO 170 i = 1,m
367 c(i,j) = beta*c(i,j)
368 170 CONTINUE
369 END IF
370 DO 190 l = 1,k
371 temp = alpha*dconjg(b(j,l))
372 DO 180 i = 1,m
373 c(i,j) = c(i,j) + temp*a(i,l)
374 180 CONTINUE
375 190 CONTINUE
376 200 CONTINUE
377 ELSE
378*
379* Form C := alpha*A*B**T + beta*C
380*
381 DO 250 j = 1,n
382 IF (beta.EQ.zero) THEN
383 DO 210 i = 1,m
384 c(i,j) = zero
385 210 CONTINUE
386 ELSE IF (beta.NE.one) THEN
387 DO 220 i = 1,m
388 c(i,j) = beta*c(i,j)
389 220 CONTINUE
390 END IF
391 DO 240 l = 1,k
392 temp = alpha*b(j,l)
393 DO 230 i = 1,m
394 c(i,j) = c(i,j) + temp*a(i,l)
395 230 CONTINUE
396 240 CONTINUE
397 250 CONTINUE
398 END IF
399 ELSE IF (conja) THEN
400 IF (conjb) THEN
401*
402* Form C := alpha*A**H*B**H + beta*C.
403*
404 DO 280 j = 1,n
405 DO 270 i = 1,m
406 temp = zero
407 DO 260 l = 1,k
408 temp = temp + dconjg(a(l,i))*dconjg(b(j,l))
409 260 CONTINUE
410 IF (beta.EQ.zero) THEN
411 c(i,j) = alpha*temp
412 ELSE
413 c(i,j) = alpha*temp + beta*c(i,j)
414 END IF
415 270 CONTINUE
416 280 CONTINUE
417 ELSE
418*
419* Form C := alpha*A**H*B**T + beta*C
420*
421 DO 310 j = 1,n
422 DO 300 i = 1,m
423 temp = zero
424 DO 290 l = 1,k
425 temp = temp + dconjg(a(l,i))*b(j,l)
426 290 CONTINUE
427 IF (beta.EQ.zero) THEN
428 c(i,j) = alpha*temp
429 ELSE
430 c(i,j) = alpha*temp + beta*c(i,j)
431 END IF
432 300 CONTINUE
433 310 CONTINUE
434 END IF
435 ELSE
436 IF (conjb) THEN
437*
438* Form C := alpha*A**T*B**H + beta*C
439*
440 DO 340 j = 1,n
441 DO 330 i = 1,m
442 temp = zero
443 DO 320 l = 1,k
444 temp = temp + a(l,i)*dconjg(b(j,l))
445 320 CONTINUE
446 IF (beta.EQ.zero) THEN
447 c(i,j) = alpha*temp
448 ELSE
449 c(i,j) = alpha*temp + beta*c(i,j)
450 END IF
451 330 CONTINUE
452 340 CONTINUE
453 ELSE
454*
455* Form C := alpha*A**T*B**T + beta*C
456*
457 DO 370 j = 1,n
458 DO 360 i = 1,m
459 temp = zero
460 DO 350 l = 1,k
461 temp = temp + a(l,i)*b(j,l)
462 350 CONTINUE
463 IF (beta.EQ.zero) THEN
464 c(i,j) = alpha*temp
465 ELSE
466 c(i,j) = alpha*temp + beta*c(i,j)
467 END IF
468 360 CONTINUE
469 370 CONTINUE
470 END IF
471 END IF
472*
473 RETURN
474*
475* End of ZGEMM
476*
#define alpha
Definition eval.h:35
subroutine xerbla(srname, info)
XERBLA
Definition xerbla.f:60
logical function lsame(ca, cb)
LSAME
Definition lsame.f:53
#define max(a, b)
Definition macros.h:21

◆ zhemm()

subroutine zhemm ( character side,
character uplo,
integer m,
integer n,
complex*16 alpha,
complex*16, dimension(lda,*) a,
integer lda,
complex*16, dimension(ldb,*) b,
integer ldb,
complex*16 beta,
complex*16, dimension(ldc,*) c,
integer ldc )

ZHEMM

Purpose:
!>
!> ZHEMM  performs one of the matrix-matrix operations
!>
!>    C := alpha*A*B + beta*C,
!>
!> or
!>
!>    C := alpha*B*A + beta*C,
!>
!> where alpha and beta are scalars, A is an hermitian matrix and  B and
!> C are m by n matrices.
!> 
Parameters
[in]SIDE
!>          SIDE is CHARACTER*1
!>           On entry,  SIDE  specifies whether  the  hermitian matrix  A
!>           appears on the  left or right  in the  operation as follows:
!>
!>              SIDE = 'L' or 'l'   C := alpha*A*B + beta*C,
!>
!>              SIDE = 'R' or 'r'   C := alpha*B*A + beta*C,
!> 
[in]UPLO
!>          UPLO is CHARACTER*1
!>           On  entry,   UPLO  specifies  whether  the  upper  or  lower
!>           triangular  part  of  the  hermitian  matrix   A  is  to  be
!>           referenced as follows:
!>
!>              UPLO = 'U' or 'u'   Only the upper triangular part of the
!>                                  hermitian matrix is to be referenced.
!>
!>              UPLO = 'L' or 'l'   Only the lower triangular part of the
!>                                  hermitian matrix is to be referenced.
!> 
[in]M
!>          M is INTEGER
!>           On entry,  M  specifies the number of rows of the matrix  C.
!>           M  must be at least zero.
!> 
[in]N
!>          N is INTEGER
!>           On entry, N specifies the number of columns of the matrix C.
!>           N  must be at least zero.
!> 
[in]ALPHA
!>          ALPHA is COMPLEX*16
!>           On entry, ALPHA specifies the scalar alpha.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension ( LDA, ka ), where ka is
!>           m  when  SIDE = 'L' or 'l'  and is n  otherwise.
!>           Before entry  with  SIDE = 'L' or 'l',  the  m by m  part of
!>           the array  A  must contain the  hermitian matrix,  such that
!>           when  UPLO = 'U' or 'u', the leading m by m upper triangular
!>           part of the array  A  must contain the upper triangular part
!>           of the  hermitian matrix and the  strictly  lower triangular
!>           part of  A  is not referenced,  and when  UPLO = 'L' or 'l',
!>           the leading  m by m  lower triangular part  of the  array  A
!>           must  contain  the  lower triangular part  of the  hermitian
!>           matrix and the  strictly upper triangular part of  A  is not
!>           referenced.
!>           Before entry  with  SIDE = 'R' or 'r',  the  n by n  part of
!>           the array  A  must contain the  hermitian matrix,  such that
!>           when  UPLO = 'U' or 'u', the leading n by n upper triangular
!>           part of the array  A  must contain the upper triangular part
!>           of the  hermitian matrix and the  strictly  lower triangular
!>           part of  A  is not referenced,  and when  UPLO = 'L' or 'l',
!>           the leading  n by n  lower triangular part  of the  array  A
!>           must  contain  the  lower triangular part  of the  hermitian
!>           matrix and the  strictly upper triangular part of  A  is not
!>           referenced.
!>           Note that the imaginary parts  of the diagonal elements need
!>           not be set, they are assumed to be zero.
!> 
[in]LDA
!>          LDA is INTEGER
!>           On entry, LDA specifies the first dimension of A as declared
!>           in the  calling (sub) program. When  SIDE = 'L' or 'l'  then
!>           LDA must be at least  max( 1, m ), otherwise  LDA must be at
!>           least max( 1, n ).
!> 
[in]B
!>          B is COMPLEX*16 array, dimension ( LDB, N )
!>           Before entry, the leading  m by n part of the array  B  must
!>           contain the matrix B.
!> 
[in]LDB
!>          LDB is INTEGER
!>           On entry, LDB specifies the first dimension of B as declared
!>           in  the  calling  (sub)  program.   LDB  must  be  at  least
!>           max( 1, m ).
!> 
[in]BETA
!>          BETA is COMPLEX*16
!>           On entry,  BETA  specifies the scalar  beta.  When  BETA  is
!>           supplied as zero then C need not be set on input.
!> 
[in,out]C
!>          C is COMPLEX*16 array, dimension ( LDC, N )
!>           Before entry, the leading  m by n  part of the array  C must
!>           contain the matrix  C,  except when  beta  is zero, in which
!>           case C need not be set on entry.
!>           On exit, the array  C  is overwritten by the  m by n updated
!>           matrix.
!> 
[in]LDC
!>          LDC is INTEGER
!>           On entry, LDC specifies the first dimension of C as declared
!>           in  the  calling  (sub)  program.   LDC  must  be  at  least
!>           max( 1, m ).
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  Level 3 Blas routine.
!>
!>  -- Written on 8-February-1989.
!>     Jack Dongarra, Argonne National Laboratory.
!>     Iain Duff, AERE Harwell.
!>     Jeremy Du Croz, Numerical Algorithms Group Ltd.
!>     Sven Hammarling, Numerical Algorithms Group Ltd.
!> 

Definition at line 190 of file zhemm.f.

191*
192* -- Reference BLAS level3 routine --
193* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
194* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
195*
196* .. Scalar Arguments ..
197 COMPLEX*16 ALPHA,BETA
198 INTEGER LDA,LDB,LDC,M,N
199 CHARACTER SIDE,UPLO
200* ..
201* .. Array Arguments ..
202 COMPLEX*16 A(LDA,*),B(LDB,*),C(LDC,*)
203* ..
204*
205* =====================================================================
206*
207* .. External Functions ..
208 LOGICAL LSAME
209 EXTERNAL lsame
210* ..
211* .. External Subroutines ..
212 EXTERNAL xerbla
213* ..
214* .. Intrinsic Functions ..
215 INTRINSIC dble,dconjg,max
216* ..
217* .. Local Scalars ..
218 COMPLEX*16 TEMP1,TEMP2
219 INTEGER I,INFO,J,K,NROWA
220 LOGICAL UPPER
221* ..
222* .. Parameters ..
223 COMPLEX*16 ONE
224 parameter(one= (1.0d+0,0.0d+0))
225 COMPLEX*16 ZERO
226 parameter(zero= (0.0d+0,0.0d+0))
227* ..
228*
229* Set NROWA as the number of rows of A.
230*
231 IF (lsame(side,'L')) THEN
232 nrowa = m
233 ELSE
234 nrowa = n
235 END IF
236 upper = lsame(uplo,'U')
237*
238* Test the input parameters.
239*
240 info = 0
241 IF ((.NOT.lsame(side,'L')) .AND. (.NOT.lsame(side,'R'))) THEN
242 info = 1
243 ELSE IF ((.NOT.upper) .AND. (.NOT.lsame(uplo,'L'))) THEN
244 info = 2
245 ELSE IF (m.LT.0) THEN
246 info = 3
247 ELSE IF (n.LT.0) THEN
248 info = 4
249 ELSE IF (lda.LT.max(1,nrowa)) THEN
250 info = 7
251 ELSE IF (ldb.LT.max(1,m)) THEN
252 info = 9
253 ELSE IF (ldc.LT.max(1,m)) THEN
254 info = 12
255 END IF
256 IF (info.NE.0) THEN
257 CALL xerbla('ZHEMM ',info)
258 RETURN
259 END IF
260*
261* Quick return if possible.
262*
263 IF ((m.EQ.0) .OR. (n.EQ.0) .OR.
264 + ((alpha.EQ.zero).AND. (beta.EQ.one))) RETURN
265*
266* And when alpha.eq.zero.
267*
268 IF (alpha.EQ.zero) THEN
269 IF (beta.EQ.zero) THEN
270 DO 20 j = 1,n
271 DO 10 i = 1,m
272 c(i,j) = zero
273 10 CONTINUE
274 20 CONTINUE
275 ELSE
276 DO 40 j = 1,n
277 DO 30 i = 1,m
278 c(i,j) = beta*c(i,j)
279 30 CONTINUE
280 40 CONTINUE
281 END IF
282 RETURN
283 END IF
284*
285* Start the operations.
286*
287 IF (lsame(side,'L')) THEN
288*
289* Form C := alpha*A*B + beta*C.
290*
291 IF (upper) THEN
292 DO 70 j = 1,n
293 DO 60 i = 1,m
294 temp1 = alpha*b(i,j)
295 temp2 = zero
296 DO 50 k = 1,i - 1
297 c(k,j) = c(k,j) + temp1*a(k,i)
298 temp2 = temp2 + b(k,j)*dconjg(a(k,i))
299 50 CONTINUE
300 IF (beta.EQ.zero) THEN
301 c(i,j) = temp1*dble(a(i,i)) + alpha*temp2
302 ELSE
303 c(i,j) = beta*c(i,j) + temp1*dble(a(i,i)) +
304 + alpha*temp2
305 END IF
306 60 CONTINUE
307 70 CONTINUE
308 ELSE
309 DO 100 j = 1,n
310 DO 90 i = m,1,-1
311 temp1 = alpha*b(i,j)
312 temp2 = zero
313 DO 80 k = i + 1,m
314 c(k,j) = c(k,j) + temp1*a(k,i)
315 temp2 = temp2 + b(k,j)*dconjg(a(k,i))
316 80 CONTINUE
317 IF (beta.EQ.zero) THEN
318 c(i,j) = temp1*dble(a(i,i)) + alpha*temp2
319 ELSE
320 c(i,j) = beta*c(i,j) + temp1*dble(a(i,i)) +
321 + alpha*temp2
322 END IF
323 90 CONTINUE
324 100 CONTINUE
325 END IF
326 ELSE
327*
328* Form C := alpha*B*A + beta*C.
329*
330 DO 170 j = 1,n
331 temp1 = alpha*dble(a(j,j))
332 IF (beta.EQ.zero) THEN
333 DO 110 i = 1,m
334 c(i,j) = temp1*b(i,j)
335 110 CONTINUE
336 ELSE
337 DO 120 i = 1,m
338 c(i,j) = beta*c(i,j) + temp1*b(i,j)
339 120 CONTINUE
340 END IF
341 DO 140 k = 1,j - 1
342 IF (upper) THEN
343 temp1 = alpha*a(k,j)
344 ELSE
345 temp1 = alpha*dconjg(a(j,k))
346 END IF
347 DO 130 i = 1,m
348 c(i,j) = c(i,j) + temp1*b(i,k)
349 130 CONTINUE
350 140 CONTINUE
351 DO 160 k = j + 1,n
352 IF (upper) THEN
353 temp1 = alpha*dconjg(a(j,k))
354 ELSE
355 temp1 = alpha*a(k,j)
356 END IF
357 DO 150 i = 1,m
358 c(i,j) = c(i,j) + temp1*b(i,k)
359 150 CONTINUE
360 160 CONTINUE
361 170 CONTINUE
362 END IF
363*
364 RETURN
365*
366* End of ZHEMM
367*

◆ zher2k()

subroutine zher2k ( character uplo,
character trans,
integer n,
integer k,
complex*16 alpha,
complex*16, dimension(lda,*) a,
integer lda,
complex*16, dimension(ldb,*) b,
integer ldb,
double precision beta,
complex*16, dimension(ldc,*) c,
integer ldc )

ZHER2K

Purpose:
!>
!> ZHER2K  performs one of the hermitian rank 2k operations
!>
!>    C := alpha*A*B**H + conjg( alpha )*B*A**H + beta*C,
!>
!> or
!>
!>    C := alpha*A**H*B + conjg( alpha )*B**H*A + beta*C,
!>
!> where  alpha and beta  are scalars with  beta  real,  C is an  n by n
!> hermitian matrix and  A and B  are  n by k matrices in the first case
!> and  k by n  matrices in the second case.
!> 
Parameters
[in]UPLO
!>          UPLO is CHARACTER*1
!>           On  entry,   UPLO  specifies  whether  the  upper  or  lower
!>           triangular  part  of the  array  C  is to be  referenced  as
!>           follows:
!>
!>              UPLO = 'U' or 'u'   Only the  upper triangular part of  C
!>                                  is to be referenced.
!>
!>              UPLO = 'L' or 'l'   Only the  lower triangular part of  C
!>                                  is to be referenced.
!> 
[in]TRANS
!>          TRANS is CHARACTER*1
!>           On entry,  TRANS  specifies the operation to be performed as
!>           follows:
!>
!>              TRANS = 'N' or 'n'    C := alpha*A*B**H          +
!>                                         conjg( alpha )*B*A**H +
!>                                         beta*C.
!>
!>              TRANS = 'C' or 'c'    C := alpha*A**H*B          +
!>                                         conjg( alpha )*B**H*A +
!>                                         beta*C.
!> 
[in]N
!>          N is INTEGER
!>           On entry,  N specifies the order of the matrix C.  N must be
!>           at least zero.
!> 
[in]K
!>          K is INTEGER
!>           On entry with  TRANS = 'N' or 'n',  K  specifies  the number
!>           of  columns  of the  matrices  A and B,  and on  entry  with
!>           TRANS = 'C' or 'c',  K  specifies  the number of rows of the
!>           matrices  A and B.  K must be at least zero.
!> 
[in]ALPHA
!>          ALPHA is COMPLEX*16 .
!>           On entry, ALPHA specifies the scalar alpha.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension ( LDA, ka ), where ka is
!>           k  when  TRANS = 'N' or 'n',  and is  n  otherwise.
!>           Before entry with  TRANS = 'N' or 'n',  the  leading  n by k
!>           part of the array  A  must contain the matrix  A,  otherwise
!>           the leading  k by n  part of the array  A  must contain  the
!>           matrix A.
!> 
[in]LDA
!>          LDA is INTEGER
!>           On entry, LDA specifies the first dimension of A as declared
!>           in  the  calling  (sub)  program.   When  TRANS = 'N' or 'n'
!>           then  LDA must be at least  max( 1, n ), otherwise  LDA must
!>           be at least  max( 1, k ).
!> 
[in]B
!>          B is COMPLEX*16 array, dimension ( LDB, kb ), where kb is
!>           k  when  TRANS = 'N' or 'n',  and is  n  otherwise.
!>           Before entry with  TRANS = 'N' or 'n',  the  leading  n by k
!>           part of the array  B  must contain the matrix  B,  otherwise
!>           the leading  k by n  part of the array  B  must contain  the
!>           matrix B.
!> 
[in]LDB
!>          LDB is INTEGER
!>           On entry, LDB specifies the first dimension of B as declared
!>           in  the  calling  (sub)  program.   When  TRANS = 'N' or 'n'
!>           then  LDB must be at least  max( 1, n ), otherwise  LDB must
!>           be at least  max( 1, k ).
!>           Unchanged on exit.
!> 
[in]BETA
!>          BETA is DOUBLE PRECISION .
!>           On entry, BETA specifies the scalar beta.
!> 
[in,out]C
!>          C is COMPLEX*16 array, dimension ( LDC, N )
!>           Before entry  with  UPLO = 'U' or 'u',  the leading  n by n
!>           upper triangular part of the array C must contain the upper
!>           triangular part  of the  hermitian matrix  and the strictly
!>           lower triangular part of C is not referenced.  On exit, the
!>           upper triangular part of the array  C is overwritten by the
!>           upper triangular part of the updated matrix.
!>           Before entry  with  UPLO = 'L' or 'l',  the leading  n by n
!>           lower triangular part of the array C must contain the lower
!>           triangular part  of the  hermitian matrix  and the strictly
!>           upper triangular part of C is not referenced.  On exit, the
!>           lower triangular part of the array  C is overwritten by the
!>           lower triangular part of the updated matrix.
!>           Note that the imaginary parts of the diagonal elements need
!>           not be set,  they are assumed to be zero,  and on exit they
!>           are set to zero.
!> 
[in]LDC
!>          LDC is INTEGER
!>           On entry, LDC specifies the first dimension of C as declared
!>           in  the  calling  (sub)  program.   LDC  must  be  at  least
!>           max( 1, n ).
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  Level 3 Blas routine.
!>
!>  -- Written on 8-February-1989.
!>     Jack Dongarra, Argonne National Laboratory.
!>     Iain Duff, AERE Harwell.
!>     Jeremy Du Croz, Numerical Algorithms Group Ltd.
!>     Sven Hammarling, Numerical Algorithms Group Ltd.
!>
!>  -- Modified 8-Nov-93 to set C(J,J) to DBLE( C(J,J) ) when BETA = 1.
!>     Ed Anderson, Cray Research Inc.
!> 

Definition at line 197 of file zher2k.f.

198*
199* -- Reference BLAS level3 routine --
200* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
201* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
202*
203* .. Scalar Arguments ..
204 COMPLEX*16 ALPHA
205 DOUBLE PRECISION BETA
206 INTEGER K,LDA,LDB,LDC,N
207 CHARACTER TRANS,UPLO
208* ..
209* .. Array Arguments ..
210 COMPLEX*16 A(LDA,*),B(LDB,*),C(LDC,*)
211* ..
212*
213* =====================================================================
214*
215* .. External Functions ..
216 LOGICAL LSAME
217 EXTERNAL lsame
218* ..
219* .. External Subroutines ..
220 EXTERNAL xerbla
221* ..
222* .. Intrinsic Functions ..
223 INTRINSIC dble,dconjg,max
224* ..
225* .. Local Scalars ..
226 COMPLEX*16 TEMP1,TEMP2
227 INTEGER I,INFO,J,L,NROWA
228 LOGICAL UPPER
229* ..
230* .. Parameters ..
231 DOUBLE PRECISION ONE
232 parameter(one=1.0d+0)
233 COMPLEX*16 ZERO
234 parameter(zero= (0.0d+0,0.0d+0))
235* ..
236*
237* Test the input parameters.
238*
239 IF (lsame(trans,'N')) THEN
240 nrowa = n
241 ELSE
242 nrowa = k
243 END IF
244 upper = lsame(uplo,'U')
245*
246 info = 0
247 IF ((.NOT.upper) .AND. (.NOT.lsame(uplo,'L'))) THEN
248 info = 1
249 ELSE IF ((.NOT.lsame(trans,'N')) .AND.
250 + (.NOT.lsame(trans,'C'))) THEN
251 info = 2
252 ELSE IF (n.LT.0) THEN
253 info = 3
254 ELSE IF (k.LT.0) THEN
255 info = 4
256 ELSE IF (lda.LT.max(1,nrowa)) THEN
257 info = 7
258 ELSE IF (ldb.LT.max(1,nrowa)) THEN
259 info = 9
260 ELSE IF (ldc.LT.max(1,n)) THEN
261 info = 12
262 END IF
263 IF (info.NE.0) THEN
264 CALL xerbla('ZHER2K',info)
265 RETURN
266 END IF
267*
268* Quick return if possible.
269*
270 IF ((n.EQ.0) .OR. (((alpha.EQ.zero).OR.
271 + (k.EQ.0)).AND. (beta.EQ.one))) RETURN
272*
273* And when alpha.eq.zero.
274*
275 IF (alpha.EQ.zero) THEN
276 IF (upper) THEN
277 IF (beta.EQ.dble(zero)) THEN
278 DO 20 j = 1,n
279 DO 10 i = 1,j
280 c(i,j) = zero
281 10 CONTINUE
282 20 CONTINUE
283 ELSE
284 DO 40 j = 1,n
285 DO 30 i = 1,j - 1
286 c(i,j) = beta*c(i,j)
287 30 CONTINUE
288 c(j,j) = beta*dble(c(j,j))
289 40 CONTINUE
290 END IF
291 ELSE
292 IF (beta.EQ.dble(zero)) THEN
293 DO 60 j = 1,n
294 DO 50 i = j,n
295 c(i,j) = zero
296 50 CONTINUE
297 60 CONTINUE
298 ELSE
299 DO 80 j = 1,n
300 c(j,j) = beta*dble(c(j,j))
301 DO 70 i = j + 1,n
302 c(i,j) = beta*c(i,j)
303 70 CONTINUE
304 80 CONTINUE
305 END IF
306 END IF
307 RETURN
308 END IF
309*
310* Start the operations.
311*
312 IF (lsame(trans,'N')) THEN
313*
314* Form C := alpha*A*B**H + conjg( alpha )*B*A**H +
315* C.
316*
317 IF (upper) THEN
318 DO 130 j = 1,n
319 IF (beta.EQ.dble(zero)) THEN
320 DO 90 i = 1,j
321 c(i,j) = zero
322 90 CONTINUE
323 ELSE IF (beta.NE.one) THEN
324 DO 100 i = 1,j - 1
325 c(i,j) = beta*c(i,j)
326 100 CONTINUE
327 c(j,j) = beta*dble(c(j,j))
328 ELSE
329 c(j,j) = dble(c(j,j))
330 END IF
331 DO 120 l = 1,k
332 IF ((a(j,l).NE.zero) .OR. (b(j,l).NE.zero)) THEN
333 temp1 = alpha*dconjg(b(j,l))
334 temp2 = dconjg(alpha*a(j,l))
335 DO 110 i = 1,j - 1
336 c(i,j) = c(i,j) + a(i,l)*temp1 +
337 + b(i,l)*temp2
338 110 CONTINUE
339 c(j,j) = dble(c(j,j)) +
340 + dble(a(j,l)*temp1+b(j,l)*temp2)
341 END IF
342 120 CONTINUE
343 130 CONTINUE
344 ELSE
345 DO 180 j = 1,n
346 IF (beta.EQ.dble(zero)) THEN
347 DO 140 i = j,n
348 c(i,j) = zero
349 140 CONTINUE
350 ELSE IF (beta.NE.one) THEN
351 DO 150 i = j + 1,n
352 c(i,j) = beta*c(i,j)
353 150 CONTINUE
354 c(j,j) = beta*dble(c(j,j))
355 ELSE
356 c(j,j) = dble(c(j,j))
357 END IF
358 DO 170 l = 1,k
359 IF ((a(j,l).NE.zero) .OR. (b(j,l).NE.zero)) THEN
360 temp1 = alpha*dconjg(b(j,l))
361 temp2 = dconjg(alpha*a(j,l))
362 DO 160 i = j + 1,n
363 c(i,j) = c(i,j) + a(i,l)*temp1 +
364 + b(i,l)*temp2
365 160 CONTINUE
366 c(j,j) = dble(c(j,j)) +
367 + dble(a(j,l)*temp1+b(j,l)*temp2)
368 END IF
369 170 CONTINUE
370 180 CONTINUE
371 END IF
372 ELSE
373*
374* Form C := alpha*A**H*B + conjg( alpha )*B**H*A +
375* C.
376*
377 IF (upper) THEN
378 DO 210 j = 1,n
379 DO 200 i = 1,j
380 temp1 = zero
381 temp2 = zero
382 DO 190 l = 1,k
383 temp1 = temp1 + dconjg(a(l,i))*b(l,j)
384 temp2 = temp2 + dconjg(b(l,i))*a(l,j)
385 190 CONTINUE
386 IF (i.EQ.j) THEN
387 IF (beta.EQ.dble(zero)) THEN
388 c(j,j) = dble(alpha*temp1+
389 + dconjg(alpha)*temp2)
390 ELSE
391 c(j,j) = beta*dble(c(j,j)) +
392 + dble(alpha*temp1+
393 + dconjg(alpha)*temp2)
394 END IF
395 ELSE
396 IF (beta.EQ.dble(zero)) THEN
397 c(i,j) = alpha*temp1 + dconjg(alpha)*temp2
398 ELSE
399 c(i,j) = beta*c(i,j) + alpha*temp1 +
400 + dconjg(alpha)*temp2
401 END IF
402 END IF
403 200 CONTINUE
404 210 CONTINUE
405 ELSE
406 DO 240 j = 1,n
407 DO 230 i = j,n
408 temp1 = zero
409 temp2 = zero
410 DO 220 l = 1,k
411 temp1 = temp1 + dconjg(a(l,i))*b(l,j)
412 temp2 = temp2 + dconjg(b(l,i))*a(l,j)
413 220 CONTINUE
414 IF (i.EQ.j) THEN
415 IF (beta.EQ.dble(zero)) THEN
416 c(j,j) = dble(alpha*temp1+
417 + dconjg(alpha)*temp2)
418 ELSE
419 c(j,j) = beta*dble(c(j,j)) +
420 + dble(alpha*temp1+
421 + dconjg(alpha)*temp2)
422 END IF
423 ELSE
424 IF (beta.EQ.dble(zero)) THEN
425 c(i,j) = alpha*temp1 + dconjg(alpha)*temp2
426 ELSE
427 c(i,j) = beta*c(i,j) + alpha*temp1 +
428 + dconjg(alpha)*temp2
429 END IF
430 END IF
431 230 CONTINUE
432 240 CONTINUE
433 END IF
434 END IF
435*
436 RETURN
437*
438* End of ZHER2K
439*

◆ zherk()

subroutine zherk ( character uplo,
character trans,
integer n,
integer k,
double precision alpha,
complex*16, dimension(lda,*) a,
integer lda,
double precision beta,
complex*16, dimension(ldc,*) c,
integer ldc )

ZHERK

Purpose:
!>
!> ZHERK  performs one of the hermitian rank k operations
!>
!>    C := alpha*A*A**H + beta*C,
!>
!> or
!>
!>    C := alpha*A**H*A + beta*C,
!>
!> where  alpha and beta  are  real scalars,  C is an  n by n  hermitian
!> matrix and  A  is an  n by k  matrix in the  first case and a  k by n
!> matrix in the second case.
!> 
Parameters
[in]UPLO
!>          UPLO is CHARACTER*1
!>           On  entry,   UPLO  specifies  whether  the  upper  or  lower
!>           triangular  part  of the  array  C  is to be  referenced  as
!>           follows:
!>
!>              UPLO = 'U' or 'u'   Only the  upper triangular part of  C
!>                                  is to be referenced.
!>
!>              UPLO = 'L' or 'l'   Only the  lower triangular part of  C
!>                                  is to be referenced.
!> 
[in]TRANS
!>          TRANS is CHARACTER*1
!>           On entry,  TRANS  specifies the operation to be performed as
!>           follows:
!>
!>              TRANS = 'N' or 'n'   C := alpha*A*A**H + beta*C.
!>
!>              TRANS = 'C' or 'c'   C := alpha*A**H*A + beta*C.
!> 
[in]N
!>          N is INTEGER
!>           On entry,  N specifies the order of the matrix C.  N must be
!>           at least zero.
!> 
[in]K
!>          K is INTEGER
!>           On entry with  TRANS = 'N' or 'n',  K  specifies  the number
!>           of  columns   of  the   matrix   A,   and  on   entry   with
!>           TRANS = 'C' or 'c',  K  specifies  the number of rows of the
!>           matrix A.  K must be at least zero.
!> 
[in]ALPHA
!>          ALPHA is DOUBLE PRECISION .
!>           On entry, ALPHA specifies the scalar alpha.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension ( LDA, ka ), where ka is
!>           k  when  TRANS = 'N' or 'n',  and is  n  otherwise.
!>           Before entry with  TRANS = 'N' or 'n',  the  leading  n by k
!>           part of the array  A  must contain the matrix  A,  otherwise
!>           the leading  k by n  part of the array  A  must contain  the
!>           matrix A.
!> 
[in]LDA
!>          LDA is INTEGER
!>           On entry, LDA specifies the first dimension of A as declared
!>           in  the  calling  (sub)  program.   When  TRANS = 'N' or 'n'
!>           then  LDA must be at least  max( 1, n ), otherwise  LDA must
!>           be at least  max( 1, k ).
!> 
[in]BETA
!>          BETA is DOUBLE PRECISION.
!>           On entry, BETA specifies the scalar beta.
!> 
[in,out]C
!>          C is COMPLEX*16 array, dimension ( LDC, N )
!>           Before entry  with  UPLO = 'U' or 'u',  the leading  n by n
!>           upper triangular part of the array C must contain the upper
!>           triangular part  of the  hermitian matrix  and the strictly
!>           lower triangular part of C is not referenced.  On exit, the
!>           upper triangular part of the array  C is overwritten by the
!>           upper triangular part of the updated matrix.
!>           Before entry  with  UPLO = 'L' or 'l',  the leading  n by n
!>           lower triangular part of the array C must contain the lower
!>           triangular part  of the  hermitian matrix  and the strictly
!>           upper triangular part of C is not referenced.  On exit, the
!>           lower triangular part of the array  C is overwritten by the
!>           lower triangular part of the updated matrix.
!>           Note that the imaginary parts of the diagonal elements need
!>           not be set,  they are assumed to be zero,  and on exit they
!>           are set to zero.
!> 
[in]LDC
!>          LDC is INTEGER
!>           On entry, LDC specifies the first dimension of C as declared
!>           in  the  calling  (sub)  program.   LDC  must  be  at  least
!>           max( 1, n ).
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  Level 3 Blas routine.
!>
!>  -- Written on 8-February-1989.
!>     Jack Dongarra, Argonne National Laboratory.
!>     Iain Duff, AERE Harwell.
!>     Jeremy Du Croz, Numerical Algorithms Group Ltd.
!>     Sven Hammarling, Numerical Algorithms Group Ltd.
!>
!>  -- Modified 8-Nov-93 to set C(J,J) to DBLE( C(J,J) ) when BETA = 1.
!>     Ed Anderson, Cray Research Inc.
!> 

Definition at line 172 of file zherk.f.

173*
174* -- Reference BLAS level3 routine --
175* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
176* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
177*
178* .. Scalar Arguments ..
179 DOUBLE PRECISION ALPHA,BETA
180 INTEGER K,LDA,LDC,N
181 CHARACTER TRANS,UPLO
182* ..
183* .. Array Arguments ..
184 COMPLEX*16 A(LDA,*),C(LDC,*)
185* ..
186*
187* =====================================================================
188*
189* .. External Functions ..
190 LOGICAL LSAME
191 EXTERNAL lsame
192* ..
193* .. External Subroutines ..
194 EXTERNAL xerbla
195* ..
196* .. Intrinsic Functions ..
197 INTRINSIC dble,dcmplx,dconjg,max
198* ..
199* .. Local Scalars ..
200 COMPLEX*16 TEMP
201 DOUBLE PRECISION RTEMP
202 INTEGER I,INFO,J,L,NROWA
203 LOGICAL UPPER
204* ..
205* .. Parameters ..
206 DOUBLE PRECISION ONE,ZERO
207 parameter(one=1.0d+0,zero=0.0d+0)
208* ..
209*
210* Test the input parameters.
211*
212 IF (lsame(trans,'N')) THEN
213 nrowa = n
214 ELSE
215 nrowa = k
216 END IF
217 upper = lsame(uplo,'U')
218*
219 info = 0
220 IF ((.NOT.upper) .AND. (.NOT.lsame(uplo,'L'))) THEN
221 info = 1
222 ELSE IF ((.NOT.lsame(trans,'N')) .AND.
223 + (.NOT.lsame(trans,'C'))) THEN
224 info = 2
225 ELSE IF (n.LT.0) THEN
226 info = 3
227 ELSE IF (k.LT.0) THEN
228 info = 4
229 ELSE IF (lda.LT.max(1,nrowa)) THEN
230 info = 7
231 ELSE IF (ldc.LT.max(1,n)) THEN
232 info = 10
233 END IF
234 IF (info.NE.0) THEN
235 CALL xerbla('ZHERK ',info)
236 RETURN
237 END IF
238*
239* Quick return if possible.
240*
241 IF ((n.EQ.0) .OR. (((alpha.EQ.zero).OR.
242 + (k.EQ.0)).AND. (beta.EQ.one))) RETURN
243*
244* And when alpha.eq.zero.
245*
246 IF (alpha.EQ.zero) THEN
247 IF (upper) THEN
248 IF (beta.EQ.zero) THEN
249 DO 20 j = 1,n
250 DO 10 i = 1,j
251 c(i,j) = zero
252 10 CONTINUE
253 20 CONTINUE
254 ELSE
255 DO 40 j = 1,n
256 DO 30 i = 1,j - 1
257 c(i,j) = beta*c(i,j)
258 30 CONTINUE
259 c(j,j) = beta*dble(c(j,j))
260 40 CONTINUE
261 END IF
262 ELSE
263 IF (beta.EQ.zero) THEN
264 DO 60 j = 1,n
265 DO 50 i = j,n
266 c(i,j) = zero
267 50 CONTINUE
268 60 CONTINUE
269 ELSE
270 DO 80 j = 1,n
271 c(j,j) = beta*dble(c(j,j))
272 DO 70 i = j + 1,n
273 c(i,j) = beta*c(i,j)
274 70 CONTINUE
275 80 CONTINUE
276 END IF
277 END IF
278 RETURN
279 END IF
280*
281* Start the operations.
282*
283 IF (lsame(trans,'N')) THEN
284*
285* Form C := alpha*A*A**H + beta*C.
286*
287 IF (upper) THEN
288 DO 130 j = 1,n
289 IF (beta.EQ.zero) THEN
290 DO 90 i = 1,j
291 c(i,j) = zero
292 90 CONTINUE
293 ELSE IF (beta.NE.one) THEN
294 DO 100 i = 1,j - 1
295 c(i,j) = beta*c(i,j)
296 100 CONTINUE
297 c(j,j) = beta*dble(c(j,j))
298 ELSE
299 c(j,j) = dble(c(j,j))
300 END IF
301 DO 120 l = 1,k
302 IF (a(j,l).NE.dcmplx(zero)) THEN
303 temp = alpha*dconjg(a(j,l))
304 DO 110 i = 1,j - 1
305 c(i,j) = c(i,j) + temp*a(i,l)
306 110 CONTINUE
307 c(j,j) = dble(c(j,j)) + dble(temp*a(i,l))
308 END IF
309 120 CONTINUE
310 130 CONTINUE
311 ELSE
312 DO 180 j = 1,n
313 IF (beta.EQ.zero) THEN
314 DO 140 i = j,n
315 c(i,j) = zero
316 140 CONTINUE
317 ELSE IF (beta.NE.one) THEN
318 c(j,j) = beta*dble(c(j,j))
319 DO 150 i = j + 1,n
320 c(i,j) = beta*c(i,j)
321 150 CONTINUE
322 ELSE
323 c(j,j) = dble(c(j,j))
324 END IF
325 DO 170 l = 1,k
326 IF (a(j,l).NE.dcmplx(zero)) THEN
327 temp = alpha*dconjg(a(j,l))
328 c(j,j) = dble(c(j,j)) + dble(temp*a(j,l))
329 DO 160 i = j + 1,n
330 c(i,j) = c(i,j) + temp*a(i,l)
331 160 CONTINUE
332 END IF
333 170 CONTINUE
334 180 CONTINUE
335 END IF
336 ELSE
337*
338* Form C := alpha*A**H*A + beta*C.
339*
340 IF (upper) THEN
341 DO 220 j = 1,n
342 DO 200 i = 1,j - 1
343 temp = zero
344 DO 190 l = 1,k
345 temp = temp + dconjg(a(l,i))*a(l,j)
346 190 CONTINUE
347 IF (beta.EQ.zero) THEN
348 c(i,j) = alpha*temp
349 ELSE
350 c(i,j) = alpha*temp + beta*c(i,j)
351 END IF
352 200 CONTINUE
353 rtemp = zero
354 DO 210 l = 1,k
355 rtemp = rtemp + dconjg(a(l,j))*a(l,j)
356 210 CONTINUE
357 IF (beta.EQ.zero) THEN
358 c(j,j) = alpha*rtemp
359 ELSE
360 c(j,j) = alpha*rtemp + beta*dble(c(j,j))
361 END IF
362 220 CONTINUE
363 ELSE
364 DO 260 j = 1,n
365 rtemp = zero
366 DO 230 l = 1,k
367 rtemp = rtemp + dconjg(a(l,j))*a(l,j)
368 230 CONTINUE
369 IF (beta.EQ.zero) THEN
370 c(j,j) = alpha*rtemp
371 ELSE
372 c(j,j) = alpha*rtemp + beta*dble(c(j,j))
373 END IF
374 DO 250 i = j + 1,n
375 temp = zero
376 DO 240 l = 1,k
377 temp = temp + dconjg(a(l,i))*a(l,j)
378 240 CONTINUE
379 IF (beta.EQ.zero) THEN
380 c(i,j) = alpha*temp
381 ELSE
382 c(i,j) = alpha*temp + beta*c(i,j)
383 END IF
384 250 CONTINUE
385 260 CONTINUE
386 END IF
387 END IF
388*
389 RETURN
390*
391* End of ZHERK
392*

◆ zsymm()

subroutine zsymm ( character side,
character uplo,
integer m,
integer n,
complex*16 alpha,
complex*16, dimension(lda,*) a,
integer lda,
complex*16, dimension(ldb,*) b,
integer ldb,
complex*16 beta,
complex*16, dimension(ldc,*) c,
integer ldc )

ZSYMM

Purpose:
!>
!> ZSYMM  performs one of the matrix-matrix operations
!>
!>    C := alpha*A*B + beta*C,
!>
!> or
!>
!>    C := alpha*B*A + beta*C,
!>
!> where  alpha and beta are scalars, A is a symmetric matrix and  B and
!> C are m by n matrices.
!> 
Parameters
[in]SIDE
!>          SIDE is CHARACTER*1
!>           On entry,  SIDE  specifies whether  the  symmetric matrix  A
!>           appears on the  left or right  in the  operation as follows:
!>
!>              SIDE = 'L' or 'l'   C := alpha*A*B + beta*C,
!>
!>              SIDE = 'R' or 'r'   C := alpha*B*A + beta*C,
!> 
[in]UPLO
!>          UPLO is CHARACTER*1
!>           On  entry,   UPLO  specifies  whether  the  upper  or  lower
!>           triangular  part  of  the  symmetric  matrix   A  is  to  be
!>           referenced as follows:
!>
!>              UPLO = 'U' or 'u'   Only the upper triangular part of the
!>                                  symmetric matrix is to be referenced.
!>
!>              UPLO = 'L' or 'l'   Only the lower triangular part of the
!>                                  symmetric matrix is to be referenced.
!> 
[in]M
!>          M is INTEGER
!>           On entry,  M  specifies the number of rows of the matrix  C.
!>           M  must be at least zero.
!> 
[in]N
!>          N is INTEGER
!>           On entry, N specifies the number of columns of the matrix C.
!>           N  must be at least zero.
!> 
[in]ALPHA
!>          ALPHA is COMPLEX*16
!>           On entry, ALPHA specifies the scalar alpha.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension ( LDA, ka ), where ka is
!>           m  when  SIDE = 'L' or 'l'  and is n  otherwise.
!>           Before entry  with  SIDE = 'L' or 'l',  the  m by m  part of
!>           the array  A  must contain the  symmetric matrix,  such that
!>           when  UPLO = 'U' or 'u', the leading m by m upper triangular
!>           part of the array  A  must contain the upper triangular part
!>           of the  symmetric matrix and the  strictly  lower triangular
!>           part of  A  is not referenced,  and when  UPLO = 'L' or 'l',
!>           the leading  m by m  lower triangular part  of the  array  A
!>           must  contain  the  lower triangular part  of the  symmetric
!>           matrix and the  strictly upper triangular part of  A  is not
!>           referenced.
!>           Before entry  with  SIDE = 'R' or 'r',  the  n by n  part of
!>           the array  A  must contain the  symmetric matrix,  such that
!>           when  UPLO = 'U' or 'u', the leading n by n upper triangular
!>           part of the array  A  must contain the upper triangular part
!>           of the  symmetric matrix and the  strictly  lower triangular
!>           part of  A  is not referenced,  and when  UPLO = 'L' or 'l',
!>           the leading  n by n  lower triangular part  of the  array  A
!>           must  contain  the  lower triangular part  of the  symmetric
!>           matrix and the  strictly upper triangular part of  A  is not
!>           referenced.
!> 
[in]LDA
!>          LDA is INTEGER
!>           On entry, LDA specifies the first dimension of A as declared
!>           in the  calling (sub) program. When  SIDE = 'L' or 'l'  then
!>           LDA must be at least  max( 1, m ), otherwise  LDA must be at
!>           least max( 1, n ).
!> 
[in]B
!>          B is COMPLEX*16 array, dimension ( LDB, N )
!>           Before entry, the leading  m by n part of the array  B  must
!>           contain the matrix B.
!> 
[in]LDB
!>          LDB is INTEGER
!>           On entry, LDB specifies the first dimension of B as declared
!>           in  the  calling  (sub)  program.   LDB  must  be  at  least
!>           max( 1, m ).
!> 
[in]BETA
!>          BETA is COMPLEX*16
!>           On entry,  BETA  specifies the scalar  beta.  When  BETA  is
!>           supplied as zero then C need not be set on input.
!> 
[in,out]C
!>          C is COMPLEX*16 array, dimension ( LDC, N )
!>           Before entry, the leading  m by n  part of the array  C must
!>           contain the matrix  C,  except when  beta  is zero, in which
!>           case C need not be set on entry.
!>           On exit, the array  C  is overwritten by the  m by n updated
!>           matrix.
!> 
[in]LDC
!>          LDC is INTEGER
!>           On entry, LDC specifies the first dimension of C as declared
!>           in  the  calling  (sub)  program.   LDC  must  be  at  least
!>           max( 1, m ).
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  Level 3 Blas routine.
!>
!>  -- Written on 8-February-1989.
!>     Jack Dongarra, Argonne National Laboratory.
!>     Iain Duff, AERE Harwell.
!>     Jeremy Du Croz, Numerical Algorithms Group Ltd.
!>     Sven Hammarling, Numerical Algorithms Group Ltd.
!> 

Definition at line 188 of file zsymm.f.

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 COMPLEX*16 ALPHA,BETA
196 INTEGER LDA,LDB,LDC,M,N
197 CHARACTER SIDE,UPLO
198* ..
199* .. Array Arguments ..
200 COMPLEX*16 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 COMPLEX*16 TEMP1,TEMP2
217 INTEGER I,INFO,J,K,NROWA
218 LOGICAL UPPER
219* ..
220* .. Parameters ..
221 COMPLEX*16 ONE
222 parameter(one= (1.0d+0,0.0d+0))
223 COMPLEX*16 ZERO
224 parameter(zero= (0.0d+0,0.0d+0))
225* ..
226*
227* Set NROWA as the number of rows of A.
228*
229 IF (lsame(side,'L')) THEN
230 nrowa = m
231 ELSE
232 nrowa = n
233 END IF
234 upper = lsame(uplo,'U')
235*
236* Test the input parameters.
237*
238 info = 0
239 IF ((.NOT.lsame(side,'L')) .AND. (.NOT.lsame(side,'R'))) THEN
240 info = 1
241 ELSE IF ((.NOT.upper) .AND. (.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('ZSYMM ',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 ZSYMM
365*

◆ zsyr2k()

subroutine zsyr2k ( character uplo,
character trans,
integer n,
integer k,
complex*16 alpha,
complex*16, dimension(lda,*) a,
integer lda,
complex*16, dimension(ldb,*) b,
integer ldb,
complex*16 beta,
complex*16, dimension(ldc,*) c,
integer ldc )

ZSYR2K

Purpose:
!>
!> ZSYR2K  performs one of the symmetric rank 2k operations
!>
!>    C := alpha*A*B**T + alpha*B*A**T + beta*C,
!>
!> or
!>
!>    C := alpha*A**T*B + alpha*B**T*A + beta*C,
!>
!> where  alpha and beta  are scalars,  C is an  n by n symmetric matrix
!> and  A and B  are  n by k  matrices  in the  first  case  and  k by n
!> matrices in the second case.
!> 
Parameters
[in]UPLO
!>          UPLO is CHARACTER*1
!>           On  entry,   UPLO  specifies  whether  the  upper  or  lower
!>           triangular  part  of the  array  C  is to be  referenced  as
!>           follows:
!>
!>              UPLO = 'U' or 'u'   Only the  upper triangular part of  C
!>                                  is to be referenced.
!>
!>              UPLO = 'L' or 'l'   Only the  lower triangular part of  C
!>                                  is to be referenced.
!> 
[in]TRANS
!>          TRANS is CHARACTER*1
!>           On entry,  TRANS  specifies the operation to be performed as
!>           follows:
!>
!>              TRANS = 'N' or 'n'    C := alpha*A*B**T + alpha*B*A**T +
!>                                         beta*C.
!>
!>              TRANS = 'T' or 't'    C := alpha*A**T*B + alpha*B**T*A +
!>                                         beta*C.
!> 
[in]N
!>          N is INTEGER
!>           On entry,  N specifies the order of the matrix C.  N must be
!>           at least zero.
!> 
[in]K
!>          K is INTEGER
!>           On entry with  TRANS = 'N' or 'n',  K  specifies  the number
!>           of  columns  of the  matrices  A and B,  and on  entry  with
!>           TRANS = 'T' or 't',  K  specifies  the number of rows of the
!>           matrices  A and B.  K must be at least zero.
!> 
[in]ALPHA
!>          ALPHA is COMPLEX*16
!>           On entry, ALPHA specifies the scalar alpha.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension ( LDA, ka ), where ka is
!>           k  when  TRANS = 'N' or 'n',  and is  n  otherwise.
!>           Before entry with  TRANS = 'N' or 'n',  the  leading  n by k
!>           part of the array  A  must contain the matrix  A,  otherwise
!>           the leading  k by n  part of the array  A  must contain  the
!>           matrix A.
!> 
[in]LDA
!>          LDA is INTEGER
!>           On entry, LDA specifies the first dimension of A as declared
!>           in  the  calling  (sub)  program.   When  TRANS = 'N' or 'n'
!>           then  LDA must be at least  max( 1, n ), otherwise  LDA must
!>           be at least  max( 1, k ).
!> 
[in]B
!>          B is COMPLEX*16 array, dimension ( LDB, kb ), where kb is
!>           k  when  TRANS = 'N' or 'n',  and is  n  otherwise.
!>           Before entry with  TRANS = 'N' or 'n',  the  leading  n by k
!>           part of the array  B  must contain the matrix  B,  otherwise
!>           the leading  k by n  part of the array  B  must contain  the
!>           matrix B.
!> 
[in]LDB
!>          LDB is INTEGER
!>           On entry, LDB specifies the first dimension of B as declared
!>           in  the  calling  (sub)  program.   When  TRANS = 'N' or 'n'
!>           then  LDB must be at least  max( 1, n ), otherwise  LDB must
!>           be at least  max( 1, k ).
!> 
[in]BETA
!>          BETA is COMPLEX*16
!>           On entry, BETA specifies the scalar beta.
!> 
[in,out]C
!>          C is COMPLEX*16 array, dimension ( LDC, N )
!>           Before entry  with  UPLO = 'U' or 'u',  the leading  n by n
!>           upper triangular part of the array C must contain the upper
!>           triangular part  of the  symmetric matrix  and the strictly
!>           lower triangular part of C is not referenced.  On exit, the
!>           upper triangular part of the array  C is overwritten by the
!>           upper triangular part of the updated matrix.
!>           Before entry  with  UPLO = 'L' or 'l',  the leading  n by n
!>           lower triangular part of the array C must contain the lower
!>           triangular part  of the  symmetric matrix  and the strictly
!>           upper triangular part of C is not referenced.  On exit, the
!>           lower triangular part of the array  C is overwritten by the
!>           lower triangular part of the updated matrix.
!> 
[in]LDC
!>          LDC is INTEGER
!>           On entry, LDC specifies the first dimension of C as declared
!>           in  the  calling  (sub)  program.   LDC  must  be  at  least
!>           max( 1, n ).
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  Level 3 Blas routine.
!>
!>  -- Written on 8-February-1989.
!>     Jack Dongarra, Argonne National Laboratory.
!>     Iain Duff, AERE Harwell.
!>     Jeremy Du Croz, Numerical Algorithms Group Ltd.
!>     Sven Hammarling, Numerical Algorithms Group Ltd.
!> 

Definition at line 187 of file zsyr2k.f.

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*16 ALPHA,BETA
195 INTEGER K,LDA,LDB,LDC,N
196 CHARACTER TRANS,UPLO
197* ..
198* .. Array Arguments ..
199 COMPLEX*16 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*16 TEMP1,TEMP2
216 INTEGER I,INFO,J,L,NROWA
217 LOGICAL UPPER
218* ..
219* .. Parameters ..
220 COMPLEX*16 ONE
221 parameter(one= (1.0d+0,0.0d+0))
222 COMPLEX*16 ZERO
223 parameter(zero= (0.0d+0,0.0d+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('ZSYR2K',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 ZSYR2K
392*

◆ zsyrk()

subroutine zsyrk ( character uplo,
character trans,
integer n,
integer k,
complex*16 alpha,
complex*16, dimension(lda,*) a,
integer lda,
complex*16 beta,
complex*16, dimension(ldc,*) c,
integer ldc )

ZSYRK

Purpose:
!>
!> ZSYRK  performs one of the symmetric rank k operations
!>
!>    C := alpha*A*A**T + beta*C,
!>
!> or
!>
!>    C := alpha*A**T*A + beta*C,
!>
!> where  alpha and beta  are scalars,  C is an  n by n symmetric matrix
!> and  A  is an  n by k  matrix in the first case and a  k by n  matrix
!> in the second case.
!> 
Parameters
[in]UPLO
!>          UPLO is CHARACTER*1
!>           On  entry,   UPLO  specifies  whether  the  upper  or  lower
!>           triangular  part  of the  array  C  is to be  referenced  as
!>           follows:
!>
!>              UPLO = 'U' or 'u'   Only the  upper triangular part of  C
!>                                  is to be referenced.
!>
!>              UPLO = 'L' or 'l'   Only the  lower triangular part of  C
!>                                  is to be referenced.
!> 
[in]TRANS
!>          TRANS is CHARACTER*1
!>           On entry,  TRANS  specifies the operation to be performed as
!>           follows:
!>
!>              TRANS = 'N' or 'n'   C := alpha*A*A**T + beta*C.
!>
!>              TRANS = 'T' or 't'   C := alpha*A**T*A + beta*C.
!> 
[in]N
!>          N is INTEGER
!>           On entry,  N specifies the order of the matrix C.  N must be
!>           at least zero.
!> 
[in]K
!>          K is INTEGER
!>           On entry with  TRANS = 'N' or 'n',  K  specifies  the number
!>           of  columns   of  the   matrix   A,   and  on   entry   with
!>           TRANS = 'T' or 't',  K  specifies  the number of rows of the
!>           matrix A.  K must be at least zero.
!> 
[in]ALPHA
!>          ALPHA is COMPLEX*16
!>           On entry, ALPHA specifies the scalar alpha.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension ( LDA, ka ), where ka is
!>           k  when  TRANS = 'N' or 'n',  and is  n  otherwise.
!>           Before entry with  TRANS = 'N' or 'n',  the  leading  n by k
!>           part of the array  A  must contain the matrix  A,  otherwise
!>           the leading  k by n  part of the array  A  must contain  the
!>           matrix A.
!> 
[in]LDA
!>          LDA is INTEGER
!>           On entry, LDA specifies the first dimension of A as declared
!>           in  the  calling  (sub)  program.   When  TRANS = 'N' or 'n'
!>           then  LDA must be at least  max( 1, n ), otherwise  LDA must
!>           be at least  max( 1, k ).
!> 
[in]BETA
!>          BETA is COMPLEX*16
!>           On entry, BETA specifies the scalar beta.
!> 
[in,out]C
!>          C is COMPLEX*16 array, dimension ( LDC, N )
!>           Before entry  with  UPLO = 'U' or 'u',  the leading  n by n
!>           upper triangular part of the array C must contain the upper
!>           triangular part  of the  symmetric matrix  and the strictly
!>           lower triangular part of C is not referenced.  On exit, the
!>           upper triangular part of the array  C is overwritten by the
!>           upper triangular part of the updated matrix.
!>           Before entry  with  UPLO = 'L' or 'l',  the leading  n by n
!>           lower triangular part of the array C must contain the lower
!>           triangular part  of the  symmetric matrix  and the strictly
!>           upper triangular part of C is not referenced.  On exit, the
!>           lower triangular part of the array  C is overwritten by the
!>           lower triangular part of the updated matrix.
!> 
[in]LDC
!>          LDC is INTEGER
!>           On entry, LDC specifies the first dimension of C as declared
!>           in  the  calling  (sub)  program.   LDC  must  be  at  least
!>           max( 1, n ).
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  Level 3 Blas routine.
!>
!>  -- Written on 8-February-1989.
!>     Jack Dongarra, Argonne National Laboratory.
!>     Iain Duff, AERE Harwell.
!>     Jeremy Du Croz, Numerical Algorithms Group Ltd.
!>     Sven Hammarling, Numerical Algorithms Group Ltd.
!> 

Definition at line 166 of file zsyrk.f.

167*
168* -- Reference BLAS level3 routine --
169* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
170* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
171*
172* .. Scalar Arguments ..
173 COMPLEX*16 ALPHA,BETA
174 INTEGER K,LDA,LDC,N
175 CHARACTER TRANS,UPLO
176* ..
177* .. Array Arguments ..
178 COMPLEX*16 A(LDA,*),C(LDC,*)
179* ..
180*
181* =====================================================================
182*
183* .. External Functions ..
184 LOGICAL LSAME
185 EXTERNAL lsame
186* ..
187* .. External Subroutines ..
188 EXTERNAL xerbla
189* ..
190* .. Intrinsic Functions ..
191 INTRINSIC max
192* ..
193* .. Local Scalars ..
194 COMPLEX*16 TEMP
195 INTEGER I,INFO,J,L,NROWA
196 LOGICAL UPPER
197* ..
198* .. Parameters ..
199 COMPLEX*16 ONE
200 parameter(one= (1.0d+0,0.0d+0))
201 COMPLEX*16 ZERO
202 parameter(zero= (0.0d+0,0.0d+0))
203* ..
204*
205* Test the input parameters.
206*
207 IF (lsame(trans,'N')) THEN
208 nrowa = n
209 ELSE
210 nrowa = k
211 END IF
212 upper = lsame(uplo,'U')
213*
214 info = 0
215 IF ((.NOT.upper) .AND. (.NOT.lsame(uplo,'L'))) THEN
216 info = 1
217 ELSE IF ((.NOT.lsame(trans,'N')) .AND.
218 + (.NOT.lsame(trans,'T'))) THEN
219 info = 2
220 ELSE IF (n.LT.0) THEN
221 info = 3
222 ELSE IF (k.LT.0) THEN
223 info = 4
224 ELSE IF (lda.LT.max(1,nrowa)) THEN
225 info = 7
226 ELSE IF (ldc.LT.max(1,n)) THEN
227 info = 10
228 END IF
229 IF (info.NE.0) THEN
230 CALL xerbla('ZSYRK ',info)
231 RETURN
232 END IF
233*
234* Quick return if possible.
235*
236 IF ((n.EQ.0) .OR. (((alpha.EQ.zero).OR.
237 + (k.EQ.0)).AND. (beta.EQ.one))) RETURN
238*
239* And when alpha.eq.zero.
240*
241 IF (alpha.EQ.zero) THEN
242 IF (upper) THEN
243 IF (beta.EQ.zero) THEN
244 DO 20 j = 1,n
245 DO 10 i = 1,j
246 c(i,j) = zero
247 10 CONTINUE
248 20 CONTINUE
249 ELSE
250 DO 40 j = 1,n
251 DO 30 i = 1,j
252 c(i,j) = beta*c(i,j)
253 30 CONTINUE
254 40 CONTINUE
255 END IF
256 ELSE
257 IF (beta.EQ.zero) THEN
258 DO 60 j = 1,n
259 DO 50 i = j,n
260 c(i,j) = zero
261 50 CONTINUE
262 60 CONTINUE
263 ELSE
264 DO 80 j = 1,n
265 DO 70 i = j,n
266 c(i,j) = beta*c(i,j)
267 70 CONTINUE
268 80 CONTINUE
269 END IF
270 END IF
271 RETURN
272 END IF
273*
274* Start the operations.
275*
276 IF (lsame(trans,'N')) THEN
277*
278* Form C := alpha*A*A**T + beta*C.
279*
280 IF (upper) THEN
281 DO 130 j = 1,n
282 IF (beta.EQ.zero) THEN
283 DO 90 i = 1,j
284 c(i,j) = zero
285 90 CONTINUE
286 ELSE IF (beta.NE.one) THEN
287 DO 100 i = 1,j
288 c(i,j) = beta*c(i,j)
289 100 CONTINUE
290 END IF
291 DO 120 l = 1,k
292 IF (a(j,l).NE.zero) THEN
293 temp = alpha*a(j,l)
294 DO 110 i = 1,j
295 c(i,j) = c(i,j) + temp*a(i,l)
296 110 CONTINUE
297 END IF
298 120 CONTINUE
299 130 CONTINUE
300 ELSE
301 DO 180 j = 1,n
302 IF (beta.EQ.zero) THEN
303 DO 140 i = j,n
304 c(i,j) = zero
305 140 CONTINUE
306 ELSE IF (beta.NE.one) THEN
307 DO 150 i = j,n
308 c(i,j) = beta*c(i,j)
309 150 CONTINUE
310 END IF
311 DO 170 l = 1,k
312 IF (a(j,l).NE.zero) THEN
313 temp = alpha*a(j,l)
314 DO 160 i = j,n
315 c(i,j) = c(i,j) + temp*a(i,l)
316 160 CONTINUE
317 END IF
318 170 CONTINUE
319 180 CONTINUE
320 END IF
321 ELSE
322*
323* Form C := alpha*A**T*A + beta*C.
324*
325 IF (upper) THEN
326 DO 210 j = 1,n
327 DO 200 i = 1,j
328 temp = zero
329 DO 190 l = 1,k
330 temp = temp + a(l,i)*a(l,j)
331 190 CONTINUE
332 IF (beta.EQ.zero) THEN
333 c(i,j) = alpha*temp
334 ELSE
335 c(i,j) = alpha*temp + beta*c(i,j)
336 END IF
337 200 CONTINUE
338 210 CONTINUE
339 ELSE
340 DO 240 j = 1,n
341 DO 230 i = j,n
342 temp = zero
343 DO 220 l = 1,k
344 temp = temp + a(l,i)*a(l,j)
345 220 CONTINUE
346 IF (beta.EQ.zero) THEN
347 c(i,j) = alpha*temp
348 ELSE
349 c(i,j) = alpha*temp + beta*c(i,j)
350 END IF
351 230 CONTINUE
352 240 CONTINUE
353 END IF
354 END IF
355*
356 RETURN
357*
358* End of ZSYRK
359*

◆ ztrmm()

subroutine ztrmm ( character side,
character uplo,
character transa,
character diag,
integer m,
integer n,
complex*16 alpha,
complex*16, dimension(lda,*) a,
integer lda,
complex*16, dimension(ldb,*) b,
integer ldb )

ZTRMM

Purpose:
!>
!> ZTRMM  performs one of the matrix-matrix operations
!>
!>    B := alpha*op( A )*B,   or   B := alpha*B*op( A )
!>
!> where  alpha  is a scalar,  B  is an m by n matrix,  A  is a unit, or
!> non-unit,  upper or lower triangular matrix  and  op( A )  is one  of
!>
!>    op( A ) = A   or   op( A ) = A**T   or   op( A ) = A**H.
!> 
Parameters
[in]SIDE
!>          SIDE is CHARACTER*1
!>           On entry,  SIDE specifies whether  op( A ) multiplies B from
!>           the left or right as follows:
!>
!>              SIDE = 'L' or 'l'   B := alpha*op( A )*B.
!>
!>              SIDE = 'R' or 'r'   B := alpha*B*op( A ).
!> 
[in]UPLO
!>          UPLO is CHARACTER*1
!>           On entry, UPLO specifies whether the matrix A is an upper or
!>           lower triangular matrix as follows:
!>
!>              UPLO = 'U' or 'u'   A is an upper triangular matrix.
!>
!>              UPLO = 'L' or 'l'   A is a lower triangular matrix.
!> 
[in]TRANSA
!>          TRANSA is CHARACTER*1
!>           On entry, TRANSA specifies the form of op( A ) to be used in
!>           the matrix multiplication as follows:
!>
!>              TRANSA = 'N' or 'n'   op( A ) = A.
!>
!>              TRANSA = 'T' or 't'   op( A ) = A**T.
!>
!>              TRANSA = 'C' or 'c'   op( A ) = A**H.
!> 
[in]DIAG
!>          DIAG is CHARACTER*1
!>           On entry, DIAG specifies whether or not A is unit triangular
!>           as follows:
!>
!>              DIAG = 'U' or 'u'   A is assumed to be unit triangular.
!>
!>              DIAG = 'N' or 'n'   A is not assumed to be unit
!>                                  triangular.
!> 
[in]M
!>          M is INTEGER
!>           On entry, M specifies the number of rows of B. M must be at
!>           least zero.
!> 
[in]N
!>          N is INTEGER
!>           On entry, N specifies the number of columns of B.  N must be
!>           at least zero.
!> 
[in]ALPHA
!>          ALPHA is COMPLEX*16
!>           On entry,  ALPHA specifies the scalar  alpha. When  alpha is
!>           zero then  A is not referenced and  B need not be set before
!>           entry.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension ( LDA, k ), where k is m
!>           when  SIDE = 'L' or 'l'  and is  n  when  SIDE = 'R' or 'r'.
!>           Before entry  with  UPLO = 'U' or 'u',  the  leading  k by k
!>           upper triangular part of the array  A must contain the upper
!>           triangular matrix  and the strictly lower triangular part of
!>           A is not referenced.
!>           Before entry  with  UPLO = 'L' or 'l',  the  leading  k by k
!>           lower triangular part of the array  A must contain the lower
!>           triangular matrix  and the strictly upper triangular part of
!>           A is not referenced.
!>           Note that when  DIAG = 'U' or 'u',  the diagonal elements of
!>           A  are not referenced either,  but are assumed to be  unity.
!> 
[in]LDA
!>          LDA is INTEGER
!>           On entry, LDA specifies the first dimension of A as declared
!>           in the calling (sub) program.  When  SIDE = 'L' or 'l'  then
!>           LDA  must be at least  max( 1, m ),  when  SIDE = 'R' or 'r'
!>           then LDA must be at least max( 1, n ).
!> 
[in,out]B
!>          B is COMPLEX*16 array, dimension ( LDB, N ).
!>           Before entry,  the leading  m by n part of the array  B must
!>           contain the matrix  B,  and  on exit  is overwritten  by the
!>           transformed matrix.
!> 
[in]LDB
!>          LDB is INTEGER
!>           On entry, LDB specifies the first dimension of B as declared
!>           in  the  calling  (sub)  program.   LDB  must  be  at  least
!>           max( 1, m ).
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  Level 3 Blas routine.
!>
!>  -- Written on 8-February-1989.
!>     Jack Dongarra, Argonne National Laboratory.
!>     Iain Duff, AERE Harwell.
!>     Jeremy Du Croz, Numerical Algorithms Group Ltd.
!>     Sven Hammarling, Numerical Algorithms Group Ltd.
!> 

Definition at line 176 of file ztrmm.f.

177*
178* -- Reference BLAS level3 routine --
179* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
180* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
181*
182* .. Scalar Arguments ..
183 COMPLEX*16 ALPHA
184 INTEGER LDA,LDB,M,N
185 CHARACTER DIAG,SIDE,TRANSA,UPLO
186* ..
187* .. Array Arguments ..
188 COMPLEX*16 A(LDA,*),B(LDB,*)
189* ..
190*
191* =====================================================================
192*
193* .. External Functions ..
194 LOGICAL LSAME
195 EXTERNAL lsame
196* ..
197* .. External Subroutines ..
198 EXTERNAL xerbla
199* ..
200* .. Intrinsic Functions ..
201 INTRINSIC dconjg,max
202* ..
203* .. Local Scalars ..
204 COMPLEX*16 TEMP
205 INTEGER I,INFO,J,K,NROWA
206 LOGICAL LSIDE,NOCONJ,NOUNIT,UPPER
207* ..
208* .. Parameters ..
209 COMPLEX*16 ONE
210 parameter(one= (1.0d+0,0.0d+0))
211 COMPLEX*16 ZERO
212 parameter(zero= (0.0d+0,0.0d+0))
213* ..
214*
215* Test the input parameters.
216*
217 lside = lsame(side,'L')
218 IF (lside) THEN
219 nrowa = m
220 ELSE
221 nrowa = n
222 END IF
223 noconj = lsame(transa,'T')
224 nounit = lsame(diag,'N')
225 upper = lsame(uplo,'U')
226*
227 info = 0
228 IF ((.NOT.lside) .AND. (.NOT.lsame(side,'R'))) THEN
229 info = 1
230 ELSE IF ((.NOT.upper) .AND. (.NOT.lsame(uplo,'L'))) THEN
231 info = 2
232 ELSE IF ((.NOT.lsame(transa,'N')) .AND.
233 + (.NOT.lsame(transa,'T')) .AND.
234 + (.NOT.lsame(transa,'C'))) THEN
235 info = 3
236 ELSE IF ((.NOT.lsame(diag,'U')) .AND. (.NOT.lsame(diag,'N'))) THEN
237 info = 4
238 ELSE IF (m.LT.0) THEN
239 info = 5
240 ELSE IF (n.LT.0) THEN
241 info = 6
242 ELSE IF (lda.LT.max(1,nrowa)) THEN
243 info = 9
244 ELSE IF (ldb.LT.max(1,m)) THEN
245 info = 11
246 END IF
247 IF (info.NE.0) THEN
248 CALL xerbla('ZTRMM ',info)
249 RETURN
250 END IF
251*
252* Quick return if possible.
253*
254 IF (m.EQ.0 .OR. n.EQ.0) RETURN
255*
256* And when alpha.eq.zero.
257*
258 IF (alpha.EQ.zero) THEN
259 DO 20 j = 1,n
260 DO 10 i = 1,m
261 b(i,j) = zero
262 10 CONTINUE
263 20 CONTINUE
264 RETURN
265 END IF
266*
267* Start the operations.
268*
269 IF (lside) THEN
270 IF (lsame(transa,'N')) THEN
271*
272* Form B := alpha*A*B.
273*
274 IF (upper) THEN
275 DO 50 j = 1,n
276 DO 40 k = 1,m
277 IF (b(k,j).NE.zero) THEN
278 temp = alpha*b(k,j)
279 DO 30 i = 1,k - 1
280 b(i,j) = b(i,j) + temp*a(i,k)
281 30 CONTINUE
282 IF (nounit) temp = temp*a(k,k)
283 b(k,j) = temp
284 END IF
285 40 CONTINUE
286 50 CONTINUE
287 ELSE
288 DO 80 j = 1,n
289 DO 70 k = m,1,-1
290 IF (b(k,j).NE.zero) THEN
291 temp = alpha*b(k,j)
292 b(k,j) = temp
293 IF (nounit) b(k,j) = b(k,j)*a(k,k)
294 DO 60 i = k + 1,m
295 b(i,j) = b(i,j) + temp*a(i,k)
296 60 CONTINUE
297 END IF
298 70 CONTINUE
299 80 CONTINUE
300 END IF
301 ELSE
302*
303* Form B := alpha*A**T*B or B := alpha*A**H*B.
304*
305 IF (upper) THEN
306 DO 120 j = 1,n
307 DO 110 i = m,1,-1
308 temp = b(i,j)
309 IF (noconj) THEN
310 IF (nounit) temp = temp*a(i,i)
311 DO 90 k = 1,i - 1
312 temp = temp + a(k,i)*b(k,j)
313 90 CONTINUE
314 ELSE
315 IF (nounit) temp = temp*dconjg(a(i,i))
316 DO 100 k = 1,i - 1
317 temp = temp + dconjg(a(k,i))*b(k,j)
318 100 CONTINUE
319 END IF
320 b(i,j) = alpha*temp
321 110 CONTINUE
322 120 CONTINUE
323 ELSE
324 DO 160 j = 1,n
325 DO 150 i = 1,m
326 temp = b(i,j)
327 IF (noconj) THEN
328 IF (nounit) temp = temp*a(i,i)
329 DO 130 k = i + 1,m
330 temp = temp + a(k,i)*b(k,j)
331 130 CONTINUE
332 ELSE
333 IF (nounit) temp = temp*dconjg(a(i,i))
334 DO 140 k = i + 1,m
335 temp = temp + dconjg(a(k,i))*b(k,j)
336 140 CONTINUE
337 END IF
338 b(i,j) = alpha*temp
339 150 CONTINUE
340 160 CONTINUE
341 END IF
342 END IF
343 ELSE
344 IF (lsame(transa,'N')) THEN
345*
346* Form B := alpha*B*A.
347*
348 IF (upper) THEN
349 DO 200 j = n,1,-1
350 temp = alpha
351 IF (nounit) temp = temp*a(j,j)
352 DO 170 i = 1,m
353 b(i,j) = temp*b(i,j)
354 170 CONTINUE
355 DO 190 k = 1,j - 1
356 IF (a(k,j).NE.zero) THEN
357 temp = alpha*a(k,j)
358 DO 180 i = 1,m
359 b(i,j) = b(i,j) + temp*b(i,k)
360 180 CONTINUE
361 END IF
362 190 CONTINUE
363 200 CONTINUE
364 ELSE
365 DO 240 j = 1,n
366 temp = alpha
367 IF (nounit) temp = temp*a(j,j)
368 DO 210 i = 1,m
369 b(i,j) = temp*b(i,j)
370 210 CONTINUE
371 DO 230 k = j + 1,n
372 IF (a(k,j).NE.zero) THEN
373 temp = alpha*a(k,j)
374 DO 220 i = 1,m
375 b(i,j) = b(i,j) + temp*b(i,k)
376 220 CONTINUE
377 END IF
378 230 CONTINUE
379 240 CONTINUE
380 END IF
381 ELSE
382*
383* Form B := alpha*B*A**T or B := alpha*B*A**H.
384*
385 IF (upper) THEN
386 DO 280 k = 1,n
387 DO 260 j = 1,k - 1
388 IF (a(j,k).NE.zero) THEN
389 IF (noconj) THEN
390 temp = alpha*a(j,k)
391 ELSE
392 temp = alpha*dconjg(a(j,k))
393 END IF
394 DO 250 i = 1,m
395 b(i,j) = b(i,j) + temp*b(i,k)
396 250 CONTINUE
397 END IF
398 260 CONTINUE
399 temp = alpha
400 IF (nounit) THEN
401 IF (noconj) THEN
402 temp = temp*a(k,k)
403 ELSE
404 temp = temp*dconjg(a(k,k))
405 END IF
406 END IF
407 IF (temp.NE.one) THEN
408 DO 270 i = 1,m
409 b(i,k) = temp*b(i,k)
410 270 CONTINUE
411 END IF
412 280 CONTINUE
413 ELSE
414 DO 320 k = n,1,-1
415 DO 300 j = k + 1,n
416 IF (a(j,k).NE.zero) THEN
417 IF (noconj) THEN
418 temp = alpha*a(j,k)
419 ELSE
420 temp = alpha*dconjg(a(j,k))
421 END IF
422 DO 290 i = 1,m
423 b(i,j) = b(i,j) + temp*b(i,k)
424 290 CONTINUE
425 END IF
426 300 CONTINUE
427 temp = alpha
428 IF (nounit) THEN
429 IF (noconj) THEN
430 temp = temp*a(k,k)
431 ELSE
432 temp = temp*dconjg(a(k,k))
433 END IF
434 END IF
435 IF (temp.NE.one) THEN
436 DO 310 i = 1,m
437 b(i,k) = temp*b(i,k)
438 310 CONTINUE
439 END IF
440 320 CONTINUE
441 END IF
442 END IF
443 END IF
444*
445 RETURN
446*
447* End of ZTRMM
448*

◆ ztrsm()

subroutine ztrsm ( character side,
character uplo,
character transa,
character diag,
integer m,
integer n,
complex*16 alpha,
complex*16, dimension(lda,*) a,
integer lda,
complex*16, dimension(ldb,*) b,
integer ldb )

ZTRSM

Purpose:
!>
!> ZTRSM  solves one of the matrix equations
!>
!>    op( A )*X = alpha*B,   or   X*op( A ) = alpha*B,
!>
!> where alpha is a scalar, X and B are m by n matrices, A is a unit, or
!> non-unit,  upper or lower triangular matrix  and  op( A )  is one  of
!>
!>    op( A ) = A   or   op( A ) = A**T   or   op( A ) = A**H.
!>
!> The matrix X is overwritten on B.
!> 
Parameters
[in]SIDE
!>          SIDE is CHARACTER*1
!>           On entry, SIDE specifies whether op( A ) appears on the left
!>           or right of X as follows:
!>
!>              SIDE = 'L' or 'l'   op( A )*X = alpha*B.
!>
!>              SIDE = 'R' or 'r'   X*op( A ) = alpha*B.
!> 
[in]UPLO
!>          UPLO is CHARACTER*1
!>           On entry, UPLO specifies whether the matrix A is an upper or
!>           lower triangular matrix as follows:
!>
!>              UPLO = 'U' or 'u'   A is an upper triangular matrix.
!>
!>              UPLO = 'L' or 'l'   A is a lower triangular matrix.
!> 
[in]TRANSA
!>          TRANSA is CHARACTER*1
!>           On entry, TRANSA specifies the form of op( A ) to be used in
!>           the matrix multiplication as follows:
!>
!>              TRANSA = 'N' or 'n'   op( A ) = A.
!>
!>              TRANSA = 'T' or 't'   op( A ) = A**T.
!>
!>              TRANSA = 'C' or 'c'   op( A ) = A**H.
!> 
[in]DIAG
!>          DIAG is CHARACTER*1
!>           On entry, DIAG specifies whether or not A is unit triangular
!>           as follows:
!>
!>              DIAG = 'U' or 'u'   A is assumed to be unit triangular.
!>
!>              DIAG = 'N' or 'n'   A is not assumed to be unit
!>                                  triangular.
!> 
[in]M
!>          M is INTEGER
!>           On entry, M specifies the number of rows of B. M must be at
!>           least zero.
!> 
[in]N
!>          N is INTEGER
!>           On entry, N specifies the number of columns of B.  N must be
!>           at least zero.
!> 
[in]ALPHA
!>          ALPHA is COMPLEX*16
!>           On entry,  ALPHA specifies the scalar  alpha. When  alpha is
!>           zero then  A is not referenced and  B need not be set before
!>           entry.
!> 
[in]A
!>          A is COMPLEX*16 array, dimension ( LDA, k ),
!>           where k is m when SIDE = 'L' or 'l'
!>             and k is n when SIDE = 'R' or 'r'.
!>           Before entry  with  UPLO = 'U' or 'u',  the  leading  k by k
!>           upper triangular part of the array  A must contain the upper
!>           triangular matrix  and the strictly lower triangular part of
!>           A is not referenced.
!>           Before entry  with  UPLO = 'L' or 'l',  the  leading  k by k
!>           lower triangular part of the array  A must contain the lower
!>           triangular matrix  and the strictly upper triangular part of
!>           A is not referenced.
!>           Note that when  DIAG = 'U' or 'u',  the diagonal elements of
!>           A  are not referenced either,  but are assumed to be  unity.
!> 
[in]LDA
!>          LDA is INTEGER
!>           On entry, LDA specifies the first dimension of A as declared
!>           in the calling (sub) program.  When  SIDE = 'L' or 'l'  then
!>           LDA  must be at least  max( 1, m ),  when  SIDE = 'R' or 'r'
!>           then LDA must be at least max( 1, n ).
!> 
[in,out]B
!>          B is COMPLEX*16 array, dimension ( LDB, N )
!>           Before entry,  the leading  m by n part of the array  B must
!>           contain  the  right-hand  side  matrix  B,  and  on exit  is
!>           overwritten by the solution matrix  X.
!> 
[in]LDB
!>          LDB is INTEGER
!>           On entry, LDB specifies the first dimension of B as declared
!>           in  the  calling  (sub)  program.   LDB  must  be  at  least
!>           max( 1, m ).
!> 
Author
Univ. of Tennessee
Univ. of California Berkeley
Univ. of Colorado Denver
NAG Ltd.
Further Details:
!>
!>  Level 3 Blas routine.
!>
!>  -- Written on 8-February-1989.
!>     Jack Dongarra, Argonne National Laboratory.
!>     Iain Duff, AERE Harwell.
!>     Jeremy Du Croz, Numerical Algorithms Group Ltd.
!>     Sven Hammarling, Numerical Algorithms Group Ltd.
!> 

Definition at line 179 of file ztrsm.f.

180*
181* -- Reference BLAS level3 routine --
182* -- Reference BLAS is a software package provided by Univ. of Tennessee, --
183* -- Univ. of California Berkeley, Univ. of Colorado Denver and NAG Ltd..--
184*
185* .. Scalar Arguments ..
186 COMPLEX*16 ALPHA
187 INTEGER LDA,LDB,M,N
188 CHARACTER DIAG,SIDE,TRANSA,UPLO
189* ..
190* .. Array Arguments ..
191 COMPLEX*16 A(LDA,*),B(LDB,*)
192* ..
193*
194* =====================================================================
195*
196* .. External Functions ..
197 LOGICAL LSAME
198 EXTERNAL lsame
199* ..
200* .. External Subroutines ..
201 EXTERNAL xerbla
202* ..
203* .. Intrinsic Functions ..
204 INTRINSIC dconjg,max
205* ..
206* .. Local Scalars ..
207 COMPLEX*16 TEMP
208 INTEGER I,INFO,J,K,NROWA
209 LOGICAL LSIDE,NOCONJ,NOUNIT,UPPER
210* ..
211* .. Parameters ..
212 COMPLEX*16 ONE
213 parameter(one= (1.0d+0,0.0d+0))
214 COMPLEX*16 ZERO
215 parameter(zero= (0.0d+0,0.0d+0))
216* ..
217*
218* Test the input parameters.
219*
220 lside = lsame(side,'L')
221 IF (lside) THEN
222 nrowa = m
223 ELSE
224 nrowa = n
225 END IF
226 noconj = lsame(transa,'T')
227 nounit = lsame(diag,'N')
228 upper = lsame(uplo,'U')
229*
230 info = 0
231 IF ((.NOT.lside) .AND. (.NOT.lsame(side,'R'))) THEN
232 info = 1
233 ELSE IF ((.NOT.upper) .AND. (.NOT.lsame(uplo,'L'))) THEN
234 info = 2
235 ELSE IF ((.NOT.lsame(transa,'N')) .AND.
236 + (.NOT.lsame(transa,'T')) .AND.
237 + (.NOT.lsame(transa,'C'))) THEN
238 info = 3
239 ELSE IF ((.NOT.lsame(diag,'U')) .AND. (.NOT.lsame(diag,'N'))) THEN
240 info = 4
241 ELSE IF (m.LT.0) THEN
242 info = 5
243 ELSE IF (n.LT.0) THEN
244 info = 6
245 ELSE IF (lda.LT.max(1,nrowa)) THEN
246 info = 9
247 ELSE IF (ldb.LT.max(1,m)) THEN
248 info = 11
249 END IF
250 IF (info.NE.0) THEN
251 CALL xerbla('ZTRSM ',info)
252 RETURN
253 END IF
254*
255* Quick return if possible.
256*
257 IF (m.EQ.0 .OR. n.EQ.0) RETURN
258*
259* And when alpha.eq.zero.
260*
261 IF (alpha.EQ.zero) THEN
262 DO 20 j = 1,n
263 DO 10 i = 1,m
264 b(i,j) = zero
265 10 CONTINUE
266 20 CONTINUE
267 RETURN
268 END IF
269*
270* Start the operations.
271*
272 IF (lside) THEN
273 IF (lsame(transa,'N')) THEN
274*
275* Form B := alpha*inv( A )*B.
276*
277 IF (upper) THEN
278 DO 60 j = 1,n
279 IF (alpha.NE.one) THEN
280 DO 30 i = 1,m
281 b(i,j) = alpha*b(i,j)
282 30 CONTINUE
283 END IF
284 DO 50 k = m,1,-1
285 IF (b(k,j).NE.zero) THEN
286 IF (nounit) b(k,j) = b(k,j)/a(k,k)
287 DO 40 i = 1,k - 1
288 b(i,j) = b(i,j) - b(k,j)*a(i,k)
289 40 CONTINUE
290 END IF
291 50 CONTINUE
292 60 CONTINUE
293 ELSE
294 DO 100 j = 1,n
295 IF (alpha.NE.one) THEN
296 DO 70 i = 1,m
297 b(i,j) = alpha*b(i,j)
298 70 CONTINUE
299 END IF
300 DO 90 k = 1,m
301 IF (b(k,j).NE.zero) THEN
302 IF (nounit) b(k,j) = b(k,j)/a(k,k)
303 DO 80 i = k + 1,m
304 b(i,j) = b(i,j) - b(k,j)*a(i,k)
305 80 CONTINUE
306 END IF
307 90 CONTINUE
308 100 CONTINUE
309 END IF
310 ELSE
311*
312* Form B := alpha*inv( A**T )*B
313* or B := alpha*inv( A**H )*B.
314*
315 IF (upper) THEN
316 DO 140 j = 1,n
317 DO 130 i = 1,m
318 temp = alpha*b(i,j)
319 IF (noconj) THEN
320 DO 110 k = 1,i - 1
321 temp = temp - a(k,i)*b(k,j)
322 110 CONTINUE
323 IF (nounit) temp = temp/a(i,i)
324 ELSE
325 DO 120 k = 1,i - 1
326 temp = temp - dconjg(a(k,i))*b(k,j)
327 120 CONTINUE
328 IF (nounit) temp = temp/dconjg(a(i,i))
329 END IF
330 b(i,j) = temp
331 130 CONTINUE
332 140 CONTINUE
333 ELSE
334 DO 180 j = 1,n
335 DO 170 i = m,1,-1
336 temp = alpha*b(i,j)
337 IF (noconj) THEN
338 DO 150 k = i + 1,m
339 temp = temp - a(k,i)*b(k,j)
340 150 CONTINUE
341 IF (nounit) temp = temp/a(i,i)
342 ELSE
343 DO 160 k = i + 1,m
344 temp = temp - dconjg(a(k,i))*b(k,j)
345 160 CONTINUE
346 IF (nounit) temp = temp/dconjg(a(i,i))
347 END IF
348 b(i,j) = temp
349 170 CONTINUE
350 180 CONTINUE
351 END IF
352 END IF
353 ELSE
354 IF (lsame(transa,'N')) THEN
355*
356* Form B := alpha*B*inv( A ).
357*
358 IF (upper) THEN
359 DO 230 j = 1,n
360 IF (alpha.NE.one) THEN
361 DO 190 i = 1,m
362 b(i,j) = alpha*b(i,j)
363 190 CONTINUE
364 END IF
365 DO 210 k = 1,j - 1
366 IF (a(k,j).NE.zero) THEN
367 DO 200 i = 1,m
368 b(i,j) = b(i,j) - a(k,j)*b(i,k)
369 200 CONTINUE
370 END IF
371 210 CONTINUE
372 IF (nounit) THEN
373 temp = one/a(j,j)
374 DO 220 i = 1,m
375 b(i,j) = temp*b(i,j)
376 220 CONTINUE
377 END IF
378 230 CONTINUE
379 ELSE
380 DO 280 j = n,1,-1
381 IF (alpha.NE.one) THEN
382 DO 240 i = 1,m
383 b(i,j) = alpha*b(i,j)
384 240 CONTINUE
385 END IF
386 DO 260 k = j + 1,n
387 IF (a(k,j).NE.zero) THEN
388 DO 250 i = 1,m
389 b(i,j) = b(i,j) - a(k,j)*b(i,k)
390 250 CONTINUE
391 END IF
392 260 CONTINUE
393 IF (nounit) THEN
394 temp = one/a(j,j)
395 DO 270 i = 1,m
396 b(i,j) = temp*b(i,j)
397 270 CONTINUE
398 END IF
399 280 CONTINUE
400 END IF
401 ELSE
402*
403* Form B := alpha*B*inv( A**T )
404* or B := alpha*B*inv( A**H ).
405*
406 IF (upper) THEN
407 DO 330 k = n,1,-1
408 IF (nounit) THEN
409 IF (noconj) THEN
410 temp = one/a(k,k)
411 ELSE
412 temp = one/dconjg(a(k,k))
413 END IF
414 DO 290 i = 1,m
415 b(i,k) = temp*b(i,k)
416 290 CONTINUE
417 END IF
418 DO 310 j = 1,k - 1
419 IF (a(j,k).NE.zero) THEN
420 IF (noconj) THEN
421 temp = a(j,k)
422 ELSE
423 temp = dconjg(a(j,k))
424 END IF
425 DO 300 i = 1,m
426 b(i,j) = b(i,j) - temp*b(i,k)
427 300 CONTINUE
428 END IF
429 310 CONTINUE
430 IF (alpha.NE.one) THEN
431 DO 320 i = 1,m
432 b(i,k) = alpha*b(i,k)
433 320 CONTINUE
434 END IF
435 330 CONTINUE
436 ELSE
437 DO 380 k = 1,n
438 IF (nounit) THEN
439 IF (noconj) THEN
440 temp = one/a(k,k)
441 ELSE
442 temp = one/dconjg(a(k,k))
443 END IF
444 DO 340 i = 1,m
445 b(i,k) = temp*b(i,k)
446 340 CONTINUE
447 END IF
448 DO 360 j = k + 1,n
449 IF (a(j,k).NE.zero) THEN
450 IF (noconj) THEN
451 temp = a(j,k)
452 ELSE
453 temp = dconjg(a(j,k))
454 END IF
455 DO 350 i = 1,m
456 b(i,j) = b(i,j) - temp*b(i,k)
457 350 CONTINUE
458 END IF
459 360 CONTINUE
460 IF (alpha.NE.one) THEN
461 DO 370 i = 1,m
462 b(i,k) = alpha*b(i,k)
463 370 CONTINUE
464 END IF
465 380 CONTINUE
466 END IF
467 END IF
468 END IF
469*
470 RETURN
471*
472* End of ZTRSM
473*