1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140
|
/*
* -- SuperLU routine (version 2.0) --
* Univ. of California Berkeley, Xerox Palo Alto Research Center,
* and Lawrence Berkeley National Lab.
* November 15, 1997
*
*/
/*
Copyright (c) 2003, The Regents of the University of California, through
Lawrence Berkeley National Laboratory (subject to receipt of any required
approvals from U.S. Dept. of Energy)
All rights reserved.
The source code is distributed under BSD license, see the file License.txt
*/
/*
Copyright (c) 1997 by Xerox Corporation. All rights reserved.
THIS MATERIAL IS PROVIDED AS IS, WITH ABSOLUTELY NO WARRANTY
EXPRESSED OR IMPLIED. ANY USE IS AT YOUR OWN RISK.
Permission is hereby granted to use or copy this program for any
purpose, provided the above notices are retained on all copies.
Permission to modify the code and to distribute modified code is
granted, provided the above notices are retained, and a notice that
the code was modified is included with the above copyright notice.
*/
/*
* File name: sp_blas3.c
* Purpose: Sparse BLAS3, using some dense BLAS3 operations.
*/
#include "slu_sdefs.h"
int
sp_sgemm(char *transa, char *transb, int m, int n, int k,
float alpha, SuperMatrix *A, float *b, int ldb,
float beta, float *c, int ldc)
{
/* Purpose
=======
sp_s 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' or op( X ) = conjg( X' ),
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
==========
TRANSA - (input) char*
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'.
TRANSA = 'C' or 'c', op( A ) = conjg( A' ).
Unchanged on exit.
TRANSB - (input) char*
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'.
TRANSB = 'C' or 'c', op( B ) = conjg( B' ).
Unchanged on exit.
M - (input) int
On entry, M specifies the number of rows of the matrix
op( A ) and of the matrix C. M must be at least zero.
Unchanged on exit.
N - (input) int
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.
Unchanged on exit.
K - (input) int
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.
Unchanged on exit.
ALPHA - (input) float
On entry, ALPHA specifies the scalar alpha.
A - (input) SuperMatrix*
Matrix A with a sparse format, of dimension (A->nrow, A->ncol).
Currently, the type of A can be:
Stype = NC or NCP; Dtype = SLU_S; Mtype = GE.
In the future, more general A can be handled.
B - FLOAT PRECISION array of 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.
Unchanged on exit.
LDB - (input) int
On entry, LDB specifies the first dimension of B as declared
in the calling (sub) program. LDB must be at least max( 1, n ).
Unchanged on exit.
BETA - (input) float
On entry, BETA specifies the scalar beta. When BETA is
supplied as zero then C need not be set on input.
C - FLOAT PRECISION array of 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 )*B + beta*C ).
LDC - (input) int
On entry, LDC specifies the first dimension of C as declared
in the calling (sub)program. LDC must be at least max(1,m).
Unchanged on exit.
==== Sparse Level 3 Blas routine.
*/
int incx = 1, incy = 1;
int j;
for (j = 0; j < n; ++j) {
sp_sgemv(transa, alpha, A, &b[ldb*j], incx, beta, &c[ldc*j], incy);
}
return 0;
}
|