Actual source code: cudavecimpl.h
petsc-3.11.4 2019-09-28
4: #if defined(__CUDACC__)
6: #include <petsccuda.h>
7: #include <petsc/private/vecimpl.h>
9: #include <cublas_v2.h>
11: #define WaitForGPU() PetscCUDASynchronize ? cudaDeviceSynchronize() : 0
13: #endif
15: typedef struct {
16: PetscScalar *GPUarray; /* this always holds the GPU data */
17: PetscScalar *GPUarray_allocated; /* if the array was allocated by PETSc this is its pointer */
18: cudaStream_t stream; /* A stream for doing asynchronous data transfers */
19: PetscBool hostDataRegisteredAsPageLocked;
20: } Vec_CUDA;
23: #include <cuda_runtime.h>
25: PETSC_INTERN PetscErrorCode VecDotNorm2_SeqCUDA(Vec,Vec,PetscScalar*, PetscScalar*);
26: PETSC_INTERN PetscErrorCode VecPointwiseDivide_SeqCUDA(Vec,Vec,Vec);
27: PETSC_INTERN PetscErrorCode VecWAXPY_SeqCUDA(Vec,PetscScalar,Vec,Vec);
28: PETSC_INTERN PetscErrorCode VecMDot_SeqCUDA(Vec,PetscInt,const Vec[],PetscScalar*);
29: PETSC_INTERN PetscErrorCode VecSet_SeqCUDA(Vec,PetscScalar);
30: PETSC_INTERN PetscErrorCode VecMAXPY_SeqCUDA(Vec,PetscInt,const PetscScalar*,Vec*);
31: PETSC_INTERN PetscErrorCode VecAXPBYPCZ_SeqCUDA(Vec,PetscScalar,PetscScalar,PetscScalar,Vec,Vec);
32: PETSC_INTERN PetscErrorCode VecPointwiseMult_SeqCUDA(Vec,Vec,Vec);
33: PETSC_INTERN PetscErrorCode VecPlaceArray_SeqCUDA(Vec,const PetscScalar*);
34: PETSC_INTERN PetscErrorCode VecResetArray_SeqCUDA(Vec);
35: PETSC_INTERN PetscErrorCode VecReplaceArray_SeqCUDA(Vec,const PetscScalar*);
36: PETSC_INTERN PetscErrorCode VecDot_SeqCUDA(Vec,Vec,PetscScalar*);
37: PETSC_INTERN PetscErrorCode VecTDot_SeqCUDA(Vec,Vec,PetscScalar*);
38: PETSC_INTERN PetscErrorCode VecScale_SeqCUDA(Vec,PetscScalar);
39: PETSC_EXTERN PetscErrorCode VecCopy_SeqCUDA(Vec,Vec);
40: PETSC_INTERN PetscErrorCode VecSwap_SeqCUDA(Vec,Vec);
41: PETSC_INTERN PetscErrorCode VecAXPY_SeqCUDA(Vec,PetscScalar,Vec);
42: PETSC_INTERN PetscErrorCode VecAXPBY_SeqCUDA(Vec,PetscScalar,PetscScalar,Vec);
43: PETSC_INTERN PetscErrorCode VecDuplicate_SeqCUDA(Vec,Vec*);
44: PETSC_INTERN PetscErrorCode VecConjugate_SeqCUDA(Vec xin);
45: PETSC_INTERN PetscErrorCode VecNorm_SeqCUDA(Vec,NormType,PetscReal*);
46: PETSC_INTERN PetscErrorCode VecCUDACopyToGPU(Vec);
47: PETSC_INTERN PetscErrorCode VecCUDAAllocateCheck(Vec);
48: PETSC_EXTERN PetscErrorCode VecCreate_SeqCUDA(Vec);
49: PETSC_INTERN PetscErrorCode VecCreate_SeqCUDA_Private(Vec,const PetscScalar*);
50: PETSC_INTERN PetscErrorCode VecCreate_MPICUDA(Vec);
51: PETSC_INTERN PetscErrorCode VecCreate_MPICUDA_Private(Vec,PetscBool,PetscInt,const PetscScalar*);
52: PETSC_INTERN PetscErrorCode VecCreate_CUDA(Vec);
53: PETSC_INTERN PetscErrorCode VecDestroy_SeqCUDA(Vec);
54: PETSC_INTERN PetscErrorCode VecDestroy_MPICUDA(Vec);
55: PETSC_INTERN PetscErrorCode VecAYPX_SeqCUDA(Vec,PetscScalar,Vec);
56: PETSC_INTERN PetscErrorCode VecSetRandom_SeqCUDA(Vec,PetscRandom);
57: PETSC_INTERN PetscErrorCode VecGetLocalVector_SeqCUDA(Vec,Vec);
58: PETSC_INTERN PetscErrorCode VecRestoreLocalVector_SeqCUDA(Vec,Vec);
59: PETSC_INTERN PetscErrorCode VecCopy_SeqCUDA_Private(Vec xin,Vec yin);
60: PETSC_INTERN PetscErrorCode VecSetRandom_SeqCUDA_Private(Vec xin,PetscRandom r);
61: PETSC_INTERN PetscErrorCode VecDestroy_SeqCUDA_Private(Vec v);
62: PETSC_INTERN PetscErrorCode VecResetArray_SeqCUDA_Private(Vec vin);
63: PETSC_INTERN PetscErrorCode VecCUDACopyToGPU_Public(Vec);
64: PETSC_INTERN PetscErrorCode VecCUDAAllocateCheck_Public(Vec);
65: PETSC_INTERN PetscErrorCode VecCUDACopyToGPUSome(Vec,PetscCUDAIndices,ScatterMode);
66: PETSC_INTERN PetscErrorCode VecCUDACopyFromGPUSome(Vec,PetscCUDAIndices,ScatterMode);
68: PETSC_INTERN PetscErrorCode VecScatterCUDAIndicesCreate_PtoP(PetscInt, PetscInt*,PetscInt, PetscInt*,PetscCUDAIndices*);
69: PETSC_INTERN PetscErrorCode VecScatterCUDAIndicesCreate_StoS(PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt*,PetscInt*,PetscCUDAIndices*);
70: PETSC_INTERN PetscErrorCode VecScatterCUDAIndicesDestroy(PetscCUDAIndices*);
71: PETSC_INTERN PetscErrorCode VecScatterCUDA_StoS(Vec,Vec,PetscCUDAIndices,InsertMode,ScatterMode);
73: typedef enum {VEC_SCATTER_CUDA_STOS, VEC_SCATTER_CUDA_PTOP} VecCUDAScatterType;
74: typedef enum {VEC_SCATTER_CUDA_GENERAL, VEC_SCATTER_CUDA_STRIDED} VecCUDASequentialScatterMode;
76: struct _p_VecScatterCUDAIndices_PtoP {
77: PetscInt ns;
78: PetscInt sendLowestIndex;
79: PetscInt nr;
80: PetscInt recvLowestIndex;
81: };
83: struct _p_VecScatterCUDAIndices_StoS {
84: /* from indices data */
85: PetscInt *fslots;
86: PetscInt fromFirst;
87: PetscInt fromStep;
88: VecCUDASequentialScatterMode fromMode;
90: /* to indices data */
91: PetscInt *tslots;
92: PetscInt toFirst;
93: PetscInt toStep;
94: VecCUDASequentialScatterMode toMode;
96: PetscInt n;
97: PetscInt MAX_BLOCKS;
98: PetscInt MAX_CORESIDENT_THREADS;
99: cudaStream_t stream;
100: };
102: struct _p_PetscCUDAIndices {
103: void * scatter;
104: VecCUDAScatterType scatterType;
105: };
107: /* complex single */
108: #if defined(PETSC_USE_COMPLEX)
109: #if defined(PETSC_USE_REAL_SINGLE)
110: #define cublasXaxpy(a,b,c,d,e,f,g) cublasCaxpy((a),(b),(cuComplex*)(c),(cuComplex*)(d),(e),(cuComplex*)(f),(g))
111: #define cublasXscal(a,b,c,d,e) cublasCscal((a),(b),(cuComplex*)(c),(cuComplex*)(d),(e))
112: #define cublasXdotu(a,b,c,d,e,f,g) cublasCdotu((a),(b),(cuComplex*)(c),(d),(cuComplex*)(e),(f),(cuComplex*)(g))
113: #define cublasXdot(a,b,c,d,e,f,g) cublasCdotc((a),(b),(cuComplex*)(c),(d),(cuComplex*)(e),(f),(cuComplex*)(g))
114: #define cublasXswap(a,b,c,d,e,f) cublasCswap((a),(b),(cuComplex*)(c),(d),(cuComplex*)(e),(f))
115: #define cublasXnrm2(a,b,c,d,e) cublasScnrm2((a),(b),(cuComplex*)(c),(d),(e))
116: #define cublasIXamax(a,b,c,d,e) cublasIcamax((a),(b),(cuComplex*)(c),(d),(e))
117: #define cublasXasum(a,b,c,d,e) cublasScasum((a),(b),(cuComplex*)(c),(d),(e))
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: #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: #else /* real double */
139: #define cublasXaxpy cublasDaxpy
140: #define cublasXscal cublasDscal
141: #define cublasXdotu cublasDdot
142: #define cublasXdot cublasDdot
143: #define cublasXswap cublasDswap
144: #define cublasXnrm2 cublasDnrm2
145: #define cublasIXamax cublasIdamax
146: #define cublasXasum cublasDasum
147: #endif
148: #endif
150: #endif