Actual source code: bvblas.c
slepc-3.6.1 2015-09-03
1: /*
2: BV private kernels that use the BLAS.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc/private/bvimpl.h>
25: #include <slepcblaslapack.h>
27: #define BLOCKSIZE 64
31: /*
32: C := alpha*A*B + beta*C
34: A is mxk (ld=m), B is kxn (ld=ldb), C is mxn (ld=m)
35: */
36: PetscErrorCode BVMult_BLAS_Private(BV bv,PetscInt m_,PetscInt n_,PetscInt k_,PetscInt ldb_,PetscScalar alpha,const PetscScalar *A,const PetscScalar *B,PetscScalar beta,PetscScalar *C)
37: {
39: PetscBLASInt m,n,k,ldb;
40: #if defined(PETSC_HAVE_FBLASLAPACK) || defined(PETSC_HAVE_F2CBLASLAPACK)
41: PetscBLASInt l,bs=BLOCKSIZE;
42: #endif
45: PetscBLASIntCast(m_,&m);
46: PetscBLASIntCast(n_,&n);
47: PetscBLASIntCast(k_,&k);
48: PetscBLASIntCast(ldb_,&ldb);
49: #if defined(PETSC_HAVE_FBLASLAPACK) || defined(PETSC_HAVE_F2CBLASLAPACK)
50: l = m % bs;
51: if (l) PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&l,&n,&k,&alpha,(PetscScalar*)A,&m,(PetscScalar*)B,&ldb,&beta,C,&m));
52: for (;l<m;l+=bs) {
53: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&bs,&n,&k,&alpha,(PetscScalar*)A+l,&m,(PetscScalar*)B,&ldb,&beta,C+l,&m));
54: }
55: #else
56: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&n,&k,&alpha,(PetscScalar*)A,&m,(PetscScalar*)B,&ldb,&beta,C,&m));
57: #endif
58: PetscLogFlops(2.0*m*n*k);
59: return(0);
60: }
64: /*
65: y := alpha*A*x + beta*y
67: A is nxk (ld=n)
68: */
69: PetscErrorCode BVMultVec_BLAS_Private(BV bv,PetscInt n_,PetscInt k_,PetscScalar alpha,const PetscScalar *A,const PetscScalar *x,PetscScalar beta,PetscScalar *y)
70: {
72: PetscBLASInt n,k,one=1;
75: PetscBLASIntCast(n_,&n);
76: PetscBLASIntCast(k_,&k);
77: if (n) PetscStackCallBLAS("BLASgemv",BLASgemv_("N",&n,&k,&alpha,A,&n,x,&one,&beta,y,&one));
78: PetscLogFlops(2.0*n*k);
79: return(0);
80: }
84: /*
85: A(:,s:e-1) := A*B(:,s:e-1)
87: A is mxk (ld=m), B is kxn (ld=ldb) n=e-s
88: */
89: PetscErrorCode BVMultInPlace_BLAS_Private(BV bv,PetscInt m_,PetscInt k_,PetscInt ldb_,PetscInt s,PetscInt e,PetscScalar *A,const PetscScalar *B,PetscBool btrans)
90: {
92: PetscScalar *pb,zero=0.0,one=1.0;
93: PetscBLASInt m,n,k,l,ldb,bs=BLOCKSIZE;
94: PetscInt j,n_=e-s;
95: const char *bt;
98: PetscBLASIntCast(m_,&m);
99: PetscBLASIntCast(n_,&n);
100: PetscBLASIntCast(k_,&k);
101: PetscBLASIntCast(ldb_,&ldb);
102: BVAllocateWork_Private(bv,BLOCKSIZE*n_);
103: if (btrans) {
104: pb = (PetscScalar*)B+s;
105: bt = "C";
106: } else {
107: pb = (PetscScalar*)B+s*ldb;
108: bt = "N";
109: }
110: l = m % bs;
111: if (l) {
112: PetscStackCallBLAS("BLASgemm",BLASgemm_("N",bt,&l,&n,&k,&one,A,&m,pb,&ldb,&zero,bv->work,&l));
113: for (j=0;j<n;j++) {
114: PetscMemcpy(A+(s+j)*m,bv->work+j*l,l*sizeof(PetscScalar));
115: }
116: }
117: for (;l<m;l+=bs) {
118: PetscStackCallBLAS("BLASgemm",BLASgemm_("N",bt,&bs,&n,&k,&one,A+l,&m,pb,&ldb,&zero,bv->work,&bs));
119: for (j=0;j<n;j++) {
120: PetscMemcpy(A+(s+j)*m+l,bv->work+j*bs,bs*sizeof(PetscScalar));
121: }
122: }
123: PetscLogFlops(2.0*m*n*k);
124: return(0);
125: }
129: /*
130: V := V*B
132: V is mxn (ld=m), B is nxn (ld=k)
133: */
134: PetscErrorCode BVMultInPlace_Vecs_Private(BV bv,PetscInt m_,PetscInt n_,PetscInt k_,Vec *V,const PetscScalar *B,PetscBool btrans)
135: {
136: PetscErrorCode ierr;
137: PetscScalar zero=0.0,one=1.0,*out,*pout;
138: const PetscScalar *pin;
139: PetscBLASInt m,n,k,l,bs=BLOCKSIZE;
140: PetscInt j;
141: const char *bt;
144: PetscBLASIntCast(m_,&m);
145: PetscBLASIntCast(n_,&n);
146: PetscBLASIntCast(k_,&k);
147: BVAllocateWork_Private(bv,2*BLOCKSIZE*n_);
148: out = bv->work+BLOCKSIZE*n_;
149: if (btrans) bt = "C";
150: else bt = "N";
151: l = m % bs;
152: if (l) {
153: for (j=0;j<n;j++) {
154: VecGetArrayRead(V[j],&pin);
155: PetscMemcpy(bv->work+j*l,pin,l*sizeof(PetscScalar));
156: VecRestoreArrayRead(V[j],&pin);
157: }
158: PetscStackCallBLAS("BLASgemm",BLASgemm_("N",bt,&l,&n,&n,&one,bv->work,&l,(PetscScalar*)B,&k,&zero,out,&l));
159: for (j=0;j<n;j++) {
160: VecGetArray(V[j],&pout);
161: PetscMemcpy(pout,out+j*l,l*sizeof(PetscScalar));
162: VecRestoreArray(V[j],&pout);
163: }
164: }
165: for (;l<m;l+=bs) {
166: for (j=0;j<n;j++) {
167: VecGetArrayRead(V[j],&pin);
168: PetscMemcpy(bv->work+j*bs,pin+l,bs*sizeof(PetscScalar));
169: VecRestoreArrayRead(V[j],&pin);
170: }
171: PetscStackCallBLAS("BLASgemm",BLASgemm_("N",bt,&bs,&n,&n,&one,bv->work,&bs,(PetscScalar*)B,&k,&zero,out,&bs));
172: for (j=0;j<n;j++) {
173: VecGetArray(V[j],&pout);
174: PetscMemcpy(pout+l,out+j*bs,bs*sizeof(PetscScalar));
175: VecRestoreArray(V[j],&pout);
176: }
177: }
178: PetscLogFlops(2.0*n*n*k);
179: return(0);
180: }
184: /*
185: B := alpha*A + B
187: A,B are nxk (ld=n)
188: */
189: PetscErrorCode BVAXPY_BLAS_Private(BV bv,PetscInt n_,PetscInt k_,PetscScalar alpha,const PetscScalar *A,PetscScalar *B)
190: {
192: PetscBLASInt m,one=1;
195: PetscBLASIntCast(n_*k_,&m);
196: PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&m,&alpha,A,&one,B,&one));
197: PetscLogFlops(2.0*n_*k_);
198: return(0);
199: }
203: /*
204: C := A'*B
206: A' is mxk (ld=k), B is kxn (ld=k), C is mxn (ld=ldc)
207: */
208: PetscErrorCode BVDot_BLAS_Private(BV bv,PetscInt m_,PetscInt n_,PetscInt k_,PetscInt ldc_,const PetscScalar *A,const PetscScalar *B,PetscScalar *C,PetscBool mpi)
209: {
211: PetscScalar zero=0.0,one=1.0,*CC;
212: PetscBLASInt m,n,k,ldc,j;
215: PetscBLASIntCast(m_,&m);
216: PetscBLASIntCast(n_,&n);
217: PetscBLASIntCast(k_,&k);
218: PetscBLASIntCast(ldc_,&ldc);
219: if (mpi) {
220: if (ldc==m) {
221: BVAllocateWork_Private(bv,m*n);
222: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&m,&n,&k,&one,(PetscScalar*)A,&k,(PetscScalar*)B,&k,&zero,bv->work,&ldc));
223: MPI_Allreduce(bv->work,C,m*n,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)bv));
224: } else {
225: BVAllocateWork_Private(bv,2*m*n);
226: CC = bv->work+m*n;
227: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&m,&n,&k,&one,(PetscScalar*)A,&k,(PetscScalar*)B,&k,&zero,bv->work,&m));
228: MPI_Allreduce(bv->work,CC,m*n,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)bv));
229: for (j=0;j<n;j++) {
230: PetscMemcpy(C+j*ldc,CC+j*m,m*sizeof(PetscScalar));
231: }
232: }
233: } else {
234: PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&m,&n,&k,&one,(PetscScalar*)A,&k,(PetscScalar*)B,&k,&zero,C,&ldc));
235: }
236: PetscLogFlops(2.0*m*n*k);
237: return(0);
238: }
242: /*
243: y := A'*x
245: A is nxk (ld=n)
246: */
247: PetscErrorCode BVDotVec_BLAS_Private(BV bv,PetscInt n_,PetscInt k_,const PetscScalar *A,const PetscScalar *x,PetscScalar *y,PetscBool mpi)
248: {
250: PetscScalar zero=0.0,done=1.0;
251: PetscBLASInt n,k,one=1;
254: PetscBLASIntCast(n_,&n);
255: PetscBLASIntCast(k_,&k);
256: if (mpi) {
257: BVAllocateWork_Private(bv,k);
258: if (n) {
259: PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&k,&done,A,&n,x,&one,&zero,bv->work,&one));
260: } else {
261: PetscMemzero(bv->work,k*sizeof(PetscScalar));
262: }
263: MPI_Allreduce(bv->work,y,k,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)bv));
264: } else {
265: if (n) PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&k,&done,A,&n,x,&one,&zero,y,&one));
266: }
267: PetscLogFlops(2.0*n*k);
268: return(0);
269: }
273: /*
274: Scale n scalars
275: */
276: PetscErrorCode BVScale_BLAS_Private(BV bv,PetscInt n_,PetscScalar *A,PetscScalar alpha)
277: {
279: PetscBLASInt n,one=1;
282: if (alpha == (PetscScalar)0.0) {
283: PetscMemzero(A,n_*sizeof(PetscScalar));
284: } else {
285: PetscBLASIntCast(n_,&n);
286: PetscStackCallBLAS("BLASscal",BLASscal_(&n,&alpha,A,&one));
287: PetscLogFlops(n);
288: }
289: return(0);
290: }
294: /*
295: Compute ||A|| for an mxn matrix
296: */
297: PetscErrorCode BVNorm_LAPACK_Private(BV bv,PetscInt m_,PetscInt n_,const PetscScalar *A,NormType type,PetscReal *nrm,PetscBool mpi)
298: {
300: PetscBLASInt m,n,i,j;
301: PetscReal lnrm,*rwork=NULL,*rwork2=NULL;
304: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
305: PetscBLASIntCast(m_,&m);
306: PetscBLASIntCast(n_,&n);
307: if (type==NORM_FROBENIUS || type==NORM_2) {
308: lnrm = LAPACKlange_("F",&m,&n,(PetscScalar*)A,&m,rwork);
309: if (mpi) {
310: lnrm = lnrm*lnrm;
311: MPI_Allreduce(&lnrm,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)bv));
312: *nrm = PetscSqrtReal(*nrm);
313: } else *nrm = lnrm;
314: PetscLogFlops(2.0*m*n);
315: } else if (type==NORM_1) {
316: if (mpi) {
317: BVAllocateWork_Private(bv,2*n_);
318: rwork = (PetscReal*)bv->work;
319: rwork2 = rwork+n_;
320: PetscMemzero(rwork,n_*sizeof(PetscReal));
321: PetscMemzero(rwork2,n_*sizeof(PetscReal));
322: for (j=0;j<n_;j++) {
323: for (i=0;i<m_;i++) {
324: rwork[j] += PetscAbsScalar(A[i+j*m_]);
325: }
326: }
327: MPI_Allreduce(rwork,rwork2,n_,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)bv));
328: for (j=0;j<n_;j++) if (rwork2[j] > *nrm) *nrm = rwork2[j];
329: } else {
330: *nrm = LAPACKlange_("O",&m,&n,(PetscScalar*)A,&m,rwork);
331: }
332: PetscLogFlops(1.0*m*n);
333: } else if (type==NORM_INFINITY) {
334: BVAllocateWork_Private(bv,m_);
335: rwork = (PetscReal*)bv->work;
336: lnrm = LAPACKlange_("I",&m,&n,(PetscScalar*)A,&m,rwork);
337: if (mpi) {
338: MPI_Allreduce(&lnrm,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)bv));
339: } else *nrm = lnrm;
340: PetscLogFlops(1.0*m*n);
341: }
342: PetscFPTrapPop();
343: return(0);
344: }
348: /*
349: QR factorization of an mxn matrix
350: */
351: PetscErrorCode BVOrthogonalize_LAPACK_Private(BV bv,PetscInt m_,PetscInt n_,PetscScalar *Q,PetscScalar *R,PetscBool mpi)
352: {
353: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(SLEPC_MISSING_LAPACK_ORGQR)
355: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/ORGQR - Lapack routines are unavailable");
356: #else
358: PetscBLASInt m,n,i,j,k,l,nb,lwork,info;
359: PetscScalar *tau,*work,*Rl=NULL,*A=NULL,*C=NULL,one=1.0,zero=0.0;
360: PetscMPIInt rank,size;
363: PetscFPTrapPush(PETSC_FP_TRAP_OFF);
364: PetscBLASIntCast(m_,&m);
365: PetscBLASIntCast(n_,&n);
366: k = PetscMin(m,n);
367: nb = 16;
368: if (mpi) {
369: MPI_Comm_rank(PetscObjectComm((PetscObject)bv),&rank);
370: MPI_Comm_size(PetscObjectComm((PetscObject)bv),&size);
371: BVAllocateWork_Private(bv,k+n*nb+n*n+n*n*size+m*n);
372: } else {
373: BVAllocateWork_Private(bv,k+n*nb);
374: }
375: tau = bv->work;
376: work = bv->work+k;
377: PetscBLASIntCast(n*nb,&lwork);
378: if (mpi) {
379: Rl = bv->work+k+n*nb;
380: A = bv->work+k+n*nb+n*n;
381: C = bv->work+k+n*nb+n*n+n*n*size;
382: }
384: /* Compute QR */
385: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&m,&n,Q,&m,tau,work,&lwork,&info));
386: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
388: /* Extract R */
389: if (R || mpi) {
390: PetscMemzero(mpi? Rl: R,n*n*sizeof(PetscScalar));
391: for (j=0;j<n;j++) {
392: for (i=0;i<=j;i++) {
393: if (mpi) Rl[i+j*n] = Q[i+j*m];
394: else R[i+j*n] = Q[i+j*m];
395: }
396: }
397: }
399: /* Compute orthogonal matrix in Q */
400: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&m,&n,&k,Q,&m,tau,work,&lwork,&info));
401: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
403: if (mpi) {
405: /* Stack triangular matrices */
406: PetscBLASIntCast(n*size,&l);
407: for (j=0;j<n;j++) {
408: MPI_Allgather(Rl+j*n,n,MPIU_SCALAR,A+j*l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)bv));
409: }
411: /* Compute QR */
412: PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&l,&n,A,&l,tau,work,&lwork,&info));
413: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);
415: /* Extract R */
416: if (R) {
417: PetscMemzero(R,n*n*sizeof(PetscScalar));
418: for (j=0;j<n;j++)
419: for (i=0;i<=j;i++)
420: R[i+j*n] = A[i+j*l];
421: }
423: /* Accumulate orthogonal matrix */
424: PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&l,&n,&n,A,&l,tau,work,&lwork,&info));
425: if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
426: PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&n,&n,&one,Q,&m,A+rank*n,&l,&zero,C,&m));
427: PetscMemcpy(Q,C,m*n*sizeof(PetscScalar));
428: }
430: PetscLogFlops(3.0*m*n*n);
431: PetscFPTrapPop();
432: return(0);
433: #endif
434: }