Actual source code: cudavecimpl.h

petsc-master 2020-05-30
Report Typos and Errors

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

  8: #include <cublas_v2.h>

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

 16: #include <cuda_runtime.h>

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

 62: PETSC_INTERN PetscErrorCode VecScatterCUDAIndicesCreate_PtoP(PetscInt, PetscInt*,PetscInt, PetscInt*,PetscCUDAIndices*);
 63: PETSC_INTERN PetscErrorCode VecScatterCUDAIndicesCreate_StoS(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*,PetscCUDAIndices*);
 64: PETSC_INTERN PetscErrorCode VecScatterCUDAIndicesDestroy(PetscCUDAIndices*);
 65: PETSC_INTERN PetscErrorCode VecScatterCUDA_StoS(Vec,Vec,PetscCUDAIndices,InsertMode,ScatterMode);

 67: typedef enum {VEC_SCATTER_CUDA_STOS, VEC_SCATTER_CUDA_PTOP} VecCUDAScatterType;
 68: typedef enum {VEC_SCATTER_CUDA_GENERAL, VEC_SCATTER_CUDA_STRIDED} VecCUDASequentialScatterMode;

 70: struct  _p_VecScatterCUDAIndices_PtoP {
 71:   PetscInt ns;
 72:   PetscInt sendLowestIndex;
 73:   PetscInt nr;
 74:   PetscInt recvLowestIndex;
 75: };

 77: struct  _p_VecScatterCUDAIndices_StoS {
 78:   /* from indices data */
 79:   PetscInt *fslots;
 80:   PetscInt fromFirst;
 81:   PetscInt fromStep;
 82:   VecCUDASequentialScatterMode fromMode;

 84:   /* to indices data */
 85:   PetscInt *tslots;
 86:   PetscInt toFirst;
 87:   PetscInt toStep;
 88:   VecCUDASequentialScatterMode toMode;

 90:   PetscInt n;
 91:   PetscInt MAX_BLOCKS;
 92:   PetscInt MAX_CORESIDENT_THREADS;
 93:   cudaStream_t stream;
 94: };

 96: struct  _p_PetscCUDAIndices {
 97:   void * scatter;
 98:   VecCUDAScatterType scatterType;
 99: };

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

156: #endif