Actual source code: bvorthogcuda.cu

slepc-3.15.0 2021-03-31
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    BV orthogonalization routines (CUDA)
 12: */

 14: #include <slepc/private/bvimpl.h>
 15: #include <slepcblaslapack.h>
 16: #include <slepccublas.h>

 18: /*
 19:    BV_CleanCoefficients_CUDA - Sets to zero all entries of column j of the bv buffer
 20: */
 21: PetscErrorCode BV_CleanCoefficients_CUDA(BV bv,PetscInt j,PetscScalar *h)
 22: {
 24:   PetscScalar    *d_hh,*d_a;
 25:   PetscInt       i;
 26:   cudaError_t    cerr;

 29:   if (!h) {
 30:     VecCUDAGetArray(bv->buffer,&d_a);
 31:     d_hh = d_a + j*(bv->nc+bv->m);
 32:     cerr = cudaMemset(d_hh,0,(bv->nc+j)*sizeof(PetscScalar));CHKERRCUDA(cerr);
 33:     cerr = WaitForCUDA();CHKERRCUDA(cerr);
 34:     VecCUDARestoreArray(bv->buffer,&d_a);
 35:   } else { /* cpu memory */
 36:     for (i=0;i<bv->nc+j;i++) h[i] = 0.0;
 37:   }
 38:   return(0);
 39: }

 41: /*
 42:    BV_AddCoefficients_CUDA - Add the contents of the scratch (0-th column) of the bv buffer
 43:    into column j of the bv buffer
 44:  */
 45: PetscErrorCode BV_AddCoefficients_CUDA(BV bv,PetscInt j,PetscScalar *h,PetscScalar *c)
 46: {
 48:   PetscScalar    *d_h,*d_c,sone=1.0;
 49:   PetscInt       i;
 50:   PetscBLASInt   idx=0,one=1;
 51:   cublasStatus_t cberr;
 52:   cublasHandle_t cublasv2handle;

 55:   if (!h) {
 56:     PetscCUBLASGetHandle(&cublasv2handle);
 57:     VecCUDAGetArray(bv->buffer,&d_c);
 58:     d_h = d_c + j*(bv->nc+bv->m);
 59:     PetscBLASIntCast(bv->nc+j,&idx);
 60:     PetscLogGpuTimeBegin();
 61:     cberr = cublasXaxpy(cublasv2handle,idx,&sone,d_c,one,d_h,one);CHKERRCUBLAS(cberr);
 62:     PetscLogGpuTimeEnd();
 63:     PetscLogGpuFlops(1.0*bv->nc+j);
 64:     VecCUDARestoreArray(bv->buffer,&d_c);
 65:   } else { /* cpu memory */
 66:     for (i=0;i<bv->nc+j;i++) h[i] += c[i];
 67:     PetscLogFlops(1.0*bv->nc+j);
 68:   }
 69:   return(0);
 70: }

 72: /*
 73:    BV_SetValue_CUDA - Sets value in row j (counted after the constraints) of column k
 74:    of the coefficients array
 75: */
 76: PetscErrorCode BV_SetValue_CUDA(BV bv,PetscInt j,PetscInt k,PetscScalar *h,PetscScalar value)
 77: {
 79:   PetscScalar    *d_h,*a;
 80:   cudaError_t    cerr;

 83:   if (!h) {
 84:     VecCUDAGetArray(bv->buffer,&a);
 85:     d_h = a + k*(bv->nc+bv->m) + bv->nc+j;
 86:     cerr = cudaMemcpy(d_h,&value,sizeof(PetscScalar),cudaMemcpyHostToDevice);CHKERRCUDA(cerr);
 87:     PetscLogCpuToGpu(sizeof(PetscScalar));
 88:     cerr = WaitForCUDA();CHKERRCUDA(cerr);
 89:     VecCUDARestoreArray(bv->buffer,&a);
 90:   } else { /* cpu memory */
 91:     h[bv->nc+j] = value;
 92:   }
 93:   return(0);
 94: }

 96: /*
 97:    BV_SquareSum_CUDA - Returns the value h'*h, where h represents the contents of the
 98:    coefficients array (up to position j)
 99: */
100: PetscErrorCode BV_SquareSum_CUDA(BV bv,PetscInt j,PetscScalar *h,PetscReal *sum)
101: {
102:   PetscErrorCode    ierr;
103:   const PetscScalar *d_h;
104:   PetscScalar       dot;
105:   PetscInt          i;
106:   PetscBLASInt      idx=0,one=1;
107:   cublasStatus_t    cberr;
108:   cublasHandle_t    cublasv2handle;

111:   if (!h) {
112:     PetscCUBLASGetHandle(&cublasv2handle);
113:     VecCUDAGetArrayRead(bv->buffer,&d_h);
114:     PetscBLASIntCast(bv->nc+j,&idx);
115:     PetscLogGpuTimeBegin();
116:     cberr = cublasXdotc(cublasv2handle,idx,d_h,one,d_h,one,&dot);CHKERRCUBLAS(cberr);
117:     PetscLogGpuTimeEnd();
118:     PetscLogGpuFlops(2.0*bv->nc+j);
119:     *sum = PetscRealPart(dot);
120:     VecCUDARestoreArrayRead(bv->buffer,&d_h);
121:   } else { /* cpu memory */
122:     *sum = 0.0;
123:     for (i=0;i<bv->nc+j;i++) *sum += PetscRealPart(h[i]*PetscConj(h[i]));
124:     PetscLogFlops(2.0*bv->nc+j);
125:   }
126:   return(0);
127: }

129: #define X_AXIS        0
130: #define BLOCK_SIZE_X 64
131: #define TILE_SIZE_X  16 /* work to be done by any thread on axis x */

133: /*
134:    Set the kernels grid dimensions
135:    xcount: number of kernel calls needed for the requested size
136:  */
137: PetscErrorCode SetGrid1D(PetscInt n, dim3 *dimGrid, dim3 *dimBlock,PetscInt *xcount)
138: {
139:   PetscInt              one=1;
140:   PetscBLASInt          card;
141:   struct cudaDeviceProp devprop;
142:   cudaError_t           cerr;

145:   *xcount = 1;
146:   if (n>BLOCK_SIZE_X) {
147:     dimBlock->x = BLOCK_SIZE_X;
148:     dimGrid->x = (n+BLOCK_SIZE_X*TILE_SIZE_X-one)/BLOCK_SIZE_X*TILE_SIZE_X;
149:   } else {
150:     dimBlock->x = (n+TILE_SIZE_X-one)/TILE_SIZE_X;
151:     dimGrid->x = one;
152:   }
153:   cerr = cudaGetDevice(&card);CHKERRCUDA(cerr);
154:   cerr = cudaGetDeviceProperties(&devprop,card);CHKERRCUDA(cerr);
155:   if (dimGrid->x>(unsigned)devprop.maxGridSize[X_AXIS]) {
156:     *xcount = (dimGrid->x+devprop.maxGridSize[X_AXIS]-one)/devprop.maxGridSize[X_AXIS];
157:     dimGrid->x = devprop.maxGridSize[X_AXIS];
158:   }
159:   return(0);
160: }

162: /* pointwise multiplication */
163: __global__ void PointwiseMult_kernel(PetscInt xcount,PetscScalar *a,const PetscScalar *b,PetscInt n)
164: {
165:   PetscInt i,x;

167:   x = xcount*gridDim.x*blockDim.x+blockIdx.x*blockDim.x*TILE_SIZE_X+threadIdx.x*TILE_SIZE_X;
168:   for (i=x;i<x+TILE_SIZE_X&&i<n;i++) {
169:     a[i] *= PetscRealPart(b[i]);
170:   }
171: }

173: /* pointwise division */
174: __global__ void PointwiseDiv_kernel(PetscInt xcount,PetscScalar *a,const PetscScalar *b,PetscInt n)
175: {
176:   PetscInt i,x;

178:   x = xcount*gridDim.x*blockDim.x+blockIdx.x*blockDim.x*TILE_SIZE_X+threadIdx.x*TILE_SIZE_X;
179:   for (i=x;i<x+TILE_SIZE_X&&i<n;i++) {
180:     a[i] /= PetscRealPart(b[i]);
181:   }
182: }

184: /*
185:    BV_ApplySignature_CUDA - Computes the pointwise product h*omega, where h represents
186:    the contents of the coefficients array (up to position j) and omega is the signature;
187:    if inverse=TRUE then the operation is h/omega
188: */
189: PetscErrorCode BV_ApplySignature_CUDA(BV bv,PetscInt j,PetscScalar *h,PetscBool inverse)
190: {
191:   PetscErrorCode    ierr;
192:   PetscScalar       *d_h;
193:   const PetscScalar *d_omega,*omega;
194:   PetscInt          i,xcount;
195:   dim3              blocks3d, threads3d;
196:   cudaError_t       cerr;

199:   if (!(bv->nc+j)) return(0);
200:   if (!h) {
201:     VecCUDAGetArray(bv->buffer,&d_h);
202:     VecCUDAGetArrayRead(bv->omega,&d_omega);
203:     SetGrid1D(bv->nc+j,&blocks3d,&threads3d,&xcount);
204:     PetscLogGpuTimeBegin();
205:     if (inverse) {
206:       for (i=0;i<xcount;i++) {
207:         PointwiseDiv_kernel<<<blocks3d,threads3d>>>(i,d_h,d_omega,bv->nc+j);
208:       }
209:     } else {
210:       for (i=0;i<xcount;i++) {
211:         PointwiseMult_kernel<<<blocks3d,threads3d>>>(i,d_h,d_omega,bv->nc+j);
212:       }
213:     }
214:     cerr = cudaGetLastError();CHKERRCUDA(cerr);
215:     PetscLogGpuTimeEnd();
216:     PetscLogGpuFlops(1.0*bv->nc+j);
217:     cerr = WaitForCUDA();CHKERRCUDA(cerr);
218:     VecCUDARestoreArrayRead(bv->omega,&d_omega);
219:     VecCUDARestoreArray(bv->buffer,&d_h);
220:   } else {
221:     VecGetArrayRead(bv->omega,&omega);
222:     if (inverse) for (i=0;i<bv->nc+j;i++) h[i] /= PetscRealPart(omega[i]);
223:     else for (i=0;i<bv->nc+j;i++) h[i] *= PetscRealPart(omega[i]);
224:     VecRestoreArrayRead(bv->omega,&omega);
225:     PetscLogFlops(1.0*bv->nc+j);
226:   }
227:   return(0);
228: }

230: /*
231:    BV_SquareRoot_CUDA - Returns the square root of position j (counted after the constraints)
232:    of the coefficients array
233: */
234: PetscErrorCode BV_SquareRoot_CUDA(BV bv,PetscInt j,PetscScalar *h,PetscReal *beta)
235: {
236:   PetscErrorCode    ierr;
237:   const PetscScalar *d_h;
238:   PetscScalar       hh;
239:   cudaError_t       cerr;

242:   if (!h) {
243:     VecCUDAGetArrayRead(bv->buffer,&d_h);
244:     cerr = cudaMemcpy(&hh,d_h+bv->nc+j,sizeof(PetscScalar),cudaMemcpyDeviceToHost);CHKERRCUDA(cerr);
245:     PetscLogGpuToCpu(sizeof(PetscScalar));
246:     cerr = WaitForCUDA();CHKERRCUDA(cerr);
247:     BV_SafeSqrt(bv,hh,beta);
248:     VecCUDARestoreArrayRead(bv->buffer,&d_h);
249:   } else {
250:     BV_SafeSqrt(bv,h[bv->nc+j],beta);
251:   }
252:   return(0);
253: }

255: /*
256:    BV_StoreCoefficients_CUDA - Copy the contents of the coefficients array to an array dest
257:    provided by the caller (only values from l to j are copied)
258: */
259: PetscErrorCode BV_StoreCoefficients_CUDA(BV bv,PetscInt j,PetscScalar *h,PetscScalar *dest)
260: {
261:   PetscErrorCode    ierr;
262:   const PetscScalar *d_h,*d_a;
263:   PetscInt          i;
264:   cudaError_t       cerr;

267:   if (!h) {
268:     VecCUDAGetArrayRead(bv->buffer,&d_a);
269:     d_h = d_a + j*(bv->nc+bv->m)+bv->nc;
270:     cerr = cudaMemcpy(dest-bv->l,d_h,(j-bv->l)*sizeof(PetscScalar),cudaMemcpyDeviceToHost);CHKERRCUDA(cerr);
271:     PetscLogGpuToCpu((j-bv->l)*sizeof(PetscScalar));
272:     cerr = WaitForCUDA();CHKERRCUDA(cerr);
273:     VecCUDARestoreArrayRead(bv->buffer,&d_a);
274:   } else {
275:     for (i=bv->l;i<j;i++) dest[i-bv->l] = h[bv->nc+i];
276:   }
277:   return(0);
278: }