Actual source code: veckokkosimpl.hpp
1: #pragma once
3: #include <../src/vec/vec/impls/mpi/pvecimpl.h>
4: #include <petsc/private/vecimpl_kokkos.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;
25: /* COO stuff */
26: PetscCountKokkosView jmap1_d; /* [m+1]: i-th entry of the vector has jmap1[i+1]-jmap1[i] repeats in COO arrays */
27: PetscCountKokkosView perm1_d; /* [tot1]: permutation array for local entries */
29: PetscCountKokkosView imap2_d; /* [nnz2]: i-th unique entry in recvbuf is imap2[i]-th entry in the vector */
30: PetscCountKokkosView jmap2_d; /* [nnz2+1] */
31: PetscCountKokkosView perm2_d; /* [recvlen] */
32: PetscCountKokkosView Cperm_d; /* [sendlen]: permutation array to fill sendbuf[]. 'C' for communication */
33: PetscScalarKokkosView sendbuf_d, recvbuf_d; /* Buffers for remote values in VecSetValuesCOO() */
35: /* Construct Vec_Kokkos with the given array(s). n is the length of the array.
36: If n != 0, host array (array_h) must not be NULL.
37: If device array (array_d) is NULL, then a proper device mirror will be allocated.
38: Otherwise, the mirror will be created using the given array_d.
39: */
40: Vec_Kokkos(PetscInt n, PetscScalar *array_h, PetscScalar *array_d = NULL)
41: {
42: PetscScalarKokkosViewHost v_h(array_h, n);
43: PetscScalarKokkosView v_d;
45: if (array_d) {
46: v_d = PetscScalarKokkosView(array_d, n); /* Use the given device array */
47: } else {
48: v_d = Kokkos::create_mirror_view(Kokkos::WithoutInitializing, DefaultMemorySpace(), v_h); /* Create a mirror in DefaultMemorySpace but do not copy values */
49: }
50: v_dual = PetscScalarKokkosDualView(v_d, v_h);
51: if (!array_d) v_dual.modify_host();
52: }
54: /* SFINAE: Update the object with an array in the given memory space,
55: assuming the given array contains the latest value for this vector.
56: */
57: template <typename MemorySpace, std::enable_if_t<std::is_same<MemorySpace, Kokkos::HostSpace>::value, bool> = true, std::enable_if_t<std::is_same<MemorySpace, DefaultMemorySpace>::value, bool> = true>
58: PetscErrorCode UpdateArray(PetscScalar *array)
59: {
60: PetscScalarKokkosViewHost v_h(array, v_dual.extent(0));
61: PetscFunctionBegin;
62: /* Kokkos said they would add error-checking so that users won't accidentally pass two different Views in this case */
63: PetscCallCXX(v_dual = PetscScalarKokkosDualView(v_h, v_h));
64: PetscFunctionReturn(PETSC_SUCCESS);
65: }
67: template <typename MemorySpace, std::enable_if_t<std::is_same<MemorySpace, Kokkos::HostSpace>::value, bool> = true, std::enable_if_t<!std::is_same<MemorySpace, DefaultMemorySpace>::value, bool> = true>
68: PetscErrorCode UpdateArray(PetscScalar *array)
69: {
70: PetscScalarKokkosViewHost v_h(array, v_dual.extent(0));
71: PetscFunctionBegin;
72: PetscCallCXX(v_dual = PetscScalarKokkosDualView(v_dual.view<DefaultMemorySpace>(), v_h));
73: PetscCallCXX(v_dual.modify_host());
74: PetscFunctionReturn(PETSC_SUCCESS);
75: }
77: template <typename MemorySpace, std::enable_if_t<!std::is_same<MemorySpace, Kokkos::HostSpace>::value, bool> = true, std::enable_if_t<std::is_same<MemorySpace, DefaultMemorySpace>::value, bool> = true>
78: PetscErrorCode UpdateArray(PetscScalar *array)
79: {
80: PetscScalarKokkosView v_d(array, v_dual.extent(0));
81: PetscFunctionBegin;
82: PetscCallCXX(v_dual = PetscScalarKokkosDualView(v_d, v_dual.view<Kokkos::HostSpace>()));
83: PetscCallCXX(v_dual.modify_device());
84: PetscFunctionReturn(PETSC_SUCCESS);
85: }
87: PetscErrorCode SetUpCOO(const Vec_Seq *vecseq, PetscInt m)
88: {
89: PetscFunctionBegin;
90: PetscCallCXX(jmap1_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecseq->jmap1, m + 1)));
91: PetscCallCXX(perm1_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecseq->perm1, vecseq->tot1)));
92: PetscFunctionReturn(PETSC_SUCCESS);
93: }
95: PetscErrorCode SetUpCOO(const Vec_MPI *vecmpi, PetscInt m)
96: {
97: PetscFunctionBegin;
98: PetscCallCXX(jmap1_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecmpi->jmap1, m + 1)));
99: PetscCallCXX(perm1_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecmpi->perm1, vecmpi->tot1)));
100: PetscCallCXX(imap2_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecmpi->imap2, vecmpi->nnz2)));
101: PetscCallCXX(jmap2_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecmpi->jmap2, vecmpi->nnz2 + 1)));
102: PetscCallCXX(perm2_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecmpi->perm2, vecmpi->recvlen)));
103: PetscCallCXX(Cperm_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscCountKokkosViewHost(vecmpi->Cperm, vecmpi->sendlen)));
104: PetscCallCXX(sendbuf_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscScalarKokkosViewHost(vecmpi->sendbuf, vecmpi->sendlen)));
105: PetscCallCXX(recvbuf_d = Kokkos::create_mirror_view_and_copy(DefaultMemorySpace(), PetscScalarKokkosViewHost(vecmpi->recvbuf, vecmpi->recvlen)));
106: PetscFunctionReturn(PETSC_SUCCESS);
107: }
108: };
110: PETSC_INTERN PetscErrorCode VecAbs_SeqKokkos(Vec);
111: PETSC_INTERN PetscErrorCode VecReciprocal_SeqKokkos(Vec);
112: PETSC_INTERN PetscErrorCode VecDotNorm2_SeqKokkos(Vec, Vec, PetscScalar *, PetscScalar *);
113: PETSC_INTERN PetscErrorCode VecPointwiseDivide_SeqKokkos(Vec, Vec, Vec);
114: PETSC_INTERN PetscErrorCode VecWAXPY_SeqKokkos(Vec, PetscScalar, Vec, Vec);
115: PETSC_INTERN PetscErrorCode VecMDot_SeqKokkos(Vec, PetscInt, const Vec[], PetscScalar *);
116: PETSC_INTERN PetscErrorCode VecMTDot_SeqKokkos(Vec, PetscInt, const Vec[], PetscScalar *);
117: PETSC_INTERN PetscErrorCode VecSet_SeqKokkos(Vec, PetscScalar);
118: PETSC_INTERN PetscErrorCode VecMAXPY_SeqKokkos(Vec, PetscInt, const PetscScalar *, Vec *);
119: PETSC_INTERN PetscErrorCode VecAXPBYPCZ_SeqKokkos(Vec, PetscScalar, PetscScalar, PetscScalar, Vec, Vec);
120: PETSC_INTERN PetscErrorCode VecPointwiseMult_SeqKokkos(Vec, Vec, Vec);
121: PETSC_INTERN PetscErrorCode VecPlaceArray_SeqKokkos(Vec, const PetscScalar *);
122: PETSC_INTERN PetscErrorCode VecResetArray_SeqKokkos(Vec);
123: PETSC_INTERN PetscErrorCode VecReplaceArray_SeqKokkos(Vec, const PetscScalar *);
124: PETSC_INTERN PetscErrorCode VecDot_SeqKokkos(Vec, Vec, PetscScalar *);
125: PETSC_INTERN PetscErrorCode VecTDot_SeqKokkos(Vec, Vec, PetscScalar *);
126: PETSC_INTERN PetscErrorCode VecScale_SeqKokkos(Vec, PetscScalar);
127: PETSC_INTERN PetscErrorCode VecCopy_SeqKokkos(Vec, Vec);
128: PETSC_INTERN PetscErrorCode VecSwap_SeqKokkos(Vec, Vec);
129: PETSC_INTERN PetscErrorCode VecAXPY_SeqKokkos(Vec, PetscScalar, Vec);
130: PETSC_INTERN PetscErrorCode VecAXPBY_SeqKokkos(Vec, PetscScalar, PetscScalar, Vec);
131: PETSC_INTERN PetscErrorCode VecConjugate_SeqKokkos(Vec xin);
132: PETSC_INTERN PetscErrorCode VecNorm_SeqKokkos(Vec, NormType, PetscReal *);
133: PETSC_INTERN PetscErrorCode VecErrorWeightedNorms_SeqKokkos(Vec, Vec, Vec, NormType, PetscReal, Vec, PetscReal, Vec, PetscReal, PetscReal *, PetscInt *, PetscReal *, PetscInt *, PetscReal *, PetscInt *);
134: PETSC_INTERN PetscErrorCode VecCreate_SeqKokkos(Vec);
135: PETSC_INTERN PetscErrorCode VecCreate_SeqKokkos_Private(Vec, const PetscScalar *);
136: PETSC_INTERN PetscErrorCode VecCreate_MPIKokkos(Vec);
137: PETSC_INTERN PetscErrorCode VecCreate_MPIKokkos_Private(Vec, PetscBool, PetscInt, const PetscScalar *);
138: PETSC_INTERN PetscErrorCode VecCreate_Kokkos(Vec);
139: PETSC_INTERN PetscErrorCode VecAYPX_SeqKokkos(Vec, PetscScalar, Vec);
140: PETSC_INTERN PetscErrorCode VecSetRandom_SeqKokkos(Vec, PetscRandom);
141: PETSC_INTERN PetscErrorCode VecGetLocalVector_SeqKokkos(Vec, Vec);
142: PETSC_INTERN PetscErrorCode VecRestoreLocalVector_SeqKokkos(Vec, Vec);
143: PETSC_INTERN PetscErrorCode VecGetArrayWrite_SeqKokkos(Vec, PetscScalar **);
144: PETSC_INTERN PetscErrorCode VecCopy_SeqKokkos_Private(Vec, Vec);
145: PETSC_INTERN PetscErrorCode VecSetRandom_SeqKokkos_Private(Vec, PetscRandom);
146: PETSC_INTERN PetscErrorCode VecResetArray_SeqKokkos_Private(Vec);
147: PETSC_INTERN PetscErrorCode VecMin_SeqKokkos(Vec, PetscInt *, PetscReal *);
148: PETSC_INTERN PetscErrorCode VecMax_SeqKokkos(Vec, PetscInt *, PetscReal *);
149: PETSC_INTERN PetscErrorCode VecSum_SeqKokkos(Vec, PetscScalar *);
150: PETSC_INTERN PetscErrorCode VecShift_SeqKokkos(Vec, PetscScalar);
151: PETSC_INTERN PetscErrorCode VecGetArray_SeqKokkos(Vec, PetscScalar **);
152: PETSC_INTERN PetscErrorCode VecRestoreArray_SeqKokkos(Vec, PetscScalar **);
154: PETSC_INTERN PetscErrorCode VecGetArrayAndMemType_SeqKokkos(Vec, PetscScalar **, PetscMemType *);
155: PETSC_INTERN PetscErrorCode VecRestoreArrayAndMemType_SeqKokkos(Vec, PetscScalar **);
156: PETSC_INTERN PetscErrorCode VecGetArrayWriteAndMemType_SeqKokkos(Vec, PetscScalar **, PetscMemType *);
157: PETSC_INTERN PetscErrorCode VecGetSubVector_Kokkos_Private(Vec, PetscBool, IS, Vec *);
158: PETSC_INTERN PetscErrorCode VecRestoreSubVector_SeqKokkos(Vec, IS, Vec *);