Actual source code: mpicuda.cu


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

  7: #include <petscconf.h>
  8: #include <../src/vec/vec/impls/mpi/pvecimpl.h>
  9: #include <petsc/private/cudavecimpl.h>

 11: /*MC
 12:    VECCUDA - VECCUDA = "cuda" - A VECSEQCUDA on a single-process communicator, and VECMPICUDA otherwise.

 14:    Options Database Keys:
 15: . -vec_type cuda - sets the vector type to VECCUDA during a call to VecSetFromOptions()

 17:   Level: beginner

 19: .seealso: VecCreate(), VecSetType(), VecSetFromOptions(), VecCreateMPIWithArray(), VECSEQCUDA, VECMPICUDA, VECSTANDARD, VecType, VecCreateMPI(), VecSetPinnedMemoryMin()
 20: M*/

 22: PetscErrorCode VecDestroy_MPICUDA(Vec v)
 23: {
 24:   Vec_MPI        *vecmpi = (Vec_MPI*)v->data;
 25:   Vec_CUDA       *veccuda;
 27:   cudaError_t    err;

 30:   if (v->spptr) {
 31:     veccuda = (Vec_CUDA*)v->spptr;
 32:     if (veccuda->GPUarray_allocated) {
 33:       err = cudaFree(((Vec_CUDA*)v->spptr)->GPUarray_allocated);CHKERRCUDA(err);
 34:       veccuda->GPUarray_allocated = NULL;
 35:     }
 36:     if (veccuda->stream) {
 37:       err = cudaStreamDestroy(((Vec_CUDA*)v->spptr)->stream);CHKERRCUDA(err);
 38:     }
 39:     if (v->pinned_memory) {
 40:       PetscMallocSetCUDAHost();
 41:       PetscFree(vecmpi->array_allocated);
 42:       PetscMallocResetCUDAHost();
 43:       v->pinned_memory = PETSC_FALSE;
 44:     }
 45:     PetscFree(v->spptr);
 46:   }
 47:   VecDestroy_MPI(v);
 48:   return(0);
 49: }

 51: PetscErrorCode VecNorm_MPICUDA(Vec xin,NormType type,PetscReal *z)
 52: {
 53:   PetscReal      sum,work = 0.0;

 57:   if (type == NORM_2 || type == NORM_FROBENIUS) {
 58:     VecNorm_SeqCUDA(xin,NORM_2,&work);
 59:     work *= work;
 60:     MPIU_Allreduce(&work,&sum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 61:     *z    = PetscSqrtReal(sum);
 62:   } else if (type == NORM_1) {
 63:     /* Find the local part */
 64:     VecNorm_SeqCUDA(xin,NORM_1,&work);
 65:     /* Find the global max */
 66:     MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 67:   } else if (type == NORM_INFINITY) {
 68:     /* Find the local max */
 69:     VecNorm_SeqCUDA(xin,NORM_INFINITY,&work);
 70:     /* Find the global max */
 71:     MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)xin));
 72:   } else if (type == NORM_1_AND_2) {
 73:     PetscReal temp[2];
 74:     VecNorm_SeqCUDA(xin,NORM_1,temp);
 75:     VecNorm_SeqCUDA(xin,NORM_2,temp+1);
 76:     temp[1] = temp[1]*temp[1];
 77:     MPIU_Allreduce(temp,z,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 78:     z[1] = PetscSqrtReal(z[1]);
 79:   }
 80:   return(0);
 81: }

 83: PetscErrorCode VecDot_MPICUDA(Vec xin,Vec yin,PetscScalar *z)
 84: {
 85:   PetscScalar    sum,work;

 89:   VecDot_SeqCUDA(xin,yin,&work);
 90:   MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
 91:   *z   = sum;
 92:   return(0);
 93: }

 95: PetscErrorCode VecTDot_MPICUDA(Vec xin,Vec yin,PetscScalar *z)
 96: {
 97:   PetscScalar    sum,work;

101:   VecTDot_SeqCUDA(xin,yin,&work);
102:   MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
103:   *z   = sum;
104:   return(0);
105: }

107: PetscErrorCode VecMDot_MPICUDA(Vec xin,PetscInt nv,const Vec y[],PetscScalar *z)
108: {
109:   PetscScalar    awork[128],*work = awork;

113:   if (nv > 128) {
114:     PetscMalloc1(nv,&work);
115:   }
116:   VecMDot_SeqCUDA(xin,nv,y,work);
117:   MPIU_Allreduce(work,z,nv,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
118:   if (nv > 128) {
119:     PetscFree(work);
120:   }
121:   return(0);
122: }

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

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

130:   Level: beginner

132: .seealso: VecCreate(), VecSetType(), VecSetFromOptions(), VecCreateMPIWithArray(), VECMPI, VecType, VecCreateMPI(), VecSetPinnedMemoryMin()
133: M*/

135: PetscErrorCode VecDuplicate_MPICUDA(Vec win,Vec *v)
136: {
138:   Vec_MPI        *vw,*w = (Vec_MPI*)win->data;
139:   PetscScalar    *array;

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

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

149:   /* save local representation of the parallel vector (and scatter) if it exists */
150:   if (w->localrep) {
151:     VecGetArray(*v,&array);
152:     VecCreateSeqWithArray(PETSC_COMM_SELF,1,win->map->n+w->nghost,array,&vw->localrep);
153:     PetscMemcpy(vw->localrep->ops,w->localrep->ops,sizeof(struct _VecOps));
154:     VecRestoreArray(*v,&array);
155:     PetscLogObjectParent((PetscObject)*v,(PetscObject)vw->localrep);
156:     vw->localupdate = w->localupdate;
157:     if (vw->localupdate) {
158:       PetscObjectReference((PetscObject)vw->localupdate);
159:     }
160:   }

162:   /* New vector should inherit stashing property of parent */
163:   (*v)->stash.donotstash   = win->stash.donotstash;
164:   (*v)->stash.ignorenegidx = win->stash.ignorenegidx;

166:   /* change type_name appropriately */
167:   VecCUDAAllocateCheck(*v);
168:   PetscObjectChangeTypeName((PetscObject)(*v),VECMPICUDA);

170:   PetscObjectListDuplicate(((PetscObject)win)->olist,&((PetscObject)(*v))->olist);
171:   PetscFunctionListDuplicate(((PetscObject)win)->qlist,&((PetscObject)(*v))->qlist);
172:   (*v)->map->bs   = PetscAbs(win->map->bs);
173:   (*v)->bstash.bs = win->bstash.bs;
174:   return(0);
175: }

177: PetscErrorCode VecDotNorm2_MPICUDA(Vec s,Vec t,PetscScalar *dp,PetscScalar *nm)
178: {
180:   PetscScalar    work[2],sum[2];

183:   VecDotNorm2_SeqCUDA(s,t,work,work+1);
184:   MPIU_Allreduce(&work,&sum,2,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)s));
185:   *dp  = sum[0];
186:   *nm  = sum[1];
187:   return(0);
188: }

190: PetscErrorCode VecCreate_MPICUDA(Vec vv)
191: {

195:   PetscCUDAInitializeCheck();
196:   PetscLayoutSetUp(vv->map);
197:   VecCUDAAllocateCheck(vv);
198:   VecCreate_MPICUDA_Private(vv,PETSC_FALSE,0,((Vec_CUDA*)vv->spptr)->GPUarray_allocated);
199:   VecCUDAAllocateCheckHost(vv);
200:   VecSet(vv,0.0);
201:   VecSet_Seq(vv,0.0);
202:   vv->offloadmask = PETSC_OFFLOAD_BOTH;
203:   return(0);
204: }

206: PetscErrorCode VecCreate_CUDA(Vec v)
207: {
209:   PetscMPIInt    size;

212:   MPI_Comm_size(PetscObjectComm((PetscObject)v),&size);
213:   if (size == 1) {
214:     VecSetType(v,VECSEQCUDA);
215:   } else {
216:     VecSetType(v,VECMPICUDA);
217:   }
218:   return(0);
219: }

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

225:    Collective

227:    Input Parameters:
228: +  comm  - the MPI communicator to use
229: .  bs    - block size, same meaning as VecSetBlockSize()
230: .  n     - local vector length, cannot be PETSC_DECIDE
231: .  N     - global vector length (or PETSC_DECIDE to have calculated)
232: -  array - the user provided GPU array to store the vector values

234:    Output Parameter:
235: .  vv - the vector

237:    Notes:
238:    Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
239:    same type as an existing vector.

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

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

247:    Level: intermediate

249: .seealso: VecCreateSeqCUDAWithArray(), VecCreateMPIWithArray(), VecCreateSeqWithArray(),
250:           VecCreate(), VecDuplicate(), VecDuplicateVecs(), VecCreateGhost(),
251:           VecCreateMPI(), VecCreateGhostWithArray(), VecPlaceArray()

253: @*/
254: PetscErrorCode  VecCreateMPICUDAWithArray(MPI_Comm comm,PetscInt bs,PetscInt n,PetscInt N,const PetscScalar array[],Vec *vv)
255: {

259:   if (n == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local size of vector");
260:   PetscCUDAInitializeCheck();
261:   VecCreate(comm,vv);
262:   VecSetSizes(*vv,n,N);
263:   VecSetBlockSize(*vv,bs);
264:   VecCreate_MPICUDA_Private(*vv,PETSC_FALSE,0,array);
265:   return(0);
266: }

268: /*@C
269:    VecCreateMPICUDAWithArrays - Creates a parallel, array-style vector,
270:    where the user provides the GPU array space to store the vector values.

272:    Collective

274:    Input Parameters:
275: +  comm  - the MPI communicator to use
276: .  bs    - block size, same meaning as VecSetBlockSize()
277: .  n     - local vector length, cannot be PETSC_DECIDE
278: .  N     - global vector length (or PETSC_DECIDE to have calculated)
279: -  cpuarray - the user provided CPU array to store the vector values
280: -  gpuarray - the user provided GPU array to store the vector values

282:    Output Parameter:
283: .  vv - the vector

285:    Notes:
286:    If both cpuarray and gpuarray are provided, the caller must ensure that
287:    the provided arrays have identical values.

289:    Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
290:    same type as an existing vector.

292:    PETSc does NOT free the provided arrays when the vector is destroyed via
293:    VecDestroy(). The user should not free the array until the vector is
294:    destroyed.

296:    Level: intermediate

298: .seealso: VecCreateSeqCUDAWithArrays(), VecCreateMPIWithArray(), VecCreateSeqWithArray(),
299:           VecCreate(), VecDuplicate(), VecDuplicateVecs(), VecCreateGhost(),
300:           VecCreateMPI(), VecCreateGhostWithArray(), VecCUDAPlaceArray(), VecPlaceArray(),
301:           VecCUDAAllocateCheckHost()
302: @*/
303: PetscErrorCode  VecCreateMPICUDAWithArrays(MPI_Comm comm,PetscInt bs,PetscInt n,PetscInt N,const PetscScalar cpuarray[],const PetscScalar gpuarray[],Vec *vv)
304: {

308:   VecCreateMPICUDAWithArray(comm,bs,n,N,gpuarray,vv);

310:   if (cpuarray && gpuarray) {
311:     Vec_MPI *s         = (Vec_MPI*)((*vv)->data);
312:     s->array           = (PetscScalar*)cpuarray;
313:     (*vv)->offloadmask = PETSC_OFFLOAD_BOTH;
314:   } else if (cpuarray) {
315:     Vec_MPI *s         = (Vec_MPI*)((*vv)->data);
316:     s->array           = (PetscScalar*)cpuarray;
317:     (*vv)->offloadmask =  PETSC_OFFLOAD_CPU;
318:   } else if (gpuarray) {
319:     (*vv)->offloadmask = PETSC_OFFLOAD_GPU;
320:   } else {
321:     (*vv)->offloadmask = PETSC_OFFLOAD_UNALLOCATED;
322:   }

324:   return(0);
325: }

327: PetscErrorCode VecMax_MPICUDA(Vec xin,PetscInt *idx,PetscReal *z)
328: {
330:   PetscReal      work;

333:   VecMax_SeqCUDA(xin,idx,&work);
334:   if (!idx) {
335:     MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)xin));
336:   } else {
337:     struct { PetscReal v; PetscInt i; } in,out;

339:     in.v  = work;
340:     in.i  = *idx + xin->map->rstart;
341:     MPIU_Allreduce(&in,&out,1,MPIU_REAL_INT,MPIU_MAXLOC,PetscObjectComm((PetscObject)xin));
342:     *z    = out.v;
343:     *idx  = out.i;
344:   }
345:   return(0);
346: }

348: PetscErrorCode VecMin_MPICUDA(Vec xin,PetscInt *idx,PetscReal *z)
349: {
351:   PetscReal      work;

354:   VecMin_SeqCUDA(xin,idx,&work);
355:   if (!idx) {
356:     MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_MIN,PetscObjectComm((PetscObject)xin));
357:   } else {
358:     struct { PetscReal v; PetscInt i; } in,out;

360:     in.v  = work;
361:     in.i  = *idx + xin->map->rstart;
362:     MPIU_Allreduce(&in,&out,1,MPIU_REAL_INT,MPIU_MINLOC,PetscObjectComm((PetscObject)xin));
363:     *z    = out.v;
364:     *idx  = out.i;
365:   }
366:   return(0);
367: }

369: PetscErrorCode VecBindToCPU_MPICUDA(Vec V,PetscBool pin)
370: {

374:   V->boundtocpu = pin;
375:   if (pin) {
376:     VecCUDACopyFromGPU(V);
377:     V->offloadmask = PETSC_OFFLOAD_CPU; /* since the CPU code will likely change values in the vector */
378:     V->ops->dotnorm2               = NULL;
379:     V->ops->waxpy                  = VecWAXPY_Seq;
380:     V->ops->dot                    = VecDot_MPI;
381:     V->ops->mdot                   = VecMDot_MPI;
382:     V->ops->tdot                   = VecTDot_MPI;
383:     V->ops->norm                   = VecNorm_MPI;
384:     V->ops->scale                  = VecScale_Seq;
385:     V->ops->copy                   = VecCopy_Seq;
386:     V->ops->set                    = VecSet_Seq;
387:     V->ops->swap                   = VecSwap_Seq;
388:     V->ops->axpy                   = VecAXPY_Seq;
389:     V->ops->axpby                  = VecAXPBY_Seq;
390:     V->ops->maxpy                  = VecMAXPY_Seq;
391:     V->ops->aypx                   = VecAYPX_Seq;
392:     V->ops->axpbypcz               = VecAXPBYPCZ_Seq;
393:     V->ops->pointwisemult          = VecPointwiseMult_Seq;
394:     V->ops->setrandom              = VecSetRandom_Seq;
395:     V->ops->placearray             = VecPlaceArray_Seq;
396:     V->ops->replacearray           = VecReplaceArray_SeqCUDA;
397:     V->ops->resetarray             = VecResetArray_Seq;
398:     V->ops->dot_local              = VecDot_Seq;
399:     V->ops->tdot_local             = VecTDot_Seq;
400:     V->ops->norm_local             = VecNorm_Seq;
401:     V->ops->mdot_local             = VecMDot_Seq;
402:     V->ops->pointwisedivide        = VecPointwiseDivide_Seq;
403:     V->ops->getlocalvector         = NULL;
404:     V->ops->restorelocalvector     = NULL;
405:     V->ops->getlocalvectorread     = NULL;
406:     V->ops->restorelocalvectorread = NULL;
407:     V->ops->getarraywrite          = NULL;
408:     V->ops->max                    = VecMax_MPI;
409:     V->ops->min                    = VecMin_MPI;
410:     V->ops->reciprocal             = VecReciprocal_Default;
411:     V->ops->sum                    = NULL;
412:     V->ops->shift                  = NULL;
413:     /* default random number generator */
414:     PetscFree(V->defaultrandtype);
415:     PetscStrallocpy(PETSCRANDER48,&V->defaultrandtype);
416:   } else {
417:     V->ops->dotnorm2               = VecDotNorm2_MPICUDA;
418:     V->ops->waxpy                  = VecWAXPY_SeqCUDA;
419:     V->ops->duplicate              = VecDuplicate_MPICUDA;
420:     V->ops->dot                    = VecDot_MPICUDA;
421:     V->ops->mdot                   = VecMDot_MPICUDA;
422:     V->ops->tdot                   = VecTDot_MPICUDA;
423:     V->ops->norm                   = VecNorm_MPICUDA;
424:     V->ops->scale                  = VecScale_SeqCUDA;
425:     V->ops->copy                   = VecCopy_SeqCUDA;
426:     V->ops->set                    = VecSet_SeqCUDA;
427:     V->ops->swap                   = VecSwap_SeqCUDA;
428:     V->ops->axpy                   = VecAXPY_SeqCUDA;
429:     V->ops->axpby                  = VecAXPBY_SeqCUDA;
430:     V->ops->maxpy                  = VecMAXPY_SeqCUDA;
431:     V->ops->aypx                   = VecAYPX_SeqCUDA;
432:     V->ops->axpbypcz               = VecAXPBYPCZ_SeqCUDA;
433:     V->ops->pointwisemult          = VecPointwiseMult_SeqCUDA;
434:     V->ops->setrandom              = VecSetRandom_SeqCUDA;
435:     V->ops->placearray             = VecPlaceArray_SeqCUDA;
436:     V->ops->replacearray           = VecReplaceArray_SeqCUDA;
437:     V->ops->resetarray             = VecResetArray_SeqCUDA;
438:     V->ops->dot_local              = VecDot_SeqCUDA;
439:     V->ops->tdot_local             = VecTDot_SeqCUDA;
440:     V->ops->norm_local             = VecNorm_SeqCUDA;
441:     V->ops->mdot_local             = VecMDot_SeqCUDA;
442:     V->ops->destroy                = VecDestroy_MPICUDA;
443:     V->ops->pointwisedivide        = VecPointwiseDivide_SeqCUDA;
444:     V->ops->getlocalvector         = VecGetLocalVector_SeqCUDA;
445:     V->ops->restorelocalvector     = VecRestoreLocalVector_SeqCUDA;
446:     V->ops->getlocalvectorread     = VecGetLocalVector_SeqCUDA;
447:     V->ops->restorelocalvectorread = VecRestoreLocalVector_SeqCUDA;
448:     V->ops->getarraywrite          = VecGetArrayWrite_SeqCUDA;
449:     V->ops->getarray               = VecGetArray_SeqCUDA;
450:     V->ops->restorearray           = VecRestoreArray_SeqCUDA;
451:     V->ops->getarrayandmemtype     = VecGetArrayAndMemType_SeqCUDA;
452:     V->ops->restorearrayandmemtype = VecRestoreArrayAndMemType_SeqCUDA;
453:     V->ops->max                    = VecMax_MPICUDA;
454:     V->ops->min                    = VecMin_MPICUDA;
455:     V->ops->reciprocal             = VecReciprocal_SeqCUDA;
456:     V->ops->sum                    = VecSum_SeqCUDA;
457:     V->ops->shift                  = VecShift_SeqCUDA;
458:     /* default random number generator */
459:     PetscFree(V->defaultrandtype);
460:     PetscStrallocpy(PETSCCURAND,&V->defaultrandtype);
461:   }
462:   return(0);
463: }

465: PetscErrorCode VecCreate_MPICUDA_Private(Vec vv,PetscBool alloc,PetscInt nghost,const PetscScalar array[])
466: {
468:   Vec_CUDA       *veccuda;

471:   VecCreate_MPI_Private(vv,PETSC_FALSE,0,0);
472:   PetscObjectChangeTypeName((PetscObject)vv,VECMPICUDA);

474:   VecBindToCPU_MPICUDA(vv,PETSC_FALSE);
475:   vv->ops->bindtocpu = VecBindToCPU_MPICUDA;

477:   /* Later, functions check for the Vec_CUDA structure existence, so do not create it without array */
478:   if (alloc && !array) {
479:     VecCUDAAllocateCheck(vv);
480:     VecCUDAAllocateCheckHost(vv);
481:     VecSet(vv,0.0);
482:     VecSet_Seq(vv,0.0);
483:     vv->offloadmask = PETSC_OFFLOAD_BOTH;
484:   }
485:   if (array) {
486:     if (!vv->spptr) {
487:       PetscReal pinned_memory_min;
488:       PetscBool flag;
489:       /* Cannot use PetscNew() here because spptr is void* */
490:       PetscCalloc(sizeof(Vec_CUDA),&vv->spptr);
491:       veccuda = (Vec_CUDA*)vv->spptr;
492:       vv->minimum_bytes_pinned_memory = 0;

494:       /* Need to parse command line for minimum size to use for pinned memory allocations on host here.
495:          Note: This same code duplicated in VecCreate_SeqCUDA_Private() and VecCUDAAllocateCheck(). Is there a good way to avoid this? */
496:       PetscOptionsBegin(PetscObjectComm((PetscObject)vv),((PetscObject)vv)->prefix,"VECCUDA Options","Vec");
497:       pinned_memory_min = vv->minimum_bytes_pinned_memory;
498:       PetscOptionsReal("-vec_pinned_memory_min","Minimum size (in bytes) for an allocation to use pinned memory on host","VecSetPinnedMemoryMin",pinned_memory_min,&pinned_memory_min,&flag);
499:       if (flag) vv->minimum_bytes_pinned_memory = pinned_memory_min;
500:       PetscOptionsEnd();
501:     }
502:     veccuda = (Vec_CUDA*)vv->spptr;
503:     veccuda->GPUarray = (PetscScalar*)array;
504:     vv->offloadmask = PETSC_OFFLOAD_GPU;
505:   }
506:   return(0);
507: }