Actual source code: veckokkosimpl.hpp
1: #pragma once
3: #include <../src/vec/vec/impls/mpi/pvecimpl.h>
4: #include <petsc/private/kokkosimpl.hpp>
6: #if defined(PETSC_USE_DEBUG)
7: #define VecErrorIfNotKokkos(v) \
8: do { \
9: PetscBool isKokkos = PETSC_FALSE; \
10: PetscCall(PetscObjectTypeCompareAny((PetscObject)(v), &isKokkos, VECSEQKOKKOS, VECMPIKOKKOS, VECKOKKOS, "")); \
11: PetscCheck(isKokkos, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Calling VECKOKKOS methods on a non-VECKOKKOS object"); \
12: } while (0)
13: #else
14: #define VecErrorIfNotKokkos(v) \
15: do { \
16: (void)(v); \
17: } while (0)
18: #endif
20: /* Stuff related to Vec_Kokkos */
22: struct Vec_Kokkos {
23: PetscScalarKokkosDualView v_dual;
24: PetscScalarKokkosView unplaced_d; // Unplaced device array in VecKokkosPlaceArray()
25: PetscScalarKokkosDualView unplaced_dual; // Unplaced v_dual, used in VecGetLocalVector_SeqKokkos()
27: /* COO stuff */
28: PetscCountKokkosView jmap1_d; /* [m+1]: i-th entry of the vector has jmap1[i+1]-jmap1[i] repeats in COO arrays */
29: PetscCountKokkosView perm1_d; /* [tot1]: permutation array for local entries */
31: PetscCountKokkosView imap2_d; /* [nnz2]: i-th unique entry in recvbuf is imap2[i]-th entry in the vector */
32: PetscCountKokkosView jmap2_d; /* [nnz2+1] */
33: PetscCountKokkosView perm2_d; /* [recvlen] */
34: PetscCountKokkosView Cperm_d; /* [sendlen]: permutation array to fill sendbuf[]. 'C' for communication */
35: PetscScalarKokkosView sendbuf_d, recvbuf_d; /* Buffers for remote values in VecSetValuesCOO() */
37: // (internal use only) sometimes we need to allocate multiple vectors from a contiguous memory block.
38: // We stash the memory in w_dual, which has the same lifespan as this vector. See VecDuplicateVecs_SeqKokkos_GEMV.
39: PetscScalarKokkosDualView w_dual;
41: /* Construct Vec_Kokkos with the given array(s). n is the length of the array.
42: If n != 0, host array (array_h) must not be NULL.
43: If device array (array_d) is NULL, then a proper device mirror will be allocated.
44: Otherwise, the mirror will be created using the given array_d.
45: If both arrays are given, we assume they contain the same value (i.e., sync'ed)
46: */
47: Vec_Kokkos(PetscInt n, PetscScalar *array_h, PetscScalar *array_d = NULL)
48: {
49: PetscScalarKokkosViewHost v_h(array_h, n);
50: PetscScalarKokkosView v_d;
52: if (array_d) {
53: v_d = PetscScalarKokkosView(array_d, n); /* Use the given device array */
54: } else {
55: v_d = Kokkos::create_mirror_view(Kokkos::WithoutInitializing, DefaultMemorySpace(), v_h); /* Create a mirror in DefaultMemorySpace but do not copy values */
56: }
57: v_dual = PetscScalarKokkosDualView(v_d, v_h);
58: if (!array_d) v_dual.modify_host();
59: }
61: // Construct Vec_Kokkos with the given DualView. Use the sync state as is. With reference counting, Kokkos manages its lifespan.
62: Vec_Kokkos(PetscScalarKokkosDualView dual) : v_dual(dual) { }
64: /* SFINAE: Update the object with an array in the given memory space,
65: assuming the given array contains the latest value for this vector.
66: */
67: template <typename MemorySpace, std::enable_if_t<std::is_same<MemorySpace, HostMirrorMemorySpace>::value, bool> = true, std::enable_if_t<std::is_same<MemorySpace, DefaultMemorySpace>::value, bool> = true>
68: PetscErrorCode UpdateArray(PetscScalar *array)
69: {
70: PetscScalarKokkosView v_d(array, v_dual.extent(0));
71: PetscScalarKokkosViewHost v_h(array, v_dual.extent(0));
73: PetscFunctionBegin;
74: /* Kokkos said they would add error-checking so that users won't accidentally pass two different Views in this case */
75: PetscCallCXX(v_dual = PetscScalarKokkosDualView(v_d, v_h));
76: PetscFunctionReturn(PETSC_SUCCESS);
77: }
79: template <typename MemorySpace, std::enable_if_t<std::is_same<MemorySpace, HostMirrorMemorySpace>::value, bool> = true, std::enable_if_t<!std::is_same<MemorySpace, DefaultMemorySpace>::value, bool> = true>
80: PetscErrorCode UpdateArray(PetscScalar *array)
81: {
82: PetscScalarKokkosViewHost v_h(array, v_dual.extent(0));
84: PetscFunctionBegin;
85: PetscCallCXX(v_dual = PetscScalarKokkosDualView(v_dual.view<DefaultMemorySpace>(), v_h));
86: PetscCallCXX(v_dual.modify_host());
87: PetscFunctionReturn(PETSC_SUCCESS);
88: }
90: template <typename MemorySpace, std::enable_if_t<!std::is_same<MemorySpace, HostMirrorMemorySpace>::value, bool> = true, std::enable_if_t<std::is_same<MemorySpace, DefaultMemorySpace>::value, bool> = true>
91: PetscErrorCode UpdateArray(PetscScalar *array)
92: {
93: PetscScalarKokkosView v_d(array, v_dual.extent(0));
95: PetscFunctionBegin;
96: PetscCallCXX(v_dual = PetscScalarKokkosDualView(v_d, v_dual.view_host()));
97: PetscCallCXX(v_dual.modify_device());
98: PetscFunctionReturn(PETSC_SUCCESS);
99: }
101: PetscErrorCode SetUpCOO(const Vec_Seq *vecseq, PetscInt m)
102: {
103: PetscFunctionBegin;
104: PetscCallCXX(jmap1_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecseq->jmap1, m + 1)));
105: PetscCallCXX(perm1_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecseq->perm1, vecseq->tot1)));
106: PetscFunctionReturn(PETSC_SUCCESS);
107: }
109: PetscErrorCode SetUpCOO(const Vec_MPI *vecmpi, PetscInt m)
110: {
111: PetscFunctionBegin;
112: PetscCallCXX(jmap1_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecmpi->jmap1, m + 1)));
113: PetscCallCXX(perm1_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecmpi->perm1, vecmpi->tot1)));
114: PetscCallCXX(imap2_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecmpi->imap2, vecmpi->nnz2)));
115: PetscCallCXX(jmap2_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecmpi->jmap2, vecmpi->nnz2 + 1)));
116: PetscCallCXX(perm2_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecmpi->perm2, vecmpi->recvlen)));
117: PetscCallCXX(Cperm_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecmpi->Cperm, vecmpi->sendlen)));
118: PetscCallCXX(sendbuf_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscScalarKokkosViewHost(vecmpi->sendbuf, vecmpi->sendlen)));
119: PetscCallCXX(recvbuf_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscScalarKokkosViewHost(vecmpi->recvbuf, vecmpi->recvlen)));
120: PetscFunctionReturn(PETSC_SUCCESS);
121: }
122: };
124: PETSC_INTERN PetscErrorCode VecAbs_SeqKokkos(Vec);
125: PETSC_INTERN PetscErrorCode VecReciprocal_SeqKokkos(Vec);
126: PETSC_INTERN PetscErrorCode VecDotNorm2_SeqKokkos(Vec, Vec, PetscScalar *, PetscScalar *);
127: PETSC_INTERN PetscErrorCode VecPointwiseDivide_SeqKokkos(Vec, Vec, Vec);
128: PETSC_INTERN PetscErrorCode VecWAXPY_SeqKokkos(Vec, PetscScalar, Vec, Vec);
129: PETSC_INTERN PetscErrorCode VecMDot_SeqKokkos(Vec, PetscInt, const Vec[], PetscScalar *);
130: PETSC_INTERN PetscErrorCode VecMTDot_SeqKokkos(Vec, PetscInt, const Vec[], PetscScalar *);
131: PETSC_INTERN PetscErrorCode VecSet_SeqKokkos(Vec, PetscScalar);
132: PETSC_INTERN PetscErrorCode VecMAXPY_SeqKokkos(Vec, PetscInt, const PetscScalar *, Vec *);
133: PETSC_INTERN PetscErrorCode VecAXPBYPCZ_SeqKokkos(Vec, PetscScalar, PetscScalar, PetscScalar, Vec, Vec);
134: PETSC_INTERN PetscErrorCode VecPointwiseMult_SeqKokkos(Vec, Vec, Vec);
135: PETSC_INTERN PetscErrorCode VecPlaceArray_SeqKokkos(Vec, const PetscScalar *);
136: PETSC_INTERN PetscErrorCode VecResetArray_SeqKokkos(Vec);
137: PETSC_INTERN PetscErrorCode VecReplaceArray_SeqKokkos(Vec, const PetscScalar *);
138: PETSC_INTERN PetscErrorCode VecDot_SeqKokkos(Vec, Vec, PetscScalar *);
139: PETSC_INTERN PetscErrorCode VecTDot_SeqKokkos(Vec, Vec, PetscScalar *);
140: PETSC_INTERN PetscErrorCode VecScale_SeqKokkos(Vec, PetscScalar);
141: PETSC_INTERN PetscErrorCode VecCopy_SeqKokkos(Vec, Vec);
142: PETSC_INTERN PetscErrorCode VecSwap_SeqKokkos(Vec, Vec);
143: PETSC_INTERN PetscErrorCode VecAXPY_SeqKokkos(Vec, PetscScalar, Vec);
144: PETSC_INTERN PetscErrorCode VecAXPBY_SeqKokkos(Vec, PetscScalar, PetscScalar, Vec);
145: PETSC_INTERN PetscErrorCode VecConjugate_SeqKokkos(Vec xin);
146: PETSC_INTERN PetscErrorCode VecNorm_SeqKokkos(Vec, NormType, PetscReal *);
147: PETSC_INTERN PetscErrorCode VecErrorWeightedNorms_SeqKokkos(Vec, Vec, Vec, NormType, PetscReal, Vec, PetscReal, Vec, PetscReal, PetscReal *, PetscInt *, PetscReal *, PetscInt *, PetscReal *, PetscInt *);
148: PETSC_INTERN PetscErrorCode VecCreate_SeqKokkos(Vec);
149: PETSC_INTERN PetscErrorCode VecCreate_SeqKokkos_Private(Vec, const PetscScalar *);
150: PETSC_INTERN PetscErrorCode VecCreate_MPIKokkos(Vec);
151: PETSC_INTERN PetscErrorCode VecCreate_MPIKokkos_Private(Vec, PetscBool, PetscInt, const PetscScalar *);
152: PETSC_INTERN PetscErrorCode VecCreate_Kokkos(Vec);
153: PETSC_INTERN PetscErrorCode VecAYPX_SeqKokkos(Vec, PetscScalar, Vec);
154: PETSC_INTERN PetscErrorCode VecSetRandom_SeqKokkos(Vec, PetscRandom);
155: PETSC_INTERN PetscErrorCode VecGetLocalVector_SeqKokkos(Vec, Vec);
156: PETSC_INTERN PetscErrorCode VecGetLocalVectorRead_SeqKokkos(Vec, Vec);
157: PETSC_INTERN PetscErrorCode VecRestoreLocalVector_SeqKokkos(Vec, Vec);
158: PETSC_INTERN PetscErrorCode VecRestoreLocalVectorRead_SeqKokkos(Vec, Vec);
159: PETSC_INTERN PetscErrorCode VecGetArrayWrite_SeqKokkos(Vec, PetscScalar **);
160: PETSC_INTERN PetscErrorCode VecCopy_SeqKokkos_Private(Vec, Vec);
161: PETSC_INTERN PetscErrorCode VecSetRandom_SeqKokkos_Private(Vec, PetscRandom);
162: PETSC_INTERN PetscErrorCode VecResetArray_SeqKokkos_Private(Vec);
163: PETSC_INTERN PetscErrorCode VecMin_SeqKokkos(Vec, PetscInt *, PetscReal *);
164: PETSC_INTERN PetscErrorCode VecMax_SeqKokkos(Vec, PetscInt *, PetscReal *);
165: PETSC_INTERN PetscErrorCode VecSum_SeqKokkos(Vec, PetscScalar *);
166: PETSC_INTERN PetscErrorCode VecShift_SeqKokkos(Vec, PetscScalar);
167: PETSC_INTERN PetscErrorCode VecGetArray_SeqKokkos(Vec, PetscScalar **);
168: PETSC_INTERN PetscErrorCode VecRestoreArray_SeqKokkos(Vec, PetscScalar **);
170: PETSC_INTERN PetscErrorCode VecGetArrayAndMemType_SeqKokkos(Vec, PetscScalar **, PetscMemType *);
171: PETSC_INTERN PetscErrorCode VecRestoreArrayAndMemType_SeqKokkos(Vec, PetscScalar **);
172: PETSC_INTERN PetscErrorCode VecGetArrayWriteAndMemType_SeqKokkos(Vec, PetscScalar **, PetscMemType *);
173: PETSC_INTERN PetscErrorCode VecGetSubVector_Kokkos_Private(Vec, PetscBool, IS, Vec *);
174: PETSC_INTERN PetscErrorCode VecRestoreSubVector_SeqKokkos(Vec, IS, Vec *);
176: PETSC_INTERN PetscErrorCode VecDuplicateVecs_Kokkos_GEMV(Vec, PetscInt, Vec **);
177: PETSC_INTERN PetscErrorCode VecMDot_SeqKokkos_GEMV(Vec, PetscInt, const Vec *, PetscScalar *);
178: PETSC_INTERN PetscErrorCode VecMTDot_SeqKokkos_GEMV(Vec, PetscInt, const Vec *, PetscScalar *);
179: PETSC_INTERN PetscErrorCode VecMAXPY_SeqKokkos_GEMV(Vec, PetscInt, const PetscScalar *, Vec *);
181: PETSC_INTERN PetscErrorCode VecCreateMPIKokkosWithLayoutAndArrays_Private(PetscLayout map, const PetscScalar *, const PetscScalar *, Vec *);
182: PETSC_INTERN PetscErrorCode VecBindToCPU_SeqKokkos(Vec, PetscBool);