Actual source code: vecscattercuda.cu

petsc-3.7.1 2016-05-15
Report Typos and Errors
  1: /*
  2:    Implements the various scatter operations on cuda vectors
  3: */

  5: #define PETSC_SKIP_SPINLOCK

  7: #include <petscconf.h>
  8: #include <petsc/private/vecimpl.h>          /*I "petscvec.h" I*/
  9: #include <../src/vec/vec/impls/dvecimpl.h>
 10: #include <../src/vec/vec/impls/seq/seqcuda/cudavecimpl.h>

 12: #include <cuda_runtime.h>

 16: PetscErrorCode VecScatterCUDAIndicesCreate_StoS(PetscInt n,PetscInt toFirst,PetscInt fromFirst,PetscInt toStep, PetscInt fromStep,PetscInt *tslots, PetscInt *fslots,PetscCUDAIndices *ci) {

 18:   PetscCUDAIndices           cci;
 19:   VecScatterCUDAIndices_StoS stos_scatter;
 20:   cudaError_t                err;
 21:   cudaStream_t               stream;
 22:   PetscInt                   *intVecGPU;
 23:   int                        device;
 24:   cudaDeviceProp             props;

 27:   cci = new struct _p_PetscCUDAIndices;
 28:   stos_scatter = new struct _p_VecScatterCUDAIndices_StoS;

 30:   /* create the "from" indices */
 31:   stos_scatter->fslots = 0;
 32:   stos_scatter->fromFirst = 0;
 33:   stos_scatter->fromStep = 0;
 34:   if (n) {
 35:     if (fslots) {
 36:       /* allocate GPU memory for the to-slots */
 37:       err = cudaMalloc((void **)&intVecGPU,n*sizeof(PetscInt));CHKERRCUDA(err);
 38:       err = cudaMemcpy(intVecGPU,fslots,n*sizeof(PetscInt),cudaMemcpyHostToDevice);CHKERRCUDA(err);

 40:       /* assign the pointer to the struct */
 41:       stos_scatter->fslots = intVecGPU;
 42:       stos_scatter->fromMode = VEC_SCATTER_CUDA_GENERAL;
 43:     } else if (fromStep) {
 44:       stos_scatter->fromFirst = fromFirst;
 45:       stos_scatter->fromStep = fromStep;
 46:       stos_scatter->fromMode = VEC_SCATTER_CUDA_STRIDED;
 47:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide fslots or fromStep.");
 48:   }

 50:   /* create the "to" indices */
 51:   stos_scatter->tslots = 0;
 52:   stos_scatter->toFirst = 0;
 53:   stos_scatter->toStep = 0;
 54:   if (n) {
 55:     if (tslots) {
 56:       /* allocate GPU memory for the to-slots */
 57:       err = cudaMalloc((void **)&intVecGPU,n*sizeof(PetscInt));CHKERRCUDA(err);
 58:       err = cudaMemcpy(intVecGPU,tslots,n*sizeof(PetscInt),cudaMemcpyHostToDevice);CHKERRCUDA(err);

 60:       /* assign the pointer to the struct */
 61:       stos_scatter->tslots = intVecGPU;
 62:       stos_scatter->toMode = VEC_SCATTER_CUDA_GENERAL;
 63:     } else if (toStep) {
 64:       stos_scatter->toFirst = toFirst;
 65:       stos_scatter->toStep = toStep;
 66:       stos_scatter->toMode = VEC_SCATTER_CUDA_STRIDED;
 67:     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must provide tslots or toStep.");
 68:   }

 70:   /* allocate the stream variable */
 71:   err = cudaStreamCreate(&stream);CHKERRCUDA(err);
 72:   stos_scatter->stream = stream;

 74:   /* the number of indices */
 75:   stos_scatter->n = n;

 77:   /* get the maximum number of coresident thread blocks */
 78:   cudaGetDevice(&device);
 79:   cudaGetDeviceProperties(&props, device);
 80:   stos_scatter->MAX_CORESIDENT_THREADS = props.maxThreadsPerMultiProcessor;
 81:   if (props.major>=3) {
 82:     stos_scatter->MAX_BLOCKS = 16*props.multiProcessorCount;
 83:   } else {
 84:     stos_scatter->MAX_BLOCKS = 8*props.multiProcessorCount;
 85:   }

 87:   /* assign the indices */
 88:   cci->scatter = (VecScatterCUDAIndices_StoS)stos_scatter;
 89:   cci->scatterType = VEC_SCATTER_CUDA_STOS;
 90:   *ci = cci;
 91:   return(0);
 92: }

 96: PetscErrorCode VecScatterCUDAIndicesCreate_PtoP(PetscInt ns,PetscInt *sendIndices,PetscInt nr,PetscInt *recvIndices,PetscCUDAIndices *ci)
 97: {
 98:   PetscCUDAIndices           cci;
 99:   VecScatterCUDAIndices_PtoP ptop_scatter;

102:   cci = new struct _p_PetscCUDAIndices;
103:   ptop_scatter = new struct _p_VecScatterCUDAIndices_PtoP;

105:   /* this calculation assumes that the input indices are sorted */
106:   ptop_scatter->ns = sendIndices[ns-1]-sendIndices[0]+1;
107:   ptop_scatter->sendLowestIndex = sendIndices[0];
108:   ptop_scatter->nr = recvIndices[nr-1]-recvIndices[0]+1;
109:   ptop_scatter->recvLowestIndex = recvIndices[0];

111:   /* assign indices */
112:   cci->scatter = (VecScatterCUDAIndices_PtoP)ptop_scatter;
113:   cci->scatterType = VEC_SCATTER_CUDA_PTOP;

115:   *ci = cci;
116:   return(0);
117: }

121: PetscErrorCode VecScatterCUDAIndicesDestroy(PetscCUDAIndices *ci)
122: {
124:   if (*ci) {
125:     if ((*ci)->scatterType == VEC_SCATTER_CUDA_PTOP) {
126:       delete (VecScatterCUDAIndices_PtoP)(*ci)->scatter;
127:       (*ci)->scatter = 0;
128:     } else {
129:       cudaError_t err;
130:       VecScatterCUDAIndices_StoS stos_scatter = (VecScatterCUDAIndices_StoS)(*ci)->scatter;
131:       if (stos_scatter->fslots) {
132:         err = cudaFree(stos_scatter->fslots);CHKERRCUDA(err);
133:         stos_scatter->fslots = 0;
134:       }

136:       /* free the GPU memory for the to-slots */
137:       if (stos_scatter->tslots) {
138:         err = cudaFree(stos_scatter->tslots);CHKERRCUDA(err);
139:         stos_scatter->tslots = 0;
140:       }

142:       /* free the stream variable */
143:       if (stos_scatter->stream) {
144:         err = cudaStreamDestroy(stos_scatter->stream);CHKERRCUDA(err);
145:         stos_scatter->stream = 0;
146:       }
147:       delete stos_scatter;
148:       (*ci)->scatter = 0;
149:     }
150:     delete *ci;
151:     *ci = 0;
152:   }
153:   return(0);
154: }

156: /* Insert operator */
157: class Insert {
158:   public:
159:     __device__ PetscScalar operator() (PetscScalar a,PetscScalar b) const {
160:       return a;
161:     }
162: };

164: /* Add operator */
165: class Add {
166:   public:
167:     __device__ PetscScalar operator() (PetscScalar a,PetscScalar b) const {
168:       return a+b;
169:     }
170: };

172: /* Add operator */
173: class Max {
174:   public:
175:     __device__ PetscScalar operator() (PetscScalar a,PetscScalar b) const {
176:       return PetscMax(PetscRealPart(a),PetscRealPart(b));
177:     }
178: };

180: /* Sequential general to sequential general GPU kernel */
181: template<class OPERATOR>
182: __global__ void VecScatterCUDA_SGtoSG_kernel(PetscInt n,PetscInt *xind,PetscScalar *x,PetscInt *yind,PetscScalar *y,OPERATOR OP) {
183:   const int tidx = blockIdx.x*blockDim.x + threadIdx.x;
184:   const int grid_size = gridDim.x * blockDim.x;
185:   for (int i = tidx; i < n; i += grid_size) {
186:     y[yind[i]] = OP(x[xind[i]],y[yind[i]]);
187:   }
188: }

190: /* Sequential general to sequential strided GPU kernel */
191: template<class OPERATOR>
192: __global__ void VecScatterCUDA_SGtoSS_kernel(PetscInt n,PetscInt *xind,PetscScalar *x,PetscInt toFirst,PetscInt toStep,PetscScalar *y,OPERATOR OP) {
193:   const int tidx = blockIdx.x*blockDim.x + threadIdx.x;
194:   const int grid_size = gridDim.x * blockDim.x;
195:   for (int i = tidx; i < n; i += grid_size) {
196:     y[toFirst+i*toStep] = OP(x[xind[i]],y[toFirst+i*toStep]);
197:   }
198: }

200: /* Sequential strided to sequential strided GPU kernel */
201: template<class OPERATOR>
202: __global__ void VecScatterCUDA_SStoSS_kernel(PetscInt n,PetscInt fromFirst,PetscInt fromStep,PetscScalar *x,PetscInt toFirst,PetscInt toStep,PetscScalar *y,OPERATOR OP) {
203:   const int tidx = blockIdx.x*blockDim.x + threadIdx.x;
204:   const int grid_size = gridDim.x * blockDim.x;
205:   for (int i = tidx; i < n; i += grid_size) {
206:     y[toFirst+i*toStep] = OP(x[fromFirst+i*fromStep],y[toFirst+i*toStep]);
207:   }
208: }

210: /* Sequential strided to sequential general GPU kernel */
211: template<class OPERATOR>
212: __global__ void VecScatterCUDA_SStoSG_kernel(PetscInt n,PetscInt fromFirst,PetscInt fromStep,PetscScalar *x,PetscInt *yind,PetscScalar *y,OPERATOR OP) {
213:   const int tidx = blockIdx.x*blockDim.x + threadIdx.x;
214:   const int grid_size = gridDim.x * blockDim.x;
215:   for (int i = tidx; i < n; i += grid_size) {
216:     y[yind[i]] = OP(x[fromFirst+i*fromStep],y[yind[i]]);
217:   }
218: }

220: template<class OPERATOR>
221: void VecScatterCUDA_StoS_Dispatcher(PetscScalar *xarray,PetscScalar *yarray,PetscCUDAIndices ci,ScatterMode mode,OPERATOR OP) {

223:   PetscInt                   nBlocks=0,nThreads=128;
224:   VecScatterCUDAIndices_StoS stos_scatter = (VecScatterCUDAIndices_StoS)ci->scatter;

226:   nBlocks=(int)ceil(((float) stos_scatter->n)/((float) nThreads))+1;
227:   if (nBlocks>stos_scatter->MAX_CORESIDENT_THREADS/nThreads) {
228:     nBlocks = stos_scatter->MAX_CORESIDENT_THREADS/nThreads;
229:   }
230:   dim3 block(nThreads,1,1);
231:   dim3 grid(nBlocks,1,1);

233:   if (mode == SCATTER_FORWARD) {
234:     if (stos_scatter->fromMode == VEC_SCATTER_CUDA_GENERAL && stos_scatter->toMode == VEC_SCATTER_CUDA_GENERAL) {
235:       VecScatterCUDA_SGtoSG_kernel<<<grid,block,0,stos_scatter->stream>>>(stos_scatter->n,stos_scatter->fslots,xarray,stos_scatter->tslots,yarray,OP);
236:     } else if (stos_scatter->fromMode == VEC_SCATTER_CUDA_GENERAL && stos_scatter->toMode == VEC_SCATTER_CUDA_STRIDED) {
237:       VecScatterCUDA_SGtoSS_kernel<<<grid,block,0,stos_scatter->stream>>>(stos_scatter->n,stos_scatter->fslots,xarray,stos_scatter->toFirst,stos_scatter->toStep,yarray,OP);
238:     } else if (stos_scatter->fromMode == VEC_SCATTER_CUDA_STRIDED && stos_scatter->toMode == VEC_SCATTER_CUDA_STRIDED) {
239:       VecScatterCUDA_SStoSS_kernel<<<grid,block,0,stos_scatter->stream>>>(stos_scatter->n,stos_scatter->fromFirst,stos_scatter->fromStep,xarray,stos_scatter->toFirst,stos_scatter->toStep,yarray,OP);
240:     } else if (stos_scatter->fromMode == VEC_SCATTER_CUDA_STRIDED && stos_scatter->toMode == VEC_SCATTER_CUDA_GENERAL) {
241:       VecScatterCUDA_SStoSG_kernel<<<grid,block,0,stos_scatter->stream>>>(stos_scatter->n,stos_scatter->fromFirst,stos_scatter->fromStep,xarray,stos_scatter->tslots,yarray,OP);
242:     }
243:   } else {
244:     if (stos_scatter->toMode == VEC_SCATTER_CUDA_GENERAL && stos_scatter->fromMode == VEC_SCATTER_CUDA_GENERAL) {
245:       VecScatterCUDA_SGtoSG_kernel<<<grid,block,0,stos_scatter->stream>>>(stos_scatter->n,stos_scatter->tslots,xarray,stos_scatter->fslots,yarray,OP);
246:     } else if (stos_scatter->toMode == VEC_SCATTER_CUDA_GENERAL && stos_scatter->fromMode == VEC_SCATTER_CUDA_STRIDED) {
247:       VecScatterCUDA_SGtoSS_kernel<<<grid,block,0,stos_scatter->stream>>>(stos_scatter->n,stos_scatter->tslots,xarray,stos_scatter->fromFirst,stos_scatter->fromStep,yarray,OP);
248:     } else if (stos_scatter->toMode == VEC_SCATTER_CUDA_STRIDED && stos_scatter->fromMode == VEC_SCATTER_CUDA_STRIDED) {
249:       VecScatterCUDA_SStoSS_kernel<<<grid,block,0,stos_scatter->stream>>>(stos_scatter->n,stos_scatter->toFirst,stos_scatter->toStep,xarray,stos_scatter->fromFirst,stos_scatter->fromStep,yarray,OP);
250:     } else if (stos_scatter->toMode == VEC_SCATTER_CUDA_STRIDED && stos_scatter->fromMode == VEC_SCATTER_CUDA_GENERAL) {
251:       VecScatterCUDA_SStoSG_kernel<<<grid,block,0,stos_scatter->stream>>>(stos_scatter->n,stos_scatter->toFirst,stos_scatter->toStep,xarray,stos_scatter->fslots,yarray,OP);
252:     }
253:   }
254: }

258: PetscErrorCode VecScatterCUDA_StoS(Vec x,Vec y,PetscCUDAIndices ci,InsertMode addv,ScatterMode mode)
259: {
260:   PetscErrorCode             ierr;
261:   PetscScalar                *xarray,*yarray;
262:   VecScatterCUDAIndices_StoS stos_scatter = (VecScatterCUDAIndices_StoS)ci->scatter;
263:   cudaError_t                err;

266:   VecCUDAAllocateCheck(x);
267:   VecCUDAAllocateCheck(y);
268:   VecCUDAGetArrayRead(x,&xarray);
269:   VecCUDAGetArrayReadWrite(y,&yarray);
270:   if (stos_scatter->n) {
271:     if (addv == INSERT_VALUES)
272:       VecScatterCUDA_StoS_Dispatcher(xarray,yarray,ci,mode,Insert());
273:     else if (addv == ADD_VALUES)
274:       VecScatterCUDA_StoS_Dispatcher(xarray,yarray,ci,mode,Add());
275:     else if (addv == MAX_VALUES)
276:       VecScatterCUDA_StoS_Dispatcher(xarray,yarray,ci,mode,Max());
277:     else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Wrong insert option");
278:     err = cudaGetLastError();CHKERRCUDA(err);
279:     err = cudaStreamSynchronize(stos_scatter->stream);CHKERRCUDA(err);
280:   }
281:   VecCUDARestoreArrayRead(x,&xarray);
282:   VecCUDARestoreArrayReadWrite(y,&yarray);
283:   return(0);
284: }