Actual source code: mpicuda.cu
petsc-3.13.3 2020-07-01
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: /*MC
13: VECCUDA - VECCUDA = "cuda" - A VECSEQCUDA on a single-process communicator, and VECMPICUDA otherwise.
15: Options Database Keys:
16: . -vec_type cuda - sets the vector type to VECCUDA during a call to VecSetFromOptions()
18: Level: beginner
20: .seealso: VecCreate(), VecSetType(), VecSetFromOptions(), VecCreateMPIWithArray(), VECSEQCUDA, VECMPICUDA, VECSTANDARD, VecType, VecCreateMPI(), VecSetPinnedMemoryMin()
21: M*/
23: PetscErrorCode VecDestroy_MPICUDA(Vec v)
24: {
25: Vec_MPI *vecmpi = (Vec_MPI*)v->data;
26: Vec_CUDA *veccuda;
28: cudaError_t err;
31: if (v->spptr) {
32: veccuda = (Vec_CUDA*)v->spptr;
33: if (veccuda->GPUarray_allocated) {
34: err = cudaFree(((Vec_CUDA*)v->spptr)->GPUarray_allocated);CHKERRCUDA(err);
35: veccuda->GPUarray_allocated = NULL;
36: }
37: if (veccuda->stream) {
38: err = cudaStreamDestroy(((Vec_CUDA*)v->spptr)->stream);CHKERRCUDA(err);
39: }
40: if (v->pinned_memory) {
41: PetscMallocSetCUDAHost();
42: PetscFree(vecmpi->array_allocated);
43: PetscMallocResetCUDAHost();
44: v->pinned_memory = PETSC_FALSE;
45: }
46: PetscFree(v->spptr);
47: }
48: VecDestroy_MPI(v);
49: return(0);
50: }
52: PetscErrorCode VecNorm_MPICUDA(Vec xin,NormType type,PetscReal *z)
53: {
54: PetscReal sum,work = 0.0;
58: if (type == NORM_2 || type == NORM_FROBENIUS) {
59: VecNorm_SeqCUDA(xin,NORM_2,&work);
60: work *= work;
61: MPIU_Allreduce(&work,&sum,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
62: *z = PetscSqrtReal(sum);
63: } else if (type == NORM_1) {
64: /* Find the local part */
65: VecNorm_SeqCUDA(xin,NORM_1,&work);
66: /* Find the global max */
67: MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
68: } else if (type == NORM_INFINITY) {
69: /* Find the local max */
70: VecNorm_SeqCUDA(xin,NORM_INFINITY,&work);
71: /* Find the global max */
72: MPIU_Allreduce(&work,z,1,MPIU_REAL,MPIU_MAX,PetscObjectComm((PetscObject)xin));
73: } else if (type == NORM_1_AND_2) {
74: PetscReal temp[2];
75: VecNorm_SeqCUDA(xin,NORM_1,temp);
76: VecNorm_SeqCUDA(xin,NORM_2,temp+1);
77: temp[1] = temp[1]*temp[1];
78: MPIU_Allreduce(temp,z,2,MPIU_REAL,MPIU_SUM,PetscObjectComm((PetscObject)xin));
79: z[1] = PetscSqrtReal(z[1]);
80: }
81: return(0);
82: }
84: PetscErrorCode VecDot_MPICUDA(Vec xin,Vec yin,PetscScalar *z)
85: {
86: PetscScalar sum,work;
90: VecDot_SeqCUDA(xin,yin,&work);
91: MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
92: *z = sum;
93: return(0);
94: }
96: PetscErrorCode VecTDot_MPICUDA(Vec xin,Vec yin,PetscScalar *z)
97: {
98: PetscScalar sum,work;
102: VecTDot_SeqCUDA(xin,yin,&work);
103: MPIU_Allreduce(&work,&sum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
104: *z = sum;
105: return(0);
106: }
108: PetscErrorCode VecMDot_MPICUDA(Vec xin,PetscInt nv,const Vec y[],PetscScalar *z)
109: {
110: PetscScalar awork[128],*work = awork;
114: if (nv > 128) {
115: PetscMalloc1(nv,&work);
116: }
117: VecMDot_SeqCUDA(xin,nv,y,work);
118: MPIU_Allreduce(work,z,nv,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)xin));
119: if (nv > 128) {
120: PetscFree(work);
121: }
122: return(0);
123: }
125: /*MC
126: VECMPICUDA - VECMPICUDA = "mpicuda" - The basic parallel vector, modified to use CUDA
128: Options Database Keys:
129: . -vec_type mpicuda - sets the vector type to VECMPICUDA during a call to VecSetFromOptions()
131: Level: beginner
133: .seealso: VecCreate(), VecSetType(), VecSetFromOptions(), VecCreateMPIWithArray(), VECMPI, VecType, VecCreateMPI(), VecSetPinnedMemoryMin()
134: M*/
137: PetscErrorCode VecDuplicate_MPICUDA(Vec win,Vec *v)
138: {
140: Vec_MPI *vw,*w = (Vec_MPI*)win->data;
141: PetscScalar *array;
144: VecCreate(PetscObjectComm((PetscObject)win),v);
145: PetscLayoutReference(win->map,&(*v)->map);
147: VecCreate_MPICUDA_Private(*v,PETSC_TRUE,w->nghost,0);
148: vw = (Vec_MPI*)(*v)->data;
149: PetscMemcpy((*v)->ops,win->ops,sizeof(struct _VecOps));
151: /* save local representation of the parallel vector (and scatter) if it exists */
152: if (w->localrep) {
153: VecGetArray(*v,&array);
154: VecCreateSeqWithArray(PETSC_COMM_SELF,1,win->map->n+w->nghost,array,&vw->localrep);
155: PetscMemcpy(vw->localrep->ops,w->localrep->ops,sizeof(struct _VecOps));
156: VecRestoreArray(*v,&array);
157: PetscLogObjectParent((PetscObject)*v,(PetscObject)vw->localrep);
158: vw->localupdate = w->localupdate;
159: if (vw->localupdate) {
160: PetscObjectReference((PetscObject)vw->localupdate);
161: }
162: }
164: /* New vector should inherit stashing property of parent */
165: (*v)->stash.donotstash = win->stash.donotstash;
166: (*v)->stash.ignorenegidx = win->stash.ignorenegidx;
168: /* change type_name appropriately */
169: VecCUDAAllocateCheck(*v);
170: PetscObjectChangeTypeName((PetscObject)(*v),VECMPICUDA);
172: PetscObjectListDuplicate(((PetscObject)win)->olist,&((PetscObject)(*v))->olist);
173: PetscFunctionListDuplicate(((PetscObject)win)->qlist,&((PetscObject)(*v))->qlist);
174: (*v)->map->bs = PetscAbs(win->map->bs);
175: (*v)->bstash.bs = win->bstash.bs;
176: return(0);
177: }
179: PetscErrorCode VecDotNorm2_MPICUDA(Vec s,Vec t,PetscScalar *dp,PetscScalar *nm)
180: {
182: PetscScalar work[2],sum[2];
185: VecDotNorm2_SeqCUDA(s,t,work,work+1);
186: MPIU_Allreduce(&work,&sum,2,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)s));
187: *dp = sum[0];
188: *nm = sum[1];
189: return(0);
190: }
192: PetscErrorCode VecCreate_MPICUDA(Vec vv)
193: {
197: PetscLayoutSetUp(vv->map);
198: VecCUDAAllocateCheck(vv);
199: VecCreate_MPICUDA_Private(vv,PETSC_FALSE,0,((Vec_CUDA*)vv->spptr)->GPUarray_allocated);
200: VecCUDAAllocateCheckHost(vv);
201: VecSet(vv,0.0);
202: VecSet_Seq(vv,0.0);
203: vv->offloadmask = PETSC_OFFLOAD_BOTH;
204: return(0);
205: }
207: PetscErrorCode VecCreate_CUDA(Vec v)
208: {
210: PetscMPIInt size;
213: MPI_Comm_size(PetscObjectComm((PetscObject)v),&size);
214: if (size == 1) {
215: VecSetType(v,VECSEQCUDA);
216: } else {
217: VecSetType(v,VECMPICUDA);
218: }
219: return(0);
220: }
222: /*@C
223: VecCreateMPICUDAWithArray - Creates a parallel, array-style vector,
224: where the user provides the GPU array space to store the vector values.
226: Collective
228: Input Parameters:
229: + comm - the MPI communicator to use
230: . bs - block size, same meaning as VecSetBlockSize()
231: . n - local vector length, cannot be PETSC_DECIDE
232: . N - global vector length (or PETSC_DECIDE to have calculated)
233: - array - the user provided GPU array to store the vector values
235: Output Parameter:
236: . vv - the vector
238: Notes:
239: Use VecDuplicate() or VecDuplicateVecs() to form additional vectors of the
240: same type as an existing vector.
242: If the user-provided array is NULL, then VecCUDAPlaceArray() can be used
243: at a later stage to SET the array for storing the vector values.
245: PETSc does NOT free the array when the vector is destroyed via VecDestroy().
246: The user should not free the array until the vector is destroyed.
248: Level: intermediate
250: .seealso: VecCreateSeqCUDAWithArray(), VecCreateMPIWithArray(), VecCreateSeqWithArray(),
251: VecCreate(), VecDuplicate(), VecDuplicateVecs(), VecCreateGhost(),
252: VecCreateMPI(), VecCreateGhostWithArray(), VecPlaceArray()
254: @*/
255: PetscErrorCode VecCreateMPICUDAWithArray(MPI_Comm comm,PetscInt bs,PetscInt n,PetscInt N,const PetscScalar array[],Vec *vv)
256: {
260: if (n == PETSC_DECIDE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must set local size of vector");
261: PetscSplitOwnership(comm,&n,&N);
262: VecCreate(comm,vv);
263: VecSetSizes(*vv,n,N);
264: VecSetBlockSize(*vv,bs);
265: VecCreate_MPICUDA_Private(*vv,PETSC_FALSE,0,array);
266: return(0);
267: }
269: PetscErrorCode VecBindToCPU_MPICUDA(Vec V,PetscBool pin)
270: {
274: V->boundtocpu = pin;
275: if (pin) {
276: VecCUDACopyFromGPU(V);
277: V->offloadmask = PETSC_OFFLOAD_CPU; /* since the CPU code will likely change values in the vector */
278: V->ops->dotnorm2 = NULL;
279: V->ops->waxpy = VecWAXPY_Seq;
280: V->ops->dot = VecDot_MPI;
281: V->ops->mdot = VecMDot_MPI;
282: V->ops->tdot = VecTDot_MPI;
283: V->ops->norm = VecNorm_MPI;
284: V->ops->scale = VecScale_Seq;
285: V->ops->copy = VecCopy_Seq;
286: V->ops->set = VecSet_Seq;
287: V->ops->swap = VecSwap_Seq;
288: V->ops->axpy = VecAXPY_Seq;
289: V->ops->axpby = VecAXPBY_Seq;
290: V->ops->maxpy = VecMAXPY_Seq;
291: V->ops->aypx = VecAYPX_Seq;
292: V->ops->axpbypcz = VecAXPBYPCZ_Seq;
293: V->ops->pointwisemult = VecPointwiseMult_Seq;
294: V->ops->setrandom = VecSetRandom_Seq;
295: V->ops->placearray = VecPlaceArray_Seq;
296: V->ops->replacearray = VecReplaceArray_Seq;
297: V->ops->resetarray = VecResetArray_Seq;
298: V->ops->dot_local = VecDot_Seq;
299: V->ops->tdot_local = VecTDot_Seq;
300: V->ops->norm_local = VecNorm_Seq;
301: V->ops->mdot_local = VecMDot_Seq;
302: V->ops->pointwisedivide = VecPointwiseDivide_Seq;
303: V->ops->getlocalvector = NULL;
304: V->ops->restorelocalvector = NULL;
305: V->ops->getlocalvectorread = NULL;
306: V->ops->restorelocalvectorread = NULL;
307: V->ops->getarraywrite = NULL;
308: } else {
309: V->ops->dotnorm2 = VecDotNorm2_MPICUDA;
310: V->ops->waxpy = VecWAXPY_SeqCUDA;
311: V->ops->duplicate = VecDuplicate_MPICUDA;
312: V->ops->dot = VecDot_MPICUDA;
313: V->ops->mdot = VecMDot_MPICUDA;
314: V->ops->tdot = VecTDot_MPICUDA;
315: V->ops->norm = VecNorm_MPICUDA;
316: V->ops->scale = VecScale_SeqCUDA;
317: V->ops->copy = VecCopy_SeqCUDA;
318: V->ops->set = VecSet_SeqCUDA;
319: V->ops->swap = VecSwap_SeqCUDA;
320: V->ops->axpy = VecAXPY_SeqCUDA;
321: V->ops->axpby = VecAXPBY_SeqCUDA;
322: V->ops->maxpy = VecMAXPY_SeqCUDA;
323: V->ops->aypx = VecAYPX_SeqCUDA;
324: V->ops->axpbypcz = VecAXPBYPCZ_SeqCUDA;
325: V->ops->pointwisemult = VecPointwiseMult_SeqCUDA;
326: V->ops->setrandom = VecSetRandom_SeqCUDA;
327: V->ops->placearray = VecPlaceArray_SeqCUDA;
328: V->ops->replacearray = VecReplaceArray_SeqCUDA;
329: V->ops->resetarray = VecResetArray_SeqCUDA;
330: V->ops->dot_local = VecDot_SeqCUDA;
331: V->ops->tdot_local = VecTDot_SeqCUDA;
332: V->ops->norm_local = VecNorm_SeqCUDA;
333: V->ops->mdot_local = VecMDot_SeqCUDA;
334: V->ops->destroy = VecDestroy_MPICUDA;
335: V->ops->pointwisedivide = VecPointwiseDivide_SeqCUDA;
336: V->ops->getlocalvector = VecGetLocalVector_SeqCUDA;
337: V->ops->restorelocalvector = VecRestoreLocalVector_SeqCUDA;
338: V->ops->getlocalvectorread = VecGetLocalVector_SeqCUDA;
339: V->ops->restorelocalvectorread = VecRestoreLocalVector_SeqCUDA;
340: V->ops->getarraywrite = VecGetArrayWrite_SeqCUDA;
341: }
342: return(0);
343: }
345: PetscErrorCode VecCreate_MPICUDA_Private(Vec vv,PetscBool alloc,PetscInt nghost,const PetscScalar array[])
346: {
348: Vec_CUDA *veccuda;
351: VecCreate_MPI_Private(vv,PETSC_FALSE,0,0);
352: PetscObjectChangeTypeName((PetscObject)vv,VECMPICUDA);
354: VecBindToCPU_MPICUDA(vv,PETSC_FALSE);
355: vv->ops->bindtocpu = VecBindToCPU_MPICUDA;
357: /* Later, functions check for the Vec_CUDA structure existence, so do not create it without array */
358: if (alloc && !array) {
359: VecCUDAAllocateCheck(vv);
360: VecCUDAAllocateCheckHost(vv);
361: VecSet(vv,0.0);
362: VecSet_Seq(vv,0.0);
363: vv->offloadmask = PETSC_OFFLOAD_BOTH;
364: }
365: if (array) {
366: if (!vv->spptr) {
367: PetscReal pinned_memory_min;
368: PetscBool flag;
369: /* Cannot use PetscNew() here because spptr is void* */
370: PetscMalloc(sizeof(Vec_CUDA),&vv->spptr);
371: veccuda = (Vec_CUDA*)vv->spptr;
372: veccuda->stream = 0; /* using default stream */
373: veccuda->GPUarray_allocated = 0;
374: vv->offloadmask = PETSC_OFFLOAD_UNALLOCATED;
375: vv->minimum_bytes_pinned_memory = 0;
377: /* Need to parse command line for minimum size to use for pinned memory allocations on host here.
378: Note: This same code duplicated in VecCreate_SeqCUDA_Private() and VecCUDAAllocateCheck(). Is there a good way to avoid this? */
379: PetscOptionsBegin(PetscObjectComm((PetscObject)vv),((PetscObject)vv)->prefix,"VECCUDA Options","Vec");
380: pinned_memory_min = vv->minimum_bytes_pinned_memory;
381: 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);
382: if (flag) vv->minimum_bytes_pinned_memory = pinned_memory_min;
383: PetscOptionsEnd();
384: }
385: veccuda = (Vec_CUDA*)vv->spptr;
386: veccuda->GPUarray = (PetscScalar*)array;
387: }
389: return(0);
390: }