Actual source code: mpicuda.cu

petsc-master 2019-08-18
Report Typos and Errors

  2: /*
  3:    This file contains routines for Parallel vector operations.
  4:  */
  5: #define PETSC_SKIP_SPINLOCK
  6: #define PETSC_SKIP_CXX_COMPLEX_FIX

  8: #include <petscconf.h>
  9:  #include <../src/vec/vec/impls/mpi/pvecimpl.h>
 10:  #include <../src/vec/vec/impls/seq/seqcuda/cudavecimpl.h>

 12: PetscErrorCode VecDestroy_MPICUDA(Vec v)
 13: {
 15:   cudaError_t    err;

 18:   if (v->spptr) {
 19:     if (((Vec_CUDA*)v->spptr)->GPUarray_allocated) {
 20:       err = cudaFree(((Vec_CUDA*)v->spptr)->GPUarray_allocated);CHKERRCUDA(err);
 21:       ((Vec_CUDA*)v->spptr)->GPUarray_allocated = NULL;
 22:     }
 23:     err = cudaStreamDestroy(((Vec_CUDA*)v->spptr)->stream);CHKERRCUDA(err);
 24:     PetscFree(v->spptr);
 25:   }
 26:   VecDestroy_MPI(v);
 27:   return(0);
 28: }

 30: PetscErrorCode VecNorm_MPICUDA(Vec xin,NormType type,PetscReal *z)
 31: {
 32:   PetscReal      sum,work = 0.0;

 36:   if (type == NORM_2 || type == NORM_FROBENIUS) {
 37:     VecNorm_SeqCUDA(xin,NORM_2,&work);
 38:     work *= work;
 39:     MPIU_Allreduce(&work,&sum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 40:     *z    = PetscSqrtReal(sum);
 41:   } else if (type == NORM_1) {
 42:     /* Find the local part */
 43:     VecNorm_SeqCUDA(xin,NORM_1,&work);
 44:     /* Find the global max */
 45:     MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 46:   } else if (type == NORM_INFINITY) {
 47:     /* Find the local max */
 48:     VecNorm_SeqCUDA(xin,NORM_INFINITY,&work);
 49:     /* Find the global max */
 50:     MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)xin));
 51:   } else if (type == NORM_1_AND_2) {
 52:     PetscReal temp[2];
 53:     VecNorm_SeqCUDA(xin,NORM_1,temp);
 54:     VecNorm_SeqCUDA(xin,NORM_2,temp+1);
 55:     temp[1] = temp[1]*temp[1];
 56:     MPIU_Allreduce(temp,z,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 57:     z[1] = PetscSqrtReal(z[1]);
 58:   }
 59:   return(0);
 60: }

 62: PetscErrorCode VecDot_MPICUDA(Vec xin,Vec yin,PetscScalar *z)
 63: {
 64:   PetscScalar    sum,work;

 68:   VecDot_SeqCUDA(xin,yin,&work);
 69:   MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 70:   *z   = sum;
 71:   return(0);
 72: }

 74: PetscErrorCode VecTDot_MPICUDA(Vec xin,Vec yin,PetscScalar *z)
 75: {
 76:   PetscScalar    sum,work;

 80:   VecTDot_SeqCUDA(xin,yin,&work);
 81:   MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 82:   *z   = sum;
 83:   return(0);
 84: }

 86: PetscErrorCode VecMDot_MPICUDA(Vec xin,PetscInt nv,const Vec y[],PetscScalar *z)
 87: {
 88:   PetscScalar    awork[128],*work = awork;

 92:   if (nv > 128) {
 93:     PetscMalloc1(nv,&work);
 94:   }
 95:   VecMDot_SeqCUDA(xin,nv,y,work);
 96:   MPIU_Allreduce(work,z,nv,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 97:   if (nv > 128) {
 98:     PetscFree(work);
 99:   }
100:   return(0);
101: }

103: /*MC
104:    VECMPICUDA - VECMPICUDA = "mpicuda" - The basic parallel vector, modified to use CUDA

106:    Options Database Keys:
107: . -vec_type mpicuda - sets the vector type to VECMPICUDA during a call to VecSetFromOptions()

109:   Level: beginner

111: .seealso: VecCreate(), VecSetType(), VecSetFromOptions(), VecCreateMPIWithArray(), VECMPI, VecType, VecCreateMPI()
112: M*/


115: PetscErrorCode VecDuplicate_MPICUDA(Vec win,Vec *v)
116: {
118:   Vec_MPI        *vw,*w = (Vec_MPI*)win->data;
119:   PetscScalar    *array;

122:   VecCreate(PetscObjectComm((PetscObject)win),v);
123:   PetscLayoutReference(win->map,&(*v)->map);

125:   VecCreate_MPICUDA_Private(*v,PETSC_TRUE,w->nghost,0);
126:   vw   = (Vec_MPI*)(*v)->data;
127:   PetscMemcpy((*v)->ops,win->ops,sizeof(struct _VecOps));

129:   /* save local representation of the parallel vector (and scatter) if it exists */
130:   if (w->localrep) {
131:     VecGetArray(*v,&array);
132:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,win->map->n+w->nghost,array,&vw->localrep);
133:     PetscMemcpy(vw->localrep->ops,w->localrep->ops,sizeof(struct _VecOps));
134:     VecRestoreArray(*v,&array);
135:     PetscLogObjectParent((PetscObject)*v,(PetscObject)vw->localrep);
136:     vw->localupdate = w->localupdate;
137:     if (vw->localupdate) {
138:       PetscObjectReference((PetscObject)vw->localupdate);
139:     }
140:   }

142:   /* New vector should inherit stashing property of parent */
143:   (*v)->stash.donotstash   = win->stash.donotstash;
144:   (*v)->stash.ignorenegidx = win->stash.ignorenegidx;

146:   /* change type_name appropriately */
147:   VecCUDAAllocateCheck(*v);
148:   PetscObjectChangeTypeName((PetscObject)(*v),VECMPICUDA);

150:   PetscObjectListDuplicate(((PetscObject)win)->olist,&((PetscObject)(*v))->olist);
151:   PetscFunctionListDuplicate(((PetscObject)win)->qlist,&((PetscObject)(*v))->qlist);
152:   (*v)->map->bs   = PetscAbs(win->map->bs);
153:   (*v)->bstash.bs = win->bstash.bs;
154:   return(0);
155: }

157: PetscErrorCode VecDotNorm2_MPICUDA(Vec s,Vec t,PetscScalar *dp,PetscScalar *nm)
158: {
160:   PetscScalar    work[2],sum[2];

163:   VecDotNorm2_SeqCUDA(s,t,work,work+1);
164:   MPIU_Allreduce(&work,&sum,2,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)s));
165:   *dp  = sum[0];
166:   *nm  = sum[1];
167:   return(0);
168: }

170: PetscErrorCode VecCreate_MPICUDA(Vec vv)
171: {

175:   PetscLayoutSetUp(vv->map);
176:   VecCUDAAllocateCheck(vv);CHKERRCUDA(ierr);
177:   vv->valid_GPU_array      = PETSC_OFFLOAD_GPU;
178:   VecCreate_MPICUDA_Private(vv,PETSC_FALSE,0,((Vec_CUDA*)vv->spptr)->GPUarray_allocated);
179:   VecSet(vv,0.0);
180:   return(0);
181: }

183: PetscErrorCode VecCreate_CUDA(Vec v)
184: {
186:   PetscMPIInt    size;

189:   MPI_Comm_size(PetscObjectComm((PetscObject)v),&size);
190:   if (size == 1) {
191:     VecSetType(v,VECSEQCUDA);
192:   } else {
193:     VecSetType(v,VECMPICUDA);
194:   }
195:   return(0);
196: }

198: /*@C
199:    VecCreateMPICUDAWithArray - Creates a parallel, array-style vector,
200:    where the user provides the GPU array space to store the vector values.

202:    Collective

204:    Input Parameters:
205: +  comm  - the MPI communicator to use
206: .  bs    - block size, same meaning as VecSetBlockSize()
207: .  n     - local vector length, cannot be PETSC_DECIDE
208: .  N     - global vector length (or PETSC_DECIDE to have calculated)
209: -  array - the user provided GPU array to store the vector values

211:    Output Parameter:
212: .  vv - the vector

214:    Notes:
215:    Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
216:    same type as an existing vector.

218:    If the user-provided array is NULL, then VecCUDAPlaceArray() can be used
219:    at a later stage to SET the array for storing the vector values.

221:    PETSc does NOT free the array when the vector is destroyed via VecDestroy().
222:    The user should not free the array until the vector is destroyed.

224:    Level: intermediate

226: .seealso: VecCreateSeqCUDAWithArray(), VecCreateMPIWithArray(), VecCreateSeqWithArray(),
227:           VecCreate(), VecDuplicate(), VecDuplicateVecs(), VecCreateGhost(),
228:           VecCreateMPI(), VecCreateGhostWithArray(), VecPlaceArray()

230: @*/
231: PetscErrorCode  VecCreateMPICUDAWithArray(MPI_Comm comm,PetscInt bs,PetscInt n,PetscInt N,const PetscScalar array[],Vec *vv)
232: {

236:   if (n == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local size of vector");
237:   PetscSplitOwnership(comm,&n,&N);
238:   VecCreate(comm,vv);
239:   VecSetSizes(*vv,n,N);
240:   VecSetBlockSize(*vv,bs);
241:   VecCreate_MPICUDA_Private(*vv,PETSC_FALSE,0,array);
242:   return(0);
243: }

245: extern "C" PetscErrorCode VecGetArrayWrite_SeqCUDA(Vec,PetscScalar**);

247: PetscErrorCode VecPinToCPU_MPICUDA(Vec V,PetscBool pin)
248: {

252:   V->pinnedtocpu = pin;
253:   if (pin) {
254:     VecCUDACopyFromGPU(V);
255:     V->valid_GPU_array = PETSC_OFFLOAD_CPU; /* since the CPU code will likely change values in the vector */
256:     V->ops->dotnorm2               = NULL;
257:     V->ops->waxpy                  = VecWAXPY_Seq;
258:     V->ops->dot                    = VecDot_MPI;
259:     V->ops->mdot                   = VecMDot_MPI;
260:     V->ops->tdot                   = VecTDot_MPI;
261:     V->ops->norm                   = VecNorm_MPI;
262:     V->ops->scale                  = VecScale_Seq;
263:     V->ops->copy                   = VecCopy_Seq;
264:     V->ops->set                    = VecSet_Seq;
265:     V->ops->swap                   = VecSwap_Seq;
266:     V->ops->axpy                   = VecAXPY_Seq;
267:     V->ops->axpby                  = VecAXPBY_Seq;
268:     V->ops->maxpy                  = VecMAXPY_Seq;
269:     V->ops->aypx                   = VecAYPX_Seq;
270:     V->ops->axpbypcz               = VecAXPBYPCZ_Seq;
271:     V->ops->pointwisemult          = VecPointwiseMult_Seq;
272:     V->ops->setrandom              = VecSetRandom_Seq;
273:     V->ops->placearray             = VecPlaceArray_Seq;
274:     V->ops->replacearray           = VecReplaceArray_Seq;
275:     V->ops->resetarray             = VecResetArray_Seq;
276:     V->ops->dot_local              = VecDot_Seq;
277:     V->ops->tdot_local             = VecTDot_Seq;
278:     V->ops->norm_local             = VecNorm_Seq;
279:     V->ops->mdot_local             = VecMDot_Seq;
280:     V->ops->pointwisedivide        = VecPointwiseDivide_Seq;
281:     V->ops->getlocalvector         = NULL;
282:     V->ops->restorelocalvector     = NULL;
283:     V->ops->getlocalvectorread     = NULL;
284:     V->ops->restorelocalvectorread = NULL;
285:     V->ops->getarraywrite          = NULL;
286:   } else {
287:     V->ops->dotnorm2               = VecDotNorm2_MPICUDA;
288:     V->ops->waxpy                  = VecWAXPY_SeqCUDA;
289:     V->ops->duplicate              = VecDuplicate_MPICUDA;
290:     V->ops->dot                    = VecDot_MPICUDA;
291:     V->ops->mdot                   = VecMDot_MPICUDA;
292:     V->ops->tdot                   = VecTDot_MPICUDA;
293:     V->ops->norm                   = VecNorm_MPICUDA;
294:     V->ops->scale                  = VecScale_SeqCUDA;
295:     V->ops->copy                   = VecCopy_SeqCUDA;
296:     V->ops->set                    = VecSet_SeqCUDA;
297:     V->ops->swap                   = VecSwap_SeqCUDA;
298:     V->ops->axpy                   = VecAXPY_SeqCUDA;
299:     V->ops->axpby                  = VecAXPBY_SeqCUDA;
300:     V->ops->maxpy                  = VecMAXPY_SeqCUDA;
301:     V->ops->aypx                   = VecAYPX_SeqCUDA;
302:     V->ops->axpbypcz               = VecAXPBYPCZ_SeqCUDA;
303:     V->ops->pointwisemult          = VecPointwiseMult_SeqCUDA;
304:     V->ops->setrandom              = VecSetRandom_SeqCUDA;
305:     V->ops->placearray             = VecPlaceArray_SeqCUDA;
306:     V->ops->replacearray           = VecReplaceArray_SeqCUDA;
307:     V->ops->resetarray             = VecResetArray_SeqCUDA;
308:     V->ops->dot_local              = VecDot_SeqCUDA;
309:     V->ops->tdot_local             = VecTDot_SeqCUDA;
310:     V->ops->norm_local             = VecNorm_SeqCUDA;
311:     V->ops->mdot_local             = VecMDot_SeqCUDA;
312:     V->ops->destroy                = VecDestroy_MPICUDA;
313:     V->ops->pointwisedivide        = VecPointwiseDivide_SeqCUDA;
314:     V->ops->getlocalvector         = VecGetLocalVector_SeqCUDA;
315:     V->ops->restorelocalvector     = VecRestoreLocalVector_SeqCUDA;
316:     V->ops->getlocalvectorread     = VecGetLocalVector_SeqCUDA;
317:     V->ops->restorelocalvectorread = VecRestoreLocalVector_SeqCUDA;
318:     V->ops->getarraywrite          = VecGetArrayWrite_SeqCUDA;
319:   }
320:   return(0);
321: }

323: PetscErrorCode VecCreate_MPICUDA_Private(Vec vv,PetscBool alloc,PetscInt nghost,const PetscScalar array[])
324: {
326:   cudaError_t    err;
327:   Vec_CUDA       *veccuda;

330:   VecCreate_MPI_Private(vv,PETSC_FALSE,0,0);
331:   PetscObjectChangeTypeName((PetscObject)vv,VECMPICUDA);

333:   VecPinToCPU_MPICUDA(vv,PETSC_FALSE);
334:   vv->ops->pintocpu = VecPinToCPU_MPICUDA;

336:   /* Later, functions check for the Vec_CUDA structure existence, so do not create it without array */
337:   if (alloc && !array) {
338:     VecCUDAAllocateCheck(vv);
339:     VecSet(vv,0.0);
340:   }
341:   if (array) {
342:     if (!vv->spptr) {
343:       /* Cannot use PetscNew() here because spptr is void* */
344:       PetscMalloc(sizeof(Vec_CUDA),&vv->spptr);
345:       veccuda = (Vec_CUDA*)vv->spptr;
346:       err = cudaStreamCreate(&veccuda->stream);CHKERRCUDA(err);
347:       veccuda->GPUarray_allocated = 0;
348:       veccuda->hostDataRegisteredAsPageLocked = PETSC_FALSE;
349:       vv->valid_GPU_array = PETSC_OFFLOAD_UNALLOCATED;
350:     }
351:     veccuda = (Vec_CUDA*)vv->spptr;
352:     veccuda->GPUarray = (PetscScalar*)array;
353:   }

355:   return(0);
356: }