Actual source code: bvblas.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    BV private kernels that use the BLAS.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2014, 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,PetscScalar *A,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,A,&m,B,&ldb,&beta,C,&m));
 52:   for (;l<m;l+=bs) {
 53:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&bs,&n,&k,&alpha,A+l,&m,B,&ldb,&beta,C+l,&m));
 54:   }
 55: #else
 56:   PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&n,&k,&alpha,A,&m,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,PetscScalar *A,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,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,PetscScalar *B,PetscBool btrans)
135: {
137:   PetscScalar    zero=0.0,one=1.0,*out,*pv;
138:   PetscBLASInt   m,n,k,l,bs=BLOCKSIZE;
139:   PetscInt       j;
140:   const char     *bt;

143:   PetscBLASIntCast(m_,&m);
144:   PetscBLASIntCast(n_,&n);
145:   PetscBLASIntCast(k_,&k);
146:   BVAllocateWork_Private(bv,2*BLOCKSIZE*n_);
147:   out = bv->work+BLOCKSIZE*n_;
148:   if (btrans) bt = "C";
149:   else bt = "N";
150:   l = m % bs;
151:   if (l) {
152:     for (j=0;j<n;j++) {
153:       VecGetArray(V[j],&pv);
154:       PetscMemcpy(bv->work+j*l,pv,l*sizeof(PetscScalar));
155:       VecRestoreArray(V[j],&pv);
156:     }
157:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N",bt,&l,&n,&n,&one,bv->work,&l,B,&k,&zero,out,&l));
158:     for (j=0;j<n;j++) {
159:       VecGetArray(V[j],&pv);
160:       PetscMemcpy(pv,out+j*l,l*sizeof(PetscScalar));
161:       VecRestoreArray(V[j],&pv);
162:     }
163:   }
164:   for (;l<m;l+=bs) {
165:     for (j=0;j<n;j++) {
166:       VecGetArray(V[j],&pv);
167:       PetscMemcpy(bv->work+j*bs,pv+l,bs*sizeof(PetscScalar));
168:       VecRestoreArray(V[j],&pv);
169:     }
170:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N",bt,&bs,&n,&n,&one,bv->work,&bs,B,&k,&zero,out,&bs));
171:     for (j=0;j<n;j++) {
172:       VecGetArray(V[j],&pv);
173:       PetscMemcpy(pv+l,out+j*bs,bs*sizeof(PetscScalar));
174:       VecRestoreArray(V[j],&pv);
175:     }
176:   }
177:   PetscLogFlops(2.0*n*n*k);
178:   return(0);
179: }

183: /*
184:     B := alpha*A + B

186:     A,B are nxk (ld=n)
187: */
188: PetscErrorCode BVAXPY_BLAS_Private(BV bv,PetscInt n_,PetscInt k_,PetscScalar alpha,PetscScalar *A,PetscScalar *B)
189: {
191:   PetscBLASInt   m,one=1;

194:   PetscBLASIntCast(n_*k_,&m);
195:   PetscStackCallBLAS("BLASaxpy",BLASaxpy_(&m,&alpha,A,&one,B,&one));
196:   PetscLogFlops(2.0*n_*k_);
197:   return(0);
198: }

202: /*
203:     C := A'*B

205:     A' is mxk (ld=k), B is kxn (ld=k), C is mxn (ld=ldc)
206: */
207: PetscErrorCode BVDot_BLAS_Private(BV bv,PetscInt m_,PetscInt n_,PetscInt k_,PetscInt ldc_,PetscScalar *A,PetscScalar *B,PetscScalar *C,PetscBool mpi)
208: {
210:   PetscScalar    zero=0.0,one=1.0,*CC;
211:   PetscBLASInt   m,n,k,ldc,j;

214:   PetscBLASIntCast(m_,&m);
215:   PetscBLASIntCast(n_,&n);
216:   PetscBLASIntCast(k_,&k);
217:   PetscBLASIntCast(ldc_,&ldc);
218:   if (mpi) {
219:     if (ldc==m) {
220:       BVAllocateWork_Private(bv,m*n);
221:       PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&m,&n,&k,&one,A,&k,B,&k,&zero,bv->work,&ldc));
222:       MPI_Allreduce(bv->work,C,m*n,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)bv));
223:     } else {
224:       BVAllocateWork_Private(bv,2*m*n);
225:       CC = bv->work+m*n;
226:       PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&m,&n,&k,&one,A,&k,B,&k,&zero,bv->work,&m));
227:       MPI_Allreduce(bv->work,CC,m*n,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)bv));
228:       for (j=0;j<n;j++) {
229:         PetscMemcpy(C+j*ldc,CC+j*m,m*sizeof(PetscScalar));
230:       }
231:     }
232:   } else {
233:     PetscStackCallBLAS("BLASgemm",BLASgemm_("C","N",&m,&n,&k,&one,A,&k,B,&k,&zero,C,&ldc));
234:   }
235:   PetscLogFlops(2.0*m*n*k);
236:   return(0);
237: }

241: /*
242:     y := A'*x

244:     A is nxk (ld=n)
245: */
246: PetscErrorCode BVDotVec_BLAS_Private(BV bv,PetscInt n_,PetscInt k_,PetscScalar *A,PetscScalar *x,PetscScalar *y,PetscBool mpi)
247: {
249:   PetscScalar    zero=0.0,done=1.0;
250:   PetscBLASInt   n,k,one=1;

253:   PetscBLASIntCast(n_,&n);
254:   PetscBLASIntCast(k_,&k);
255:   if (mpi) {
256:     BVAllocateWork_Private(bv,k);
257:     if (n) PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&k,&done,A,&n,x,&one,&zero,bv->work,&one));
258:     MPI_Allreduce(bv->work,y,k,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)bv));
259:   } else {
260:     if (n) PetscStackCallBLAS("BLASgemv",BLASgemv_("C",&n,&k,&done,A,&n,x,&one,&zero,y,&one));
261:   }
262:   PetscLogFlops(2.0*n*k);
263:   return(0);
264: }

268: /*
269:     Scale n scalars
270: */
271: PetscErrorCode BVScale_BLAS_Private(BV bv,PetscInt n_,PetscScalar *A,PetscScalar alpha)
272: {
274:   PetscBLASInt   n,one=1;

277:   if (alpha == (PetscScalar)0.0) {
278:     PetscMemzero(A,n_*sizeof(PetscScalar));
279:   } else {
280:     PetscBLASIntCast(n_,&n);
281:     PetscStackCallBLAS("BLASscal",BLASscal_(&n,&alpha,A,&one));
282:     PetscLogFlops(n);
283:   }
284:   return(0);
285: }

289: /*
290:     Compute ||A|| for an mxn matrix
291: */
292: PetscErrorCode BVNorm_LAPACK_Private(BV bv,PetscInt m_,PetscInt n_,PetscScalar *A,NormType type,PetscReal *nrm,PetscBool mpi)
293: {
295:   PetscBLASInt   m,n,i,j;
296:   PetscReal      lnrm,*rwork=NULL,*rwork2=NULL;

299:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
300:   PetscBLASIntCast(m_,&m);
301:   PetscBLASIntCast(n_,&n);
302:   if (type==NORM_FROBENIUS || type==NORM_2) {
303:     lnrm = LAPACKlange_("F",&m,&n,A,&m,rwork);
304:     if (mpi) {
305:       lnrm = lnrm*lnrm;
306:       MPI_Allreduce(&lnrm,nrm,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)bv));
307:       *nrm = PetscSqrtReal(*nrm);
308:     } else *nrm = lnrm;
309:     PetscLogFlops(2.0*m*n);
310:   } else if (type==NORM_1) {
311:     if (mpi) {
312:       BVAllocateWork_Private(bv,2*n_);
313:       rwork = (PetscReal*)bv->work;
314:       rwork2 = rwork+n_;
315:       PetscMemzero(rwork,n_*sizeof(PetscReal));
316:       PetscMemzero(rwork2,n_*sizeof(PetscReal));
317:       for (j=0;j<n_;j++) {
318:         for (i=0;i<m_;i++) {
319:           rwork[j] += PetscAbsScalar(A[i+j*m_]);
320:         }
321:       }
322:       MPI_Allreduce(rwork,rwork2,n_,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)bv));
323:       for (j=0;j<n_;j++) if (rwork2[j] > *nrm) *nrm = rwork2[j];
324:     } else {
325:       *nrm = LAPACKlange_("O",&m,&n,A,&m,rwork);
326:     }
327:     PetscLogFlops(1.0*m*n);
328:   } else if (type==NORM_INFINITY) {
329:     BVAllocateWork_Private(bv,m_);
330:     rwork = (PetscReal*)bv->work;
331:     lnrm = LAPACKlange_("I",&m,&n,A,&m,rwork);
332:     if (mpi) {
333:       MPI_Allreduce(&lnrm,nrm,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)bv));
334:     } else *nrm = lnrm;
335:     PetscLogFlops(1.0*m*n);
336:   }
337:   PetscFPTrapPop();
338:   return(0);
339: }

343: /*
344:     QR factorization of an mxn matrix
345: */
346: PetscErrorCode BVOrthogonalize_LAPACK_Private(BV bv,PetscInt m_,PetscInt n_,PetscScalar *Q,PetscScalar *R,PetscBool mpi)
347: {
348: #if defined(PETSC_MISSING_LAPACK_GEQRF) || defined(SLEPC_MISSING_LAPACK_ORGQR)
350:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"GEQRF/ORGQR - Lapack routines are unavailable");
351: #else
353:   PetscBLASInt   m,n,i,j,k,l,nb,lwork,info;
354:   PetscScalar    *tau,*work,*Rl=NULL,*A=NULL,*C=NULL,one=1.0,zero=0.0;
355:   PetscMPIInt    rank,size;

358:   PetscFPTrapPush(PETSC_FP_TRAP_OFF);
359:   PetscBLASIntCast(m_,&m);
360:   PetscBLASIntCast(n_,&n);
361:   k = PetscMin(m,n);
362:   nb = 16;
363:   if (mpi) {
364:     MPI_Comm_rank(PetscObjectComm((PetscObject)bv),&rank);
365:     MPI_Comm_size(PetscObjectComm((PetscObject)bv),&size);
366:     BVAllocateWork_Private(bv,k+n*nb+n*n+n*n*size+m*n);
367:   } else {
368:     BVAllocateWork_Private(bv,k+n*nb);
369:    }
370:   tau = bv->work;
371:   work = bv->work+k;
372:   PetscBLASIntCast(n*nb,&lwork);
373:   if (mpi) {
374:     Rl = bv->work+k+n*nb;
375:     A  = bv->work+k+n*nb+n*n;
376:     C  = bv->work+k+n*nb+n*n+n*n*size;
377:   }

379:   /* Compute QR */
380:   PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&m,&n,Q,&m,tau,work,&lwork,&info));
381:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);

383:   /* Extract R */
384:   if (R || mpi) {
385:     PetscMemzero(mpi? Rl: R,n*n*sizeof(PetscScalar));
386:     for (j=0;j<n;j++) {
387:       for (i=0;i<=j;i++) {
388:         if (mpi) Rl[i+j*n] = Q[i+j*m];
389:         else R[i+j*n] = Q[i+j*m];
390:       }
391:     }
392:   }

394:   /* Compute orthogonal matrix in Q */
395:   PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&m,&n,&k,Q,&m,tau,work,&lwork,&info));
396:   if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);

398:   if (mpi) {

400:     /* Stack triangular matrices */
401:     PetscBLASIntCast(n*size,&l);
402:     for (j=0;j<n;j++) {
403:       MPI_Allgather(Rl+j*n,n,MPIU_SCALAR,A+j*l,n,MPIU_SCALAR,PetscObjectComm((PetscObject)bv));
404:     }

406:     /* Compute QR */
407:     PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&l,&n,A,&l,tau,work,&lwork,&info));
408:     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xGEQRF %d",info);

410:     /* Extract R */
411:     if (R) {
412:       PetscMemzero(R,n*n*sizeof(PetscScalar));
413:       for (j=0;j<n;j++)
414:         for (i=0;i<=j;i++)
415:           R[i+j*n] = A[i+j*l];
416:     }

418:     /* Accumulate orthogonal matrix */
419:     PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&l,&n,&n,A,&l,tau,work,&lwork,&info));
420:     if (info) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in Lapack xORGQR %d",info);
421:     PetscStackCallBLAS("BLASgemm",BLASgemm_("N","N",&m,&n,&n,&one,Q,&m,A+rank*n,&l,&zero,C,&m));
422:     PetscMemcpy(Q,C,m*n*sizeof(PetscScalar));
423:   }

425:   PetscLogFlops(3.0*m*n*n);
426:   PetscFPTrapPop();
427:   return(0);
428: #endif
429: }