Actual source code: cudavecimpl.h

petsc-3.12.5 2020-03-29
Report Typos and Errors

  4:  #include <petscvec.h>
  5:  #include <petsccublas.h>
  6:  #include <petsc/private/vecimpl.h>

  8: #include <cublas_v2.h>

 10: #define WaitForGPU() PetscCUDASynchronize ? cudaDeviceSynchronize() : 0

 12: typedef struct {
 13:   PetscScalar  *GPUarray;           /* this always holds the GPU data */
 14:   PetscScalar  *GPUarray_allocated; /* if the array was allocated by PETSc this is its pointer */
 15:   cudaStream_t stream;              /* A stream for doing asynchronous data transfers */
 16:   PetscBool    hostDataRegisteredAsPageLocked;
 17: } Vec_CUDA;


 20: #include <cuda_runtime.h>

 22: PETSC_INTERN PetscErrorCode VecDotNorm2_SeqCUDA(Vec,Vec,PetscScalar*, PetscScalar*);
 23: PETSC_INTERN PetscErrorCode VecPointwiseDivide_SeqCUDA(Vec,Vec,Vec);
 24: PETSC_INTERN PetscErrorCode VecWAXPY_SeqCUDA(Vec,PetscScalar,Vec,Vec);
 25: PETSC_INTERN PetscErrorCode VecMDot_SeqCUDA(Vec,PetscInt,const Vec[],PetscScalar*);
 26: PETSC_EXTERN PetscErrorCode VecSet_SeqCUDA(Vec,PetscScalar);
 27: PETSC_INTERN PetscErrorCode VecMAXPY_SeqCUDA(Vec,PetscInt,const PetscScalar*,Vec*);
 28: PETSC_INTERN PetscErrorCode VecAXPBYPCZ_SeqCUDA(Vec,PetscScalar,PetscScalar,PetscScalar,Vec,Vec);
 29: PETSC_INTERN PetscErrorCode VecPointwiseMult_SeqCUDA(Vec,Vec,Vec);
 30: PETSC_INTERN PetscErrorCode VecPlaceArray_SeqCUDA(Vec,const PetscScalar*);
 31: PETSC_INTERN PetscErrorCode VecResetArray_SeqCUDA(Vec);
 32: PETSC_INTERN PetscErrorCode VecReplaceArray_SeqCUDA(Vec,const PetscScalar*);
 33: PETSC_INTERN PetscErrorCode VecDot_SeqCUDA(Vec,Vec,PetscScalar*);
 34: PETSC_INTERN PetscErrorCode VecTDot_SeqCUDA(Vec,Vec,PetscScalar*);
 35: PETSC_INTERN PetscErrorCode VecScale_SeqCUDA(Vec,PetscScalar);
 36: PETSC_EXTERN PetscErrorCode VecCopy_SeqCUDA(Vec,Vec);
 37: PETSC_INTERN PetscErrorCode VecSwap_SeqCUDA(Vec,Vec);
 38: PETSC_EXTERN PetscErrorCode VecAXPY_SeqCUDA(Vec,PetscScalar,Vec);
 39: PETSC_INTERN PetscErrorCode VecAXPBY_SeqCUDA(Vec,PetscScalar,PetscScalar,Vec);
 40: PETSC_INTERN PetscErrorCode VecDuplicate_SeqCUDA(Vec,Vec*);
 41: PETSC_INTERN PetscErrorCode VecConjugate_SeqCUDA(Vec xin);
 42: PETSC_INTERN PetscErrorCode VecNorm_SeqCUDA(Vec,NormType,PetscReal*);
 43: PETSC_INTERN PetscErrorCode VecCUDACopyToGPU(Vec);
 44: PETSC_INTERN PetscErrorCode VecCUDAAllocateCheck(Vec);
 45: PETSC_EXTERN PetscErrorCode VecCreate_SeqCUDA(Vec);
 46: PETSC_INTERN PetscErrorCode VecCreate_SeqCUDA_Private(Vec,const PetscScalar*);
 47: PETSC_INTERN PetscErrorCode VecCreate_MPICUDA(Vec);
 48: PETSC_INTERN PetscErrorCode VecCreate_MPICUDA_Private(Vec,PetscBool,PetscInt,const PetscScalar*);
 49: PETSC_INTERN PetscErrorCode VecCreate_CUDA(Vec);
 50: PETSC_INTERN PetscErrorCode VecDestroy_SeqCUDA(Vec);
 51: PETSC_INTERN PetscErrorCode VecDestroy_MPICUDA(Vec);
 52: PETSC_INTERN PetscErrorCode VecAYPX_SeqCUDA(Vec,PetscScalar,Vec);
 53: PETSC_INTERN PetscErrorCode VecSetRandom_SeqCUDA(Vec,PetscRandom);
 54: PETSC_INTERN PetscErrorCode VecGetLocalVector_SeqCUDA(Vec,Vec);
 55: PETSC_INTERN PetscErrorCode VecRestoreLocalVector_SeqCUDA(Vec,Vec);
 56: PETSC_INTERN PetscErrorCode VecGetArrayWrite_SeqCUDA(Vec,PetscScalar**);
 57: PETSC_INTERN PetscErrorCode VecCopy_SeqCUDA_Private(Vec xin,Vec yin);
 58: PETSC_INTERN PetscErrorCode VecSetRandom_SeqCUDA_Private(Vec xin,PetscRandom r);
 59: PETSC_INTERN PetscErrorCode VecDestroy_SeqCUDA_Private(Vec v);
 60: PETSC_INTERN PetscErrorCode VecResetArray_SeqCUDA_Private(Vec vin);
 61: PETSC_INTERN PetscErrorCode VecCUDACopyToGPU_Public(Vec);
 62: PETSC_INTERN PetscErrorCode VecCUDAAllocateCheck_Public(Vec);
 63: PETSC_INTERN PetscErrorCode VecCUDACopyToGPUSome(Vec,PetscCUDAIndices,ScatterMode);
 64: PETSC_INTERN PetscErrorCode VecCUDACopyFromGPUSome(Vec,PetscCUDAIndices,ScatterMode);

 66: PETSC_INTERN PetscErrorCode VecScatterCUDAIndicesCreate_PtoP(PetscInt, PetscInt*,PetscInt, PetscInt*,PetscCUDAIndices*);
 67: PETSC_INTERN PetscErrorCode VecScatterCUDAIndicesCreate_StoS(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*,PetscCUDAIndices*);
 68: PETSC_INTERN PetscErrorCode VecScatterCUDAIndicesDestroy(PetscCUDAIndices*);
 69: PETSC_INTERN PetscErrorCode VecScatterCUDA_StoS(Vec,Vec,PetscCUDAIndices,InsertMode,ScatterMode);

 71: typedef enum {VEC_SCATTER_CUDA_STOS, VEC_SCATTER_CUDA_PTOP} VecCUDAScatterType;
 72: typedef enum {VEC_SCATTER_CUDA_GENERAL, VEC_SCATTER_CUDA_STRIDED} VecCUDASequentialScatterMode;

 74: struct  _p_VecScatterCUDAIndices_PtoP {
 75:   PetscInt ns;
 76:   PetscInt sendLowestIndex;
 77:   PetscInt nr;
 78:   PetscInt recvLowestIndex;
 79: };

 81: struct  _p_VecScatterCUDAIndices_StoS {
 82:   /* from indices data */
 83:   PetscInt *fslots;
 84:   PetscInt fromFirst;
 85:   PetscInt fromStep;
 86:   VecCUDASequentialScatterMode fromMode;

 88:   /* to indices data */
 89:   PetscInt *tslots;
 90:   PetscInt toFirst;
 91:   PetscInt toStep;
 92:   VecCUDASequentialScatterMode toMode;

 94:   PetscInt n;
 95:   PetscInt MAX_BLOCKS;
 96:   PetscInt MAX_CORESIDENT_THREADS;
 97:   cudaStream_t stream;
 98: };

100: struct  _p_PetscCUDAIndices {
101:   void * scatter;
102:   VecCUDAScatterType scatterType;
103: };

105: /* complex single */
106: #if defined(PETSC_USE_COMPLEX)
107: #if defined(PETSC_USE_REAL_SINGLE)
108: #define cublasXaxpy(a,b,c,d,e,f,g) cublasCaxpy((a),(b),(cuComplex*)(c),(cuComplex*)(d),(e),(cuComplex*)(f),(g))
109: #define cublasXscal(a,b,c,d,e)     cublasCscal((a),(b),(cuComplex*)(c),(cuComplex*)(d),(e))
110: #define cublasXdotu(a,b,c,d,e,f,g) cublasCdotu((a),(b),(cuComplex*)(c),(d),(cuComplex*)(e),(f),(cuComplex*)(g))
111: #define cublasXdot(a,b,c,d,e,f,g)  cublasCdotc((a),(b),(cuComplex*)(c),(d),(cuComplex*)(e),(f),(cuComplex*)(g))
112: #define cublasXswap(a,b,c,d,e,f)   cublasCswap((a),(b),(cuComplex*)(c),(d),(cuComplex*)(e),(f))
113: #define cublasXnrm2(a,b,c,d,e)     cublasScnrm2((a),(b),(cuComplex*)(c),(d),(e))
114: #define cublasIXamax(a,b,c,d,e)    cublasIcamax((a),(b),(cuComplex*)(c),(d),(e))
115: #define cublasXasum(a,b,c,d,e)     cublasScasum((a),(b),(cuComplex*)(c),(d),(e))
116: #define cublasXgemv(a,b,c,d,e,f,g,h,i,j,k,l) cublasCgemv((a),(b),(c),(d),(cuComplex*)(e),(cuComplex*)(f),(g),(cuComplex*)(h),(i),(cuComplex*)(j),(cuComplex*)(k),(l))
117: #define cublasXgemm(a,b,c,d,e,f,g,h,i,j,k,l,m,n) cublasCgemm((a),(b),(c),(d),(e),(f),(cuComplex*)(g),(cuComplex*)(h),(i),(cuComplex*)(j),(k),(cuComplex*)(l),(cuComplex*)(m),(n))
118: #else /* complex double */
119: #define cublasXaxpy(a,b,c,d,e,f,g) cublasZaxpy((a),(b),(cuDoubleComplex*)(c),(cuDoubleComplex*)(d),(e),(cuDoubleComplex*)(f),(g))
120: #define cublasXscal(a,b,c,d,e)     cublasZscal((a),(b),(cuDoubleComplex*)(c),(cuDoubleComplex*)(d),(e))
121: #define cublasXdotu(a,b,c,d,e,f,g) cublasZdotu((a),(b),(cuDoubleComplex*)(c),(d),(cuDoubleComplex*)(e),(f),(cuDoubleComplex*)(g))
122: #define cublasXdot(a,b,c,d,e,f,g)  cublasZdotc((a),(b),(cuDoubleComplex*)(c),(d),(cuDoubleComplex*)(e),(f),(cuDoubleComplex*)(g))
123: #define cublasXswap(a,b,c,d,e,f)   cublasZswap((a),(b),(cuDoubleComplex*)(c),(d),(cuDoubleComplex*)(e),(f))
124: #define cublasXnrm2(a,b,c,d,e)     cublasDznrm2((a),(b),(cuDoubleComplex*)(c),(d),(e))
125: #define cublasIXamax(a,b,c,d,e)    cublasIzamax((a),(b),(cuDoubleComplex*)(c),(d),(e))
126: #define cublasXasum(a,b,c,d,e)     cublasDzasum((a),(b),(cuDoubleComplex*)(c),(d),(e))
127: #define cublasXgemv(a,b,c,d,e,f,g,h,i,j,k,l) cublasZgemv((a),(b),(c),(d),(cuDoubleComplex*)(e),(cuDoubleComplex*)(f),(g),(cuDoubleComplex*)(h),(i),(cuDoubleComplex*)(j),(cuDoubleComplex*)(k),(l))
128: #define cublasXgemm(a,b,c,d,e,f,g,h,i,j,k,l,m,n) cublasZgemm((a),(b),(c),(d),(e),(f),(cuDoubleComplex*)(g),(cuDoubleComplex*)(h),(i),(cuDoubleComplex*)(j),(k),(cuDoubleComplex*)(l),(cuDoubleComplex*)(m),(n))
129: #endif
130: #else /* real single */
131: #if defined(PETSC_USE_REAL_SINGLE)
132: #define cublasXaxpy  cublasSaxpy
133: #define cublasXscal  cublasSscal
134: #define cublasXdotu  cublasSdot
135: #define cublasXdot   cublasSdot
136: #define cublasXswap  cublasSswap
137: #define cublasXnrm2  cublasSnrm2
138: #define cublasIXamax cublasIsamax
139: #define cublasXasum  cublasSasum
140: #define cublasXgemv  cublasSgemv
141: #define cublasXgemm  cublasSgemm
142: #else /* real double */
143: #define cublasXaxpy  cublasDaxpy
144: #define cublasXscal  cublasDscal
145: #define cublasXdotu  cublasDdot
146: #define cublasXdot   cublasDdot
147: #define cublasXswap  cublasDswap
148: #define cublasXnrm2  cublasDnrm2
149: #define cublasIXamax cublasIdamax
150: #define cublasXasum  cublasDasum
151: #define cublasXgemv  cublasDgemv
152: #define cublasXgemm  cublasDgemm
153: #endif
154: #endif

156: #endif