Ref #285: added axpby kernels

This commit is contained in:
wernsaar 2014-06-08 11:54:24 +02:00
parent a40116de25
commit faf3ac0aad
18 changed files with 512 additions and 11 deletions

10
cblas.h
View File

@ -305,6 +305,16 @@ void cblas_zher2k(OPENBLAS_CONST enum CBLAS_ORDER Order, OPENBLAS_CONST enum CBL
void cblas_xerbla(blasint p, char *rout, char *form, ...);
/*** BLAS extensions ***/
void cblas_saxpby(OPENBLAS_CONST blasint n, OPENBLAS_CONST float alpha, OPENBLAS_CONST float *x, OPENBLAS_CONST blasint incx,OPENBLAS_CONST float beta, float *y, OPENBLAS_CONST blasint incy);
void cblas_daxpby(OPENBLAS_CONST blasint n, OPENBLAS_CONST double alpha, OPENBLAS_CONST double *x, OPENBLAS_CONST blasint incx,OPENBLAS_CONST double beta, double *y, OPENBLAS_CONST blasint incy);
void cblas_caxpby(OPENBLAS_CONST blasint n, OPENBLAS_CONST float *alpha, OPENBLAS_CONST float *x, OPENBLAS_CONST blasint incx,OPENBLAS_CONST float *beta, float *y, OPENBLAS_CONST blasint incy);
void cblas_zaxpby(OPENBLAS_CONST blasint n, OPENBLAS_CONST double *alpha, OPENBLAS_CONST double *x, OPENBLAS_CONST blasint incx,OPENBLAS_CONST double *beta, double *y, OPENBLAS_CONST blasint incy);
#ifdef __cplusplus
}
#endif /* __cplusplus */

View File

@ -296,6 +296,17 @@ void cblas_zher2k(enum CBLAS_ORDER Order, enum CBLAS_UPLO Uplo, enum CBLAS_TRANS
void cblas_xerbla(blasint p, char *rout, char *form, ...);
/*** BLAS extensions ***/
void cblas_saxpby(blasint n, float alpha, float *x, blasint incx,float beta, float *y, blasint incy);
void cblas_daxpby(blasint n, double alpha, double *x, blasint incx,double beta, double *y, blasint incy);
void cblas_caxpby(blasint n, float *alpha, float *x, blasint incx,float *beta, float *y, blasint incy);
void cblas_zaxpby(blasint n, double *alpha, double *x, blasint incx,double *beta, double *y, blasint incy);
#ifdef __cplusplus
}
#endif /* __cplusplus */

View File

@ -209,6 +209,8 @@
#define CNEG_TCOPY cneg_tcopy
#define CLASWP_NCOPY claswp_ncopy
#define CAXPBY_K caxpby_k
#else
#define CAMAX_K gotoblas -> camax_k
@ -380,6 +382,7 @@
#define CNEG_TCOPY gotoblas -> cneg_tcopy
#define CLASWP_NCOPY gotoblas -> claswp_ncopy
#define CAXPBY_K gotoblas -> caxpby_k
#endif
#define CGEMM_NN cgemm_nn

View File

@ -144,6 +144,8 @@
#define DNEG_TCOPY dneg_tcopy
#define DLASWP_NCOPY dlaswp_ncopy
#define DAXPBY_K daxpby_k
#else
#define DAMAX_K gotoblas -> damax_k
@ -255,6 +257,7 @@
#define DNEG_TCOPY gotoblas -> dneg_tcopy
#define DLASWP_NCOPY gotoblas -> dlaswp_ncopy
#define DAXPBY_K gotoblas -> daxpby_k
#endif
#define DGEMM_NN dgemm_nn

View File

@ -757,6 +757,14 @@ FLOATRET BLASFUNC(slamc3)(float *, float *);
double BLASFUNC(dlamc3)(double *, double *);
xdouble BLASFUNC(qlamc3)(xdouble *, xdouble *);
/* BLAS extensions */
void BLASFUNC(saxpby) (blasint *, float *, float *, blasint *, float *, float *, blasint *);
void BLASFUNC(daxpby) (blasint *, double *, double *, blasint *, double *, double *, blasint *);
void BLASFUNC(caxpby) (blasint *, float *, float *, blasint *, float *, float *, blasint *);
void BLASFUNC(zaxpby) (blasint *, double *, double *, blasint *, double *, double *, blasint *);
#ifdef __cplusplus
}

View File

@ -204,6 +204,13 @@ int srotm_k (BLASLONG, float, BLASLONG, float, BLASLONG, float);
int drotm_k (BLASLONG, double, BLASLONG, double, BLASLONG, double);
int qrotm_k (BLASLONG, xdouble, BLASLONG, xdouble, BLASLONG, xdouble);
int saxpby_k (BLASLONG, float, float *, BLASLONG, float, float *, BLASLONG);
int daxpby_k (BLASLONG, double, double *, BLASLONG, double, double *, BLASLONG);
int caxpby_k (BLASLONG, float, float, float *, BLASLONG, float, float, float *, BLASLONG);
int zaxpby_k (BLASLONG, double, double, double *, BLASLONG, double, double, double *, BLASLONG);
#ifdef __CUDACC__
}
#endif

View File

@ -628,6 +628,8 @@
#define HERK_THREAD_LR DSYRK_THREAD_LN
#define HERK_THREAD_LC DSYRK_THREAD_LT
#define AXPBY_K DAXPBY_K
#else
#define AMAX_K SAMAX_K
@ -918,6 +920,8 @@
#define HERK_THREAD_LR SSYRK_THREAD_LN
#define HERK_THREAD_LC SSYRK_THREAD_LT
#define AXPBY_K SAXPBY_K
#endif
#else
#ifdef XDOUBLE
@ -1722,6 +1726,8 @@
#define SYMM_OUTCOPY ZSYMM_OUTCOPY
#define SYMM_OLTCOPY ZSYMM_OLTCOPY
#define AXPBY_K ZAXPBY_K
#else
#define AMAX_K CAMAX_K
@ -2123,6 +2129,8 @@
#define SYMM_OUTCOPY CSYMM_OUTCOPY
#define SYMM_OLTCOPY CSYMM_OLTCOPY
#define AXPBY_K CAXPBY_K
#endif
#endif

View File

@ -806,10 +806,16 @@ BLASLONG (*ixamin_k)(BLASLONG, xdouble *, BLASLONG);
#endif
void (*init)(void);
int snum_opt, dnum_opt, qnum_opt;
int (*saxpby_k) (BLASLONG, float, float*, BLASLONG,float, float*, BLASLONG);
int (*daxpby_k) (BLASLONG, double, double*, BLASLONG,double, double*, BLASLONG);
int (*caxpby_k) (BLASLONG, float, float, float*, BLASLONG,float,float, float*, BLASLONG);
int (*zaxpby_k) (BLASLONG, double, double, double*, BLASLONG,double,double, double*, BLASLONG);
} gotoblas_t;
extern gotoblas_t *gotoblas;

View File

@ -146,6 +146,8 @@
#define SNEG_TCOPY sneg_tcopy
#define SLASWP_NCOPY slaswp_ncopy
#define SAXPBY_K saxpby_k
#else
#define SAMAX_K gotoblas -> samax_k
@ -259,6 +261,8 @@
#define SNEG_TCOPY gotoblas -> sneg_tcopy
#define SLASWP_NCOPY gotoblas -> slaswp_ncopy
#define SAXPBY_K gotoblas -> saxpby_k
#endif
#define SGEMM_NN sgemm_nn

View File

@ -209,6 +209,8 @@
#define ZNEG_TCOPY zneg_tcopy
#define ZLASWP_NCOPY zlaswp_ncopy
#define ZAXPBY_K zaxpby_k
#else
#define ZAMAX_K gotoblas -> zamax_k
@ -380,6 +382,7 @@
#define ZNEG_TCOPY gotoblas -> zneg_tcopy
#define ZLASWP_NCOPY gotoblas -> zlaswp_ncopy
#define ZAXPBY_K gotoblas -> zaxpby_k
#endif
#define ZGEMM_NN zgemm_nn

View File

@ -22,7 +22,9 @@
zhbmv,zhemm,zhemv,zher2,zher2k,zher,zherk,zhpmv,zhpr2,
zhpr,zrotg,zscal,zswap,zsymm,zsyr2k,zsyrk,ztbmv,
ztbsv,ztpmv,ztpsv,ztrmm,ztrmv,ztrsm,ztrsv, zsymv,
xerbla);
xerbla,
saxpby,daxpby,caxpby,zaxpby
);
@cblasobjs = (
cblas_caxpy, cblas_ccopy, cblas_cdotc, cblas_cdotu, cblas_cgbmv, cblas_cgemm, cblas_cgemv,
@ -49,7 +51,9 @@
cblas_zhemv, cblas_zher2, cblas_zher2k, cblas_zher, cblas_zherk, cblas_zhpmv, cblas_zhpr2,
cblas_zhpr, cblas_zscal, cblas_zswap, cblas_zsymm, cblas_zsyr2k, cblas_zsyrk,
cblas_ztbmv, cblas_ztbsv, cblas_ztpmv, cblas_ztpsv, cblas_ztrmm, cblas_ztrmv, cblas_ztrsm,
cblas_ztrsv, cblas_cdotc_sub, cblas_cdotu_sub, cblas_zdotc_sub, cblas_zdotu_sub );
cblas_ztrsv, cblas_cdotc_sub, cblas_cdotu_sub, cblas_zdotc_sub, cblas_zdotu_sub,
cblas_saxpby,cblas_daxpby,cblas_caxpby,cblas_zaxpby
);
@exblasobjs = (
qamax,qamin,qasum,qaxpy,qcabs1,qcopy,qdot,qgbmv,qgemm,

View File

@ -27,6 +27,7 @@ SBLAS1OBJS = \
smax.$(SUFFIX) samax.$(SUFFIX) ismax.$(SUFFIX) isamax.$(SUFFIX) \
smin.$(SUFFIX) samin.$(SUFFIX) ismin.$(SUFFIX) isamin.$(SUFFIX) \
srot.$(SUFFIX) srotg.$(SUFFIX) srotm.$(SUFFIX) srotmg.$(SUFFIX) \
saxpby.$(SUFFIX)
SBLAS2OBJS = \
sgemv.$(SUFFIX) sger.$(SUFFIX) \
@ -44,11 +45,12 @@ SBLAS3OBJS = \
DBLAS1OBJS = \
daxpy.$(SUFFIX) dswap.$(SUFFIX) \
dcopy.$(SUFFIX) dscal.$(SUFFIX) \
ddot.$(SUFFIX) \
ddot.$(SUFFIX) \
dasum.$(SUFFIX) dnrm2.$(SUFFIX) \
dmax.$(SUFFIX) damax.$(SUFFIX) idmax.$(SUFFIX) idamax.$(SUFFIX) \
dmin.$(SUFFIX) damin.$(SUFFIX) idmin.$(SUFFIX) idamin.$(SUFFIX) \
drot.$(SUFFIX) drotg.$(SUFFIX) drotm.$(SUFFIX) drotmg.$(SUFFIX) \
daxpby.$(SUFFIX)
DBLAS2OBJS = \
dgemv.$(SUFFIX) dger.$(SUFFIX) \
@ -71,6 +73,7 @@ CBLAS1OBJS = \
scamax.$(SUFFIX) icamax.$(SUFFIX) \
scamin.$(SUFFIX) icamin.$(SUFFIX) \
csrot.$(SUFFIX) crotg.$(SUFFIX) \
caxpby.$(SUFFIX)
CBLAS2OBJS = \
cgemv.$(SUFFIX) cgeru.$(SUFFIX) cgerc.$(SUFFIX) \
@ -97,6 +100,7 @@ ZBLAS1OBJS = \
dzamax.$(SUFFIX) izamax.$(SUFFIX) \
dzamin.$(SUFFIX) izamin.$(SUFFIX) \
zdrot.$(SUFFIX) zrotg.$(SUFFIX) \
zaxpby.$(SUFFIX)
ZBLAS2OBJS = \
zgemv.$(SUFFIX) zgeru.$(SUFFIX) zgerc.$(SUFFIX) \
@ -246,7 +250,7 @@ CSBLAS1OBJS = \
cblas_isamax.$(SUFFIX) cblas_sasum.$(SUFFIX) cblas_saxpy.$(SUFFIX) \
cblas_scopy.$(SUFFIX) cblas_sdot.$(SUFFIX) cblas_sdsdot.$(SUFFIX) cblas_dsdot.$(SUFFIX) \
cblas_srot.$(SUFFIX) cblas_srotg.$(SUFFIX) cblas_srotm.$(SUFFIX) cblas_srotmg.$(SUFFIX) \
cblas_sscal.$(SUFFIX) cblas_sswap.$(SUFFIX) cblas_snrm2.$(SUFFIX)
cblas_sscal.$(SUFFIX) cblas_sswap.$(SUFFIX) cblas_snrm2.$(SUFFIX) cblas_saxpby.$(SUFFIX)
CSBLAS2OBJS = \
cblas_sgemv.$(SUFFIX) cblas_sger.$(SUFFIX) cblas_ssymv.$(SUFFIX) cblas_strmv.$(SUFFIX) \
@ -262,7 +266,7 @@ CDBLAS1OBJS = \
cblas_idamax.$(SUFFIX) cblas_dasum.$(SUFFIX) cblas_daxpy.$(SUFFIX) \
cblas_dcopy.$(SUFFIX) cblas_ddot.$(SUFFIX) \
cblas_drot.$(SUFFIX) cblas_drotg.$(SUFFIX) cblas_drotm.$(SUFFIX) cblas_drotmg.$(SUFFIX) \
cblas_dscal.$(SUFFIX) cblas_dswap.$(SUFFIX) cblas_dnrm2.$(SUFFIX)
cblas_dscal.$(SUFFIX) cblas_dswap.$(SUFFIX) cblas_dnrm2.$(SUFFIX) cblas_daxpby.$(SUFFIX)
CDBLAS2OBJS = \
cblas_dgemv.$(SUFFIX) cblas_dger.$(SUFFIX) cblas_dsymv.$(SUFFIX) cblas_dtrmv.$(SUFFIX) \
@ -280,7 +284,8 @@ CCBLAS1OBJS = \
cblas_cdotc.$(SUFFIX) cblas_cdotu.$(SUFFIX) \
cblas_cdotc_sub.$(SUFFIX) cblas_cdotu_sub.$(SUFFIX) \
cblas_cscal.$(SUFFIX) cblas_csscal.$(SUFFIX) \
cblas_cswap.$(SUFFIX) cblas_scnrm2.$(SUFFIX)
cblas_cswap.$(SUFFIX) cblas_scnrm2.$(SUFFIX) \
cblas_caxpby.$(SUFFIX)
CCBLAS2OBJS = \
cblas_cgemv.$(SUFFIX) cblas_cgerc.$(SUFFIX) cblas_cgeru.$(SUFFIX) \
@ -301,7 +306,8 @@ CZBLAS1OBJS = \
cblas_zdotc.$(SUFFIX) cblas_zdotu.$(SUFFIX) \
cblas_zdotc_sub.$(SUFFIX) cblas_zdotu_sub.$(SUFFIX) \
cblas_zscal.$(SUFFIX) cblas_zdscal.$(SUFFIX) \
cblas_zswap.$(SUFFIX) cblas_dznrm2.$(SUFFIX)
cblas_zswap.$(SUFFIX) cblas_dznrm2.$(SUFFIX) \
cblas_zaxpby.$(SUFFIX)
CZBLAS2OBJS = \
cblas_zgemv.$(SUFFIX) cblas_zgerc.$(SUFFIX) cblas_zgeru.$(SUFFIX) \
@ -1991,3 +1997,32 @@ zlarf.$(SUFFIX) zlarf.$(PSUFFIX) : larf.c
xlarf.$(SUFFIX) xlarf.$(PSUFFIX) : larf.c
$(CC) -c $(CFLAGS) $< -o $(@F)
############# BLAS EXTENSIONS #####################################
daxpby.$(SUFFIX) daxpby.$(PSUFFIX) : axpby.c
$(CC) $(CFLAGS) -c $< -o $(@F)
cblas_daxpby.$(SUFFIX) cblas_daxpby.$(PSUFFIX) : axpby.c
$(CC) $(CFLAGS) -DCBLAS -c $< -o $(@F)
saxpby.$(SUFFIX) saxpby.$(PSUFFIX) : axpby.c
$(CC) $(CFLAGS) -c $< -o $(@F)
cblas_saxpby.$(SUFFIX) cblas_saxpby.$(PSUFFIX) : axpby.c
$(CC) $(CFLAGS) -DCBLAS -c $< -o $(@F)
zaxpby.$(SUFFIX) zaxpby.$(PSUFFIX) : zaxpby.c
$(CC) $(CFLAGS) -c $< -o $(@F)
cblas_zaxpby.$(SUFFIX) cblas_zaxpby.$(PSUFFIX) : zaxpby.c
$(CC) $(CFLAGS) -DCBLAS -c $< -o $(@F)
caxpby.$(SUFFIX) caxpby.$(PSUFFIX) : zaxpby.c
$(CC) $(CFLAGS) -c $< -o $(@F)
cblas_caxpby.$(SUFFIX) cblas_caxpby.$(PSUFFIX) : zaxpby.c
$(CC) $(CFLAGS) -DCBLAS -c $< -o $(@F)

72
interface/axpby.c Normal file
View File

@ -0,0 +1,72 @@
/***************************************************************************
Copyright (c) 2013, The OpenBLAS Project
All rights reserved.
Redistribution and use in source and binary forms, with or without
modification, are permitted provided that the following conditions are
met:
1. Redistributions of source code must retain the above copyright
notice, this list of conditions and the following disclaimer.
2. Redistributions in binary form must reproduce the above copyright
notice, this list of conditions and the following disclaimer in
the documentation and/or other materials provided with the
distribution.
3. Neither the name of the OpenBLAS project nor the names of
its contributors may be used to endorse or promote products
derived from this software without specific prior written permission.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
ARE DISCLAIMED. IN NO EVENT SHALL THE OPENBLAS PROJECT OR CONTRIBUTORS BE
LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR
SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY,
OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE
USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
*****************************************************************************/
/******************************************************************
2014/06/07 Saar
******************************************************************/
#include <stdio.h>
#include "common.h"
#ifdef FUNCTION_PROFILE
#include "functable.h"
#endif
#ifndef CBLAS
void NAME(blasint *N, FLOAT *ALPHA, FLOAT *x, blasint *INCX, FLOAT *BETA, FLOAT *y, blasint *INCY)
{
BLASLONG n = *N;
BLASLONG incx = *INCX;
BLASLONG incy = *INCY;
FLOAT alpha = *ALPHA;
FLOAT beta = *BETA;
#else
void CNAME(blasint n, FLOAT alpha, FLOAT *x, blasint incx, FLOAT beta, FLOAT *y, blasint incy)
{
#endif
if (n <= 0) return;
FUNCTION_PROFILE_START();
if (incx < 0) x -= (n - 1) * incx;
if (incy < 0) y -= (n - 1) * incy;
AXPBY_K(n, alpha, x, incx, beta, y, incy);
FUNCTION_PROFILE_END(1, 2 * n, 2 * n);
return;
}

74
interface/zaxpby.c Normal file
View File

@ -0,0 +1,74 @@
/***************************************************************************
Copyright (c) 2013, The OpenBLAS Project
All rights reserved.
Redistribution and use in source and binary forms, with or without
modification, are permitted provided that the following conditions are
met:
1. Redistributions of source code must retain the above copyright
notice, this list of conditions and the following disclaimer.
2. Redistributions in binary form must reproduce the above copyright
notice, this list of conditions and the following disclaimer in
the documentation and/or other materials provided with the
distribution.
3. Neither the name of the OpenBLAS project nor the names of
its contributors may be used to endorse or promote products
derived from this software without specific prior written permission.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
ARE DISCLAIMED. IN NO EVENT SHALL THE OPENBLAS PROJECT OR CONTRIBUTORS BE
LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR
SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY,
OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE
USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
*****************************************************************************/
/**********************************************************************
2014/06/07 Saar
**********************************************************************/
#include <stdio.h>
#include "common.h"
#ifdef FUNCTION_PROFILE
#include "functable.h"
#endif
#ifndef CBLAS
void NAME(blasint *N, FLOAT *ALPHA, FLOAT *x, blasint *INCX, FLOAT *BETA, FLOAT *y, blasint *INCY)
{
blasint n = *N;
blasint incx = *INCX;
blasint incy = *INCY;
#else
void CNAME(blasint n, FLOAT *ALPHA, FLOAT *x, blasint incx, FLOAT *BETA, FLOAT *y, blasint incy)
{
#endif
if (n <= 0) return;
FLOAT alpha_r = *(ALPHA + 0);
FLOAT alpha_i = *(ALPHA + 1);
FLOAT beta_r = *(BETA + 0);
FLOAT beta_i = *(BETA + 1);
FUNCTION_PROFILE_START();
if (incx < 0) x -= (n - 1) * incx * 2;
if (incy < 0) y -= (n - 1) * incy * 2;
AXPBY_K (n, alpha_r, alpha_i, x, incx, beta_r, beta_i, y, incy);
FUNCTION_PROFILE_END(4, 2 * n, 2 * n);
return;
}

View File

@ -432,18 +432,38 @@ ifndef LSAME_KERNEL
LSAME_KERNEL = lsame.S
endif
### AXPBY ###
ifndef SAXPBYKERNEL
SAXPBYKERNEL = ../arm/axpby.c
endif
ifndef DAXPBYKERNEL
DAXPBYKERNEL = ../arm/axpby.c
endif
ifndef CAXPBYKERNEL
CAXPBYKERNEL = ../arm/zaxpby.c
endif
ifndef ZAXPBYKERNEL
ZAXPBYKERNEL = ../arm/zaxpby.c
endif
SBLASOBJS += \
samax_k$(TSUFFIX).$(SUFFIX) samin_k$(TSUFFIX).$(SUFFIX) smax_k$(TSUFFIX).$(SUFFIX) smin_k$(TSUFFIX).$(SUFFIX) \
isamax_k$(TSUFFIX).$(SUFFIX) isamin_k$(TSUFFIX).$(SUFFIX) ismax_k$(TSUFFIX).$(SUFFIX) ismin_k$(TSUFFIX).$(SUFFIX) \
sasum_k$(TSUFFIX).$(SUFFIX) saxpy_k$(TSUFFIX).$(SUFFIX) scopy_k$(TSUFFIX).$(SUFFIX) \
sdot_k$(TSUFFIX).$(SUFFIX) sdsdot_k$(TSUFFIX).$(SUFFIX) dsdot_k$(TSUFFIX).$(SUFFIX) \
snrm2_k$(TSUFFIX).$(SUFFIX) srot_k$(TSUFFIX).$(SUFFIX) sscal_k$(TSUFFIX).$(SUFFIX) sswap_k$(TSUFFIX).$(SUFFIX)
snrm2_k$(TSUFFIX).$(SUFFIX) srot_k$(TSUFFIX).$(SUFFIX) sscal_k$(TSUFFIX).$(SUFFIX) sswap_k$(TSUFFIX).$(SUFFIX) \
saxpby_k$(TSUFFIX).$(SUFFIX)
DBLASOBJS += \
damax_k$(TSUFFIX).$(SUFFIX) damin_k$(TSUFFIX).$(SUFFIX) dmax_k$(TSUFFIX).$(SUFFIX) dmin_k$(TSUFFIX).$(SUFFIX) \
idamax_k$(TSUFFIX).$(SUFFIX) idamin_k$(TSUFFIX).$(SUFFIX) idmax_k$(TSUFFIX).$(SUFFIX) idmin_k$(TSUFFIX).$(SUFFIX) \
dasum_k$(TSUFFIX).$(SUFFIX) daxpy_k$(TSUFFIX).$(SUFFIX) dcopy_k$(TSUFFIX).$(SUFFIX) ddot_k$(TSUFFIX).$(SUFFIX) \
dnrm2_k$(TSUFFIX).$(SUFFIX) drot_k$(TSUFFIX).$(SUFFIX) dscal_k$(TSUFFIX).$(SUFFIX) dswap_k$(TSUFFIX).$(SUFFIX)
dnrm2_k$(TSUFFIX).$(SUFFIX) drot_k$(TSUFFIX).$(SUFFIX) dscal_k$(TSUFFIX).$(SUFFIX) dswap_k$(TSUFFIX).$(SUFFIX) \
daxpby_k$(TSUFFIX).$(SUFFIX)
QBLASOBJS += \
qamax_k$(TSUFFIX).$(SUFFIX) qamin_k$(TSUFFIX).$(SUFFIX) qmax_k$(TSUFFIX).$(SUFFIX) qmin_k$(TSUFFIX).$(SUFFIX) \
@ -455,13 +475,13 @@ CBLASOBJS += \
camax_k$(TSUFFIX).$(SUFFIX) camin_k$(TSUFFIX).$(SUFFIX) icamax_k$(TSUFFIX).$(SUFFIX) icamin_k$(TSUFFIX).$(SUFFIX) \
casum_k$(TSUFFIX).$(SUFFIX) caxpy_k$(TSUFFIX).$(SUFFIX) caxpyc_k$(TSUFFIX).$(SUFFIX) ccopy_k$(TSUFFIX).$(SUFFIX) \
cdotc_k$(TSUFFIX).$(SUFFIX) cdotu_k$(TSUFFIX).$(SUFFIX) cnrm2_k$(TSUFFIX).$(SUFFIX) csrot_k$(TSUFFIX).$(SUFFIX) \
cscal_k$(TSUFFIX).$(SUFFIX) cswap_k$(TSUFFIX).$(SUFFIX)
cscal_k$(TSUFFIX).$(SUFFIX) cswap_k$(TSUFFIX).$(SUFFIX) caxpby_k$(TSUFFIX).$(SUFFIX)
ZBLASOBJS += \
zamax_k$(TSUFFIX).$(SUFFIX) zamin_k$(TSUFFIX).$(SUFFIX) izamax_k$(TSUFFIX).$(SUFFIX) izamin_k$(TSUFFIX).$(SUFFIX) \
zasum_k$(TSUFFIX).$(SUFFIX) zaxpy_k$(TSUFFIX).$(SUFFIX) zaxpyc_k$(TSUFFIX).$(SUFFIX) zcopy_k$(TSUFFIX).$(SUFFIX) \
zdotc_k$(TSUFFIX).$(SUFFIX) zdotu_k$(TSUFFIX).$(SUFFIX) znrm2_k$(TSUFFIX).$(SUFFIX) zdrot_k$(TSUFFIX).$(SUFFIX) \
zscal_k$(TSUFFIX).$(SUFFIX) zswap_k$(TSUFFIX).$(SUFFIX)
zscal_k$(TSUFFIX).$(SUFFIX) zswap_k$(TSUFFIX).$(SUFFIX) zaxpby_k$(TSUFFIX).$(SUFFIX)
XBLASOBJS += \
xamax_k$(TSUFFIX).$(SUFFIX) xamin_k$(TSUFFIX).$(SUFFIX) ixamax_k$(TSUFFIX).$(SUFFIX) ixamin_k$(TSUFFIX).$(SUFFIX) \
@ -765,3 +785,17 @@ $(KDIR)zswap_k$(TSUFFIX).$(SUFFIX) $(KDIR)zswap_k$(TPSUFFIX).$(PSUFFIX) : $(KE
$(KDIR)xswap_k$(TSUFFIX).$(SUFFIX) $(KDIR)xswap_k$(TPSUFFIX).$(PSUFFIX) : $(KERNELDIR)/$(XSWAPKERNEL)
$(CC) -c $(CFLAGS) -DCOMPLEX -DXDOUBLE $< -o $@
$(KDIR)saxpby_k$(TSUFFIX).$(SUFFIX) $(KDIR)saxpby_k$(TPSUFFIX).$(PSUFFIX) : $(KERNELDIR)/$(SAXPBYKERNEL)
$(CC) -c $(CFLAGS) -UCOMPLEX -UCOMPLEX -UDOUBLE $< -o $@
$(KDIR)daxpby_k$(TSUFFIX).$(SUFFIX) $(KDIR)daxpby_k$(TPSUFFIX).$(PSUFFIX) : $(KERNELDIR)/$(DAXPBYKERNEL)
$(CC) -c $(CFLAGS) -UCOMPLEX -UCOMPLEX -DDOUBLE $< -o $@
$(KDIR)caxpby_k$(TSUFFIX).$(SUFFIX) $(KDIR)caxpby_k$(TPSUFFIX).$(PSUFFIX) : $(KERNELDIR)/$(CAXPBYKERNEL)
$(CC) -c $(CFLAGS) -DCOMPLEX -DCOMPLEX -UCONJ -UDOUBLE $< -o $@
$(KDIR)zaxpby_k$(TSUFFIX).$(SUFFIX) $(KDIR)zaxpby_k$(TPSUFFIX).$(PSUFFIX) : $(KERNELDIR)/$(ZAXPBYKERNEL)
$(CC) -c $(CFLAGS) -DCOMPLEX -DCOMPLEX -UCONJ -DDOUBLE $< -o $@

96
kernel/arm/axpby.c Normal file
View File

@ -0,0 +1,96 @@
/***************************************************************************
Copyright (c) 2013, The OpenBLAS Project
All rights reserved.
Redistribution and use in source and binary forms, with or without
modification, are permitted provided that the following conditions are
met:
1. Redistributions of source code must retain the above copyright
notice, this list of conditions and the following disclaimer.
2. Redistributions in binary form must reproduce the above copyright
notice, this list of conditions and the following disclaimer in
the documentation and/or other materials provided with the
distribution.
3. Neither the name of the OpenBLAS project nor the names of
its contributors may be used to endorse or promote products
derived from this software without specific prior written permission.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
ARE DISCLAIMED. IN NO EVENT SHALL THE OPENBLAS PROJECT OR CONTRIBUTORS BE
LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR
SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY,
OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE
USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
*****************************************************************************/
#include "common.h"
int CNAME(BLASLONG n, FLOAT alpha, FLOAT *x, BLASLONG inc_x, FLOAT beta, FLOAT *y, BLASLONG inc_y)
{
BLASLONG i=0;
BLASLONG ix,iy;
if ( n < 0 ) return(0);
ix = 0;
iy = 0;
if ( beta == 0.0 )
{
if ( alpha == 0.0 )
{
while(i < n)
{
y[iy] = 0.0 ;
iy += inc_y ;
i++ ;
}
}
else
{
while(i < n)
{
y[iy] = alpha * x[ix] ;
ix += inc_x ;
iy += inc_y ;
i++ ;
}
}
}
else
{
if ( alpha == 0.0 )
{
while(i < n)
{
y[iy] = beta * y[iy] ;
iy += inc_y ;
i++ ;
}
}
else
{
while(i < n)
{
y[iy] = alpha * x[ix] + beta * y[iy] ;
ix += inc_x ;
iy += inc_y ;
i++ ;
}
}
}
return(0);
}

117
kernel/arm/zaxpby.c Normal file
View File

@ -0,0 +1,117 @@
/***************************************************************************
Copyright (c) 2013, The OpenBLAS Project
All rights reserved.
Redistribution and use in source and binary forms, with or without
modification, are permitted provided that the following conditions are
met:
1. Redistributions of source code must retain the above copyright
notice, this list of conditions and the following disclaimer.
2. Redistributions in binary form must reproduce the above copyright
notice, this list of conditions and the following disclaimer in
the documentation and/or other materials provided with the
distribution.
3. Neither the name of the OpenBLAS project nor the names of
its contributors may be used to endorse or promote products
derived from this software without specific prior written permission.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
ARE DISCLAIMED. IN NO EVENT SHALL THE OPENBLAS PROJECT OR CONTRIBUTORS BE
LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR
SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY,
OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE
USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
*****************************************************************************/
/***************************************************************************
* 2014/06/07 Saar
*
***************************************************************************/
#include "common.h"
int CNAME(BLASLONG n, FLOAT alpha_r, FLOAT alpha_i, FLOAT *x, BLASLONG inc_x, FLOAT beta_r, FLOAT beta_i,FLOAT *y, BLASLONG inc_y)
{
BLASLONG i=0;
BLASLONG ix,iy;
FLOAT temp;
if ( n < 0 ) return(0);
ix = 0;
iy = 0;
BLASLONG inc_x2 = 2 * inc_x;
BLASLONG inc_y2 = 2 * inc_y;
if ( beta_r == 0.0 && beta_i == 0.0)
{
if ( alpha_r == 0.0 && alpha_i == 0.0 )
{
while(i < n)
{
y[iy] = 0.0 ;
y[iy+1] = 0.0 ;
iy += inc_y2 ;
i++ ;
}
}
else
{
while(i < n)
{
y[iy] = ( alpha_r * x[ix] - alpha_i * x[ix+1] ) ;
y[iy+1] = ( alpha_r * x[ix+1] + alpha_i * x[ix] ) ;
ix += inc_x2 ;
iy += inc_y2 ;
i++ ;
}
}
}
else
{
if ( alpha_r == 0.0 && alpha_i == 0.0 )
{
while(i < n)
{
temp = ( beta_r * y[iy] - beta_i * y[iy+1] ) ;
y[iy+1] = ( beta_r * y[iy+1] + beta_i * y[iy] ) ;
y[iy] = temp;
iy += inc_y2 ;
i++ ;
}
}
else
{
while(i < n)
{
temp = ( alpha_r * x[ix] - alpha_i * x[ix+1] ) + ( beta_r * y[iy] - beta_i * y[iy+1] ) ;
y[iy+1] = ( alpha_r * x[ix+1] + alpha_i * x[ix] ) + ( beta_r * y[iy+1] + beta_i * y[iy] ) ;
y[iy] = temp;
ix += inc_x2 ;
iy += inc_y2 ;
i++ ;
}
}
}
return(0);
}

View File

@ -500,6 +500,12 @@ gotoblas_t TABLE_NAME = {
SNUMOPT, DNUMOPT, QNUMOPT,
saxpby_kTS,
daxpby_kTS,
caxpby_kTS,
zaxpby_kTS
};
#ifdef ARCH_X86