Actual source code: cudavecimpl.h

petsc-3.13.6 2020-09-29
Report Typos and Errors

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

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

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

 58: PETSC_INTERN PetscErrorCode VecScatterCUDAIndicesCreate_PtoP(PetscInt, PetscInt*,PetscInt, PetscInt*,PetscCUDAIndices*);
 59: PETSC_INTERN PetscErrorCode VecScatterCUDAIndicesCreate_StoS(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*,PetscCUDAIndices*);
 60: PETSC_INTERN PetscErrorCode VecScatterCUDAIndicesDestroy(PetscCUDAIndices*);
 61: PETSC_INTERN PetscErrorCode VecScatterCUDA_StoS(Vec,Vec,PetscCUDAIndices,InsertMode,ScatterMode);

 63: typedef enum {VEC_SCATTER_CUDA_STOS, VEC_SCATTER_CUDA_PTOP} VecCUDAScatterType;
 64: typedef enum {VEC_SCATTER_CUDA_GENERAL, VEC_SCATTER_CUDA_STRIDED} VecCUDASequentialScatterMode;

 66: struct  _p_VecScatterCUDAIndices_PtoP {
 67:   PetscInt ns;
 68:   PetscInt sendLowestIndex;
 69:   PetscInt nr;
 70:   PetscInt recvLowestIndex;
 71: };

 73: struct  _p_VecScatterCUDAIndices_StoS {
 74:   /* from indices data */
 75:   PetscInt *fslots;
 76:   PetscInt fromFirst;
 77:   PetscInt fromStep;
 78:   VecCUDASequentialScatterMode fromMode;

 80:   /* to indices data */
 81:   PetscInt *tslots;
 82:   PetscInt toFirst;
 83:   PetscInt toStep;
 84:   VecCUDASequentialScatterMode toMode;

 86:   PetscInt n;
 87:   PetscInt MAX_BLOCKS;
 88:   PetscInt MAX_CORESIDENT_THREADS;
 89:   cudaStream_t stream;
 90: };

 92: struct  _p_PetscCUDAIndices {
 93:   void * scatter;
 94:   VecCUDAScatterType scatterType;
 95: };

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

148: #endif